| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | /* | ||
| 2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
| 4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
| 5 | |||
| 6 | This file is part of SLEPc. | ||
| 7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
| 8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 9 | */ | ||
| 10 | /* | ||
| 11 | Basic DS routines | ||
| 12 | */ | ||
| 13 | |||
| 14 | #include <slepc/private/dsimpl.h> /*I "slepcds.h" I*/ | ||
| 15 | |||
| 16 | PetscFunctionList DSList = NULL; | ||
| 17 | PetscBool DSRegisterAllCalled = PETSC_FALSE; | ||
| 18 | PetscClassId DS_CLASSID = 0; | ||
| 19 | PetscLogEvent DS_Solve = 0,DS_Vectors = 0,DS_Synchronize = 0,DS_Other = 0; | ||
| 20 | static PetscBool DSPackageInitialized = PETSC_FALSE; | ||
| 21 | |||
| 22 | const char *DSStateTypes[] = {"RAW","INTERMEDIATE","CONDENSED","TRUNCATED","DSStateType","DS_STATE_",NULL}; | ||
| 23 | const char *DSParallelTypes[] = {"REDUNDANT","SYNCHRONIZED","DISTRIBUTED","DSParallelType","DS_PARALLEL_",NULL}; | ||
| 24 | const char *DSMatName[DS_NUM_MAT] = {"A","B","C","T","D","Q","Z","X","Y","U","V","W","E0","E1","E2","E3","E4","E5","E6","E7","E8","E9"}; | ||
| 25 | DSMatType DSMatExtra[DS_NUM_EXTRA] = {DS_MAT_E0,DS_MAT_E1,DS_MAT_E2,DS_MAT_E3,DS_MAT_E4,DS_MAT_E5,DS_MAT_E6,DS_MAT_E7,DS_MAT_E8,DS_MAT_E9}; | ||
| 26 | |||
| 27 | /*@C | ||
| 28 | DSFinalizePackage - This function destroys everything in the SLEPc interface | ||
| 29 | to the `DS` package. It is called from `SlepcFinalize()`. | ||
| 30 | |||
| 31 | Level: developer | ||
| 32 | |||
| 33 | .seealso: [](sec:ds), `SlepcFinalize()`, `DSInitializePackage()` | ||
| 34 | @*/ | ||
| 35 | 11074 | PetscErrorCode DSFinalizePackage(void) | |
| 36 | { | ||
| 37 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
11074 | PetscFunctionBegin; |
| 38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscFunctionListDestroy(&DSList)); |
| 39 | 11074 | DSPackageInitialized = PETSC_FALSE; | |
| 40 | 11074 | DSRegisterAllCalled = PETSC_FALSE; | |
| 41 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
11074 | PetscFunctionReturn(PETSC_SUCCESS); |
| 42 | } | ||
| 43 | |||
| 44 | /*@C | ||
| 45 | DSInitializePackage - This function initializes everything in the `DS` package. | ||
| 46 | It is called from `PetscDLLibraryRegister_slepc()` when using dynamic libraries, and | ||
| 47 | on the first call to `DSCreate()` when using shared or static libraries. | ||
| 48 | |||
| 49 | Note: | ||
| 50 | This function never needs to be called by SLEPc users. | ||
| 51 | |||
| 52 | Level: developer | ||
| 53 | |||
| 54 | .seealso: [](sec:ds), `DS`, `SlepcInitialize()`, `DSFinalizePackage()` | ||
| 55 | @*/ | ||
| 56 | 135203 | PetscErrorCode DSInitializePackage(void) | |
| 57 | { | ||
| 58 | 135203 | char logList[256]; | |
| 59 | 135203 | PetscBool opt,pkg; | |
| 60 | 135203 | PetscClassId classids[1]; | |
| 61 | |||
| 62 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
135203 | PetscFunctionBegin; |
| 63 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
135203 | if (DSPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); |
| 64 | 11074 | DSPackageInitialized = PETSC_TRUE; | |
| 65 | /* Register Classes */ | ||
| 66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscClassIdRegister("Direct Solver",&DS_CLASSID)); |
| 67 | /* Register Constructors */ | ||
| 68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegisterAll()); |
| 69 | /* Register Events */ | ||
| 70 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscLogEventRegister("DSSolve",DS_CLASSID,&DS_Solve)); |
| 71 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscLogEventRegister("DSVectors",DS_CLASSID,&DS_Vectors)); |
| 72 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscLogEventRegister("DSSynchronize",DS_CLASSID,&DS_Synchronize)); |
| 73 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscLogEventRegister("DSOther",DS_CLASSID,&DS_Other)); |
| 74 | /* Process Info */ | ||
| 75 | 11074 | classids[0] = DS_CLASSID; | |
| 76 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscInfoProcessClass("ds",1,&classids[0])); |
| 77 | /* Process summary exclusions */ | ||
| 78 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt)); |
| 79 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
11074 | if (opt) { |
| 80 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
84 | PetscCall(PetscStrInList("ds",logList,',',&pkg)); |
| 81 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
84 | if (pkg) PetscCall(PetscLogEventDeactivateClass(DS_CLASSID)); |
| 82 | } | ||
| 83 | /* Register package finalizer */ | ||
| 84 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(PetscRegisterFinalize(DSFinalizePackage)); |
| 85 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2144 | PetscFunctionReturn(PETSC_SUCCESS); |
| 86 | } | ||
| 87 | |||
| 88 | /*@ | ||
| 89 | DSCreate - Creates a `DS` context. | ||
| 90 | |||
| 91 | Collective | ||
| 92 | |||
| 93 | Input Parameter: | ||
| 94 | . comm - MPI communicator | ||
| 95 | |||
| 96 | Output Parameter: | ||
| 97 | . newds - location to put the `DS` context | ||
| 98 | |||
| 99 | Level: beginner | ||
| 100 | |||
| 101 | Note: | ||
| 102 | `DS` objects are not intended for normal users but only for | ||
| 103 | advanced user that for instance implement their own solvers. | ||
| 104 | |||
| 105 | .seealso: [](sec:ds), `DSDestroy()`, `DS` | ||
| 106 | @*/ | ||
| 107 | 13342 | PetscErrorCode DSCreate(MPI_Comm comm,DS *newds) | |
| 108 | { | ||
| 109 | 13342 | DS ds; | |
| 110 | 13342 | PetscInt i; | |
| 111 | |||
| 112 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
13342 | PetscFunctionBegin; |
| 113 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
13342 | PetscAssertPointer(newds,2); |
| 114 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13342 | PetscCall(DSInitializePackage()); |
| 115 |
7/12✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
13342 | PetscCall(SlepcHeaderCreate(ds,DS_CLASSID,"DS","Direct Solver (or Dense System)","DS",comm,DSDestroy,DSView)); |
| 116 | |||
| 117 | 13342 | ds->state = DS_STATE_RAW; | |
| 118 | 13342 | ds->method = 0; | |
| 119 | 13342 | ds->compact = PETSC_FALSE; | |
| 120 | 13342 | ds->refined = PETSC_FALSE; | |
| 121 | 13342 | ds->extrarow = PETSC_FALSE; | |
| 122 | 13342 | ds->ld = 0; | |
| 123 | 13342 | ds->l = 0; | |
| 124 | 13342 | ds->n = 0; | |
| 125 | 13342 | ds->k = 0; | |
| 126 | 13342 | ds->t = 0; | |
| 127 | 13342 | ds->bs = 1; | |
| 128 | 13342 | ds->sc = NULL; | |
| 129 | 13342 | ds->pmode = DS_PARALLEL_REDUNDANT; | |
| 130 | |||
| 131 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
306866 | for (i=0;i<DS_NUM_MAT;i++) ds->omat[i] = NULL; |
| 132 | 13342 | ds->perm = NULL; | |
| 133 | 13342 | ds->data = NULL; | |
| 134 | 13342 | ds->scset = PETSC_FALSE; | |
| 135 | 13342 | ds->work = NULL; | |
| 136 | 13342 | ds->rwork = NULL; | |
| 137 | 13342 | ds->iwork = NULL; | |
| 138 | 13342 | ds->lwork = 0; | |
| 139 | 13342 | ds->lrwork = 0; | |
| 140 | 13342 | ds->liwork = 0; | |
| 141 | |||
| 142 | 13342 | *newds = ds; | |
| 143 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
13342 | PetscFunctionReturn(PETSC_SUCCESS); |
| 144 | } | ||
| 145 | |||
| 146 | /*@ | ||
| 147 | DSSetOptionsPrefix - Sets the prefix used for searching for all | ||
| 148 | `DS` options in the database. | ||
| 149 | |||
| 150 | Logically Collective | ||
| 151 | |||
| 152 | Input Parameters: | ||
| 153 | + ds - the direct solver context | ||
| 154 | - prefix - the prefix string to prepend to all `DS` option requests | ||
| 155 | |||
| 156 | Notes: | ||
| 157 | A hyphen (-) must NOT be given at the beginning of the prefix name. | ||
| 158 | The first character of all runtime options is AUTOMATICALLY the | ||
| 159 | hyphen. | ||
| 160 | |||
| 161 | Level: advanced | ||
| 162 | |||
| 163 | .seealso: [](sec:ds), `DSAppendOptionsPrefix()` | ||
| 164 | @*/ | ||
| 165 | 2485 | PetscErrorCode DSSetOptionsPrefix(DS ds,const char prefix[]) | |
| 166 | { | ||
| 167 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
2485 | PetscFunctionBegin; |
| 168 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2485 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 169 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2485 | PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ds,prefix)); |
| 170 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
478 | PetscFunctionReturn(PETSC_SUCCESS); |
| 171 | } | ||
| 172 | |||
| 173 | /*@ | ||
| 174 | DSAppendOptionsPrefix - Appends to the prefix used for searching for all | ||
| 175 | `DS` options in the database. | ||
| 176 | |||
| 177 | Logically Collective | ||
| 178 | |||
| 179 | Input Parameters: | ||
| 180 | + ds - the direct solver context | ||
| 181 | - prefix - the prefix string to prepend to all DS option requests | ||
| 182 | |||
| 183 | Notes: | ||
| 184 | A hyphen (-) must NOT be given at the beginning of the prefix name. | ||
| 185 | The first character of all runtime options is AUTOMATICALLY the hyphen. | ||
| 186 | |||
| 187 | Level: advanced | ||
| 188 | |||
| 189 | .seealso: [](sec:ds), `DSSetOptionsPrefix()` | ||
| 190 | @*/ | ||
| 191 | 1992 | PetscErrorCode DSAppendOptionsPrefix(DS ds,const char prefix[]) | |
| 192 | { | ||
| 193 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1992 | PetscFunctionBegin; |
| 194 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
1992 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 195 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1992 | PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)ds,prefix)); |
| 196 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
396 | PetscFunctionReturn(PETSC_SUCCESS); |
| 197 | } | ||
| 198 | |||
| 199 | /*@ | ||
| 200 | DSGetOptionsPrefix - Gets the prefix used for searching for all | ||
| 201 | `DS` options in the database. | ||
| 202 | |||
| 203 | Not Collective | ||
| 204 | |||
| 205 | Input Parameter: | ||
| 206 | . ds - the direct solver context | ||
| 207 | |||
| 208 | Output Parameter: | ||
| 209 | . prefix - pointer to the prefix string used is returned | ||
| 210 | |||
| 211 | Level: advanced | ||
| 212 | |||
| 213 | .seealso: [](sec:ds), `DSSetOptionsPrefix()`, `DSAppendOptionsPrefix()` | ||
| 214 | @*/ | ||
| 215 | ✗ | PetscErrorCode DSGetOptionsPrefix(DS ds,const char *prefix[]) | |
| 216 | { | ||
| 217 | ✗ | PetscFunctionBegin; | |
| 218 | ✗ | PetscValidHeaderSpecific(ds,DS_CLASSID,1); | |
| 219 | ✗ | PetscAssertPointer(prefix,2); | |
| 220 | ✗ | PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ds,prefix)); | |
| 221 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
| 222 | } | ||
| 223 | |||
| 224 | /*@ | ||
| 225 | DSSetType - Selects the type for the DS object. | ||
| 226 | |||
| 227 | Logically Collective | ||
| 228 | |||
| 229 | Input Parameters: | ||
| 230 | + ds - the direct solver context | ||
| 231 | - type - a known type | ||
| 232 | |||
| 233 | Level: intermediate | ||
| 234 | |||
| 235 | .seealso: [](sec:ds), `DSGetType()` | ||
| 236 | @*/ | ||
| 237 | 26832 | PetscErrorCode DSSetType(DS ds,DSType type) | |
| 238 | { | ||
| 239 | 26832 | PetscErrorCode (*r)(DS); | |
| 240 | 26832 | PetscBool match; | |
| 241 | |||
| 242 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
26832 | PetscFunctionBegin; |
| 243 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
26832 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 244 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
26832 | PetscAssertPointer(type,2); |
| 245 | |||
| 246 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
26832 | PetscCall(PetscObjectTypeCompare((PetscObject)ds,type,&match)); |
| 247 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
26832 | if (match) PetscFunctionReturn(PETSC_SUCCESS); |
| 248 | |||
| 249 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18003 | PetscCall(PetscFunctionListFind(DSList,type,&r)); |
| 250 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
18003 | PetscCheck(r,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested DS type %s",type); |
| 251 | |||
| 252 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
18003 | PetscTryTypeMethod(ds,destroy); |
| 253 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18003 | PetscCall(DSReset(ds)); |
| 254 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18003 | PetscCall(PetscMemzero(ds->ops,sizeof(struct _DSOps))); |
| 255 | |||
| 256 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18003 | PetscCall(PetscObjectChangeTypeName((PetscObject)ds,type)); |
| 257 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
18003 | PetscCall((*r)(ds)); |
| 258 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3465 | PetscFunctionReturn(PETSC_SUCCESS); |
| 259 | } | ||
| 260 | |||
| 261 | /*@ | ||
| 262 | DSGetType - Gets the `DS` type name (as a string) from the `DS` context. | ||
| 263 | |||
| 264 | Not Collective | ||
| 265 | |||
| 266 | Input Parameter: | ||
| 267 | . ds - the direct solver context | ||
| 268 | |||
| 269 | Output Parameter: | ||
| 270 | . type - name of the direct solver | ||
| 271 | |||
| 272 | Note: | ||
| 273 | `type` should not be retained for later use as it will be an invalid pointer | ||
| 274 | if the `DSType` of `ds` is changed. | ||
| 275 | |||
| 276 | Level: intermediate | ||
| 277 | |||
| 278 | .seealso: [](sec:ds), `DSSetType()`, `PetscObjectTypeCompare()`, `PetscObjectTypeCompareAny()` | ||
| 279 | @*/ | ||
| 280 | 20 | PetscErrorCode DSGetType(DS ds,DSType *type) | |
| 281 | { | ||
| 282 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
| 283 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
20 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 284 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
20 | PetscAssertPointer(type,2); |
| 285 | 20 | *type = ((PetscObject)ds)->type_name; | |
| 286 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
20 | PetscFunctionReturn(PETSC_SUCCESS); |
| 287 | } | ||
| 288 | |||
| 289 | /*@ | ||
| 290 | DSDuplicate - Creates a new direct solver object with the same options as | ||
| 291 | an existing one. | ||
| 292 | |||
| 293 | Collective | ||
| 294 | |||
| 295 | Input Parameter: | ||
| 296 | . ds - direct solver context | ||
| 297 | |||
| 298 | Output Parameter: | ||
| 299 | . dsnew - location to put the new `DS` | ||
| 300 | |||
| 301 | Notes: | ||
| 302 | `DSDuplicate()` DOES NOT COPY the matrices, and the new `DS` does not even have | ||
| 303 | internal arrays allocated. Use `DSAllocate()` to use the new `DS`. | ||
| 304 | |||
| 305 | The sorting criterion options are not copied, see `DSSetSlepcSC()`. | ||
| 306 | |||
| 307 | Level: intermediate | ||
| 308 | |||
| 309 | .seealso: [](sec:ds), `DSCreate()`, `DSAllocate()`, `DSSetSlepcSC()` | ||
| 310 | @*/ | ||
| 311 | ✗ | PetscErrorCode DSDuplicate(DS ds,DS *dsnew) | |
| 312 | { | ||
| 313 | ✗ | PetscFunctionBegin; | |
| 314 | ✗ | PetscValidHeaderSpecific(ds,DS_CLASSID,1); | |
| 315 | ✗ | PetscAssertPointer(dsnew,2); | |
| 316 | ✗ | PetscCall(DSCreate(PetscObjectComm((PetscObject)ds),dsnew)); | |
| 317 | ✗ | if (((PetscObject)ds)->type_name) PetscCall(DSSetType(*dsnew,((PetscObject)ds)->type_name)); | |
| 318 | ✗ | (*dsnew)->method = ds->method; | |
| 319 | ✗ | (*dsnew)->compact = ds->compact; | |
| 320 | ✗ | (*dsnew)->refined = ds->refined; | |
| 321 | ✗ | (*dsnew)->extrarow = ds->extrarow; | |
| 322 | ✗ | (*dsnew)->bs = ds->bs; | |
| 323 | ✗ | (*dsnew)->pmode = ds->pmode; | |
| 324 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
| 325 | } | ||
| 326 | |||
| 327 | /*@ | ||
| 328 | DSSetMethod - Selects the method to be used to solve the problem. | ||
| 329 | |||
| 330 | Logically Collective | ||
| 331 | |||
| 332 | Input Parameters: | ||
| 333 | + ds - the direct solver context | ||
| 334 | - meth - an index identifying the method | ||
| 335 | |||
| 336 | Options Database Key: | ||
| 337 | . -ds_method meth - sets the method | ||
| 338 | |||
| 339 | Level: intermediate | ||
| 340 | |||
| 341 | .seealso: [](sec:ds), `DSGetMethod()` | ||
| 342 | @*/ | ||
| 343 | 679 | PetscErrorCode DSSetMethod(DS ds,PetscInt meth) | |
| 344 | { | ||
| 345 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
679 | PetscFunctionBegin; |
| 346 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
679 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 347 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
679 | PetscValidLogicalCollectiveInt(ds,meth,2); |
| 348 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
679 | PetscCheck(meth>=0,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The method must be a non-negative integer"); |
| 349 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
679 | PetscCheck(meth<=DS_MAX_SOLVE,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Too large value for the method"); |
| 350 | 679 | ds->method = meth; | |
| 351 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
679 | PetscFunctionReturn(PETSC_SUCCESS); |
| 352 | } | ||
| 353 | |||
| 354 | /*@ | ||
| 355 | DSGetMethod - Gets the method currently used in the `DS`. | ||
| 356 | |||
| 357 | Not Collective | ||
| 358 | |||
| 359 | Input Parameter: | ||
| 360 | . ds - the direct solver context | ||
| 361 | |||
| 362 | Output Parameter: | ||
| 363 | . meth - identifier of the method | ||
| 364 | |||
| 365 | Level: intermediate | ||
| 366 | |||
| 367 | .seealso: [](sec:ds), `DSSetMethod()` | ||
| 368 | @*/ | ||
| 369 | 227 | PetscErrorCode DSGetMethod(DS ds,PetscInt *meth) | |
| 370 | { | ||
| 371 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
227 | PetscFunctionBegin; |
| 372 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
227 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 373 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
227 | PetscAssertPointer(meth,2); |
| 374 | 227 | *meth = ds->method; | |
| 375 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
227 | PetscFunctionReturn(PETSC_SUCCESS); |
| 376 | } | ||
| 377 | |||
| 378 | /*@ | ||
| 379 | DSSetParallel - Selects the mode of operation in parallel runs. | ||
| 380 | |||
| 381 | Logically Collective | ||
| 382 | |||
| 383 | Input Parameters: | ||
| 384 | + ds - the direct solver context | ||
| 385 | - pmode - the parallel mode | ||
| 386 | |||
| 387 | Options Database Key: | ||
| 388 | . -ds_parallel (redundant|synchronized|distributed) - sets the parallel mode | ||
| 389 | |||
| 390 | Notes: | ||
| 391 | In the `redundant` parallel mode, all processes will make the computation | ||
| 392 | redundantly, starting from the same data, and producing the same result. | ||
| 393 | This result may be slightly different in the different processes if using a | ||
| 394 | multithreaded BLAS library, which may cause issues in ill-conditioned problems. | ||
| 395 | |||
| 396 | In the `synchronized` parallel mode, only the first MPI process performs the | ||
| 397 | computation and then the computed quantities are broadcast to the other | ||
| 398 | processes in the communicator. This communication is not done automatically, | ||
| 399 | an explicit call to `DSSynchronize()` is required. | ||
| 400 | |||
| 401 | The `distributed` parallel mode can be used in some `DS` types only, such | ||
| 402 | as the contour integral method of `DSNEP`. In this case, every MPI process | ||
| 403 | will be in charge of part of the computation. | ||
| 404 | |||
| 405 | Level: advanced | ||
| 406 | |||
| 407 | .seealso: [](sec:ds), `DSSynchronize()`, `DSGetParallel()` | ||
| 408 | @*/ | ||
| 409 | 837 | PetscErrorCode DSSetParallel(DS ds,DSParallelType pmode) | |
| 410 | { | ||
| 411 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
837 | PetscFunctionBegin; |
| 412 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
837 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 413 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
837 | PetscValidLogicalCollectiveEnum(ds,pmode,2); |
| 414 | 837 | ds->pmode = pmode; | |
| 415 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
837 | PetscFunctionReturn(PETSC_SUCCESS); |
| 416 | } | ||
| 417 | |||
| 418 | /*@ | ||
| 419 | DSGetParallel - Gets the mode of operation in parallel runs. | ||
| 420 | |||
| 421 | Not Collective | ||
| 422 | |||
| 423 | Input Parameter: | ||
| 424 | . ds - the direct solver context | ||
| 425 | |||
| 426 | Output Parameter: | ||
| 427 | . pmode - the parallel mode | ||
| 428 | |||
| 429 | Level: advanced | ||
| 430 | |||
| 431 | .seealso: [](sec:ds), `DSSetParallel()` | ||
| 432 | @*/ | ||
| 433 | 197 | PetscErrorCode DSGetParallel(DS ds,DSParallelType *pmode) | |
| 434 | { | ||
| 435 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
197 | PetscFunctionBegin; |
| 436 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
197 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 437 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
197 | PetscAssertPointer(pmode,2); |
| 438 | 197 | *pmode = ds->pmode; | |
| 439 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
197 | PetscFunctionReturn(PETSC_SUCCESS); |
| 440 | } | ||
| 441 | |||
| 442 | /*@ | ||
| 443 | DSSetCompact - Switch to compact storage of matrices. | ||
| 444 | |||
| 445 | Logically Collective | ||
| 446 | |||
| 447 | Input Parameters: | ||
| 448 | + ds - the direct solver context | ||
| 449 | - comp - a boolean flag | ||
| 450 | |||
| 451 | Notes: | ||
| 452 | Compact storage is used in some `DS` types such as `DSHEP` when the matrix | ||
| 453 | is tridiagonal. This flag can be used to indicate whether the user | ||
| 454 | provides the matrix entries via the compact form (the tridiagonal `DS_MAT_T`) | ||
| 455 | or the non-compact one (`DS_MAT_A`). | ||
| 456 | |||
| 457 | The default is `PETSC_FALSE`. | ||
| 458 | |||
| 459 | Level: advanced | ||
| 460 | |||
| 461 | .seealso: [](sec:ds), `DSGetCompact()` | ||
| 462 | @*/ | ||
| 463 | 6569 | PetscErrorCode DSSetCompact(DS ds,PetscBool comp) | |
| 464 | { | ||
| 465 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
6569 | PetscFunctionBegin; |
| 466 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
6569 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 467 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
6569 | PetscValidLogicalCollectiveBool(ds,comp,2); |
| 468 |
10/12✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 8 times.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
6569 | if (ds->compact != comp && ds->ld) PetscTryTypeMethod(ds,setcompact,comp); |
| 469 | 6569 | ds->compact = comp; | |
| 470 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
6569 | PetscFunctionReturn(PETSC_SUCCESS); |
| 471 | } | ||
| 472 | |||
| 473 | /*@ | ||
| 474 | DSGetCompact - Gets the compact storage flag. | ||
| 475 | |||
| 476 | Not Collective | ||
| 477 | |||
| 478 | Input Parameter: | ||
| 479 | . ds - the direct solver context | ||
| 480 | |||
| 481 | Output Parameter: | ||
| 482 | . comp - the flag | ||
| 483 | |||
| 484 | Level: advanced | ||
| 485 | |||
| 486 | .seealso: [](sec:ds), `DSSetCompact()` | ||
| 487 | @*/ | ||
| 488 | 216 | PetscErrorCode DSGetCompact(DS ds,PetscBool *comp) | |
| 489 | { | ||
| 490 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
216 | PetscFunctionBegin; |
| 491 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
216 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 492 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
216 | PetscAssertPointer(comp,2); |
| 493 | 216 | *comp = ds->compact; | |
| 494 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
216 | PetscFunctionReturn(PETSC_SUCCESS); |
| 495 | } | ||
| 496 | |||
| 497 | /*@ | ||
| 498 | DSSetExtraRow - Sets a flag to indicate that the matrix has one extra | ||
| 499 | row. | ||
| 500 | |||
| 501 | Logically Collective | ||
| 502 | |||
| 503 | Input Parameters: | ||
| 504 | + ds - the direct solver context | ||
| 505 | - ext - a boolean flag | ||
| 506 | |||
| 507 | Notes: | ||
| 508 | In Krylov methods it is useful that the matrix representing the direct solver | ||
| 509 | has one extra row, i.e., has dimension $(n+1) \times n$. If this flag is activated, all | ||
| 510 | transformations applied to the right of the matrix also affect this additional | ||
| 511 | row. In that case, $(n+1)$ must be less or equal than the leading dimension. | ||
| 512 | |||
| 513 | The default is `PETSC_FALSE`. | ||
| 514 | |||
| 515 | Level: advanced | ||
| 516 | |||
| 517 | .seealso: [](sec:ds), `DSSolve()`, `DSAllocate()`, `DSGetExtraRow()` | ||
| 518 | @*/ | ||
| 519 | 9898 | PetscErrorCode DSSetExtraRow(DS ds,PetscBool ext) | |
| 520 | { | ||
| 521 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
9898 | PetscFunctionBegin; |
| 522 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
9898 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 523 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
9898 | PetscValidLogicalCollectiveBool(ds,ext,2); |
| 524 |
5/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 10 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
9898 | PetscCheck(!ext || ds->n==0 || ds->n!=ds->ld,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"Cannot set extra row after setting n=ld"); |
| 525 | 9898 | ds->extrarow = ext; | |
| 526 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
9898 | PetscFunctionReturn(PETSC_SUCCESS); |
| 527 | } | ||
| 528 | |||
| 529 | /*@ | ||
| 530 | DSGetExtraRow - Gets the extra row flag. | ||
| 531 | |||
| 532 | Not Collective | ||
| 533 | |||
| 534 | Input Parameter: | ||
| 535 | . ds - the direct solver context | ||
| 536 | |||
| 537 | Output Parameter: | ||
| 538 | . ext - the flag | ||
| 539 | |||
| 540 | Level: advanced | ||
| 541 | |||
| 542 | .seealso: [](sec:ds), `DSSetExtraRow()` | ||
| 543 | @*/ | ||
| 544 | 29260 | PetscErrorCode DSGetExtraRow(DS ds,PetscBool *ext) | |
| 545 | { | ||
| 546 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
29260 | PetscFunctionBegin; |
| 547 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
29260 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 548 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
29260 | PetscAssertPointer(ext,2); |
| 549 | 29260 | *ext = ds->extrarow; | |
| 550 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
29260 | PetscFunctionReturn(PETSC_SUCCESS); |
| 551 | } | ||
| 552 | |||
| 553 | /*@ | ||
| 554 | DSSetRefined - Sets a flag to indicate that refined vectors must be | ||
| 555 | computed. | ||
| 556 | |||
| 557 | Logically Collective | ||
| 558 | |||
| 559 | Input Parameters: | ||
| 560 | + ds - the direct solver context | ||
| 561 | - ref - a boolean flag | ||
| 562 | |||
| 563 | Notes: | ||
| 564 | Normally the vectors returned in `DS_MAT_X` are eigenvectors of the | ||
| 565 | projected matrix. With this flag activated, `DSVectors()` will return | ||
| 566 | the right singular vector of the smallest singular value of matrix | ||
| 567 | $\tilde{A}-\theta I$, where $\tilde{A}$ is the extended $(n+1)\times n$ | ||
| 568 | matrix and $\theta$ is the Ritz value. This is used in the refined Ritz | ||
| 569 | approximation {cite:p}`Jia97`. | ||
| 570 | |||
| 571 | The default is `PETSC_FALSE`. | ||
| 572 | |||
| 573 | Level: advanced | ||
| 574 | |||
| 575 | .seealso: [](sec:ds), `DSVectors()`, `DSGetRefined()` | ||
| 576 | @*/ | ||
| 577 | 10 | PetscErrorCode DSSetRefined(DS ds,PetscBool ref) | |
| 578 | { | ||
| 579 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 580 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 581 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
10 | PetscValidLogicalCollectiveBool(ds,ref,2); |
| 582 | 10 | ds->refined = ref; | |
| 583 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 584 | } | ||
| 585 | |||
| 586 | /*@ | ||
| 587 | DSGetRefined - Gets the refined vectors flag. | ||
| 588 | |||
| 589 | Not Collective | ||
| 590 | |||
| 591 | Input Parameter: | ||
| 592 | . ds - the direct solver context | ||
| 593 | |||
| 594 | Output Parameter: | ||
| 595 | . ref - the flag | ||
| 596 | |||
| 597 | Level: advanced | ||
| 598 | |||
| 599 | .seealso: [](sec:ds), `DSSetRefined()` | ||
| 600 | @*/ | ||
| 601 | 56646 | PetscErrorCode DSGetRefined(DS ds,PetscBool *ref) | |
| 602 | { | ||
| 603 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
56646 | PetscFunctionBegin; |
| 604 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
56646 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 605 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
56646 | PetscAssertPointer(ref,2); |
| 606 | 56646 | *ref = ds->refined; | |
| 607 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
56646 | PetscFunctionReturn(PETSC_SUCCESS); |
| 608 | } | ||
| 609 | |||
| 610 | /*@ | ||
| 611 | DSSetBlockSize - Sets the block size. | ||
| 612 | |||
| 613 | Logically Collective | ||
| 614 | |||
| 615 | Input Parameters: | ||
| 616 | + ds - the direct solver context | ||
| 617 | - bs - the block size | ||
| 618 | |||
| 619 | Options Database Key: | ||
| 620 | . -ds_block_size bs - sets the block size | ||
| 621 | |||
| 622 | Level: intermediate | ||
| 623 | |||
| 624 | .seealso: [](sec:ds), `DSGetBlockSize()` | ||
| 625 | @*/ | ||
| 626 | 10 | PetscErrorCode DSSetBlockSize(DS ds,PetscInt bs) | |
| 627 | { | ||
| 628 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 629 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 630 |
27/62✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 1 times.
✓ Branch 16 taken 1 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 1 times.
✓ Branch 20 taken 1 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 1 times.
✓ Branch 24 taken 1 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 1 times.
✓ Branch 28 taken 1 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 1 times.
✓ Branch 32 taken 1 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 1 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 1 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 1 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 1 times.
✓ Branch 46 taken 1 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 1 times.
✓ Branch 50 taken 1 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 1 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 1 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 1 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
10 | PetscValidLogicalCollectiveInt(ds,bs,2); |
| 631 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
10 | PetscCheck(bs>0,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"The block size must be at least one"); |
| 632 | 10 | ds->bs = bs; | |
| 633 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 634 | } | ||
| 635 | |||
| 636 | /*@ | ||
| 637 | DSGetBlockSize - Gets the block size. | ||
| 638 | |||
| 639 | Not Collective | ||
| 640 | |||
| 641 | Input Parameter: | ||
| 642 | . ds - the direct solver context | ||
| 643 | |||
| 644 | Output Parameter: | ||
| 645 | . bs - block size | ||
| 646 | |||
| 647 | Level: intermediate | ||
| 648 | |||
| 649 | .seealso: [](sec:ds), `DSSetBlockSize()` | ||
| 650 | @*/ | ||
| 651 | ✗ | PetscErrorCode DSGetBlockSize(DS ds,PetscInt *bs) | |
| 652 | { | ||
| 653 | ✗ | PetscFunctionBegin; | |
| 654 | ✗ | PetscValidHeaderSpecific(ds,DS_CLASSID,1); | |
| 655 | ✗ | PetscAssertPointer(bs,2); | |
| 656 | ✗ | *bs = ds->bs; | |
| 657 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
| 658 | } | ||
| 659 | |||
| 660 | /*@C | ||
| 661 | DSSetSlepcSC - Sets the sorting criterion context. | ||
| 662 | |||
| 663 | Logically Collective | ||
| 664 | |||
| 665 | Input Parameters: | ||
| 666 | + ds - the direct solver context | ||
| 667 | - sc - a pointer to the sorting criterion context | ||
| 668 | |||
| 669 | Note: | ||
| 670 | Not available in Fortran. | ||
| 671 | |||
| 672 | Level: developer | ||
| 673 | |||
| 674 | .seealso: [](sec:ds), `DSGetSlepcSC()`, `DSSort()` | ||
| 675 | @*/ | ||
| 676 | 10 | PetscErrorCode DSSetSlepcSC(DS ds,SlepcSC sc) | |
| 677 | { | ||
| 678 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBegin; |
| 679 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
10 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 680 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
10 | PetscAssertPointer(sc,2); |
| 681 |
1/12✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
10 | if (ds->sc && !ds->scset) PetscCall(PetscFree(ds->sc)); |
| 682 | 10 | ds->sc = sc; | |
| 683 | 10 | ds->scset = PETSC_TRUE; | |
| 684 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
10 | PetscFunctionReturn(PETSC_SUCCESS); |
| 685 | } | ||
| 686 | |||
| 687 | /*@C | ||
| 688 | DSGetSlepcSC - Gets the sorting criterion context. | ||
| 689 | |||
| 690 | Not Collective | ||
| 691 | |||
| 692 | Input Parameter: | ||
| 693 | . ds - the direct solver context | ||
| 694 | |||
| 695 | Output Parameter: | ||
| 696 | . sc - a pointer to the sorting criterion context | ||
| 697 | |||
| 698 | Note: | ||
| 699 | Not available in Fortran. | ||
| 700 | |||
| 701 | Level: developer | ||
| 702 | |||
| 703 | .seealso: [](sec:ds), `DSSetSlepcSC()`, `DSSort()` | ||
| 704 | @*/ | ||
| 705 | 17575 | PetscErrorCode DSGetSlepcSC(DS ds,SlepcSC *sc) | |
| 706 | { | ||
| 707 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
17575 | PetscFunctionBegin; |
| 708 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
17575 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 709 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
17575 | PetscAssertPointer(sc,2); |
| 710 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
17575 | if (!ds->sc) PetscCall(PetscNew(&ds->sc)); |
| 711 | 17575 | *sc = ds->sc; | |
| 712 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
17575 | PetscFunctionReturn(PETSC_SUCCESS); |
| 713 | } | ||
| 714 | |||
| 715 | /*@ | ||
| 716 | DSSetFromOptions - Sets `DS` options from the options database. | ||
| 717 | |||
| 718 | Collective | ||
| 719 | |||
| 720 | Input Parameter: | ||
| 721 | . ds - the direct solver context | ||
| 722 | |||
| 723 | Notes: | ||
| 724 | To see all options, run your program with the `-help` option. | ||
| 725 | |||
| 726 | Level: beginner | ||
| 727 | |||
| 728 | .seealso: [](sec:ds), `DSSetOptionsPrefix()` | ||
| 729 | @*/ | ||
| 730 | 12501 | PetscErrorCode DSSetFromOptions(DS ds) | |
| 731 | { | ||
| 732 | 12501 | PetscInt bs,meth; | |
| 733 | 12501 | PetscBool flag; | |
| 734 | 12501 | DSParallelType pmode; | |
| 735 | |||
| 736 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
12501 | PetscFunctionBegin; |
| 737 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
12501 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 738 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12501 | PetscCall(DSRegisterAll()); |
| 739 | /* Set default type (we do not allow changing it with -ds_type) */ | ||
| 740 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
12501 | if (!((PetscObject)ds)->type_name) PetscCall(DSSetType(ds,DSNHEP)); |
| 741 |
8/10✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
|
37503 | PetscObjectOptionsBegin((PetscObject)ds); |
| 742 | |||
| 743 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12501 | PetscCall(PetscOptionsInt("-ds_block_size","Block size for the dense system solver","DSSetBlockSize",ds->bs,&bs,&flag)); |
| 744 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
12501 | if (flag) PetscCall(DSSetBlockSize(ds,bs)); |
| 745 | |||
| 746 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12501 | PetscCall(PetscOptionsInt("-ds_method","Method to be used for the dense system","DSSetMethod",ds->method,&meth,&flag)); |
| 747 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
12501 | if (flag) PetscCall(DSSetMethod(ds,meth)); |
| 748 | |||
| 749 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12501 | PetscCall(PetscOptionsEnum("-ds_parallel","Operation mode in parallel runs","DSSetParallel",DSParallelTypes,(PetscEnum)ds->pmode,(PetscEnum*)&pmode,&flag)); |
| 750 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
12501 | if (flag) PetscCall(DSSetParallel(ds,pmode)); |
| 751 | |||
| 752 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
12501 | PetscTryTypeMethod(ds,setfromoptions,PetscOptionsObject); |
| 753 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12501 | PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)ds,PetscOptionsObject)); |
| 754 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12501 | PetscOptionsEnd(); |
| 755 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
12501 | PetscFunctionReturn(PETSC_SUCCESS); |
| 756 | } | ||
| 757 | |||
| 758 | /*@ | ||
| 759 | DSView - Prints the `DS` data structure. | ||
| 760 | |||
| 761 | Collective | ||
| 762 | |||
| 763 | Input Parameters: | ||
| 764 | + ds - the direct solver context | ||
| 765 | - viewer - optional visualization context | ||
| 766 | |||
| 767 | Notes: | ||
| 768 | The available visualization contexts include | ||
| 769 | + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) | ||
| 770 | - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard output where only the | ||
| 771 | first process opens the file; all other processes send their data to the | ||
| 772 | first one to print | ||
| 773 | |||
| 774 | The user can open an alternative visualization context with `PetscViewerASCIIOpen()` | ||
| 775 | to output to a specified file. | ||
| 776 | |||
| 777 | Use `DSViewFromOptions()` to allow the user to select many different `PetscViewerType` | ||
| 778 | and formats from the options database. | ||
| 779 | |||
| 780 | Level: beginner | ||
| 781 | |||
| 782 | .seealso: [](sec:ds), `DSCreate()`, `DSViewFromOptions()`, `DSViewMat()` | ||
| 783 | @*/ | ||
| 784 | 1243 | PetscErrorCode DSView(DS ds,PetscViewer viewer) | |
| 785 | { | ||
| 786 | 1243 | PetscBool isascii; | |
| 787 | 1243 | PetscViewerFormat format; | |
| 788 | 1243 | PetscMPIInt size; | |
| 789 | |||
| 790 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1243 | PetscFunctionBegin; |
| 791 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
1243 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 792 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1243 | if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ds),&viewer)); |
| 793 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
1243 | PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); |
| 794 |
13/32✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
|
1243 | PetscCheckSameComm(ds,1,viewer,2); |
| 795 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1243 | PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); |
| 796 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
1243 | if (isascii) { |
| 797 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1243 | PetscCall(PetscViewerGetFormat(viewer,&format)); |
| 798 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1243 | PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ds,viewer)); |
| 799 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
1243 | PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ds),&size)); |
| 800 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1243 | if (size>1) PetscCall(PetscViewerASCIIPrintf(viewer," parallel operation mode: %s\n",DSParallelTypes[ds->pmode])); |
| 801 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1243 | if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { |
| 802 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
655 | PetscCall(PetscViewerASCIIPrintf(viewer," current state: %s\n",DSStateTypes[ds->state])); |
| 803 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
655 | PetscCall(PetscViewerASCIIPrintf(viewer," dimensions: ld=%" PetscInt_FMT ", n=%" PetscInt_FMT ", l=%" PetscInt_FMT ", k=%" PetscInt_FMT,ds->ld,ds->n,ds->l,ds->k)); |
| 804 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
655 | if (ds->state==DS_STATE_TRUNCATED) PetscCall(PetscViewerASCIIPrintf(viewer,", t=%" PetscInt_FMT "\n",ds->t)); |
| 805 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
655 | else PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); |
| 806 |
10/12✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 10 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2170 | PetscCall(PetscViewerASCIIPrintf(viewer," flags:%s%s%s\n",ds->compact?" compact":"",ds->extrarow?" extrarow":"",ds->refined?" refined":"")); |
| 807 | } | ||
| 808 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1243 | PetscCall(PetscViewerASCIIPushTab(viewer)); |
| 809 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
1243 | PetscTryTypeMethod(ds,view,viewer); |
| 810 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1243 | PetscCall(PetscViewerASCIIPopTab(viewer)); |
| 811 | } | ||
| 812 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
245 | PetscFunctionReturn(PETSC_SUCCESS); |
| 813 | } | ||
| 814 | |||
| 815 | /*@ | ||
| 816 | DSViewFromOptions - View (print) a `DS` object based on values in the options database. | ||
| 817 | |||
| 818 | Collective | ||
| 819 | |||
| 820 | Input Parameters: | ||
| 821 | + ds - the direct solver context | ||
| 822 | . obj - optional object that provides the options prefix used to query the options database | ||
| 823 | - name - command line option | ||
| 824 | |||
| 825 | Level: intermediate | ||
| 826 | |||
| 827 | .seealso: [](sec:ds), `DSView()`, `DSCreate()`, `PetscObjectViewFromOptions()` | ||
| 828 | @*/ | ||
| 829 | 15 | PetscErrorCode DSViewFromOptions(DS ds,PetscObject obj,const char name[]) | |
| 830 | { | ||
| 831 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
15 | PetscFunctionBegin; |
| 832 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
15 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 833 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
15 | PetscCall(PetscObjectViewFromOptions((PetscObject)ds,obj,name)); |
| 834 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3 | PetscFunctionReturn(PETSC_SUCCESS); |
| 835 | } | ||
| 836 | |||
| 837 | /*@ | ||
| 838 | DSAllocate - Allocates memory for internal storage or matrices in `DS`. | ||
| 839 | |||
| 840 | Logically Collective | ||
| 841 | |||
| 842 | Input Parameters: | ||
| 843 | + ds - the direct solver context | ||
| 844 | - ld - leading dimension (maximum allowed dimension for the matrices, including | ||
| 845 | the extra row if present) | ||
| 846 | |||
| 847 | Note: | ||
| 848 | If the leading dimension is different from a previously set value, then | ||
| 849 | all matrices are destroyed with `DSReset()`. | ||
| 850 | |||
| 851 | Level: intermediate | ||
| 852 | |||
| 853 | .seealso: [](sec:ds), `DSGetLeadingDimension()`, `DSSetDimensions()`, `DSSetExtraRow()`, `DSReset()`, `DSReallocate()` | ||
| 854 | @*/ | ||
| 855 | 15652 | PetscErrorCode DSAllocate(DS ds,PetscInt ld) | |
| 856 | { | ||
| 857 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
15652 | PetscFunctionBegin; |
| 858 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
15652 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 859 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
15652 | PetscValidLogicalCollectiveInt(ds,ld,2); |
| 860 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
15652 | PetscValidType(ds,1); |
| 861 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
15652 | PetscCheck(ld>0,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"Leading dimension should be at least one"); |
| 862 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
15652 | if (ld!=ds->ld) { |
| 863 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12200 | PetscCall(PetscInfo(ds,"Allocating memory with leading dimension=%" PetscInt_FMT "\n",ld)); |
| 864 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
12200 | PetscCall(DSReset(ds)); |
| 865 | 12200 | ds->ld = ld; | |
| 866 |
5/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
12200 | PetscUseTypeMethod(ds,allocate,ld); |
| 867 | } | ||
| 868 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
3006 | PetscFunctionReturn(PETSC_SUCCESS); |
| 869 | } | ||
| 870 | |||
| 871 | /*@ | ||
| 872 | DSReallocate - Reallocates memory for internal storage or matrices in `DS`, | ||
| 873 | keeping the previously set data. | ||
| 874 | |||
| 875 | Logically Collective | ||
| 876 | |||
| 877 | Input Parameters: | ||
| 878 | + ds - the direct solver context | ||
| 879 | - ld - new leading dimension | ||
| 880 | |||
| 881 | Notes: | ||
| 882 | The new leading dimension must be larger than the previous one. The relevant | ||
| 883 | data previously set is copied over to the new data structures. | ||
| 884 | |||
| 885 | This operation is not available in all `DS` types. | ||
| 886 | |||
| 887 | Level: developer | ||
| 888 | |||
| 889 | .seealso: [](sec:ds), `DSAllocate()` | ||
| 890 | @*/ | ||
| 891 | 133 | PetscErrorCode DSReallocate(DS ds,PetscInt ld) | |
| 892 | { | ||
| 893 | 133 | PetscInt i; | |
| 894 | |||
| 895 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
133 | PetscFunctionBegin; |
| 896 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
133 | PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 897 |
27/62✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
|
133 | PetscValidLogicalCollectiveInt(ds,ld,2); |
| 898 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
133 | PetscValidType(ds,1); |
| 899 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
133 | PetscCheck(ds->ld,PetscObjectComm((PetscObject)ds),PETSC_ERR_ORDER,"DSAllocate() must be called first"); |
| 900 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
133 | PetscCheck(ld>ds->ld,PetscObjectComm((PetscObject)ds),PETSC_ERR_ARG_OUTOFRANGE,"New leading dimension %" PetscInt_FMT " must be larger than the previous one %" PetscInt_FMT,ld,ds->ld); |
| 901 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
133 | PetscCall(PetscInfo(ds,"Reallocating memory with new leading dimension=%" PetscInt_FMT "\n",ld)); |
| 902 |
5/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
133 | PetscUseTypeMethod(ds,reallocate,ld); |
| 903 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1041 | for (i=ds->ld;i<ld;i++) ds->perm[i] = i; |
| 904 | 133 | ds->ld = ld; | |
| 905 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
133 | PetscFunctionReturn(PETSC_SUCCESS); |
| 906 | } | ||
| 907 | |||
| 908 | /*@ | ||
| 909 | DSReset - Resets the `DS` context to the initial state. | ||
| 910 | |||
| 911 | Collective | ||
| 912 | |||
| 913 | Input Parameter: | ||
| 914 | . ds - the direct solver context | ||
| 915 | |||
| 916 | Note: | ||
| 917 | All data structures with size depending on the leading dimension | ||
| 918 | of `DSAllocate()` are released. | ||
| 919 | |||
| 920 | Level: advanced | ||
| 921 | |||
| 922 | .seealso: [](sec:ds), `DSDestroy()`, `DSAllocate()` | ||
| 923 | @*/ | ||
| 924 | 43797 | PetscErrorCode DSReset(DS ds) | |
| 925 | { | ||
| 926 | 43797 | PetscInt i; | |
| 927 | |||
| 928 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
43797 | PetscFunctionBegin; |
| 929 |
3/14✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
43797 | if (ds) PetscValidHeaderSpecific(ds,DS_CLASSID,1); |
| 930 |
1/12✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
35342 | if (!ds) PetscFunctionReturn(PETSC_SUCCESS); |
| 931 | 43797 | ds->state = DS_STATE_RAW; | |
| 932 | 43797 | ds->ld = 0; | |
| 933 | 43797 | ds->l = 0; | |
| 934 | 43797 | ds->n = 0; | |
| 935 | 43797 | ds->k = 0; | |
| 936 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
1007331 | for (i=0;i<DS_NUM_MAT;i++) PetscCall(MatDestroy(&ds->omat[i])); |
| 937 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
43797 | PetscCall(PetscFree(ds->perm)); |
| 938 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
8455 | PetscFunctionReturn(PETSC_SUCCESS); |
| 939 | } | ||
| 940 | |||
| 941 | /*@ | ||
| 942 | DSDestroy - Destroys `DS` context that was created with `DSCreate()`. | ||
| 943 | |||
| 944 | Collective | ||
| 945 | |||
| 946 | Input Parameter: | ||
| 947 | . ds - the direct solver context | ||
| 948 | |||
| 949 | Level: beginner | ||
| 950 | |||
| 951 | .seealso: [](sec:ds), `DSCreate()` | ||
| 952 | @*/ | ||
| 953 | 14341 | PetscErrorCode DSDestroy(DS *ds) | |
| 954 | { | ||
| 955 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
14341 | PetscFunctionBegin; |
| 956 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
14341 | if (!*ds) PetscFunctionReturn(PETSC_SUCCESS); |
| 957 |
2/12✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
13342 | PetscValidHeaderSpecific(*ds,DS_CLASSID,1); |
| 958 |
1/14✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
|
13342 | if (--((PetscObject)*ds)->refct > 0) { *ds = NULL; PetscFunctionReturn(PETSC_SUCCESS); } |
| 959 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13342 | PetscCall(DSReset(*ds)); |
| 960 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
13342 | PetscTryTypeMethod(*ds,destroy); |
| 961 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
13342 | PetscCall(PetscFree((*ds)->work)); |
| 962 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
13342 | PetscCall(PetscFree((*ds)->rwork)); |
| 963 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
13342 | PetscCall(PetscFree((*ds)->iwork)); |
| 964 |
7/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
13342 | if (!(*ds)->scset) PetscCall(PetscFree((*ds)->sc)); |
| 965 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
13342 | PetscCall(PetscHeaderDestroy(ds)); |
| 966 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2579 | PetscFunctionReturn(PETSC_SUCCESS); |
| 967 | } | ||
| 968 | |||
| 969 | /*@C | ||
| 970 | DSRegister - Adds a direct solver to the `DS` package. | ||
| 971 | |||
| 972 | Not Collective | ||
| 973 | |||
| 974 | Input Parameters: | ||
| 975 | + name - name of a new user-defined `DS` | ||
| 976 | - function - routine to create context | ||
| 977 | |||
| 978 | Note: | ||
| 979 | `DSRegister()` may be called multiple times to add several user-defined | ||
| 980 | direct solvers. | ||
| 981 | |||
| 982 | Level: advanced | ||
| 983 | |||
| 984 | .seealso: [](sec:ds), `DSRegisterAll()` | ||
| 985 | @*/ | ||
| 986 | 121814 | PetscErrorCode DSRegister(const char *name,PetscErrorCode (*function)(DS)) | |
| 987 | { | ||
| 988 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
121814 | PetscFunctionBegin; |
| 989 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
121814 | PetscCall(DSInitializePackage()); |
| 990 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
121814 | PetscCall(PetscFunctionListAdd(&DSList,name,function)); |
| 991 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
23584 | PetscFunctionReturn(PETSC_SUCCESS); |
| 992 | } | ||
| 993 | |||
| 994 | SLEPC_EXTERN PetscErrorCode DSCreate_HEP(DS); | ||
| 995 | SLEPC_EXTERN PetscErrorCode DSCreate_NHEP(DS); | ||
| 996 | SLEPC_EXTERN PetscErrorCode DSCreate_GHEP(DS); | ||
| 997 | SLEPC_EXTERN PetscErrorCode DSCreate_GHIEP(DS); | ||
| 998 | SLEPC_EXTERN PetscErrorCode DSCreate_GNHEP(DS); | ||
| 999 | SLEPC_EXTERN PetscErrorCode DSCreate_NHEPTS(DS); | ||
| 1000 | SLEPC_EXTERN PetscErrorCode DSCreate_SVD(DS); | ||
| 1001 | SLEPC_EXTERN PetscErrorCode DSCreate_HSVD(DS); | ||
| 1002 | SLEPC_EXTERN PetscErrorCode DSCreate_GSVD(DS); | ||
| 1003 | SLEPC_EXTERN PetscErrorCode DSCreate_PEP(DS); | ||
| 1004 | SLEPC_EXTERN PetscErrorCode DSCreate_NEP(DS); | ||
| 1005 | |||
| 1006 | /*@C | ||
| 1007 | DSRegisterAll - Registers all of the direct solvers in the `DS` package. | ||
| 1008 | |||
| 1009 | Not Collective | ||
| 1010 | |||
| 1011 | Level: advanced | ||
| 1012 | |||
| 1013 | .seealso: [](sec:ds), `DSRegister()` | ||
| 1014 | @*/ | ||
| 1015 | 23575 | PetscErrorCode DSRegisterAll(void) | |
| 1016 | { | ||
| 1017 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
23575 | PetscFunctionBegin; |
| 1018 |
8/14✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
|
23575 | if (DSRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); |
| 1019 | 11074 | DSRegisterAllCalled = PETSC_TRUE; | |
| 1020 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSHEP,DSCreate_HEP)); |
| 1021 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSNHEP,DSCreate_NHEP)); |
| 1022 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSGHEP,DSCreate_GHEP)); |
| 1023 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSGHIEP,DSCreate_GHIEP)); |
| 1024 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSGNHEP,DSCreate_GNHEP)); |
| 1025 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSNHEPTS,DSCreate_NHEPTS)); |
| 1026 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSSVD,DSCreate_SVD)); |
| 1027 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSHSVD,DSCreate_HSVD)); |
| 1028 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSGSVD,DSCreate_GSVD)); |
| 1029 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSPEP,DSCreate_PEP)); |
| 1030 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
11074 | PetscCall(DSRegister(DSNEP,DSCreate_NEP)); |
| 1031 |
6/12✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
|
2144 | PetscFunctionReturn(PETSC_SUCCESS); |
| 1032 | } | ||
| 1033 |