xref: /petsc/src/ts/interface/tsrhssplit.c (revision 487a658c8b32ba712a1dc8280daad2fd70c1dcd9)
1 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2 
3 static PetscErrorCode TSRHSSplitGetRHSSplit(TS ts,const char splitname[],TS_RHSSplitLink *isplit)
4 {
5   PetscBool       found = PETSC_FALSE;
6   PetscErrorCode  ierr;
7 
8   PetscFunctionBegin;
9   *isplit = ts->tsrhssplit;
10   /* look up the split */
11   while (*isplit) {
12     ierr = PetscStrcmp((*isplit)->splitname,splitname,&found);CHKERRQ(ierr);
13     if (found) break;
14     *isplit = (*isplit)->next;
15   }
16   PetscFunctionReturn(0);
17 }
18 
19 /*@C
20    TSRHSSplitSetIS - Set the index set for the specified split
21 
22    Logically Collective on TS
23 
24    Input Parameters:
25 +  ts        - the TS context obtained from TSCreate()
26 .  splitname - name of this split, if NULL the number of the split is used
27 -  is        - the index set for part of the solution vector
28 
29    Level: intermediate
30 
31 .seealso: TSRHSSplitGetIS()
32 
33 .keywords: TS, TSRHSSplit
34 @*/
35 PetscErrorCode TSRHSSplitSetIS(TS ts,const char splitname[],IS is)
36 {
37   TS_RHSSplitLink newsplit,next = ts->tsrhssplit;
38   char            prefix[128];
39   PetscErrorCode  ierr;
40 
41   PetscFunctionBegin;
42   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
43   PetscValidHeaderSpecific(is,IS_CLASSID,3);
44 
45   ierr = PetscNew(&newsplit);CHKERRQ(ierr);
46   if (splitname) {
47     ierr = PetscStrallocpy(splitname,&newsplit->splitname);CHKERRQ(ierr);
48   } else {
49     ierr = PetscMalloc1(8,&newsplit->splitname);CHKERRQ(ierr);
50     ierr = PetscSNPrintf(newsplit->splitname,7,"%D",ts->num_rhs_splits);CHKERRQ(ierr);
51   }
52   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
53   newsplit->is = is;
54   ierr = TSCreate(PetscObjectComm((PetscObject)ts),&newsplit->ts);CHKERRQ(ierr);
55   ierr = PetscObjectIncrementTabLevel((PetscObject)newsplit->ts,(PetscObject)ts,1);CHKERRQ(ierr);
56   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)newsplit->ts);CHKERRQ(ierr);
57   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%srhsplit_%s_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "",newsplit->splitname);
58   ierr = TSSetOptionsPrefix(newsplit->ts,prefix);CHKERRQ(ierr);
59   if (!next) ts->tsrhssplit = newsplit;
60   else {
61     while (next->next) next = next->next;
62     next->next = newsplit;
63   }
64   ts->num_rhs_splits++;
65   PetscFunctionReturn(0);
66 }
67 
68 /*@C
69    TSRHSSplitGetIS - Retrieves the elements for a split as an IS
70 
71    Logically Collective on TS
72 
73    Input Parameters:
74 +  ts        - the TS context obtained from TSCreate()
75 -  splitname - name of this split
76 
77    Output Parameters:
78 -  is        - the index set for part of the solution vector
79 
80    Level: intermediate
81 
82 .seealso: TSRHSSplitSetIS()
83 
84 .keywords: TS, TSRHSSplit
85 @*/
86 PetscErrorCode TSRHSSplitGetIS(TS ts,const char splitname[],IS *is)
87 {
88   TS_RHSSplitLink isplit;
89   PetscErrorCode  ierr;
90 
91   PetscFunctionBegin;
92   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
93   *is = NULL;
94   /* look up the split */
95   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
96   if (isplit) *is = isplit->is;
97   PetscFunctionReturn(0);
98 }
99 
100 /*@C
101    TSRHSSplitSetRHSFunction - Set the split right-hand-side functions.
102 
103    Logically Collective on TS
104 
105    Input Parameters:
106 +  ts        - the TS context obtained from TSCreate()
107 .  splitname - name of this split
108 .  r         - vector to hold the residual (or NULL to have it created internally)
109 .  rhsfunc   - the RHS function evaluation routine
110 -  ctx       - user-defined context for private data for the split function evaluation routine (may be NULL)
111 
112  Calling sequence of fun:
113 $  rhsfunc(TS ts,PetscReal t,Vec u,Vec f,ctx);
114 
115 +  t    - time at step/stage being solved
116 .  u    - state vector
117 .  f    - function vector
118 -  ctx  - [optional] user-defined context for matrix evaluation routine (may be NULL)
119 
120  Level: beginner
121 
122 .keywords: TS, timestep, set, ODE, Hamiltonian, Function
123 @*/
124 PetscErrorCode TSRHSSplitSetRHSFunction(TS ts,const char splitname[],Vec r,TSRHSFunction rhsfunc,void *ctx)
125 {
126   TS_RHSSplitLink isplit;
127   Vec             subvec,ralloc = NULL;
128   PetscErrorCode  ierr;
129 
130   PetscFunctionBegin;
131   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
132   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
133 
134   /* look up the split */
135   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
136   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);
137 
138   if (!r && ts->vec_sol) {
139     ierr = VecGetSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr);
140     ierr = VecDuplicate(subvec,&ralloc);CHKERRQ(ierr);
141     r    = ralloc;
142     ierr = VecRestoreSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr);
143   }
144   ierr = TSSetRHSFunction(isplit->ts,r,rhsfunc,ctx);CHKERRQ(ierr);
145   ierr = VecDestroy(&ralloc);CHKERRQ(ierr);
146   PetscFunctionReturn(0);
147 }
148 
149 /*@C
150    TSRHSSplitGetSubTS - Get the sub-TS by split name.
151 
152    Logically Collective on TS
153 
154    Output Parameters:
155 +  splitname - the number of the split
156 -  subts - the array of TS contexts
157 
158    Level: advanced
159 
160 .seealso: TSGetRHSSplitFunction()
161 @*/
162 PetscErrorCode TSRHSSplitGetSubTS(TS ts,const char splitname[],TS *subts)
163 {
164   TS_RHSSplitLink isplit;
165   PetscErrorCode  ierr;
166 
167   PetscFunctionBegin;
168   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
169   PetscValidPointer(subts,3);
170   *subts = NULL;
171   /* look up the split */
172   ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr);
173   if (isplit) *subts = isplit->ts;
174   PetscFunctionReturn(0);
175 }
176 
177 /*@C
178    TSRHSSplitGetSubTSs - Get an array of all sub-TS contexts.
179 
180    Logically Collective on TS
181 
182    Output Parameters:
183 +  n - the number of splits
184 -  subksp - the array of TS contexts
185 
186    Note:
187    After TSRHSSplitGetSubTS() the array of TSs is to be freed by the user with PetscFree()
188    (not the TS just the array that contains them).
189 
190    Level: advanced
191 
192 .seealso: TSGetRHSSplitFunction()
193 @*/
194 PetscErrorCode TSRHSSplitGetSubTSs(TS ts,PetscInt *n,TS *subts[])
195 {
196   TS_RHSSplitLink ilink = ts->tsrhssplit;
197   PetscInt        i = 0;
198   PetscErrorCode  ierr;
199 
200   PetscFunctionBegin;
201   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
202   if (subts) {
203     ierr = PetscMalloc1(ts->num_rhs_splits,subts);CHKERRQ(ierr);
204     while (ilink) {
205       (*subts)[i++] = ilink->ts;
206       ilink = ilink->next;
207     }
208   }
209   if (n) *n = ts->num_rhs_splits;
210   PetscFunctionReturn(0);
211 }
212