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 18cc4c1da9SBarry 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 64cc4c1da9SBarry 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 1463a2a065fSHong Zhang /*@C 1473a2a065fSHong Zhang TSRHSSplitSetIFunction - Set the split implicit function for `TSARKIMEX` 1483a2a065fSHong Zhang 1493a2a065fSHong Zhang Logically Collective 1503a2a065fSHong Zhang 1513a2a065fSHong Zhang Input Parameters: 1523a2a065fSHong Zhang + ts - the `TS` context obtained from `TSCreate()` 1533a2a065fSHong Zhang . splitname - name of this split 1543a2a065fSHong Zhang . r - vector to hold the residual (or `NULL` to have it created internally) 1553a2a065fSHong Zhang . ifunc - the implicit function evaluation routine 1563a2a065fSHong Zhang - ctx - user-defined context for private data for the split function evaluation routine (may be `NULL`) 1573a2a065fSHong Zhang 1583a2a065fSHong Zhang Level: intermediate 1593a2a065fSHong Zhang 1603a2a065fSHong Zhang .seealso: [](ch_ts), `TS`, `TSIFunctionFn`, `IS`, `TSRHSSplitSetIS()`, `TSARKIMEX` 1613a2a065fSHong Zhang @*/ 1623a2a065fSHong Zhang PetscErrorCode TSRHSSplitSetIFunction(TS ts, const char splitname[], Vec r, TSIFunctionFn *ifunc, void *ctx) 1633a2a065fSHong Zhang { 1643a2a065fSHong Zhang TS_RHSSplitLink isplit; 1653a2a065fSHong Zhang DM dmc; 1663a2a065fSHong Zhang Vec subvec, ralloc = NULL; 1673a2a065fSHong Zhang 1683a2a065fSHong Zhang PetscFunctionBegin; 1693a2a065fSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1703a2a065fSHong Zhang if (r) PetscValidHeaderSpecific(r, VEC_CLASSID, 3); 1713a2a065fSHong Zhang 1723a2a065fSHong Zhang /* look up the split */ 1733a2a065fSHong Zhang PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 1743a2a065fSHong Zhang 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); 1753a2a065fSHong Zhang 1763a2a065fSHong Zhang if (!r && ts->vec_sol) { 1773a2a065fSHong Zhang PetscCall(VecGetSubVector(ts->vec_sol, isplit->is, &subvec)); 1783a2a065fSHong Zhang PetscCall(VecDuplicate(subvec, &ralloc)); 1793a2a065fSHong Zhang r = ralloc; 1803a2a065fSHong Zhang PetscCall(VecRestoreSubVector(ts->vec_sol, isplit->is, &subvec)); 1813a2a065fSHong Zhang } 1823a2a065fSHong Zhang 1833a2a065fSHong Zhang if (ts->dm) { 1843a2a065fSHong Zhang PetscInt dim; 1853a2a065fSHong Zhang 1863a2a065fSHong Zhang PetscCall(DMGetDimension(ts->dm, &dim)); 1873a2a065fSHong Zhang if (dim != -1) { 1883a2a065fSHong Zhang PetscCall(DMClone(ts->dm, &dmc)); 1893a2a065fSHong Zhang PetscCall(TSSetDM(isplit->ts, dmc)); 1903a2a065fSHong Zhang PetscCall(DMDestroy(&dmc)); 1913a2a065fSHong Zhang } 1923a2a065fSHong Zhang } 1933a2a065fSHong Zhang 1943a2a065fSHong Zhang PetscCall(TSSetIFunction(isplit->ts, r, ifunc, ctx)); 1953a2a065fSHong Zhang PetscCall(VecDestroy(&ralloc)); 1963a2a065fSHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 1973a2a065fSHong Zhang } 1983a2a065fSHong Zhang 1993a2a065fSHong Zhang /*@C 2003a2a065fSHong Zhang TSRHSSplitSetIJacobian - Set the Jacobian for the split implicit function with `TSARKIMEX` 2013a2a065fSHong Zhang 2023a2a065fSHong Zhang Logically Collective 2033a2a065fSHong Zhang 2043a2a065fSHong Zhang Input Parameters: 2053a2a065fSHong Zhang + ts - the `TS` context obtained from `TSCreate()` 2063a2a065fSHong Zhang . splitname - name of this split 2073a2a065fSHong Zhang . Amat - (approximate) matrix to store Jacobian entries computed by `f` 2083a2a065fSHong Zhang . Pmat - matrix used to compute preconditioner (usually the same as `Amat`) 2093a2a065fSHong Zhang . ijac - the Jacobian evaluation routine 2103a2a065fSHong Zhang - ctx - user-defined context for private data for the split function evaluation routine (may be `NULL`) 2113a2a065fSHong Zhang 2123a2a065fSHong Zhang Level: intermediate 2133a2a065fSHong Zhang 2143a2a065fSHong Zhang .seealso: [](ch_ts), `TS`, `TSRHSSplitSetIFunction`, `TSIJacobianFn`, `IS`, `TSRHSSplitSetIS()` 2153a2a065fSHong Zhang @*/ 2163a2a065fSHong Zhang PetscErrorCode TSRHSSplitSetIJacobian(TS ts, const char splitname[], Mat Amat, Mat Pmat, TSIJacobianFn *ijac, void *ctx) 2173a2a065fSHong Zhang { 2183a2a065fSHong Zhang TS_RHSSplitLink isplit; 2193a2a065fSHong Zhang DM dmc; 2203a2a065fSHong Zhang 2213a2a065fSHong Zhang PetscFunctionBegin; 2223a2a065fSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 2233a2a065fSHong Zhang if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 3); 2243a2a065fSHong Zhang if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 4); 2253a2a065fSHong Zhang if (Amat) PetscCheckSameComm(ts, 1, Amat, 3); 2263a2a065fSHong Zhang if (Pmat) PetscCheckSameComm(ts, 1, Pmat, 4); 2273a2a065fSHong Zhang 2283a2a065fSHong Zhang /* look up the split */ 2293a2a065fSHong Zhang PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 2303a2a065fSHong Zhang 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); 2313a2a065fSHong Zhang 2323a2a065fSHong Zhang if (ts->dm) { 2333a2a065fSHong Zhang PetscInt dim; 2343a2a065fSHong Zhang 2353a2a065fSHong Zhang PetscCall(DMGetDimension(ts->dm, &dim)); 2363a2a065fSHong Zhang if (dim != -1) { 2373a2a065fSHong Zhang PetscCall(DMClone(ts->dm, &dmc)); 2383a2a065fSHong Zhang PetscCall(TSSetDM(isplit->ts, dmc)); 2393a2a065fSHong Zhang PetscCall(DMDestroy(&dmc)); 2403a2a065fSHong Zhang } 2413a2a065fSHong Zhang } 2423a2a065fSHong Zhang 2433a2a065fSHong Zhang PetscCall(TSSetIJacobian(isplit->ts, Amat, Pmat, ijac, ctx)); 2443a2a065fSHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 2453a2a065fSHong Zhang } 2463a2a065fSHong Zhang 2473a2a065fSHong Zhang /*@C 248bcf0153eSBarry Smith TSRHSSplitGetSubTS - Get the sub-`TS` by split name. 2491d06f6b3SHong Zhang 250c3339decSBarry Smith Logically Collective 2511d06f6b3SHong Zhang 2526b867d5aSJose E. Roman Input Parameter: 253bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 2546b867d5aSJose E. Roman 2551d06f6b3SHong Zhang Output Parameters: 2561d06f6b3SHong Zhang + splitname - the number of the split 25720f4b53cSBarry Smith - subts - the sub-`TS` 2581d06f6b3SHong Zhang 2591d06f6b3SHong Zhang Level: advanced 2601d06f6b3SHong Zhang 2611cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `IS`, `TSGetRHSSplitFunction()` 2621d06f6b3SHong Zhang @*/ 263d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetSubTS(TS ts, const char splitname[], TS *subts) 264d71ae5a4SJacob Faibussowitsch { 2651d06f6b3SHong Zhang TS_RHSSplitLink isplit; 2661d06f6b3SHong Zhang 2671d06f6b3SHong Zhang PetscFunctionBegin; 2681d06f6b3SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 2694f572ea9SToby Isaac PetscAssertPointer(subts, 3); 2701d06f6b3SHong Zhang *subts = NULL; 2711d06f6b3SHong Zhang /* look up the split */ 2729566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit)); 2731d06f6b3SHong Zhang if (isplit) *subts = isplit->ts; 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2751d06f6b3SHong Zhang } 2761d06f6b3SHong Zhang 2771d06f6b3SHong Zhang /*@C 278bcf0153eSBarry Smith TSRHSSplitGetSubTSs - Get an array of all sub-`TS` contexts. 279545aaa6fSHong Zhang 280c3339decSBarry Smith Logically Collective 281545aaa6fSHong Zhang 2826b867d5aSJose E. Roman Input Parameter: 283bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 2846b867d5aSJose E. Roman 285545aaa6fSHong Zhang Output Parameters: 286545aaa6fSHong Zhang + n - the number of splits 287b43aa488SJacob Faibussowitsch - subts - the array of `TS` contexts 288545aaa6fSHong Zhang 289545aaa6fSHong Zhang Level: advanced 290545aaa6fSHong Zhang 291bcf0153eSBarry Smith Note: 292bcf0153eSBarry Smith After `TSRHSSplitGetSubTS()` the array of `TS`s is to be freed by the user with `PetscFree()` 29320f4b53cSBarry Smith (not the `TS` in the array just the array that contains them). 294bcf0153eSBarry Smith 2951cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `IS`, `TSGetRHSSplitFunction()` 296545aaa6fSHong Zhang @*/ 297d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetSubTSs(TS ts, PetscInt *n, TS *subts[]) 298d71ae5a4SJacob Faibussowitsch { 2991d06f6b3SHong Zhang TS_RHSSplitLink ilink = ts->tsrhssplit; 300545aaa6fSHong Zhang PetscInt i = 0; 301545aaa6fSHong Zhang 302545aaa6fSHong Zhang PetscFunctionBegin; 303545aaa6fSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 304545aaa6fSHong Zhang if (subts) { 3059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ts->num_rhs_splits, subts)); 3061d06f6b3SHong Zhang while (ilink) { 3071d06f6b3SHong Zhang (*subts)[i++] = ilink->ts; 3081d06f6b3SHong Zhang ilink = ilink->next; 309545aaa6fSHong Zhang } 310545aaa6fSHong Zhang } 311545aaa6fSHong Zhang if (n) *n = ts->num_rhs_splits; 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313545aaa6fSHong Zhang } 314*4bd3aaa3SHong Zhang 315*4bd3aaa3SHong Zhang /*@ 316*4bd3aaa3SHong Zhang TSRHSSplitGetSNES - Returns the `SNES` (nonlinear solver) associated with 317*4bd3aaa3SHong Zhang a `TS` (timestepper) context when RHS splits are used. 318*4bd3aaa3SHong Zhang 319*4bd3aaa3SHong Zhang Not Collective, but snes is parallel if ts is parallel 320*4bd3aaa3SHong Zhang 321*4bd3aaa3SHong Zhang Input Parameter: 322*4bd3aaa3SHong Zhang . ts - the `TS` context obtained from `TSCreate()` 323*4bd3aaa3SHong Zhang 324*4bd3aaa3SHong Zhang Output Parameter: 325*4bd3aaa3SHong Zhang . snes - the nonlinear solver context 326*4bd3aaa3SHong Zhang 327*4bd3aaa3SHong Zhang Level: intermediate 328*4bd3aaa3SHong Zhang 329*4bd3aaa3SHong Zhang Note: 330*4bd3aaa3SHong Zhang The returned `SNES` may have a different `DM` with the `TS` `DM`. 331*4bd3aaa3SHong Zhang 332*4bd3aaa3SHong Zhang .seealso: [](ch_ts), `TS`, `SNES`, `TSCreate()`, `TSRHSSplitSetSNES()` 333*4bd3aaa3SHong Zhang @*/ 334*4bd3aaa3SHong Zhang PetscErrorCode TSRHSSplitGetSNES(TS ts, SNES *snes) 335*4bd3aaa3SHong Zhang { 336*4bd3aaa3SHong Zhang PetscFunctionBegin; 337*4bd3aaa3SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 338*4bd3aaa3SHong Zhang PetscAssertPointer(snes, 2); 339*4bd3aaa3SHong Zhang if (!ts->snesrhssplit) { 340*4bd3aaa3SHong Zhang PetscCall(SNESCreate(PetscObjectComm((PetscObject)ts), &ts->snesrhssplit)); 341*4bd3aaa3SHong Zhang PetscCall(PetscObjectSetOptions((PetscObject)ts->snesrhssplit, ((PetscObject)ts)->options)); 342*4bd3aaa3SHong Zhang PetscCall(SNESSetFunction(ts->snesrhssplit, NULL, SNESTSFormFunction, ts)); 343*4bd3aaa3SHong Zhang PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->snesrhssplit, (PetscObject)ts, 1)); 344*4bd3aaa3SHong Zhang if (ts->problem_type == TS_LINEAR) PetscCall(SNESSetType(ts->snesrhssplit, SNESKSPONLY)); 345*4bd3aaa3SHong Zhang } 346*4bd3aaa3SHong Zhang *snes = ts->snesrhssplit; 347*4bd3aaa3SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 348*4bd3aaa3SHong Zhang } 349*4bd3aaa3SHong Zhang 350*4bd3aaa3SHong Zhang /*@ 351*4bd3aaa3SHong Zhang TSRHSSplitSetSNES - Set the `SNES` (nonlinear solver) to be used by the 352*4bd3aaa3SHong Zhang timestepping context when RHS splits are used. 353*4bd3aaa3SHong Zhang 354*4bd3aaa3SHong Zhang Collective 355*4bd3aaa3SHong Zhang 356*4bd3aaa3SHong Zhang Input Parameters: 357*4bd3aaa3SHong Zhang + ts - the `TS` context obtained from `TSCreate()` 358*4bd3aaa3SHong Zhang - snes - the nonlinear solver context 359*4bd3aaa3SHong Zhang 360*4bd3aaa3SHong Zhang Level: intermediate 361*4bd3aaa3SHong Zhang 362*4bd3aaa3SHong Zhang Note: 363*4bd3aaa3SHong Zhang Most users should have the `TS` created by calling `TSRHSSplitGetSNES()` 364*4bd3aaa3SHong Zhang 365*4bd3aaa3SHong Zhang .seealso: [](ch_ts), `TS`, `SNES`, `TSCreate()`, `TSRHSSplitGetSNES()` 366*4bd3aaa3SHong Zhang @*/ 367*4bd3aaa3SHong Zhang PetscErrorCode TSRHSSplitSetSNES(TS ts, SNES snes) 368*4bd3aaa3SHong Zhang { 369*4bd3aaa3SHong Zhang PetscErrorCode (*func)(SNES, Vec, Mat, Mat, void *); 370*4bd3aaa3SHong Zhang 371*4bd3aaa3SHong Zhang PetscFunctionBegin; 372*4bd3aaa3SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 373*4bd3aaa3SHong Zhang PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 374*4bd3aaa3SHong Zhang PetscCall(PetscObjectReference((PetscObject)snes)); 375*4bd3aaa3SHong Zhang PetscCall(SNESDestroy(&ts->snesrhssplit)); 376*4bd3aaa3SHong Zhang 377*4bd3aaa3SHong Zhang ts->snesrhssplit = snes; 378*4bd3aaa3SHong Zhang 379*4bd3aaa3SHong Zhang PetscCall(SNESSetFunction(ts->snesrhssplit, NULL, SNESTSFormFunction, ts)); 380*4bd3aaa3SHong Zhang PetscCall(SNESGetJacobian(ts->snesrhssplit, NULL, NULL, &func, NULL)); 381*4bd3aaa3SHong Zhang if (func == SNESTSFormJacobian) PetscCall(SNESSetJacobian(ts->snesrhssplit, NULL, NULL, SNESTSFormJacobian, ts)); 382*4bd3aaa3SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 383*4bd3aaa3SHong Zhang } 384