xref: /petsc/src/ts/utils/dmplexts.c (revision ef68eab947152bda6b3e24da72c19f6c7311ff5f)
1 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2 #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
3 #include <petsc/private/snesimpl.h>
4 #include <petscds.h>
5 #include <petscfv.h>
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "DMTSConvertPlex"
9 static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy)
10 {
11   PetscBool      isPlex;
12   PetscErrorCode ierr;
13 
14   PetscFunctionBegin;
15   ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);CHKERRQ(ierr);
16   if (isPlex) {
17     *plex = dm;
18     ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
19   }
20   else {
21     ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
22     if (copy) {
23       PetscInt    i;
24       PetscObject obj;
25       const char *comps[3] = {"A","dmAux","dmCh"};
26 
27       ierr = DMCopyDMTS(dm,*plex);CHKERRQ(ierr);
28       ierr = DMCopyDMSNES(dm,*plex);CHKERRQ(ierr);
29       for (i = 0; i < 3; i++) {
30         ierr = PetscObjectQuery((PetscObject)dm,comps[i],&obj);CHKERRQ(ierr);
31         ierr = PetscObjectCompose((PetscObject)*plex,comps[i],obj);CHKERRQ(ierr);
32       }
33     }
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "DMPlexTSGetGeometryFVM"
41 /*@
42   DMPlexTSGetGeometryFVM - Return precomputed geometric data
43 
44   Input Parameter:
45 . dm - The DM
46 
47   Output Parameters:
48 + facegeom - The values precomputed from face geometry
49 . cellgeom - The values precomputed from cell geometry
50 - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
51 
52   Level: developer
53 
54 .seealso: DMPlexTSSetRHSFunctionLocal()
55 @*/
56 PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
57 {
58   DMTS           dmts;
59   PetscObject    obj;
60   DM             plex;
61   PetscErrorCode ierr;
62 
63   PetscFunctionBegin;
64   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
65   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
66   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
67   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
68   if (!obj) {
69     Vec cellgeom, facegeom;
70 
71     ierr = DMPlexComputeGeometryFVM(plex, &cellgeom, &facegeom);CHKERRQ(ierr);
72     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
73     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
74     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
75     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
76   }
77   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
78   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
79   if (minRadius) {ierr = DMPlexGetMinRadius(plex, minRadius);CHKERRQ(ierr);}
80   ierr = DMDestroy(&plex);CHKERRQ(ierr);
81   PetscFunctionReturn(0);
82 }
83 
84 #undef __FUNCT__
85 #define __FUNCT__ "DMPlexTSGetGradientDM"
86 /*@C
87   DMPlexTSGetGradientDM - Return gradient data layout
88 
89   Input Parameters:
90 + dm - The DM
91 - fv - The PetscFV
92 
93   Output Parameter:
94 . dmGrad - The layout for gradient values
95 
96   Level: developer
97 
98 .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
99 @*/
100 PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
101 {
102   DMTS           dmts;
103   PetscObject    obj;
104   PetscBool      computeGradients;
105   DM             plex;
106   PetscErrorCode ierr;
107 
108   PetscFunctionBegin;
109   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
110   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
111   PetscValidPointer(dmGrad,3);
112   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
113   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
114   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
115   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
116   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
117   if (!obj) {
118     DM  dmGrad;
119     Vec faceGeometry, cellGeometry;
120 
121     ierr = DMPlexTSGetGeometryFVM(plex, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
122     ierr = DMPlexComputeGradientFVM(plex, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
123     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
124     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
125   }
126   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
127   ierr = DMDestroy(&plex);CHKERRQ(ierr);
128   PetscFunctionReturn(0);
129 }
130 
131 #undef __FUNCT__
132 #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
133 /*@
134   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
135 
136   Input Parameters:
137 + dm - The mesh
138 . t - The time
139 . locX  - Local solution
140 - user - The user context
141 
142   Output Parameter:
143 . F  - Global output vector
144 
145   Level: developer
146 
147 .seealso: DMPlexComputeJacobianActionFEM()
148 @*/
149 PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
150 {
151   Vec            locF;
152   PetscInt       cStart, cEnd, cEndInterior;
153   DM             plex;
154   PetscErrorCode ierr;
155 
156   PetscFunctionBegin;
157   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
158   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
159   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
160   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
161   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
162   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
163   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, NULL, locF, user);CHKERRQ(ierr);
164   ierr = DMLocalToGlobalBegin(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
165   ierr = DMLocalToGlobalEnd(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
166   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
167   ierr = DMDestroy(&plex);CHKERRQ(ierr);
168   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "DMPlexTSComputeBoundary"
173 /*@
174   DMPlexTSComputeBoundary - Insert the essential boundary values for the local input X and/or its time derivative X_t using pointwise functions specified by the user
175 
176   Input Parameters:
177 + dm - The mesh
178 . t - The time
179 . locX  - Local solution
180 . locX_t - Local solution time derivative, or NULL
181 - user - The user context
182 
183   Level: developer
184 
185 .seealso: DMPlexComputeJacobianActionFEM()
186 @*/
187 PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
188 {
189   DM             plex;
190   Vec            faceGeometryFVM = NULL;
191   PetscInt       Nf, f;
192   PetscErrorCode ierr;
193 
194   PetscFunctionBegin;
195   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
196   ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr);
197   for (f = 0; f < Nf; f++) {
198     PetscObject  obj;
199     PetscClassId id;
200 
201     ierr = DMGetField(plex, f, &obj);CHKERRQ(ierr);
202     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
203     if (id == PETSCFV_CLASSID) {
204       ierr = DMPlexSNESGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
205       break;
206     }
207   }
208   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, NULL, NULL, NULL);CHKERRQ(ierr);
209   /* TODO: locX_t */
210   ierr = DMDestroy(&plex);CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
216 /*@
217   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
218 
219   Input Parameters:
220 + dm - The mesh
221 . t - The time
222 . locX  - Local solution
223 . locX_t - Local solution time derivative, or NULL
224 - user - The user context
225 
226   Output Parameter:
227 . locF  - Local output vector
228 
229   Level: developer
230 
231 .seealso: DMPlexComputeJacobianActionFEM()
232 @*/
233 PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
234 {
235   PetscInt       cStart, cEnd, cEndInterior;
236   DM             plex;
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
241   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
242   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
243   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
244   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, locX_t, locF, user);CHKERRQ(ierr);
245   ierr = DMDestroy(&plex);CHKERRQ(ierr);
246   PetscFunctionReturn(0);
247 }
248 
249 #undef __FUNCT__
250 #define __FUNCT__ "DMPlexTSComputeIJacobianFEM"
251 /*@
252   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
253 
254   Input Parameters:
255 + dm - The mesh
256 . t - The time
257 . locX  - Local solution
258 . locX_t - Local solution time derivative, or NULL
259 . X_tshift - The multiplicative parameter for dF/du_t
260 - user - The user context
261 
262   Output Parameter:
263 . locF  - Local output vector
264 
265   Level: developer
266 
267 .seealso: DMPlexComputeJacobianActionFEM()
268 @*/
269 PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
270 {
271   PetscInt       cStart, cEnd, cEndInterior;
272   DM             plex;
273   PetscErrorCode ierr;
274 
275   PetscFunctionBegin;
276   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
277   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
278   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
279   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
280   ierr = DMPlexComputeJacobian_Internal(plex, cStart, cEnd, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
281   ierr = DMDestroy(&plex);CHKERRQ(ierr);
282   PetscFunctionReturn(0);
283 }
284 
285 #undef __FUNCT__
286 #define __FUNCT__ "DMTSCheckFromOptions"
287 PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
288 {
289   DM             dm;
290   SNES           snes;
291   Vec            sol;
292   PetscBool      check;
293   PetscErrorCode ierr;
294 
295   PetscFunctionBegin;
296   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
297   if (!check) PetscFunctionReturn(0);
298   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
299   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
300   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
301   ierr = TSSetUp(ts);CHKERRQ(ierr);
302   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
303   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
304   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
305   ierr = VecDestroy(&sol);CHKERRQ(ierr);
306   PetscFunctionReturn(0);
307 }
308