xref: /petsc/src/ts/interface/tsrhssplit.c (revision 424fc49db2618b7ab2212aa82f9e1492b9f632f0)
1545aaa6fSHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2*424fc49dSJoe Pusztay #include <petscdm.h>
31d06f6b3SHong Zhang static PetscErrorCode TSRHSSplitGetRHSSplit(TS ts,const char splitname[],TS_RHSSplitLink *isplit)
41d06f6b3SHong Zhang {
51d06f6b3SHong Zhang   PetscBool       found = PETSC_FALSE;
61d06f6b3SHong Zhang   PetscErrorCode  ierr;
71d06f6b3SHong Zhang 
81d06f6b3SHong Zhang   PetscFunctionBegin;
91d06f6b3SHong Zhang   *isplit = ts->tsrhssplit;
101d06f6b3SHong Zhang   /* look up the split */
111d06f6b3SHong Zhang   while (*isplit) {
121d06f6b3SHong Zhang     ierr = PetscStrcmp((*isplit)->splitname,splitname,&found);CHKERRQ(ierr);
131d06f6b3SHong Zhang     if (found) break;
141d06f6b3SHong Zhang     *isplit = (*isplit)->next;
151d06f6b3SHong Zhang   }
161d06f6b3SHong Zhang   PetscFunctionReturn(0);
171d06f6b3SHong Zhang }
181d06f6b3SHong Zhang 
19545aaa6fSHong Zhang /*@C
20545aaa6fSHong Zhang    TSRHSSplitSetIS - Set the index set for the specified split
21545aaa6fSHong Zhang 
22545aaa6fSHong Zhang    Logically Collective on TS
23545aaa6fSHong Zhang 
24545aaa6fSHong Zhang    Input Parameters:
25545aaa6fSHong Zhang +  ts        - the TS context obtained from TSCreate()
26545aaa6fSHong Zhang .  splitname - name of this split, if NULL the number of the split is used
27545aaa6fSHong Zhang -  is        - the index set for part of the solution vector
28545aaa6fSHong Zhang 
29545aaa6fSHong Zhang    Level: intermediate
30545aaa6fSHong Zhang 
311d06f6b3SHong Zhang .seealso: TSRHSSplitGetIS()
32545aaa6fSHong Zhang 
33545aaa6fSHong Zhang @*/
34545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetIS(TS ts,const char splitname[],IS is)
35545aaa6fSHong Zhang {
361d06f6b3SHong Zhang   TS_RHSSplitLink newsplit,next = ts->tsrhssplit;
37545aaa6fSHong Zhang   char            prefix[128];
38545aaa6fSHong Zhang   PetscErrorCode  ierr;
39545aaa6fSHong Zhang 
40545aaa6fSHong Zhang   PetscFunctionBegin;
41545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
42545aaa6fSHong Zhang   PetscValidHeaderSpecific(is,IS_CLASSID,3);
431d06f6b3SHong Zhang 
44545aaa6fSHong Zhang   ierr = PetscNew(&newsplit);CHKERRQ(ierr);
45545aaa6fSHong Zhang   if (splitname) {
46545aaa6fSHong Zhang     ierr = PetscStrallocpy(splitname,&newsplit->splitname);CHKERRQ(ierr);
47545aaa6fSHong Zhang   } else {
48545aaa6fSHong Zhang     ierr = PetscMalloc1(8,&newsplit->splitname);CHKERRQ(ierr);
49545aaa6fSHong Zhang     ierr = PetscSNPrintf(newsplit->splitname,7,"%D",ts->num_rhs_splits);CHKERRQ(ierr);
50545aaa6fSHong Zhang   }
51545aaa6fSHong Zhang   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
52545aaa6fSHong Zhang   newsplit->is = is;
53545aaa6fSHong Zhang   ierr = TSCreate(PetscObjectComm((PetscObject)ts),&newsplit->ts);CHKERRQ(ierr);
54109dc152SHong Zhang 
55545aaa6fSHong Zhang   ierr = PetscObjectIncrementTabLevel((PetscObject)newsplit->ts,(PetscObject)ts,1);CHKERRQ(ierr);
56545aaa6fSHong Zhang   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)newsplit->ts);CHKERRQ(ierr);
57928bb9adSStefano Zampini   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%srhsplit_%s_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "",newsplit->splitname);CHKERRQ(ierr);
58545aaa6fSHong Zhang   ierr = TSSetOptionsPrefix(newsplit->ts,prefix);CHKERRQ(ierr);
591d06f6b3SHong Zhang   if (!next) ts->tsrhssplit = newsplit;
601d06f6b3SHong Zhang   else {
611d06f6b3SHong Zhang     while (next->next) next = next->next;
621d06f6b3SHong Zhang     next->next = newsplit;
631d06f6b3SHong Zhang   }
641d06f6b3SHong Zhang   ts->num_rhs_splits++;
65545aaa6fSHong Zhang   PetscFunctionReturn(0);
66545aaa6fSHong Zhang }
67545aaa6fSHong Zhang 
68545aaa6fSHong Zhang /*@C
69545aaa6fSHong Zhang    TSRHSSplitGetIS - Retrieves the elements for a split as an IS
70545aaa6fSHong Zhang 
71545aaa6fSHong Zhang    Logically Collective on TS
72545aaa6fSHong Zhang 
73545aaa6fSHong Zhang    Input Parameters:
74545aaa6fSHong Zhang +  ts        - the TS context obtained from TSCreate()
75545aaa6fSHong Zhang -  splitname - name of this split
76545aaa6fSHong Zhang 
77545aaa6fSHong Zhang    Output Parameters:
78545aaa6fSHong Zhang -  is        - the index set for part of the solution vector
79545aaa6fSHong Zhang 
80545aaa6fSHong Zhang    Level: intermediate
81545aaa6fSHong Zhang 
821d06f6b3SHong Zhang .seealso: TSRHSSplitSetIS()
83545aaa6fSHong Zhang 
84545aaa6fSHong Zhang @*/
85545aaa6fSHong Zhang PetscErrorCode TSRHSSplitGetIS(TS ts,const char splitname[],IS *is)
86545aaa6fSHong Zhang {
871d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
88545aaa6fSHong Zhang   PetscErrorCode  ierr;
89545aaa6fSHong Zhang 
90545aaa6fSHong Zhang   PetscFunctionBegin;
91545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
921d06f6b3SHong Zhang   *is = NULL;
93545aaa6fSHong Zhang   /* look up the split */
941d06f6b3SHong Zhang   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
951d06f6b3SHong Zhang   if (isplit) *is = isplit->is;
96545aaa6fSHong Zhang   PetscFunctionReturn(0);
97545aaa6fSHong Zhang }
98545aaa6fSHong Zhang 
99545aaa6fSHong Zhang /*@C
100545aaa6fSHong Zhang    TSRHSSplitSetRHSFunction - Set the split right-hand-side functions.
101545aaa6fSHong Zhang 
102545aaa6fSHong Zhang    Logically Collective on TS
103545aaa6fSHong Zhang 
104545aaa6fSHong Zhang    Input Parameters:
105545aaa6fSHong Zhang +  ts        - the TS context obtained from TSCreate()
106545aaa6fSHong Zhang .  splitname - name of this split
107545aaa6fSHong Zhang .  r         - vector to hold the residual (or NULL to have it created internally)
108545aaa6fSHong Zhang .  rhsfunc   - the RHS function evaluation routine
109545aaa6fSHong Zhang -  ctx       - user-defined context for private data for the split function evaluation routine (may be NULL)
110545aaa6fSHong Zhang 
111545aaa6fSHong Zhang  Calling sequence of fun:
112545aaa6fSHong Zhang $  rhsfunc(TS ts,PetscReal t,Vec u,Vec f,ctx);
113545aaa6fSHong Zhang 
114545aaa6fSHong Zhang +  t    - time at step/stage being solved
115545aaa6fSHong Zhang .  u    - state vector
116545aaa6fSHong Zhang .  f    - function vector
117545aaa6fSHong Zhang -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
118545aaa6fSHong Zhang 
119545aaa6fSHong Zhang  Level: beginner
120545aaa6fSHong Zhang 
121545aaa6fSHong Zhang @*/
122545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetRHSFunction(TS ts,const char splitname[],Vec r,TSRHSFunction rhsfunc,void *ctx)
123545aaa6fSHong Zhang {
1241d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
125*424fc49dSJoe Pusztay   DM              dmc;
126545aaa6fSHong Zhang   Vec             subvec,ralloc = NULL;
127545aaa6fSHong Zhang   PetscErrorCode  ierr;
128545aaa6fSHong Zhang 
129545aaa6fSHong Zhang   PetscFunctionBegin;
130545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
131064a246eSJacob Faibussowitsch   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,3);
132545aaa6fSHong Zhang 
133545aaa6fSHong Zhang   /* look up the split */
1341d06f6b3SHong Zhang   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
1351d06f6b3SHong Zhang   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);
136545aaa6fSHong Zhang 
1371d06f6b3SHong Zhang   if (!r && ts->vec_sol) {
1381d06f6b3SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr);
139545aaa6fSHong Zhang     ierr = VecDuplicate(subvec,&ralloc);CHKERRQ(ierr);
140545aaa6fSHong Zhang     r    = ralloc;
1411d06f6b3SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr);
142545aaa6fSHong Zhang   }
143*424fc49dSJoe Pusztay 
144*424fc49dSJoe Pusztay   if (ts->dm) {
145*424fc49dSJoe Pusztay     PetscInt dim;
146*424fc49dSJoe Pusztay 
147*424fc49dSJoe Pusztay     ierr = DMGetDimension(ts->dm, &dim);CHKERRQ(ierr);
148*424fc49dSJoe Pusztay     if (dim != -1) {
149*424fc49dSJoe Pusztay       ierr = DMClone(ts->dm, &dmc);CHKERRQ(ierr);
150*424fc49dSJoe Pusztay       ierr = TSSetDM(isplit->ts, dmc);CHKERRQ(ierr);
151*424fc49dSJoe Pusztay       ierr = DMDestroy(&dmc);CHKERRQ(ierr);
152*424fc49dSJoe Pusztay     }
153*424fc49dSJoe Pusztay   }
154*424fc49dSJoe Pusztay 
1551d06f6b3SHong Zhang   ierr = TSSetRHSFunction(isplit->ts,r,rhsfunc,ctx);CHKERRQ(ierr);
156545aaa6fSHong Zhang   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
157545aaa6fSHong Zhang   PetscFunctionReturn(0);
158545aaa6fSHong Zhang }
159545aaa6fSHong Zhang 
160545aaa6fSHong Zhang /*@C
1611d06f6b3SHong Zhang    TSRHSSplitGetSubTS - Get the sub-TS by split name.
1621d06f6b3SHong Zhang 
1631d06f6b3SHong Zhang    Logically Collective on TS
1641d06f6b3SHong Zhang 
1656b867d5aSJose E. Roman    Input Parameter:
1666b867d5aSJose E. Roman .  ts - the TS context obtained from TSCreate()
1676b867d5aSJose E. Roman 
1681d06f6b3SHong Zhang    Output Parameters:
1691d06f6b3SHong Zhang +  splitname - the number of the split
1701d06f6b3SHong Zhang -  subts - the array of TS contexts
1711d06f6b3SHong Zhang 
1721d06f6b3SHong Zhang    Level: advanced
1731d06f6b3SHong Zhang 
1741d06f6b3SHong Zhang .seealso: TSGetRHSSplitFunction()
1751d06f6b3SHong Zhang @*/
1761d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTS(TS ts,const char splitname[],TS *subts)
1771d06f6b3SHong Zhang {
1781d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
1791d06f6b3SHong Zhang   PetscErrorCode  ierr;
1801d06f6b3SHong Zhang 
1811d06f6b3SHong Zhang   PetscFunctionBegin;
1821d06f6b3SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1831d06f6b3SHong Zhang   PetscValidPointer(subts,3);
1841d06f6b3SHong Zhang   *subts = NULL;
1851d06f6b3SHong Zhang   /* look up the split */
1861d06f6b3SHong Zhang   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
1871d06f6b3SHong Zhang   if (isplit) *subts = isplit->ts;
1881d06f6b3SHong Zhang   PetscFunctionReturn(0);
1891d06f6b3SHong Zhang }
1901d06f6b3SHong Zhang 
1911d06f6b3SHong Zhang /*@C
1921d06f6b3SHong Zhang    TSRHSSplitGetSubTSs - Get an array of all sub-TS contexts.
193545aaa6fSHong Zhang 
194545aaa6fSHong Zhang    Logically Collective on TS
195545aaa6fSHong Zhang 
1966b867d5aSJose E. Roman    Input Parameter:
1976b867d5aSJose E. Roman .  ts - the TS context obtained from TSCreate()
1986b867d5aSJose E. Roman 
199545aaa6fSHong Zhang    Output Parameters:
200545aaa6fSHong Zhang +  n - the number of splits
201545aaa6fSHong Zhang -  subksp - the array of TS contexts
202545aaa6fSHong Zhang 
203545aaa6fSHong Zhang    Note:
204545aaa6fSHong Zhang    After TSRHSSplitGetSubTS() the array of TSs is to be freed by the user with PetscFree()
205545aaa6fSHong Zhang    (not the TS just the array that contains them).
206545aaa6fSHong Zhang 
207545aaa6fSHong Zhang    Level: advanced
208545aaa6fSHong Zhang 
209545aaa6fSHong Zhang .seealso: TSGetRHSSplitFunction()
210545aaa6fSHong Zhang @*/
2111d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTSs(TS ts,PetscInt *n,TS *subts[])
212545aaa6fSHong Zhang {
2131d06f6b3SHong Zhang   TS_RHSSplitLink ilink = ts->tsrhssplit;
214545aaa6fSHong Zhang   PetscInt        i = 0;
215545aaa6fSHong Zhang   PetscErrorCode  ierr;
216545aaa6fSHong Zhang 
217545aaa6fSHong Zhang   PetscFunctionBegin;
218545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
219545aaa6fSHong Zhang   if (subts) {
220545aaa6fSHong Zhang     ierr = PetscMalloc1(ts->num_rhs_splits,subts);CHKERRQ(ierr);
2211d06f6b3SHong Zhang     while (ilink) {
2221d06f6b3SHong Zhang       (*subts)[i++] = ilink->ts;
2231d06f6b3SHong Zhang       ilink = ilink->next;
224545aaa6fSHong Zhang     }
225545aaa6fSHong Zhang   }
226545aaa6fSHong Zhang   if (n) *n = ts->num_rhs_splits;
227545aaa6fSHong Zhang   PetscFunctionReturn(0);
228545aaa6fSHong Zhang }
229