xref: /petsc/src/ts/utils/dmplexts.c (revision a40652d484231dadb820fc3a9703f41e31810d43)
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__ "DMPlexTSComputeRHSBoundaryFVM"
133 /*@
134   DMPlexTSComputeRHSBoundaryFVM - Insert the essential boundary conditions in 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   Level: developer
143 
144 .seealso: DMPlexComputeJacobianActionFEM()
145 @*/
146 PetscErrorCode DMPlexTSComputeRHSBoundaryFVM(DM dm, PetscReal time, Vec locX, void *user)
147 {
148   DM             plex;
149   PetscErrorCode ierr;
150 
151   PetscFunctionBegin;
152   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
153   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, NULL, NULL, NULL);CHKERRQ(ierr);
154   ierr = DMDestroy(&plex);CHKERRQ(ierr);
155   PetscFunctionReturn(0);
156 }
157 
158 #undef __FUNCT__
159 #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
160 /*@
161   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
162 
163   Input Parameters:
164 + dm - The mesh
165 . t - The time
166 . locX  - Local solution
167 - user - The user context
168 
169   Output Parameter:
170 . F  - Global output vector
171 
172   Level: developer
173 
174 .seealso: DMPlexComputeJacobianActionFEM()
175 @*/
176 PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
177 {
178   Vec            locF;
179   PetscInt       cStart, cEnd, cEndInterior;
180   DM             plex;
181   PetscErrorCode ierr;
182 
183   PetscFunctionBegin;
184   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
185   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
186   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
187   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
188   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
189   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
190   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, NULL, locF, user);CHKERRQ(ierr);
191   ierr = DMLocalToGlobalBegin(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
192   ierr = DMLocalToGlobalEnd(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
193   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
194   ierr = DMDestroy(&plex);CHKERRQ(ierr);
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "DMPlexTSComputeIBoundaryFEM"
200 /*@
201   DMPlexTSComputeIBoundaryFEM - Insert the essential boundary values for the local input X and/or its time derivative X_t using pointwise functions specified by the user
202 
203   Input Parameters:
204 + dm - The mesh
205 . t - The time
206 . locX  - Local solution
207 . locX_t - Local solution time derivative, or NULL
208 - user - The user context
209 
210   Level: developer
211 
212 .seealso: DMPlexComputeJacobianActionFEM()
213 @*/
214 PetscErrorCode DMPlexTSComputeIBoundaryFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
215 {
216   DM             plex;
217   PetscErrorCode ierr;
218 
219   PetscFunctionBegin;
220   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
221   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, NULL, NULL, NULL);CHKERRQ(ierr);
222   /* TODO: locX_t */
223   ierr = DMDestroy(&plex);CHKERRQ(ierr);
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
229 /*@
230   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
231 
232   Input Parameters:
233 + dm - The mesh
234 . t - The time
235 . locX  - Local solution
236 . locX_t - Local solution time derivative, or NULL
237 - user - The user context
238 
239   Output Parameter:
240 . locF  - Local output vector
241 
242   Level: developer
243 
244 .seealso: DMPlexComputeJacobianActionFEM()
245 @*/
246 PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
247 {
248   PetscInt       cStart, cEnd, cEndInterior;
249   DM             plex;
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
254   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
255   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
256   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
257   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, locX_t, locF, user);CHKERRQ(ierr);
258   ierr = DMDestroy(&plex);CHKERRQ(ierr);
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "DMPlexTSComputeIJacobianFEM"
264 /*@
265   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
266 
267   Input Parameters:
268 + dm - The mesh
269 . t - The time
270 . locX  - Local solution
271 . locX_t - Local solution time derivative, or NULL
272 . X_tshift - The multiplicative parameter for dF/du_t
273 - user - The user context
274 
275   Output Parameter:
276 . locF  - Local output vector
277 
278   Level: developer
279 
280 .seealso: DMPlexComputeJacobianActionFEM()
281 @*/
282 PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
283 {
284   PetscInt       cStart, cEnd, cEndInterior;
285   DM             plex;
286   PetscErrorCode ierr;
287 
288   PetscFunctionBegin;
289   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
290   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
291   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
292   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
293   ierr = DMPlexComputeJacobian_Internal(plex, cStart, cEnd, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
294   ierr = DMDestroy(&plex);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "DMTSCheckFromOptions"
300 PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
301 {
302   DM             dm;
303   SNES           snes;
304   Vec            sol;
305   PetscBool      check;
306   PetscErrorCode ierr;
307 
308   PetscFunctionBegin;
309   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
310   if (!check) PetscFunctionReturn(0);
311   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
312   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
313   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
314   ierr = TSSetUp(ts);CHKERRQ(ierr);
315   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
316   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
317   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
318   ierr = VecDestroy(&sol);CHKERRQ(ierr);
319   PetscFunctionReturn(0);
320 }
321