xref: /petsc/src/ts/utils/dmplexts.c (revision 3e9753d60b406ac2fec1249fa8fb9cfa31a5ff24)
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 
76da023fcSToby Isaac static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy)
86da023fcSToby Isaac {
96da023fcSToby Isaac   PetscBool      isPlex;
106da023fcSToby Isaac   PetscErrorCode ierr;
116da023fcSToby Isaac 
126da023fcSToby Isaac   PetscFunctionBegin;
136da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
146da023fcSToby Isaac   if (isPlex) {
156da023fcSToby Isaac     *plex = dm;
166da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
17f7148743SMatthew G. Knepley   } else {
18f7148743SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
19f7148743SMatthew G. Knepley     if (!*plex) {
206da023fcSToby Isaac       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
21f7148743SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
226da023fcSToby Isaac       if (copy) {
236da023fcSToby Isaac         PetscInt    i;
246da023fcSToby Isaac         PetscObject obj;
256da023fcSToby Isaac         const char *comps[3] = {"A","dmAux","dmCh"};
266da023fcSToby Isaac 
276da023fcSToby Isaac         ierr = DMCopyDMTS(dm, *plex);CHKERRQ(ierr);
286da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
296da023fcSToby Isaac         for (i = 0; i < 3; i++) {
306da023fcSToby Isaac           ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr);
316da023fcSToby Isaac           ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr);
326da023fcSToby Isaac         }
336da023fcSToby Isaac       }
34f7148743SMatthew G. Knepley     } else {
35f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
36f7148743SMatthew G. Knepley     }
376da023fcSToby Isaac   }
386da023fcSToby Isaac   PetscFunctionReturn(0);
396da023fcSToby Isaac }
406da023fcSToby Isaac 
4108449791SMatthew G. Knepley /*@
4208449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
4308449791SMatthew G. Knepley 
4408449791SMatthew G. Knepley   Input Parameters:
4508449791SMatthew G. Knepley + dm - The mesh
4608449791SMatthew G. Knepley . t - The time
4708449791SMatthew G. Knepley . locX  - Local solution
4808449791SMatthew G. Knepley - user - The user context
4908449791SMatthew G. Knepley 
5008449791SMatthew G. Knepley   Output Parameter:
513b16df42SMatthew G. Knepley . F  - Global output vector
5208449791SMatthew G. Knepley 
5308449791SMatthew G. Knepley   Level: developer
5408449791SMatthew G. Knepley 
5508449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
5608449791SMatthew G. Knepley @*/
573b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
58254c1ad2SMatthew G. Knepley {
593b16df42SMatthew G. Knepley   Vec            locF;
604a3e9fdbSToby Isaac   IS             cellIS;
616da023fcSToby Isaac   DM             plex;
624a3e9fdbSToby Isaac   PetscInt       depth;
63254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
64254c1ad2SMatthew G. Knepley 
65254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
666da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
674a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
68aeadca18SToby Isaac   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
694a3e9fdbSToby Isaac   if (!cellIS) {
704a3e9fdbSToby Isaac     ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);
714a3e9fdbSToby Isaac   }
726da023fcSToby Isaac   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
733b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
744a3e9fdbSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cellIS, time, locX, NULL, time, locF, user);CHKERRQ(ierr);
759bda831aSToby Isaac   ierr = DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
769bda831aSToby Isaac   ierr = DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
776da023fcSToby Isaac   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
784a3e9fdbSToby Isaac   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
79b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
80254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
81254c1ad2SMatthew G. Knepley }
82254c1ad2SMatthew G. Knepley 
83c5d70e09SMatthew G. Knepley /*@
84ef68eab9SMatthew G. Knepley   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
85c5d70e09SMatthew G. Knepley 
86c5d70e09SMatthew G. Knepley   Input Parameters:
87c5d70e09SMatthew G. Knepley + dm - The mesh
88a40652d4SToby Isaac . t - The time
89a40652d4SToby Isaac . locX  - Local solution
90a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL
91c5d70e09SMatthew G. Knepley - user - The user context
92c5d70e09SMatthew G. Knepley 
93c5d70e09SMatthew G. Knepley   Level: developer
94c5d70e09SMatthew G. Knepley 
95a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM()
96c5d70e09SMatthew G. Knepley @*/
97ef68eab9SMatthew G. Knepley PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
98c5d70e09SMatthew G. Knepley {
99c5d70e09SMatthew G. Knepley   DM             plex;
100ef68eab9SMatthew G. Knepley   Vec            faceGeometryFVM = NULL;
101ef68eab9SMatthew G. Knepley   PetscInt       Nf, f;
102c5d70e09SMatthew G. Knepley   PetscErrorCode ierr;
103c5d70e09SMatthew G. Knepley 
104c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
105c5d70e09SMatthew G. Knepley   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
106ef68eab9SMatthew G. Knepley   ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr);
1079586001cSMatthew G. Knepley   if (!locX_t) {
1089586001cSMatthew G. Knepley     /* This is the RHS part */
109ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
110ef68eab9SMatthew G. Knepley       PetscObject  obj;
111ef68eab9SMatthew G. Knepley       PetscClassId id;
112ef68eab9SMatthew G. Knepley 
11344a7f3ddSMatthew G. Knepley       ierr = DMGetField(plex, f, NULL, &obj);CHKERRQ(ierr);
114ef68eab9SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
115ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
116*3e9753d6SMatthew G. Knepley         ierr = DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
117ef68eab9SMatthew G. Knepley         break;
118ef68eab9SMatthew G. Knepley       }
119ef68eab9SMatthew G. Knepley     }
1209586001cSMatthew G. Knepley   }
1219586001cSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
122a40652d4SToby Isaac   /* TODO: locX_t */
123c5d70e09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
124c5d70e09SMatthew G. Knepley   PetscFunctionReturn(0);
125c5d70e09SMatthew G. Knepley }
126c5d70e09SMatthew G. Knepley 
12724cdb843SMatthew G. Knepley /*@
12824cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
12924cdb843SMatthew G. Knepley 
13024cdb843SMatthew G. Knepley   Input Parameters:
13124cdb843SMatthew G. Knepley + dm - The mesh
13224cdb843SMatthew G. Knepley . t - The time
13308449791SMatthew G. Knepley . locX  - Local solution
13408449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
13524cdb843SMatthew G. Knepley - user - The user context
13624cdb843SMatthew G. Knepley 
13724cdb843SMatthew G. Knepley   Output Parameter:
13808449791SMatthew G. Knepley . locF  - Local output vector
13924cdb843SMatthew G. Knepley 
14024cdb843SMatthew G. Knepley   Level: developer
14124cdb843SMatthew G. Knepley 
14224cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
14324cdb843SMatthew G. Knepley @*/
14408449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
14524cdb843SMatthew G. Knepley {
1466da023fcSToby Isaac   DM             plex;
1474a3e9fdbSToby Isaac   IS             cellIS;
1484a3e9fdbSToby Isaac   PetscInt       depth;
14924cdb843SMatthew G. Knepley   PetscErrorCode ierr;
15024cdb843SMatthew G. Knepley 
15124cdb843SMatthew G. Knepley   PetscFunctionBegin;
1526da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1534a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
154aeadca18SToby Isaac   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
1554a3e9fdbSToby Isaac   if (!cellIS) {
1564a3e9fdbSToby Isaac     ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);
1574a3e9fdbSToby Isaac   }
1584a3e9fdbSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cellIS, time, locX, locX_t, time, locF, user);CHKERRQ(ierr);
1594a3e9fdbSToby Isaac   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
160b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
16124cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
16224cdb843SMatthew G. Knepley }
1637cdb2a12SMatthew G. Knepley 
164756a1f44SMatthew G. Knepley /*@
165756a1f44SMatthew G. Knepley   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
166756a1f44SMatthew G. Knepley 
167756a1f44SMatthew G. Knepley   Input Parameters:
168756a1f44SMatthew G. Knepley + dm - The mesh
169756a1f44SMatthew G. Knepley . t - The time
170756a1f44SMatthew G. Knepley . locX  - Local solution
171756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
172756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t
173756a1f44SMatthew G. Knepley - user - The user context
174756a1f44SMatthew G. Knepley 
175756a1f44SMatthew G. Knepley   Output Parameter:
176756a1f44SMatthew G. Knepley . locF  - Local output vector
177756a1f44SMatthew G. Knepley 
178756a1f44SMatthew G. Knepley   Level: developer
179756a1f44SMatthew G. Knepley 
180756a1f44SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
181756a1f44SMatthew G. Knepley @*/
182756a1f44SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
183756a1f44SMatthew G. Knepley {
184756a1f44SMatthew G. Knepley   DM             plex;
185f04eb4edSMatthew G. Knepley   PetscDS        prob;
186f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
1879044fa66SMatthew G. Knepley   IS             cellIS;
1884a3e9fdbSToby Isaac   PetscInt       depth;
189756a1f44SMatthew G. Knepley   PetscErrorCode ierr;
190756a1f44SMatthew G. Knepley 
191756a1f44SMatthew G. Knepley   PetscFunctionBegin;
192756a1f44SMatthew G. Knepley   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
193f04eb4edSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
194f04eb4edSMatthew G. Knepley   ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr);
195f04eb4edSMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr);
196f04eb4edSMatthew G. Knepley   if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
197f04eb4edSMatthew G. Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1984a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex,&depth);CHKERRQ(ierr);
199aeadca18SToby Isaac   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
2009044fa66SMatthew G. Knepley   if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);}
2014a3e9fdbSToby Isaac   ierr = DMPlexComputeJacobian_Internal(plex, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
2024a3e9fdbSToby Isaac   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
203756a1f44SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
204756a1f44SMatthew G. Knepley   PetscFunctionReturn(0);
205756a1f44SMatthew G. Knepley }
206756a1f44SMatthew G. Knepley 
207bee9a294SMatthew G. Knepley /*@C
208bee9a294SMatthew G. Knepley   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
209bee9a294SMatthew G. Knepley 
210bee9a294SMatthew G. Knepley   Input Parameters:
211bee9a294SMatthew G. Knepley + ts - the TS object
2127f96f943SMatthew G. Knepley - u  - representative TS vector
2137f96f943SMatthew G. Knepley 
2147f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
215bee9a294SMatthew G. Knepley 
216bee9a294SMatthew G. Knepley   Level: developer
217bee9a294SMatthew G. Knepley @*/
2187f96f943SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u)
2197cdb2a12SMatthew G. Knepley {
2207cdb2a12SMatthew G. Knepley   DM             dm;
2217cdb2a12SMatthew G. Knepley   SNES           snes;
2221878804aSMatthew G. Knepley   Vec            sol;
2237cdb2a12SMatthew G. Knepley   PetscBool      check;
2247cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
2257cdb2a12SMatthew G. Knepley 
2267cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
227c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
2287cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
22972fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
2306fb91e28SMatthew G. Knepley   ierr = VecCopy(u, sol);CHKERRQ(ierr);
2316fb91e28SMatthew G. Knepley   ierr = TSSetSolution(ts, u);CHKERRQ(ierr);
2327cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
2337cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
2347cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
2356fb91e28SMatthew G. Knepley   ierr = SNESSetSolution(snes, u);CHKERRQ(ierr);
2367f96f943SMatthew G. Knepley   ierr = DMSNESCheck_Internal(snes, dm, sol);CHKERRQ(ierr);
23772fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
2387cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
2397cdb2a12SMatthew G. Knepley }
240