xref: /petsc/src/ts/interface/tsrhssplit.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1545aaa6fSHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h"  I*/
2424fc49dSJoe Pusztay #include <petscdm.h>
3*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRHSSplitGetRHSSplit(TS ts, const char splitname[], TS_RHSSplitLink *isplit)
4*d71ae5a4SJacob 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   }
151d06f6b3SHong Zhang   PetscFunctionReturn(0);
161d06f6b3SHong Zhang }
171d06f6b3SHong Zhang 
18545aaa6fSHong Zhang /*@C
19545aaa6fSHong Zhang    TSRHSSplitSetIS - Set the index set for the specified split
20545aaa6fSHong Zhang 
21545aaa6fSHong Zhang    Logically Collective on TS
22545aaa6fSHong Zhang 
23545aaa6fSHong Zhang    Input Parameters:
24545aaa6fSHong Zhang +  ts        - the TS context obtained from TSCreate()
25545aaa6fSHong Zhang .  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 
30db781477SPatrick Sanan .seealso: `TSRHSSplitGetIS()`
31545aaa6fSHong Zhang 
32545aaa6fSHong Zhang @*/
33*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitSetIS(TS ts, const char splitname[], IS is)
34*d71ae5a4SJacob Faibussowitsch {
351d06f6b3SHong Zhang   TS_RHSSplitLink newsplit, next = ts->tsrhssplit;
36545aaa6fSHong Zhang   char            prefix[128];
37545aaa6fSHong Zhang 
38545aaa6fSHong Zhang   PetscFunctionBegin;
39545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
40545aaa6fSHong Zhang   PetscValidHeaderSpecific(is, IS_CLASSID, 3);
411d06f6b3SHong Zhang 
429566063dSJacob Faibussowitsch   PetscCall(PetscNew(&newsplit));
43545aaa6fSHong Zhang   if (splitname) {
449566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(splitname, &newsplit->splitname));
45545aaa6fSHong Zhang   } else {
469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(8, &newsplit->splitname));
4763a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(newsplit->splitname, 7, "%" PetscInt_FMT, ts->num_rhs_splits));
48545aaa6fSHong Zhang   }
499566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)is));
50545aaa6fSHong Zhang   newsplit->is = is;
519566063dSJacob Faibussowitsch   PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &newsplit->ts));
52109dc152SHong Zhang 
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)newsplit->ts, (PetscObject)ts, 1));
549566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%srhsplit_%s_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "", newsplit->splitname));
559566063dSJacob Faibussowitsch   PetscCall(TSSetOptionsPrefix(newsplit->ts, prefix));
561d06f6b3SHong Zhang   if (!next) ts->tsrhssplit = newsplit;
571d06f6b3SHong Zhang   else {
581d06f6b3SHong Zhang     while (next->next) next = next->next;
591d06f6b3SHong Zhang     next->next = newsplit;
601d06f6b3SHong Zhang   }
611d06f6b3SHong Zhang   ts->num_rhs_splits++;
62545aaa6fSHong Zhang   PetscFunctionReturn(0);
63545aaa6fSHong Zhang }
64545aaa6fSHong Zhang 
65545aaa6fSHong Zhang /*@C
66545aaa6fSHong Zhang    TSRHSSplitGetIS - Retrieves the elements for a split as an IS
67545aaa6fSHong Zhang 
68545aaa6fSHong Zhang    Logically Collective on TS
69545aaa6fSHong Zhang 
70545aaa6fSHong Zhang    Input Parameters:
71545aaa6fSHong Zhang +  ts        - the TS context obtained from TSCreate()
72545aaa6fSHong Zhang -  splitname - name of this split
73545aaa6fSHong Zhang 
74545aaa6fSHong Zhang    Output Parameters:
75545aaa6fSHong Zhang -  is        - the index set for part of the solution vector
76545aaa6fSHong Zhang 
77545aaa6fSHong Zhang    Level: intermediate
78545aaa6fSHong Zhang 
79db781477SPatrick Sanan .seealso: `TSRHSSplitSetIS()`
80545aaa6fSHong Zhang 
81545aaa6fSHong Zhang @*/
82*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetIS(TS ts, const char splitname[], IS *is)
83*d71ae5a4SJacob Faibussowitsch {
841d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
85545aaa6fSHong Zhang 
86545aaa6fSHong Zhang   PetscFunctionBegin;
87545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
881d06f6b3SHong Zhang   *is = NULL;
89545aaa6fSHong Zhang   /* look up the split */
909566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit));
911d06f6b3SHong Zhang   if (isplit) *is = isplit->is;
92545aaa6fSHong Zhang   PetscFunctionReturn(0);
93545aaa6fSHong Zhang }
94545aaa6fSHong Zhang 
95545aaa6fSHong Zhang /*@C
96545aaa6fSHong Zhang    TSRHSSplitSetRHSFunction - Set the split right-hand-side functions.
97545aaa6fSHong Zhang 
98545aaa6fSHong Zhang    Logically Collective on TS
99545aaa6fSHong Zhang 
100545aaa6fSHong Zhang    Input Parameters:
101545aaa6fSHong Zhang +  ts        - the TS context obtained from TSCreate()
102545aaa6fSHong Zhang .  splitname - name of this split
103545aaa6fSHong Zhang .  r         - vector to hold the residual (or NULL to have it created internally)
104545aaa6fSHong Zhang .  rhsfunc   - the RHS function evaluation routine
105545aaa6fSHong Zhang -  ctx       - user-defined context for private data for the split function evaluation routine (may be NULL)
106545aaa6fSHong Zhang 
107545aaa6fSHong Zhang  Calling sequence of fun:
108545aaa6fSHong Zhang $  rhsfunc(TS ts,PetscReal t,Vec u,Vec f,ctx);
109545aaa6fSHong Zhang 
110545aaa6fSHong Zhang +  t    - time at step/stage being solved
111545aaa6fSHong Zhang .  u    - state vector
112545aaa6fSHong Zhang .  f    - function vector
113545aaa6fSHong Zhang -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
114545aaa6fSHong Zhang 
115545aaa6fSHong Zhang  Level: beginner
116545aaa6fSHong Zhang 
117545aaa6fSHong Zhang @*/
118*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitSetRHSFunction(TS ts, const char splitname[], Vec r, TSRHSFunction rhsfunc, void *ctx)
119*d71ae5a4SJacob Faibussowitsch {
1201d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
121424fc49dSJoe Pusztay   DM              dmc;
122545aaa6fSHong Zhang   Vec             subvec, ralloc = NULL;
123545aaa6fSHong Zhang 
124545aaa6fSHong Zhang   PetscFunctionBegin;
125545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
126064a246eSJacob Faibussowitsch   if (r) PetscValidHeaderSpecific(r, VEC_CLASSID, 3);
127545aaa6fSHong Zhang 
128545aaa6fSHong Zhang   /* look up the split */
1299566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit));
1303c633725SBarry 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);
131545aaa6fSHong Zhang 
1321d06f6b3SHong Zhang   if (!r && ts->vec_sol) {
1339566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(ts->vec_sol, isplit->is, &subvec));
1349566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(subvec, &ralloc));
135545aaa6fSHong Zhang     r = ralloc;
1369566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(ts->vec_sol, isplit->is, &subvec));
137545aaa6fSHong Zhang   }
138424fc49dSJoe Pusztay 
139424fc49dSJoe Pusztay   if (ts->dm) {
140424fc49dSJoe Pusztay     PetscInt dim;
141424fc49dSJoe Pusztay 
1429566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(ts->dm, &dim));
143424fc49dSJoe Pusztay     if (dim != -1) {
1449566063dSJacob Faibussowitsch       PetscCall(DMClone(ts->dm, &dmc));
1459566063dSJacob Faibussowitsch       PetscCall(TSSetDM(isplit->ts, dmc));
1469566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dmc));
147424fc49dSJoe Pusztay     }
148424fc49dSJoe Pusztay   }
149424fc49dSJoe Pusztay 
1509566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(isplit->ts, r, rhsfunc, ctx));
1519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ralloc));
152545aaa6fSHong Zhang   PetscFunctionReturn(0);
153545aaa6fSHong Zhang }
154545aaa6fSHong Zhang 
155545aaa6fSHong Zhang /*@C
1561d06f6b3SHong Zhang    TSRHSSplitGetSubTS - Get the sub-TS by split name.
1571d06f6b3SHong Zhang 
1581d06f6b3SHong Zhang    Logically Collective on TS
1591d06f6b3SHong Zhang 
1606b867d5aSJose E. Roman    Input Parameter:
1616b867d5aSJose E. Roman .  ts - the TS context obtained from TSCreate()
1626b867d5aSJose E. Roman 
1631d06f6b3SHong Zhang    Output Parameters:
1641d06f6b3SHong Zhang +  splitname - the number of the split
1651d06f6b3SHong Zhang -  subts - the array of TS contexts
1661d06f6b3SHong Zhang 
1671d06f6b3SHong Zhang    Level: advanced
1681d06f6b3SHong Zhang 
169db781477SPatrick Sanan .seealso: `TSGetRHSSplitFunction()`
1701d06f6b3SHong Zhang @*/
171*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetSubTS(TS ts, const char splitname[], TS *subts)
172*d71ae5a4SJacob Faibussowitsch {
1731d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
1741d06f6b3SHong Zhang 
1751d06f6b3SHong Zhang   PetscFunctionBegin;
1761d06f6b3SHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
1771d06f6b3SHong Zhang   PetscValidPointer(subts, 3);
1781d06f6b3SHong Zhang   *subts = NULL;
1791d06f6b3SHong Zhang   /* look up the split */
1809566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitGetRHSSplit(ts, splitname, &isplit));
1811d06f6b3SHong Zhang   if (isplit) *subts = isplit->ts;
1821d06f6b3SHong Zhang   PetscFunctionReturn(0);
1831d06f6b3SHong Zhang }
1841d06f6b3SHong Zhang 
1851d06f6b3SHong Zhang /*@C
1861d06f6b3SHong Zhang    TSRHSSplitGetSubTSs - Get an array of all sub-TS contexts.
187545aaa6fSHong Zhang 
188545aaa6fSHong Zhang    Logically Collective on TS
189545aaa6fSHong Zhang 
1906b867d5aSJose E. Roman    Input Parameter:
1916b867d5aSJose E. Roman .  ts - the TS context obtained from TSCreate()
1926b867d5aSJose E. Roman 
193545aaa6fSHong Zhang    Output Parameters:
194545aaa6fSHong Zhang +  n - the number of splits
195545aaa6fSHong Zhang -  subksp - the array of TS contexts
196545aaa6fSHong Zhang 
197545aaa6fSHong Zhang    Note:
198545aaa6fSHong Zhang    After TSRHSSplitGetSubTS() the array of TSs is to be freed by the user with PetscFree()
199545aaa6fSHong Zhang    (not the TS just the array that contains them).
200545aaa6fSHong Zhang 
201545aaa6fSHong Zhang    Level: advanced
202545aaa6fSHong Zhang 
203db781477SPatrick Sanan .seealso: `TSGetRHSSplitFunction()`
204545aaa6fSHong Zhang @*/
205*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSRHSSplitGetSubTSs(TS ts, PetscInt *n, TS *subts[])
206*d71ae5a4SJacob Faibussowitsch {
2071d06f6b3SHong Zhang   TS_RHSSplitLink ilink = ts->tsrhssplit;
208545aaa6fSHong Zhang   PetscInt        i     = 0;
209545aaa6fSHong Zhang 
210545aaa6fSHong Zhang   PetscFunctionBegin;
211545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
212545aaa6fSHong Zhang   if (subts) {
2139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ts->num_rhs_splits, subts));
2141d06f6b3SHong Zhang     while (ilink) {
2151d06f6b3SHong Zhang       (*subts)[i++] = ilink->ts;
2161d06f6b3SHong Zhang       ilink         = ilink->next;
217545aaa6fSHong Zhang     }
218545aaa6fSHong Zhang   }
219545aaa6fSHong Zhang   if (n) *n = ts->num_rhs_splits;
220545aaa6fSHong Zhang   PetscFunctionReturn(0);
221545aaa6fSHong Zhang }
222