xref: /petsc/src/ts/utils/dmplexts.c (revision 6da023fc85477100b0ee824c70dea76712871f2a)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
337c2070cSMatthew G. Knepley #include <petsc/private/snesimpl.h>
4924a1b8fSMatthew G. Knepley #include <petscds.h>
56dbbd306SMatthew G. Knepley #include <petscfv.h>
66dbbd306SMatthew G. Knepley 
7254c1ad2SMatthew G. Knepley #undef __FUNCT__
8*6da023fcSToby Isaac #define __FUNCT__ "DMTSConvertPlex"
9*6da023fcSToby Isaac static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy)
10*6da023fcSToby Isaac {
11*6da023fcSToby Isaac   PetscBool      isPlex;
12*6da023fcSToby Isaac   PetscErrorCode ierr;
13*6da023fcSToby Isaac 
14*6da023fcSToby Isaac   PetscFunctionBegin;
15*6da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);CHKERRQ(ierr);
16*6da023fcSToby Isaac   if (isPlex) {
17*6da023fcSToby Isaac     *plex = dm;
18*6da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
19*6da023fcSToby Isaac   }
20*6da023fcSToby Isaac   else {
21*6da023fcSToby Isaac     ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
22*6da023fcSToby Isaac     if (copy) {
23*6da023fcSToby Isaac       PetscInt    i;
24*6da023fcSToby Isaac       PetscObject obj;
25*6da023fcSToby Isaac       const char *comps[3] = {"A","dmAux","dmCh"};
26*6da023fcSToby Isaac 
27*6da023fcSToby Isaac       ierr = DMCopyDMTS(dm,*plex);CHKERRQ(ierr);
28*6da023fcSToby Isaac       ierr = DMCopyDMSNES(dm,*plex);CHKERRQ(ierr);
29*6da023fcSToby Isaac       for (i = 0; i < 3; i++) {
30*6da023fcSToby Isaac         ierr = PetscObjectQuery((PetscObject)dm,comps[i],&obj);CHKERRQ(ierr);
31*6da023fcSToby Isaac         ierr = PetscObjectCompose((PetscObject)*plex,comps[i],obj);CHKERRQ(ierr);
32*6da023fcSToby Isaac       }
33*6da023fcSToby Isaac     }
34*6da023fcSToby Isaac   }
35*6da023fcSToby Isaac   PetscFunctionReturn(0);
36*6da023fcSToby Isaac }
37*6da023fcSToby Isaac 
38*6da023fcSToby Isaac 
39*6da023fcSToby Isaac #undef __FUNCT__
40b6aca0f9SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometryFVM"
41c510411aSMatthew G. Knepley /*@
42b6aca0f9SMatthew G. Knepley   DMPlexTSGetGeometryFVM - Return precomputed geometric data
43c510411aSMatthew G. Knepley 
44c510411aSMatthew G. Knepley   Input Parameter:
45c510411aSMatthew G. Knepley . dm - The DM
46c510411aSMatthew G. Knepley 
47c510411aSMatthew G. Knepley   Output Parameters:
48c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry
49c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry
50c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
51c510411aSMatthew G. Knepley 
52c510411aSMatthew G. Knepley   Level: developer
53c510411aSMatthew G. Knepley 
54c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal()
55c510411aSMatthew G. Knepley @*/
56b6aca0f9SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
57a0ac79e7SMatthew G. Knepley {
58a0ac79e7SMatthew G. Knepley   DMTS           dmts;
591faf85eaSMatthew G. Knepley   PetscObject    obj;
60*6da023fcSToby Isaac   DM             plex;
61a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
62a0ac79e7SMatthew G. Knepley 
63a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
64924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
65*6da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
66*6da023fcSToby Isaac   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
67b6aca0f9SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
681faf85eaSMatthew G. Knepley   if (!obj) {
691faf85eaSMatthew G. Knepley     Vec cellgeom, facegeom;
701faf85eaSMatthew G. Knepley 
71*6da023fcSToby Isaac     ierr = DMPlexComputeGeometryFVM(plex, &cellgeom, &facegeom);CHKERRQ(ierr);
72b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
73b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
741faf85eaSMatthew G. Knepley     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
751faf85eaSMatthew G. Knepley     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
761faf85eaSMatthew G. Knepley   }
77b6aca0f9SMatthew G. Knepley   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
78b6aca0f9SMatthew G. Knepley   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
79*6da023fcSToby Isaac   if (minRadius) {ierr = DMPlexGetMinRadius(plex, minRadius);CHKERRQ(ierr);}
80b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
81924a1b8fSMatthew G. Knepley   PetscFunctionReturn(0);
82924a1b8fSMatthew G. Knepley }
83924a1b8fSMatthew G. Knepley 
84924a1b8fSMatthew G. Knepley #undef __FUNCT__
85c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM"
86c49ccbb3SMatthew G. Knepley /*@C
87c49ccbb3SMatthew G. Knepley   DMPlexTSGetGradientDM - Return gradient data layout
88c49ccbb3SMatthew G. Knepley 
89c49ccbb3SMatthew G. Knepley   Input Parameters:
90c49ccbb3SMatthew G. Knepley + dm - The DM
91c49ccbb3SMatthew G. Knepley - fv - The PetscFV
92c49ccbb3SMatthew G. Knepley 
93c49ccbb3SMatthew G. Knepley   Output Parameter:
94c49ccbb3SMatthew G. Knepley . dmGrad - The layout for gradient values
95c49ccbb3SMatthew G. Knepley 
96c49ccbb3SMatthew G. Knepley   Level: developer
97c49ccbb3SMatthew G. Knepley 
98c49ccbb3SMatthew G. Knepley .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
99c49ccbb3SMatthew G. Knepley @*/
100c49ccbb3SMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
101c49ccbb3SMatthew G. Knepley {
102c49ccbb3SMatthew G. Knepley   DMTS           dmts;
103c49ccbb3SMatthew G. Knepley   PetscObject    obj;
104c49ccbb3SMatthew G. Knepley   PetscBool      computeGradients;
105*6da023fcSToby Isaac   DM             plex;
106c49ccbb3SMatthew G. Knepley   PetscErrorCode ierr;
107c49ccbb3SMatthew G. Knepley 
108c49ccbb3SMatthew G. Knepley   PetscFunctionBegin;
109c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
110c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
111c49ccbb3SMatthew G. Knepley   PetscValidPointer(dmGrad,3);
112c49ccbb3SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
113c49ccbb3SMatthew G. Knepley   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
114*6da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
115*6da023fcSToby Isaac   ierr = DMGetDMTS(plex, &dmts);CHKERRQ(ierr);
116c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
117c49ccbb3SMatthew G. Knepley   if (!obj) {
118c49ccbb3SMatthew G. Knepley     DM  dmGrad;
11908449791SMatthew G. Knepley     Vec faceGeometry, cellGeometry;
120c49ccbb3SMatthew G. Knepley 
121*6da023fcSToby Isaac     ierr = DMPlexTSGetGeometryFVM(plex, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
122*6da023fcSToby Isaac     ierr = DMPlexComputeGradientFVM(plex, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
123c49ccbb3SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
124924a1b8fSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
125c49ccbb3SMatthew G. Knepley   }
126c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
127b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
128254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
129254c1ad2SMatthew G. Knepley }
130254c1ad2SMatthew G. Knepley 
131254c1ad2SMatthew G. Knepley #undef __FUNCT__
132924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
13308449791SMatthew G. Knepley /*@
13408449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
13508449791SMatthew G. Knepley 
13608449791SMatthew G. Knepley   Input Parameters:
13708449791SMatthew G. Knepley + dm - The mesh
13808449791SMatthew G. Knepley . t - The time
13908449791SMatthew G. Knepley . locX  - Local solution
14008449791SMatthew G. Knepley - user - The user context
14108449791SMatthew G. Knepley 
14208449791SMatthew G. Knepley   Output Parameter:
1433b16df42SMatthew G. Knepley . F  - Global output vector
14408449791SMatthew G. Knepley 
14508449791SMatthew G. Knepley   Level: developer
14608449791SMatthew G. Knepley 
14708449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
14808449791SMatthew G. Knepley @*/
1493b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
150254c1ad2SMatthew G. Knepley {
1513b16df42SMatthew G. Knepley   Vec            locF;
152c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
153*6da023fcSToby Isaac   DM             plex;
154254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
155254c1ad2SMatthew G. Knepley 
156254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
157*6da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
158*6da023fcSToby Isaac   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
159*6da023fcSToby Isaac   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
160c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
161*6da023fcSToby Isaac   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
1623b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
163*6da023fcSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, NULL, locF, user);CHKERRQ(ierr);
164*6da023fcSToby Isaac   ierr = DMLocalToGlobalBegin(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
165*6da023fcSToby Isaac   ierr = DMLocalToGlobalEnd(plex, locF, INSERT_VALUES, F);CHKERRQ(ierr);
166*6da023fcSToby Isaac   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
167b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
168254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
169254c1ad2SMatthew G. Knepley }
170254c1ad2SMatthew G. Knepley 
171254c1ad2SMatthew G. Knepley #undef __FUNCT__
17224cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
17324cdb843SMatthew G. Knepley /*@
17424cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
17524cdb843SMatthew G. Knepley 
17624cdb843SMatthew G. Knepley   Input Parameters:
17724cdb843SMatthew G. Knepley + dm - The mesh
17824cdb843SMatthew G. Knepley . t - The time
17908449791SMatthew G. Knepley . locX  - Local solution
18008449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
18124cdb843SMatthew G. Knepley - user - The user context
18224cdb843SMatthew G. Knepley 
18324cdb843SMatthew G. Knepley   Output Parameter:
18408449791SMatthew G. Knepley . locF  - Local output vector
18524cdb843SMatthew G. Knepley 
18624cdb843SMatthew G. Knepley   Level: developer
18724cdb843SMatthew G. Knepley 
18824cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
18924cdb843SMatthew G. Knepley @*/
19008449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
19124cdb843SMatthew G. Knepley {
192c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
193*6da023fcSToby Isaac   DM             plex;
19424cdb843SMatthew G. Knepley   PetscErrorCode ierr;
19524cdb843SMatthew G. Knepley 
19624cdb843SMatthew G. Knepley   PetscFunctionBegin;
197*6da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
198*6da023fcSToby Isaac   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
199*6da023fcSToby Isaac   ierr = DMPlexGetHybridBounds(plex, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
200c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
201*6da023fcSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cStart, cEnd, time, locX, locX_t, locF, user);CHKERRQ(ierr);
202b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
20324cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
20424cdb843SMatthew G. Knepley }
2057cdb2a12SMatthew G. Knepley 
2067cdb2a12SMatthew G. Knepley #undef __FUNCT__
2077cdb2a12SMatthew G. Knepley #define __FUNCT__ "DMTSCheckFromOptions"
2080163fd50SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2097cdb2a12SMatthew G. Knepley {
2107cdb2a12SMatthew G. Knepley   DM             dm;
2117cdb2a12SMatthew G. Knepley   SNES           snes;
2121878804aSMatthew G. Knepley   Vec            sol;
2137cdb2a12SMatthew G. Knepley   PetscBool      check;
2147cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
2157cdb2a12SMatthew G. Knepley 
2167cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
217c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
2187cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
21972fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
22072fd7fc8SMatthew G. Knepley   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
2217cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
2227cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
2237cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
2247cdb2a12SMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
2251878804aSMatthew G. Knepley   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
22672fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
2277cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
2287cdb2a12SMatthew G. Knepley }
229