1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 #include <petscdm.h> 3 static PetscErrorCode TSRHSSplitGetRHSSplit(TS ts, const char splitname[], TS_RHSSplitLink *isplit) 4 { 5 PetscBool found = PETSC_FALSE; 6 7 PetscFunctionBegin; 8 *isplit = ts->tsrhssplit; 9 /* look up the split */ 10 while (*isplit) { 11 PetscCall(PetscStrcmp((*isplit)->splitname, splitname, &found)); 12 if (found) break; 13 *isplit = (*isplit)->next; 14 } 15 PetscFunctionReturn(PETSC_SUCCESS); 16 } 17 18 /*@C 19 TSRHSSplitSetIS - Set the index set for the specified split 20 21 Logically Collective 22 23 Input Parameters: 24 + ts - the `TS` context obtained from `TSCreate()` 25 . splitname - name of this split, if NULL the number of the split is used 26 - is - the index set for part of the solution vector 27 28 Level: intermediate 29 30 .seealso: [](chapter_ts), `TS`, `IS`, `TSRHSSplitGetIS()` 31 @*/ 32 PetscErrorCode TSRHSSplitSetIS(TS ts, const char splitname[], IS is) 33 { 34 TS_RHSSplitLink newsplit, next = ts->tsrhssplit; 35 char prefix[128]; 36 37 PetscFunctionBegin; 38 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 39 PetscValidHeaderSpecific(is, IS_CLASSID, 3); 40 41 PetscCall(PetscNew(&newsplit)); 42 if (splitname) { 43 PetscCall(PetscStrallocpy(splitname, &newsplit->splitname)); 44 } else { 45 PetscCall(PetscMalloc1(8, &newsplit->splitname)); 46 PetscCall(PetscSNPrintf(newsplit->splitname, 7, "%" PetscInt_FMT, ts->num_rhs_splits)); 47 } 48 PetscCall(PetscObjectReference((PetscObject)is)); 49 newsplit->is = is; 50 PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &newsplit->ts)); 51 52 PetscCall(PetscObjectIncrementTabLevel((PetscObject)newsplit->ts, (PetscObject)ts, 1)); 53 PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%srhsplit_%s_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "", newsplit->splitname)); 54 PetscCall(TSSetOptionsPrefix(newsplit->ts, prefix)); 55 if (!next) ts->tsrhssplit = newsplit; 56 else { 57 while (next->next) next = next->next; 58 next->next = newsplit; 59 } 60 ts->num_rhs_splits++; 61 PetscFunctionReturn(PETSC_SUCCESS); 62 } 63 64 /*@C 65 TSRHSSplitGetIS - Retrieves the elements for a split as an `IS` 66 67 Logically Collective 68 69 Input Parameters: 70 + ts - the `TS` context obtained from `TSCreate()` 71 - splitname - name of this split 72 73 Output Parameters: 74 - is - the index set for part of the solution vector 75 76 Level: intermediate 77 78 .seealso: [](chapter_ts), `TS`, `IS`, `TSRHSSplitSetIS()` 79 @*/ 80 PetscErrorCode TSRHSSplitGetIS(TS ts, const char splitname[], IS *is) 81 { 82 TS_RHSSplitLink isplit; 83 84 PetscFunctionBegin; 85 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 86 *is = NULL; 87 /* look up the split */ 88 PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 89 if (isplit) *is = isplit->is; 90 PetscFunctionReturn(PETSC_SUCCESS); 91 } 92 93 /*@C 94 TSRHSSplitSetRHSFunction - Set the split right-hand-side functions. 95 96 Logically Collective 97 98 Input Parameters: 99 + ts - the `TS` context obtained from `TSCreate()` 100 . splitname - name of this split 101 . r - vector to hold the residual (or NULL to have it created internally) 102 . rhsfunc - the RHS function evaluation routine 103 - ctx - user-defined context for private data for the split function evaluation routine (may be NULL) 104 105 Calling sequence of fun: 106 $ rhsfunc(TS ts,PetscReal t,Vec u,Vec f,ctx); 107 108 + t - time at step/stage being solved 109 . u - state vector 110 . f - function vector 111 - ctx - [optional] user-defined context for matrix evaluation routine (may be NULL) 112 113 Level: beginner 114 115 .seealso: [](chapter_ts), `TS`, `IS`, `TSRHSSplitSetIS()` 116 @*/ 117 PetscErrorCode TSRHSSplitSetRHSFunction(TS ts, const char splitname[], Vec r, TSRHSFunction rhsfunc, void *ctx) 118 { 119 TS_RHSSplitLink isplit; 120 DM dmc; 121 Vec subvec, ralloc = NULL; 122 123 PetscFunctionBegin; 124 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 125 if (r) PetscValidHeaderSpecific(r, VEC_CLASSID, 3); 126 127 /* look up the split */ 128 PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 129 PetscCheck(isplit, PETSC_COMM_SELF, PETSC_ERR_USER, "The split %s is not created, check the split name or call TSRHSSplitSetIS() to create one", splitname); 130 131 if (!r && ts->vec_sol) { 132 PetscCall(VecGetSubVector(ts->vec_sol, isplit->is, &subvec)); 133 PetscCall(VecDuplicate(subvec, &ralloc)); 134 r = ralloc; 135 PetscCall(VecRestoreSubVector(ts->vec_sol, isplit->is, &subvec)); 136 } 137 138 if (ts->dm) { 139 PetscInt dim; 140 141 PetscCall(DMGetDimension(ts->dm, &dim)); 142 if (dim != -1) { 143 PetscCall(DMClone(ts->dm, &dmc)); 144 PetscCall(TSSetDM(isplit->ts, dmc)); 145 PetscCall(DMDestroy(&dmc)); 146 } 147 } 148 149 PetscCall(TSSetRHSFunction(isplit->ts, r, rhsfunc, ctx)); 150 PetscCall(VecDestroy(&ralloc)); 151 PetscFunctionReturn(PETSC_SUCCESS); 152 } 153 154 /*@C 155 TSRHSSplitGetSubTS - Get the sub-`TS` by split name. 156 157 Logically Collective 158 159 Input Parameter: 160 . ts - the `TS` context obtained from `TSCreate()` 161 162 Output Parameters: 163 + splitname - the number of the split 164 - subts - the array of `TS` contexts 165 166 Level: advanced 167 168 .seealso: [](chapter_ts), `TS`, `IS`, `TSGetRHSSplitFunction()` 169 @*/ 170 PetscErrorCode TSRHSSplitGetSubTS(TS ts, const char splitname[], TS *subts) 171 { 172 TS_RHSSplitLink isplit; 173 174 PetscFunctionBegin; 175 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 176 PetscValidPointer(subts, 3); 177 *subts = NULL; 178 /* look up the split */ 179 PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 180 if (isplit) *subts = isplit->ts; 181 PetscFunctionReturn(PETSC_SUCCESS); 182 } 183 184 /*@C 185 TSRHSSplitGetSubTSs - Get an array of all sub-`TS` contexts. 186 187 Logically Collective 188 189 Input Parameter: 190 . ts - the `TS` context obtained from `TSCreate()` 191 192 Output Parameters: 193 + n - the number of splits 194 - subksp - the array of `TS` contexts 195 196 Level: advanced 197 198 Note: 199 After `TSRHSSplitGetSubTS()` the array of `TS`s is to be freed by the user with `PetscFree()` 200 (not the `TS` just the array that contains them). 201 202 .seealso: [](chapter_ts), `TS`, `IS`, `TSGetRHSSplitFunction()` 203 @*/ 204 PetscErrorCode TSRHSSplitGetSubTSs(TS ts, PetscInt *n, TS *subts[]) 205 { 206 TS_RHSSplitLink ilink = ts->tsrhssplit; 207 PetscInt i = 0; 208 209 PetscFunctionBegin; 210 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 211 if (subts) { 212 PetscCall(PetscMalloc1(ts->num_rhs_splits, subts)); 213 while (ilink) { 214 (*subts)[i++] = ilink->ts; 215 ilink = ilink->next; 216 } 217 } 218 if (n) *n = ts->num_rhs_splits; 219 PetscFunctionReturn(PETSC_SUCCESS); 220 } 221