1545aaa6fSHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2424fc49dSJoe Pusztay #include <petscdm.h> 3d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRHSSplitGetRHSSplit(TS ts, const char splitname[], TS_RHSSplitLink *isplit) 4d71ae5a4SJacob Faibussowitsch { 51d06f6b3SHong Zhang PetscBool found = PETSC_FALSE; 61d06f6b3SHong Zhang 71d06f6b3SHong Zhang PetscFunctionBegin; 81d06f6b3SHong Zhang *isplit = ts->tsrhssplit; 91d06f6b3SHong Zhang /* look up the split */ 101d06f6b3SHong Zhang while (*isplit) { 119566063dSJacob Faibussowitsch PetscCall(PetscStrcmp((*isplit)->splitname, splitname, &found)); 121d06f6b3SHong Zhang if (found) break; 131d06f6b3SHong Zhang *isplit = (*isplit)->next; 141d06f6b3SHong Zhang } 153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161d06f6b3SHong Zhang } 171d06f6b3SHong Zhang 18*cc4c1da9SBarry Smith /*@ 19545aaa6fSHong Zhang TSRHSSplitSetIS - Set the index set for the specified split 20545aaa6fSHong Zhang 21c3339decSBarry Smith Logically Collective 22545aaa6fSHong Zhang 23545aaa6fSHong Zhang Input Parameters: 24bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 2520f4b53cSBarry Smith . splitname - name of this split, if `NULL` the number of the split is used 26545aaa6fSHong Zhang - is - the index set for part of the solution vector 27545aaa6fSHong Zhang 28545aaa6fSHong Zhang Level: intermediate 29545aaa6fSHong Zhang 301cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `IS`, `TSRHSSplitGetIS()` 31545aaa6fSHong Zhang @*/ 32d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitSetIS(TS ts, const char splitname[], IS is) 33d71ae5a4SJacob Faibussowitsch { 341d06f6b3SHong Zhang TS_RHSSplitLink newsplit, next = ts->tsrhssplit; 35545aaa6fSHong Zhang char prefix[128]; 36545aaa6fSHong Zhang 37545aaa6fSHong Zhang PetscFunctionBegin; 38545aaa6fSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 39545aaa6fSHong Zhang PetscValidHeaderSpecific(is, IS_CLASSID, 3); 401d06f6b3SHong Zhang 419566063dSJacob Faibussowitsch PetscCall(PetscNew(&newsplit)); 42545aaa6fSHong Zhang if (splitname) { 439566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname, &newsplit->splitname)); 44545aaa6fSHong Zhang } else { 459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(8, &newsplit->splitname)); 4663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(newsplit->splitname, 7, "%" PetscInt_FMT, ts->num_rhs_splits)); 47545aaa6fSHong Zhang } 489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is)); 49545aaa6fSHong Zhang newsplit->is = is; 509566063dSJacob Faibussowitsch PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &newsplit->ts)); 51109dc152SHong Zhang 529566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)newsplit->ts, (PetscObject)ts, 1)); 539566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%srhsplit_%s_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "", newsplit->splitname)); 549566063dSJacob Faibussowitsch PetscCall(TSSetOptionsPrefix(newsplit->ts, prefix)); 551d06f6b3SHong Zhang if (!next) ts->tsrhssplit = newsplit; 561d06f6b3SHong Zhang else { 571d06f6b3SHong Zhang while (next->next) next = next->next; 581d06f6b3SHong Zhang next->next = newsplit; 591d06f6b3SHong Zhang } 601d06f6b3SHong Zhang ts->num_rhs_splits++; 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62545aaa6fSHong Zhang } 63545aaa6fSHong Zhang 64*cc4c1da9SBarry Smith /*@ 65bcf0153eSBarry Smith TSRHSSplitGetIS - Retrieves the elements for a split as an `IS` 66545aaa6fSHong Zhang 67c3339decSBarry Smith Logically Collective 68545aaa6fSHong Zhang 69545aaa6fSHong Zhang Input Parameters: 70bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 71545aaa6fSHong Zhang - splitname - name of this split 72545aaa6fSHong Zhang 732fe279fdSBarry Smith Output Parameter: 742fe279fdSBarry Smith . is - the index set for part of the solution vector 75545aaa6fSHong Zhang 76545aaa6fSHong Zhang Level: intermediate 77545aaa6fSHong Zhang 781cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `IS`, `TSRHSSplitSetIS()` 79545aaa6fSHong Zhang @*/ 80d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetIS(TS ts, const char splitname[], IS *is) 81d71ae5a4SJacob Faibussowitsch { 821d06f6b3SHong Zhang TS_RHSSplitLink isplit; 83545aaa6fSHong Zhang 84545aaa6fSHong Zhang PetscFunctionBegin; 85545aaa6fSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 861d06f6b3SHong Zhang *is = NULL; 87545aaa6fSHong Zhang /* look up the split */ 889566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 891d06f6b3SHong Zhang if (isplit) *is = isplit->is; 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91545aaa6fSHong Zhang } 92545aaa6fSHong Zhang 93545aaa6fSHong Zhang /*@C 94545aaa6fSHong Zhang TSRHSSplitSetRHSFunction - Set the split right-hand-side functions. 95545aaa6fSHong Zhang 96c3339decSBarry Smith Logically Collective 97545aaa6fSHong Zhang 98545aaa6fSHong Zhang Input Parameters: 99bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 100545aaa6fSHong Zhang . splitname - name of this split 10120f4b53cSBarry Smith . r - vector to hold the residual (or `NULL` to have it created internally) 102545aaa6fSHong Zhang . rhsfunc - the RHS function evaluation routine 10320f4b53cSBarry Smith - ctx - user-defined context for private data for the split function evaluation routine (may be `NULL`) 104545aaa6fSHong Zhang 10520f4b53cSBarry Smith Level: intermediate 106545aaa6fSHong Zhang 1078434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSRHSFunctionFn`, `IS`, `TSRHSSplitSetIS()` 108545aaa6fSHong Zhang @*/ 1098434afd1SBarry Smith PetscErrorCode TSRHSSplitSetRHSFunction(TS ts, const char splitname[], Vec r, TSRHSFunctionFn *rhsfunc, void *ctx) 110d71ae5a4SJacob Faibussowitsch { 1111d06f6b3SHong Zhang TS_RHSSplitLink isplit; 112424fc49dSJoe Pusztay DM dmc; 113545aaa6fSHong Zhang Vec subvec, ralloc = NULL; 114545aaa6fSHong Zhang 115545aaa6fSHong Zhang PetscFunctionBegin; 116545aaa6fSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 117064a246eSJacob Faibussowitsch if (r) PetscValidHeaderSpecific(r, VEC_CLASSID, 3); 118545aaa6fSHong Zhang 119545aaa6fSHong Zhang /* look up the split */ 1209566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 1213c633725SBarry Smith 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); 122545aaa6fSHong Zhang 1231d06f6b3SHong Zhang if (!r && ts->vec_sol) { 1249566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, isplit->is, &subvec)); 1259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(subvec, &ralloc)); 126545aaa6fSHong Zhang r = ralloc; 1279566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, isplit->is, &subvec)); 128545aaa6fSHong Zhang } 129424fc49dSJoe Pusztay 130424fc49dSJoe Pusztay if (ts->dm) { 131424fc49dSJoe Pusztay PetscInt dim; 132424fc49dSJoe Pusztay 1339566063dSJacob Faibussowitsch PetscCall(DMGetDimension(ts->dm, &dim)); 134424fc49dSJoe Pusztay if (dim != -1) { 1359566063dSJacob Faibussowitsch PetscCall(DMClone(ts->dm, &dmc)); 1369566063dSJacob Faibussowitsch PetscCall(TSSetDM(isplit->ts, dmc)); 1379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmc)); 138424fc49dSJoe Pusztay } 139424fc49dSJoe Pusztay } 140424fc49dSJoe Pusztay 1419566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(isplit->ts, r, rhsfunc, ctx)); 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ralloc)); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144545aaa6fSHong Zhang } 145545aaa6fSHong Zhang 146*cc4c1da9SBarry Smith /*@ 147bcf0153eSBarry Smith TSRHSSplitGetSubTS - Get the sub-`TS` by split name. 1481d06f6b3SHong Zhang 149c3339decSBarry Smith Logically Collective 1501d06f6b3SHong Zhang 1516b867d5aSJose E. Roman Input Parameter: 152bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1536b867d5aSJose E. Roman 1541d06f6b3SHong Zhang Output Parameters: 1551d06f6b3SHong Zhang + splitname - the number of the split 15620f4b53cSBarry Smith - subts - the sub-`TS` 1571d06f6b3SHong Zhang 1581d06f6b3SHong Zhang Level: advanced 1591d06f6b3SHong Zhang 1601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `IS`, `TSGetRHSSplitFunction()` 1611d06f6b3SHong Zhang @*/ 162d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetSubTS(TS ts, const char splitname[], TS *subts) 163d71ae5a4SJacob Faibussowitsch { 1641d06f6b3SHong Zhang TS_RHSSplitLink isplit; 1651d06f6b3SHong Zhang 1661d06f6b3SHong Zhang PetscFunctionBegin; 1671d06f6b3SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1684f572ea9SToby Isaac PetscAssertPointer(subts, 3); 1691d06f6b3SHong Zhang *subts = NULL; 1701d06f6b3SHong Zhang /* look up the split */ 1719566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 1721d06f6b3SHong Zhang if (isplit) *subts = isplit->ts; 1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1741d06f6b3SHong Zhang } 1751d06f6b3SHong Zhang 1761d06f6b3SHong Zhang /*@C 177bcf0153eSBarry Smith TSRHSSplitGetSubTSs - Get an array of all sub-`TS` contexts. 178545aaa6fSHong Zhang 179c3339decSBarry Smith Logically Collective 180545aaa6fSHong Zhang 1816b867d5aSJose E. Roman Input Parameter: 182bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1836b867d5aSJose E. Roman 184545aaa6fSHong Zhang Output Parameters: 185545aaa6fSHong Zhang + n - the number of splits 186b43aa488SJacob Faibussowitsch - subts - the array of `TS` contexts 187545aaa6fSHong Zhang 188545aaa6fSHong Zhang Level: advanced 189545aaa6fSHong Zhang 190bcf0153eSBarry Smith Note: 191bcf0153eSBarry Smith After `TSRHSSplitGetSubTS()` the array of `TS`s is to be freed by the user with `PetscFree()` 19220f4b53cSBarry Smith (not the `TS` in the array just the array that contains them). 193bcf0153eSBarry Smith 1941cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `IS`, `TSGetRHSSplitFunction()` 195545aaa6fSHong Zhang @*/ 196d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetSubTSs(TS ts, PetscInt *n, TS *subts[]) 197d71ae5a4SJacob Faibussowitsch { 1981d06f6b3SHong Zhang TS_RHSSplitLink ilink = ts->tsrhssplit; 199545aaa6fSHong Zhang PetscInt i = 0; 200545aaa6fSHong Zhang 201545aaa6fSHong Zhang PetscFunctionBegin; 202545aaa6fSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 203545aaa6fSHong Zhang if (subts) { 2049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ts->num_rhs_splits, subts)); 2051d06f6b3SHong Zhang while (ilink) { 2061d06f6b3SHong Zhang (*subts)[i++] = ilink->ts; 2071d06f6b3SHong Zhang ilink = ilink->next; 208545aaa6fSHong Zhang } 209545aaa6fSHong Zhang } 210545aaa6fSHong Zhang if (n) *n = ts->num_rhs_splits; 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 212545aaa6fSHong Zhang } 213