xref: /petsc/src/ts/interface/tsrhssplit.c (revision 928bb9ad770004178bd1c162529fb34e5630f279)
1545aaa6fSHong Zhang #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2545aaa6fSHong Zhang 
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 .keywords: TS, TSRHSSplit
34545aaa6fSHong Zhang @*/
35545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetIS(TS ts,const char splitname[],IS is)
36545aaa6fSHong Zhang {
371d06f6b3SHong Zhang   TS_RHSSplitLink newsplit,next = ts->tsrhssplit;
38545aaa6fSHong Zhang   char            prefix[128];
39545aaa6fSHong Zhang   PetscErrorCode  ierr;
40545aaa6fSHong Zhang 
41545aaa6fSHong Zhang   PetscFunctionBegin;
42545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
43545aaa6fSHong Zhang   PetscValidHeaderSpecific(is,IS_CLASSID,3);
441d06f6b3SHong Zhang 
45545aaa6fSHong Zhang   ierr = PetscNew(&newsplit);CHKERRQ(ierr);
46545aaa6fSHong Zhang   if (splitname) {
47545aaa6fSHong Zhang     ierr = PetscStrallocpy(splitname,&newsplit->splitname);CHKERRQ(ierr);
48545aaa6fSHong Zhang   } else {
49545aaa6fSHong Zhang     ierr = PetscMalloc1(8,&newsplit->splitname);CHKERRQ(ierr);
50545aaa6fSHong Zhang     ierr = PetscSNPrintf(newsplit->splitname,7,"%D",ts->num_rhs_splits);CHKERRQ(ierr);
51545aaa6fSHong Zhang   }
52545aaa6fSHong Zhang   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
53545aaa6fSHong Zhang   newsplit->is = is;
54545aaa6fSHong Zhang   ierr = TSCreate(PetscObjectComm((PetscObject)ts),&newsplit->ts);CHKERRQ(ierr);
55545aaa6fSHong Zhang   ierr = PetscObjectIncrementTabLevel((PetscObject)newsplit->ts,(PetscObject)ts,1);CHKERRQ(ierr);
56545aaa6fSHong Zhang   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)newsplit->ts);CHKERRQ(ierr);
57*928bb9adSStefano 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 .keywords: TS, TSRHSSplit
85545aaa6fSHong Zhang @*/
86545aaa6fSHong Zhang PetscErrorCode TSRHSSplitGetIS(TS ts,const char splitname[],IS *is)
87545aaa6fSHong Zhang {
881d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
89545aaa6fSHong Zhang   PetscErrorCode  ierr;
90545aaa6fSHong Zhang 
91545aaa6fSHong Zhang   PetscFunctionBegin;
92545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
931d06f6b3SHong Zhang   *is = NULL;
94545aaa6fSHong Zhang   /* look up the split */
951d06f6b3SHong Zhang   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
961d06f6b3SHong Zhang   if (isplit) *is = isplit->is;
97545aaa6fSHong Zhang   PetscFunctionReturn(0);
98545aaa6fSHong Zhang }
99545aaa6fSHong Zhang 
100545aaa6fSHong Zhang /*@C
101545aaa6fSHong Zhang    TSRHSSplitSetRHSFunction - Set the split right-hand-side functions.
102545aaa6fSHong Zhang 
103545aaa6fSHong Zhang    Logically Collective on TS
104545aaa6fSHong Zhang 
105545aaa6fSHong Zhang    Input Parameters:
106545aaa6fSHong Zhang +  ts        - the TS context obtained from TSCreate()
107545aaa6fSHong Zhang .  splitname - name of this split
108545aaa6fSHong Zhang .  r         - vector to hold the residual (or NULL to have it created internally)
109545aaa6fSHong Zhang .  rhsfunc   - the RHS function evaluation routine
110545aaa6fSHong Zhang -  ctx       - user-defined context for private data for the split function evaluation routine (may be NULL)
111545aaa6fSHong Zhang 
112545aaa6fSHong Zhang  Calling sequence of fun:
113545aaa6fSHong Zhang $  rhsfunc(TS ts,PetscReal t,Vec u,Vec f,ctx);
114545aaa6fSHong Zhang 
115545aaa6fSHong Zhang +  t    - time at step/stage being solved
116545aaa6fSHong Zhang .  u    - state vector
117545aaa6fSHong Zhang .  f    - function vector
118545aaa6fSHong Zhang -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
119545aaa6fSHong Zhang 
120545aaa6fSHong Zhang  Level: beginner
121545aaa6fSHong Zhang 
122545aaa6fSHong Zhang .keywords: TS, timestep, set, ODE, Hamiltonian, Function
123545aaa6fSHong Zhang @*/
124545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetRHSFunction(TS ts,const char splitname[],Vec r,TSRHSFunction rhsfunc,void *ctx)
125545aaa6fSHong Zhang {
1261d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
127545aaa6fSHong Zhang   Vec             subvec,ralloc = NULL;
128545aaa6fSHong Zhang   PetscErrorCode  ierr;
129545aaa6fSHong Zhang 
130545aaa6fSHong Zhang   PetscFunctionBegin;
131545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
132545aaa6fSHong Zhang   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
133545aaa6fSHong Zhang 
134545aaa6fSHong Zhang   /* look up the split */
1351d06f6b3SHong Zhang   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
1361d06f6b3SHong 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);
137545aaa6fSHong Zhang 
1381d06f6b3SHong Zhang   if (!r && ts->vec_sol) {
1391d06f6b3SHong Zhang     ierr = VecGetSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr);
140545aaa6fSHong Zhang     ierr = VecDuplicate(subvec,&ralloc);CHKERRQ(ierr);
141545aaa6fSHong Zhang     r    = ralloc;
1421d06f6b3SHong Zhang     ierr = VecRestoreSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr);
143545aaa6fSHong Zhang   }
1441d06f6b3SHong Zhang   ierr = TSSetRHSFunction(isplit->ts,r,rhsfunc,ctx);CHKERRQ(ierr);
145545aaa6fSHong Zhang   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
146545aaa6fSHong Zhang   PetscFunctionReturn(0);
147545aaa6fSHong Zhang }
148545aaa6fSHong Zhang 
149545aaa6fSHong Zhang /*@C
1501d06f6b3SHong Zhang    TSRHSSplitGetSubTS - Get the sub-TS by split name.
1511d06f6b3SHong Zhang 
1521d06f6b3SHong Zhang    Logically Collective on TS
1531d06f6b3SHong Zhang 
1541d06f6b3SHong Zhang    Output Parameters:
1551d06f6b3SHong Zhang +  splitname - the number of the split
1561d06f6b3SHong Zhang -  subts - the array of TS contexts
1571d06f6b3SHong Zhang 
1581d06f6b3SHong Zhang    Level: advanced
1591d06f6b3SHong Zhang 
1601d06f6b3SHong Zhang .seealso: TSGetRHSSplitFunction()
1611d06f6b3SHong Zhang @*/
1621d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTS(TS ts,const char splitname[],TS *subts)
1631d06f6b3SHong Zhang {
1641d06f6b3SHong Zhang   TS_RHSSplitLink isplit;
1651d06f6b3SHong Zhang   PetscErrorCode  ierr;
1661d06f6b3SHong Zhang 
1671d06f6b3SHong Zhang   PetscFunctionBegin;
1681d06f6b3SHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1691d06f6b3SHong Zhang   PetscValidPointer(subts,3);
1701d06f6b3SHong Zhang   *subts = NULL;
1711d06f6b3SHong Zhang   /* look up the split */
1721d06f6b3SHong Zhang   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
1731d06f6b3SHong Zhang   if (isplit) *subts = isplit->ts;
1741d06f6b3SHong Zhang   PetscFunctionReturn(0);
1751d06f6b3SHong Zhang }
1761d06f6b3SHong Zhang 
1771d06f6b3SHong Zhang /*@C
1781d06f6b3SHong Zhang    TSRHSSplitGetSubTSs - Get an array of all sub-TS contexts.
179545aaa6fSHong Zhang 
180545aaa6fSHong Zhang    Logically Collective on TS
181545aaa6fSHong Zhang 
182545aaa6fSHong Zhang    Output Parameters:
183545aaa6fSHong Zhang +  n - the number of splits
184545aaa6fSHong Zhang -  subksp - the array of TS contexts
185545aaa6fSHong Zhang 
186545aaa6fSHong Zhang    Note:
187545aaa6fSHong Zhang    After TSRHSSplitGetSubTS() the array of TSs is to be freed by the user with PetscFree()
188545aaa6fSHong Zhang    (not the TS just the array that contains them).
189545aaa6fSHong Zhang 
190545aaa6fSHong Zhang    Level: advanced
191545aaa6fSHong Zhang 
192545aaa6fSHong Zhang .seealso: TSGetRHSSplitFunction()
193545aaa6fSHong Zhang @*/
1941d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTSs(TS ts,PetscInt *n,TS *subts[])
195545aaa6fSHong Zhang {
1961d06f6b3SHong Zhang   TS_RHSSplitLink ilink = ts->tsrhssplit;
197545aaa6fSHong Zhang   PetscInt        i = 0;
198545aaa6fSHong Zhang   PetscErrorCode  ierr;
199545aaa6fSHong Zhang 
200545aaa6fSHong Zhang   PetscFunctionBegin;
201545aaa6fSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
202545aaa6fSHong Zhang   if (subts) {
203545aaa6fSHong Zhang     ierr = PetscMalloc1(ts->num_rhs_splits,subts);CHKERRQ(ierr);
2041d06f6b3SHong Zhang     while (ilink) {
2051d06f6b3SHong Zhang       (*subts)[i++] = ilink->ts;
2061d06f6b3SHong Zhang       ilink = ilink->next;
207545aaa6fSHong Zhang     }
208545aaa6fSHong Zhang   }
209545aaa6fSHong Zhang   if (n) *n = ts->num_rhs_splits;
210545aaa6fSHong Zhang   PetscFunctionReturn(0);
211545aaa6fSHong Zhang }
212