xref: /petsc/src/ts/utils/dmplexts.c (revision b2338dbc0ae157d1defa802b8d9a5ac446629e46)
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__
8b6aca0f9SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGeometryFVM"
9c510411aSMatthew G. Knepley /*@
10b6aca0f9SMatthew G. Knepley   DMPlexTSGetGeometryFVM - Return precomputed geometric data
11c510411aSMatthew G. Knepley 
12c510411aSMatthew G. Knepley   Input Parameter:
13c510411aSMatthew G. Knepley . dm - The DM
14c510411aSMatthew G. Knepley 
15c510411aSMatthew G. Knepley   Output Parameters:
16c510411aSMatthew G. Knepley + facegeom - The values precomputed from face geometry
17c510411aSMatthew G. Knepley . cellgeom - The values precomputed from cell geometry
18c510411aSMatthew G. Knepley - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell
19c510411aSMatthew G. Knepley 
20c510411aSMatthew G. Knepley   Level: developer
21c510411aSMatthew G. Knepley 
22c510411aSMatthew G. Knepley .seealso: DMPlexTSSetRHSFunctionLocal()
23c510411aSMatthew G. Knepley @*/
24b6aca0f9SMatthew G. Knepley PetscErrorCode DMPlexTSGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
25a0ac79e7SMatthew G. Knepley {
26a0ac79e7SMatthew G. Knepley   DMTS           dmts;
271faf85eaSMatthew G. Knepley   PetscObject    obj;
28*b2338dbcSToby Isaac   PetscBool      isPlex;
29*b2338dbcSToby Isaac   DM             plex = NULL;
30a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
31a0ac79e7SMatthew G. Knepley 
32a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
33924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
34*b2338dbcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);CHKERRQ(ierr);
35*b2338dbcSToby Isaac   if (!isPlex) {
36*b2338dbcSToby Isaac     ierr = DMConvert(dm,DMPLEX,&plex);CHKERRQ(ierr);
37*b2338dbcSToby Isaac     if (!plex) {SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot convert to DMPlex");}
38*b2338dbcSToby Isaac     ierr = DMCopyDMTS(dm,plex);CHKERRQ(ierr);
39*b2338dbcSToby Isaac     dm = plex;
40*b2338dbcSToby Isaac   }
41a0ac79e7SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
42b6aca0f9SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
431faf85eaSMatthew G. Knepley   if (!obj) {
441faf85eaSMatthew G. Knepley     Vec cellgeom, facegeom;
451faf85eaSMatthew G. Knepley 
461faf85eaSMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr);
47b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
48b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
491faf85eaSMatthew G. Knepley     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
501faf85eaSMatthew G. Knepley     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
511faf85eaSMatthew G. Knepley   }
52b6aca0f9SMatthew G. Knepley   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
53b6aca0f9SMatthew G. Knepley   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
54113c68e6SMatthew G. Knepley   if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);}
55*b2338dbcSToby Isaac   if (!isPlex) {
56*b2338dbcSToby Isaac     ierr = DMDestroy(&plex);CHKERRQ(ierr);
57*b2338dbcSToby Isaac   }
58924a1b8fSMatthew G. Knepley   PetscFunctionReturn(0);
59924a1b8fSMatthew G. Knepley }
60924a1b8fSMatthew G. Knepley 
61924a1b8fSMatthew G. Knepley #undef __FUNCT__
62c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM"
63c49ccbb3SMatthew G. Knepley /*@C
64c49ccbb3SMatthew G. Knepley   DMPlexTSGetGradientDM - Return gradient data layout
65c49ccbb3SMatthew G. Knepley 
66c49ccbb3SMatthew G. Knepley   Input Parameters:
67c49ccbb3SMatthew G. Knepley + dm - The DM
68c49ccbb3SMatthew G. Knepley - fv - The PetscFV
69c49ccbb3SMatthew G. Knepley 
70c49ccbb3SMatthew G. Knepley   Output Parameter:
71c49ccbb3SMatthew G. Knepley . dmGrad - The layout for gradient values
72c49ccbb3SMatthew G. Knepley 
73c49ccbb3SMatthew G. Knepley   Level: developer
74c49ccbb3SMatthew G. Knepley 
75c49ccbb3SMatthew G. Knepley .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
76c49ccbb3SMatthew G. Knepley @*/
77c49ccbb3SMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
78c49ccbb3SMatthew G. Knepley {
79c49ccbb3SMatthew G. Knepley   DMTS           dmts;
80c49ccbb3SMatthew G. Knepley   PetscObject    obj;
81c49ccbb3SMatthew G. Knepley   PetscBool      computeGradients;
82*b2338dbcSToby Isaac   PetscBool      isPlex;
83*b2338dbcSToby Isaac   DM             plex = NULL;
84c49ccbb3SMatthew G. Knepley   PetscErrorCode ierr;
85c49ccbb3SMatthew G. Knepley 
86c49ccbb3SMatthew G. Knepley   PetscFunctionBegin;
87c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
88c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
89c49ccbb3SMatthew G. Knepley   PetscValidPointer(dmGrad,3);
90c49ccbb3SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
91c49ccbb3SMatthew G. Knepley   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
92*b2338dbcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);CHKERRQ(ierr);
93*b2338dbcSToby Isaac   if (!isPlex) {
94*b2338dbcSToby Isaac     ierr = DMConvert(dm,DMPLEX,&plex);CHKERRQ(ierr);
95*b2338dbcSToby Isaac     if (!plex) {SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot convert to DMPlex");}
96*b2338dbcSToby Isaac     ierr = DMCopyDMTS(dm,plex);CHKERRQ(ierr);
97*b2338dbcSToby Isaac     dm = plex;
98*b2338dbcSToby Isaac   }
99c49ccbb3SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
100c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
101c49ccbb3SMatthew G. Knepley   if (!obj) {
102c49ccbb3SMatthew G. Knepley     DM  dmGrad;
10308449791SMatthew G. Knepley     Vec faceGeometry, cellGeometry;
104c49ccbb3SMatthew G. Knepley 
10508449791SMatthew G. Knepley     ierr = DMPlexTSGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
10608449791SMatthew G. Knepley     ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
107c49ccbb3SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
108924a1b8fSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
109c49ccbb3SMatthew G. Knepley   }
110c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
111*b2338dbcSToby Isaac   if (!isPlex) {
112*b2338dbcSToby Isaac     ierr = DMDestroy(&plex);CHKERRQ(ierr);
113*b2338dbcSToby Isaac   }
114254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
115254c1ad2SMatthew G. Knepley }
116254c1ad2SMatthew G. Knepley 
117254c1ad2SMatthew G. Knepley #undef __FUNCT__
118924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
11908449791SMatthew G. Knepley /*@
12008449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
12108449791SMatthew G. Knepley 
12208449791SMatthew G. Knepley   Input Parameters:
12308449791SMatthew G. Knepley + dm - The mesh
12408449791SMatthew G. Knepley . t - The time
12508449791SMatthew G. Knepley . locX  - Local solution
12608449791SMatthew G. Knepley - user - The user context
12708449791SMatthew G. Knepley 
12808449791SMatthew G. Knepley   Output Parameter:
1293b16df42SMatthew G. Knepley . F  - Global output vector
13008449791SMatthew G. Knepley 
13108449791SMatthew G. Knepley   Level: developer
13208449791SMatthew G. Knepley 
13308449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
13408449791SMatthew G. Knepley @*/
1353b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
136254c1ad2SMatthew G. Knepley {
1373b16df42SMatthew G. Knepley   Vec            locF;
138c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
139*b2338dbcSToby Isaac   PetscBool      isPlex;
140*b2338dbcSToby Isaac   DM             plex = NULL;
141254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
142254c1ad2SMatthew G. Knepley 
143254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
144*b2338dbcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);CHKERRQ(ierr);
145*b2338dbcSToby Isaac   if (!isPlex) {
146*b2338dbcSToby Isaac     ierr = DMConvert(dm,DMPLEX,&plex);CHKERRQ(ierr);
147*b2338dbcSToby Isaac     if (!plex) {SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot convert to DMPlex");}
148*b2338dbcSToby Isaac     ierr = DMCopyDMTS(dm,plex);CHKERRQ(ierr);
149*b2338dbcSToby Isaac     dm = plex;
150*b2338dbcSToby Isaac   }
151c4d4a4f8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
152c4d4a4f8SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
153c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
1543b16df42SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locF);CHKERRQ(ierr);
1553b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
156c4d4a4f8SMatthew G. Knepley   ierr = DMPlexComputeResidual_Internal(dm, cStart, cEnd, time, locX, NULL, locF, user);CHKERRQ(ierr);
1573b16df42SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1583b16df42SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1593b16df42SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locF);CHKERRQ(ierr);
160*b2338dbcSToby Isaac   if (!isPlex) {
161*b2338dbcSToby Isaac     ierr = DMDestroy(&plex);CHKERRQ(ierr);
162*b2338dbcSToby Isaac   }
163254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
164254c1ad2SMatthew G. Knepley }
165254c1ad2SMatthew G. Knepley 
166254c1ad2SMatthew G. Knepley #undef __FUNCT__
16724cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
16824cdb843SMatthew G. Knepley /*@
16924cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
17024cdb843SMatthew G. Knepley 
17124cdb843SMatthew G. Knepley   Input Parameters:
17224cdb843SMatthew G. Knepley + dm - The mesh
17324cdb843SMatthew G. Knepley . t - The time
17408449791SMatthew G. Knepley . locX  - Local solution
17508449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
17624cdb843SMatthew G. Knepley - user - The user context
17724cdb843SMatthew G. Knepley 
17824cdb843SMatthew G. Knepley   Output Parameter:
17908449791SMatthew G. Knepley . locF  - Local output vector
18024cdb843SMatthew G. Knepley 
18124cdb843SMatthew G. Knepley   Level: developer
18224cdb843SMatthew G. Knepley 
18324cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
18424cdb843SMatthew G. Knepley @*/
18508449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
18624cdb843SMatthew G. Knepley {
187c4d4a4f8SMatthew G. Knepley   PetscInt       cStart, cEnd, cEndInterior;
188*b2338dbcSToby Isaac   PetscBool      isPlex;
189*b2338dbcSToby Isaac   DM             plex = NULL;
19024cdb843SMatthew G. Knepley   PetscErrorCode ierr;
19124cdb843SMatthew G. Knepley 
19224cdb843SMatthew G. Knepley   PetscFunctionBegin;
193*b2338dbcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);CHKERRQ(ierr);
194*b2338dbcSToby Isaac   if (!isPlex) {
195*b2338dbcSToby Isaac     ierr = DMConvert(dm,DMPLEX,&plex);CHKERRQ(ierr);
196*b2338dbcSToby Isaac     if (!plex) {SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot convert to DMPlex");}
197*b2338dbcSToby Isaac     ierr = DMCopyDMTS(dm,plex);CHKERRQ(ierr);
198*b2338dbcSToby Isaac     dm = plex;
199*b2338dbcSToby Isaac   }
200c4d4a4f8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
201c4d4a4f8SMatthew G. Knepley   ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr);
202c4d4a4f8SMatthew G. Knepley   cEnd = cEndInterior < 0 ? cEnd : cEndInterior;
203c4d4a4f8SMatthew G. Knepley   ierr = DMPlexComputeResidual_Internal(dm, cStart, cEnd, time, locX, locX_t, locF, user);CHKERRQ(ierr);
204*b2338dbcSToby Isaac   if (!isPlex) {
205*b2338dbcSToby Isaac     ierr = DMDestroy(&plex);CHKERRQ(ierr);
206*b2338dbcSToby Isaac   }
20724cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
20824cdb843SMatthew G. Knepley }
2097cdb2a12SMatthew G. Knepley 
2107cdb2a12SMatthew G. Knepley #undef __FUNCT__
2117cdb2a12SMatthew G. Knepley #define __FUNCT__ "DMTSCheckFromOptions"
2120163fd50SMatthew 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)
2137cdb2a12SMatthew G. Knepley {
2147cdb2a12SMatthew G. Knepley   DM             dm;
2157cdb2a12SMatthew G. Knepley   SNES           snes;
2161878804aSMatthew G. Knepley   Vec            sol;
2177cdb2a12SMatthew G. Knepley   PetscBool      check;
2187cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
2197cdb2a12SMatthew G. Knepley 
2207cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
221c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
2227cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
22372fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
22472fd7fc8SMatthew G. Knepley   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
2257cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
2267cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
2277cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
2287cdb2a12SMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
2291878804aSMatthew G. Knepley   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
23072fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
2317cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
2327cdb2a12SMatthew G. Knepley }
233