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