xref: /petsc/src/ts/utils/dmplexts.c (revision 37c2070c1994685cd908891d4d4f0bed1f8a27a1)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
3*37c2070cSMatthew 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;
28a0ac79e7SMatthew G. Knepley   PetscErrorCode ierr;
29a0ac79e7SMatthew G. Knepley 
30a0ac79e7SMatthew G. Knepley   PetscFunctionBegin;
31924a1b8fSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
32a0ac79e7SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
33b6aca0f9SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", &obj);CHKERRQ(ierr);
341faf85eaSMatthew G. Knepley   if (!obj) {
351faf85eaSMatthew G. Knepley     Vec cellgeom, facegeom;
361faf85eaSMatthew G. Knepley 
371faf85eaSMatthew G. Knepley     ierr = DMPlexComputeGeometryFVM(dm, &cellgeom, &facegeom);CHKERRQ(ierr);
38b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject) facegeom);CHKERRQ(ierr);
39b6aca0f9SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject) cellgeom);CHKERRQ(ierr);
401faf85eaSMatthew G. Knepley     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
411faf85eaSMatthew G. Knepley     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
421faf85eaSMatthew G. Knepley   }
43b6aca0f9SMatthew G. Knepley   if (facegeom) {PetscValidPointer(facegeom, 2); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_facegeom_fvm", (PetscObject *) facegeom);CHKERRQ(ierr);}
44b6aca0f9SMatthew G. Knepley   if (cellgeom) {PetscValidPointer(cellgeom, 3); ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_cellgeom_fvm", (PetscObject *) cellgeom);CHKERRQ(ierr);}
45113c68e6SMatthew G. Knepley   if (minRadius) {ierr = DMPlexGetMinRadius(dm, minRadius);CHKERRQ(ierr);}
46924a1b8fSMatthew G. Knepley   PetscFunctionReturn(0);
47924a1b8fSMatthew G. Knepley }
48924a1b8fSMatthew G. Knepley 
49924a1b8fSMatthew G. Knepley #undef __FUNCT__
50c49ccbb3SMatthew G. Knepley #define __FUNCT__ "DMPlexTSGetGradientDM"
51c49ccbb3SMatthew G. Knepley /*@C
52c49ccbb3SMatthew G. Knepley   DMPlexTSGetGradientDM - Return gradient data layout
53c49ccbb3SMatthew G. Knepley 
54c49ccbb3SMatthew G. Knepley   Input Parameters:
55c49ccbb3SMatthew G. Knepley + dm - The DM
56c49ccbb3SMatthew G. Knepley - fv - The PetscFV
57c49ccbb3SMatthew G. Knepley 
58c49ccbb3SMatthew G. Knepley   Output Parameter:
59c49ccbb3SMatthew G. Knepley . dmGrad - The layout for gradient values
60c49ccbb3SMatthew G. Knepley 
61c49ccbb3SMatthew G. Knepley   Level: developer
62c49ccbb3SMatthew G. Knepley 
63c49ccbb3SMatthew G. Knepley .seealso: DMPlexTSGetGeometryFVM(), DMPlexTSSetRHSFunctionLocal()
64c49ccbb3SMatthew G. Knepley @*/
65c49ccbb3SMatthew G. Knepley PetscErrorCode DMPlexTSGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
66c49ccbb3SMatthew G. Knepley {
67c49ccbb3SMatthew G. Knepley   DMTS           dmts;
68c49ccbb3SMatthew G. Knepley   PetscObject    obj;
69c49ccbb3SMatthew G. Knepley   PetscBool      computeGradients;
70c49ccbb3SMatthew G. Knepley   PetscErrorCode ierr;
71c49ccbb3SMatthew G. Knepley 
72c49ccbb3SMatthew G. Knepley   PetscFunctionBegin;
73c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
74c49ccbb3SMatthew G. Knepley   PetscValidHeaderSpecific(fv,PETSCFV_CLASSID,2);
75c49ccbb3SMatthew G. Knepley   PetscValidPointer(dmGrad,3);
76c49ccbb3SMatthew G. Knepley   ierr = PetscFVGetComputeGradients(fv, &computeGradients);CHKERRQ(ierr);
77c49ccbb3SMatthew G. Knepley   if (!computeGradients) {*dmGrad = NULL; PetscFunctionReturn(0);}
78c49ccbb3SMatthew G. Knepley   ierr = DMGetDMTS(dm, &dmts);CHKERRQ(ierr);
79c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", &obj);CHKERRQ(ierr);
80c49ccbb3SMatthew G. Knepley   if (!obj) {
81c49ccbb3SMatthew G. Knepley     DM  dmGrad;
8208449791SMatthew G. Knepley     Vec faceGeometry, cellGeometry;
83c49ccbb3SMatthew G. Knepley 
8408449791SMatthew G. Knepley     ierr = DMPlexTSGetGeometryFVM(dm, &faceGeometry, &cellGeometry, NULL);CHKERRQ(ierr);
8508449791SMatthew G. Knepley     ierr = DMPlexComputeGradientFVM(dm, fv, faceGeometry, cellGeometry, &dmGrad);CHKERRQ(ierr);
86c49ccbb3SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject) dmGrad);CHKERRQ(ierr);
87924a1b8fSMatthew G. Knepley     ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
88c49ccbb3SMatthew G. Knepley   }
89c49ccbb3SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dmts, "DMPlexTS_dmgrad_fvm", (PetscObject *) dmGrad);CHKERRQ(ierr);
90254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
91254c1ad2SMatthew G. Knepley }
92254c1ad2SMatthew G. Knepley 
93254c1ad2SMatthew G. Knepley #undef __FUNCT__
94924a1b8fSMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeRHSFunctionFVM"
9508449791SMatthew G. Knepley /*@
9608449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
9708449791SMatthew G. Knepley 
9808449791SMatthew G. Knepley   Input Parameters:
9908449791SMatthew G. Knepley + dm - The mesh
10008449791SMatthew G. Knepley . t - The time
10108449791SMatthew G. Knepley . locX  - Local solution
10208449791SMatthew G. Knepley - user - The user context
10308449791SMatthew G. Knepley 
10408449791SMatthew G. Knepley   Output Parameter:
1053b16df42SMatthew G. Knepley . F  - Global output vector
10608449791SMatthew G. Knepley 
10708449791SMatthew G. Knepley   Level: developer
10808449791SMatthew G. Knepley 
10908449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
11008449791SMatthew G. Knepley @*/
1113b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
112254c1ad2SMatthew G. Knepley {
1133b16df42SMatthew G. Knepley   Vec            locF;
114254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
115254c1ad2SMatthew G. Knepley 
116254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
1173b16df42SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &locF);CHKERRQ(ierr);
1183b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
11908449791SMatthew G. Knepley   ierr = DMPlexComputeResidual_Internal(dm, time, locX, NULL, locF, user);CHKERRQ(ierr);
1203b16df42SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1213b16df42SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, locF, INSERT_VALUES, F);CHKERRQ(ierr);
1223b16df42SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &locF);CHKERRQ(ierr);
123254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
124254c1ad2SMatthew G. Knepley }
125254c1ad2SMatthew G. Knepley 
126254c1ad2SMatthew G. Knepley #undef __FUNCT__
12724cdb843SMatthew G. Knepley #define __FUNCT__ "DMPlexTSComputeIFunctionFEM"
12824cdb843SMatthew G. Knepley /*@
12924cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
13024cdb843SMatthew G. Knepley 
13124cdb843SMatthew G. Knepley   Input Parameters:
13224cdb843SMatthew G. Knepley + dm - The mesh
13324cdb843SMatthew G. Knepley . t - The time
13408449791SMatthew G. Knepley . locX  - Local solution
13508449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
13624cdb843SMatthew G. Knepley - user - The user context
13724cdb843SMatthew G. Knepley 
13824cdb843SMatthew G. Knepley   Output Parameter:
13908449791SMatthew G. Knepley . locF  - Local output vector
14024cdb843SMatthew G. Knepley 
14124cdb843SMatthew G. Knepley   Level: developer
14224cdb843SMatthew G. Knepley 
14324cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
14424cdb843SMatthew G. Knepley @*/
14508449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
14624cdb843SMatthew G. Knepley {
14724cdb843SMatthew G. Knepley   PetscErrorCode ierr;
14824cdb843SMatthew G. Knepley 
14924cdb843SMatthew G. Knepley   PetscFunctionBegin;
15008449791SMatthew G. Knepley   ierr = DMPlexComputeResidual_Internal(dm, time, locX, locX_t, locF, user);CHKERRQ(ierr);
15124cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
15224cdb843SMatthew G. Knepley }
1537cdb2a12SMatthew G. Knepley 
1547cdb2a12SMatthew G. Knepley #undef __FUNCT__
1557cdb2a12SMatthew G. Knepley #define __FUNCT__ "DMTSCheckFromOptions"
156ad917190SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
1577cdb2a12SMatthew G. Knepley {
1587cdb2a12SMatthew G. Knepley   DM             dm;
1597cdb2a12SMatthew G. Knepley   SNES           snes;
1601878804aSMatthew G. Knepley   Vec            sol;
1617cdb2a12SMatthew G. Knepley   PetscBool      check;
1627cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
1637cdb2a12SMatthew G. Knepley 
1647cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
1657cdb2a12SMatthew G. Knepley   ierr = PetscOptionsHasName(ts->hdr.prefix, "-dmts_check", &check);CHKERRQ(ierr);
1667cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
16772fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
16872fd7fc8SMatthew G. Knepley   ierr = TSSetSolution(ts, sol);CHKERRQ(ierr);
1697cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
1707cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
1717cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
1727cdb2a12SMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
1731878804aSMatthew G. Knepley   ierr = DMSNESCheckFromOptions_Internal(snes, dm, u, sol, exactFuncs, ctxs);CHKERRQ(ierr);
17472fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
1757cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
1767cdb2a12SMatthew G. Knepley }
177