xref: /petsc/src/ts/utils/dmplexts.c (revision 5962854d720b7b8d98c62edd758f00bbb980e600)
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 
7d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy)
8d71ae5a4SJacob Faibussowitsch {
96da023fcSToby Isaac   PetscBool isPlex;
106da023fcSToby Isaac 
116da023fcSToby Isaac   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
136da023fcSToby Isaac   if (isPlex) {
146da023fcSToby Isaac     *plex = dm;
159566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
16f7148743SMatthew G. Knepley   } else {
179566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
18f7148743SMatthew G. Knepley     if (!*plex) {
199566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm, DMPLEX, plex));
209566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
216da023fcSToby Isaac       if (copy) {
229566063dSJacob Faibussowitsch         PetscCall(DMCopyDMTS(dm, *plex));
239566063dSJacob Faibussowitsch         PetscCall(DMCopyDMSNES(dm, *plex));
249566063dSJacob Faibussowitsch         PetscCall(DMCopyAuxiliaryVec(dm, *plex));
256da023fcSToby Isaac       }
26f7148743SMatthew G. Knepley     } else {
279566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)*plex));
28f7148743SMatthew G. Knepley     }
296da023fcSToby Isaac   }
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
316da023fcSToby Isaac }
326da023fcSToby Isaac 
33*5962854dSMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVMCEED(DM dm, PetscReal time, Vec locX, Vec F, void *user)
34*5962854dSMatthew G. Knepley {
35*5962854dSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
36*5962854dSMatthew G. Knepley   PetscFV    fv;
37*5962854dSMatthew G. Knepley   Vec        locF;
38*5962854dSMatthew G. Knepley   Ceed       ceed;
39*5962854dSMatthew G. Knepley   DMCeed     sd = dm->dmceed;
40*5962854dSMatthew G. Knepley   CeedVector clocX, clocF;
41*5962854dSMatthew G. Knepley #endif
42*5962854dSMatthew G. Knepley 
43*5962854dSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
44*5962854dSMatthew G. Knepley   PetscFunctionBegin;
45*5962854dSMatthew G. Knepley   PetscCall(DMGetCeed(dm, &ceed));
46*5962854dSMatthew G. Knepley   PetscCheck(sd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This DM has no CEED data. Call DMCeedCreate() before computing the residual.");
47*5962854dSMatthew G. Knepley   if (time == 0.) PetscCall(DMCeedComputeGeometry(dm, sd));
48*5962854dSMatthew G. Knepley   PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fv));
49*5962854dSMatthew G. Knepley   PetscCall(DMPlexInsertBoundaryValuesFVM(dm, fv, locX, time, NULL));
50*5962854dSMatthew G. Knepley   PetscCall(DMGetLocalVector(dm, &locF));
51*5962854dSMatthew G. Knepley   PetscCall(VecZeroEntries(locF));
52*5962854dSMatthew G. Knepley   PetscCall(VecGetCeedVectorRead(locX, ceed, &clocX));
53*5962854dSMatthew G. Knepley   PetscCall(VecGetCeedVector(locF, ceed, &clocF));
54*5962854dSMatthew G. Knepley   PetscCallCEED(CeedOperatorApplyAdd(sd->op, clocX, clocF, CEED_REQUEST_IMMEDIATE));
55*5962854dSMatthew G. Knepley   PetscCall(VecRestoreCeedVectorRead(locX, &clocX));
56*5962854dSMatthew G. Knepley   PetscCall(VecRestoreCeedVector(locF, &clocF));
57*5962854dSMatthew G. Knepley   PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F));
58*5962854dSMatthew G. Knepley   PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F));
59*5962854dSMatthew G. Knepley   PetscCall(DMRestoreLocalVector(dm, &locF));
60*5962854dSMatthew G. Knepley   PetscCall(VecViewFromOptions(F, NULL, "-fv_rhs_view"));
61*5962854dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
62*5962854dSMatthew G. Knepley #else
63*5962854dSMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "This requires libCEED. Reconfigure using --download-libceed");
64*5962854dSMatthew G. Knepley #endif
65*5962854dSMatthew G. Knepley }
66*5962854dSMatthew G. Knepley 
6708449791SMatthew G. Knepley /*@
6820f4b53cSBarry Smith   DMPlexTSComputeRHSFunctionFVM - Form the forcing `F` from the local input `locX` using pointwise functions specified by the user
6908449791SMatthew G. Knepley 
7008449791SMatthew G. Knepley   Input Parameters:
7108449791SMatthew G. Knepley + dm   - The mesh
72b43aa488SJacob Faibussowitsch . time - The time
7308449791SMatthew G. Knepley . locX - Local solution
7408449791SMatthew G. Knepley - user - The user context
7508449791SMatthew G. Knepley 
7608449791SMatthew G. Knepley   Output Parameter:
773b16df42SMatthew G. Knepley . F - Global output vector
7808449791SMatthew G. Knepley 
7908449791SMatthew G. Knepley   Level: developer
8008449791SMatthew G. Knepley 
811cc06b55SBarry Smith .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexComputeJacobianActionFEM()`
8208449791SMatthew G. Knepley @*/
83d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
84d71ae5a4SJacob Faibussowitsch {
853b16df42SMatthew G. Knepley   Vec          locF;
864a3e9fdbSToby Isaac   IS           cellIS;
876da023fcSToby Isaac   DM           plex;
884a3e9fdbSToby Isaac   PetscInt     depth;
89d70f29a3SPierre Jolivet   PetscFormKey key = {NULL, 0, 0, 0};
90254c1ad2SMatthew G. Knepley 
91254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
929566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
949566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(plex, "dim", depth, &cellIS));
959566063dSJacob Faibussowitsch   if (!cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &cellIS));
969566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(plex, &locF));
979566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
989566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locF, user));
999566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F));
1009566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F));
1019566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(plex, &locF));
1029566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
1039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
104*5962854dSMatthew G. Knepley   PetscCall(VecViewFromOptions(F, NULL, "-fv_rhs_view"));
1053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106254c1ad2SMatthew G. Knepley }
107254c1ad2SMatthew G. Knepley 
108c5d70e09SMatthew G. Knepley /*@
10920f4b53cSBarry Smith   DMPlexTSComputeBoundary - Insert the essential boundary values into the local input `locX` and/or its time derivative `locX_t` using pointwise functions specified by the user
110c5d70e09SMatthew G. Knepley 
111c5d70e09SMatthew G. Knepley   Input Parameters:
112c5d70e09SMatthew G. Knepley + dm     - The mesh
113b43aa488SJacob Faibussowitsch . time   - The time
114a40652d4SToby Isaac . locX   - Local solution
11520f4b53cSBarry Smith . locX_t - Local solution time derivative, or `NULL`
116c5d70e09SMatthew G. Knepley - user   - The user context
117c5d70e09SMatthew G. Knepley 
118c5d70e09SMatthew G. Knepley   Level: developer
119c5d70e09SMatthew G. Knepley 
1201cc06b55SBarry Smith .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexComputeJacobianActionFEM()`
121c5d70e09SMatthew G. Knepley @*/
122d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
123d71ae5a4SJacob Faibussowitsch {
124c5d70e09SMatthew G. Knepley   DM       plex;
125ef68eab9SMatthew G. Knepley   Vec      faceGeometryFVM = NULL;
126ef68eab9SMatthew G. Knepley   PetscInt Nf, f;
127c5d70e09SMatthew G. Knepley 
128c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
1309566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(plex, &Nf));
1319586001cSMatthew G. Knepley   if (!locX_t) {
1329586001cSMatthew G. Knepley     /* This is the RHS part */
133ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
134ef68eab9SMatthew G. Knepley       PetscObject  obj;
135ef68eab9SMatthew G. Knepley       PetscClassId id;
136ef68eab9SMatthew G. Knepley 
1379566063dSJacob Faibussowitsch       PetscCall(DMGetField(plex, f, NULL, &obj));
1389566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
139ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
1409566063dSJacob Faibussowitsch         PetscCall(DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL));
141ef68eab9SMatthew G. Knepley         break;
142ef68eab9SMatthew G. Knepley       }
143ef68eab9SMatthew G. Knepley     }
1449586001cSMatthew G. Knepley   }
1459566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL));
1469566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL));
1479566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
149c5d70e09SMatthew G. Knepley }
150c5d70e09SMatthew G. Knepley 
15124cdb843SMatthew G. Knepley /*@
15220f4b53cSBarry Smith   DMPlexTSComputeIFunctionFEM - Form the local residual `locF` from the local input `locX` using pointwise functions specified by the user
15324cdb843SMatthew G. Knepley 
15424cdb843SMatthew G. Knepley   Input Parameters:
15524cdb843SMatthew G. Knepley + dm     - The mesh
156b43aa488SJacob Faibussowitsch . time   - The time
15708449791SMatthew G. Knepley . locX   - Local solution
15820f4b53cSBarry Smith . locX_t - Local solution time derivative, or `NULL`
15924cdb843SMatthew G. Knepley - user   - The user context
16024cdb843SMatthew G. Knepley 
16124cdb843SMatthew G. Knepley   Output Parameter:
16208449791SMatthew G. Knepley . locF - Local output vector
16324cdb843SMatthew G. Knepley 
16424cdb843SMatthew G. Knepley   Level: developer
16524cdb843SMatthew G. Knepley 
16642747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexTSComputeRHSFunctionFEM()`
16724cdb843SMatthew G. Knepley @*/
168d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
169d71ae5a4SJacob Faibussowitsch {
1706da023fcSToby Isaac   DM       plex;
1716528b96dSMatthew G. Knepley   IS       allcellIS;
1726528b96dSMatthew G. Knepley   PetscInt Nds, s;
17324cdb843SMatthew G. Knepley 
17424cdb843SMatthew G. Knepley   PetscFunctionBegin;
1759566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
1769566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
1779566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1786528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1796528b96dSMatthew G. Knepley     PetscDS      ds;
1806528b96dSMatthew G. Knepley     IS           cellIS;
18106ad1575SMatthew G. Knepley     PetscFormKey key;
1826528b96dSMatthew G. Knepley 
18307218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
1846528b96dSMatthew G. Knepley     key.value = 0;
1856528b96dSMatthew G. Knepley     key.field = 0;
18606ad1575SMatthew G. Knepley     key.part  = 0;
1876528b96dSMatthew G. Knepley     if (!key.label) {
1889566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
1896528b96dSMatthew G. Knepley       cellIS = allcellIS;
1906528b96dSMatthew G. Knepley     } else {
1916528b96dSMatthew G. Knepley       IS pointIS;
1926528b96dSMatthew G. Knepley 
1936528b96dSMatthew G. Knepley       key.value = 1;
1949566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
1959566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
1969566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1974a3e9fdbSToby Isaac     }
1989566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, locX_t, time, locF, user));
1999566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
2006528b96dSMatthew G. Knepley   }
2019566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2029566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20424cdb843SMatthew G. Knepley }
2057cdb2a12SMatthew G. Knepley 
206756a1f44SMatthew G. Knepley /*@
20720f4b53cSBarry Smith   DMPlexTSComputeIJacobianFEM - Form the Jacobian `Jac` from the local input `locX` using pointwise functions specified by the user
208756a1f44SMatthew G. Knepley 
209756a1f44SMatthew G. Knepley   Input Parameters:
210756a1f44SMatthew G. Knepley + dm       - The mesh
211b43aa488SJacob Faibussowitsch . time     - The time
212756a1f44SMatthew G. Knepley . locX     - Local solution
21320f4b53cSBarry Smith . locX_t   - Local solution time derivative, or `NULL`
21420f4b53cSBarry Smith . X_tShift - The multiplicative parameter for dF/du_t
215756a1f44SMatthew G. Knepley - user     - The user context
216756a1f44SMatthew G. Knepley 
21720f4b53cSBarry Smith   Output Parameters:
21820f4b53cSBarry Smith + Jac  - the Jacobian
21920f4b53cSBarry Smith - JacP - an additional approximation for the Jacobian to be used to compute the preconditioner (often is `Jac`)
220756a1f44SMatthew G. Knepley 
221756a1f44SMatthew G. Knepley   Level: developer
222756a1f44SMatthew G. Knepley 
2231cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeRHSFunctionFEM()`
224756a1f44SMatthew G. Knepley @*/
225d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
226d71ae5a4SJacob Faibussowitsch {
227756a1f44SMatthew G. Knepley   DM        plex;
2286528b96dSMatthew G. Knepley   IS        allcellIS;
229f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
2306528b96dSMatthew G. Knepley   PetscInt  Nds, s;
231756a1f44SMatthew G. Knepley 
232756a1f44SMatthew G. Knepley   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
2349566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2359566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
2366528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
2376528b96dSMatthew G. Knepley     PetscDS      ds;
2386528b96dSMatthew G. Knepley     IS           cellIS;
23906ad1575SMatthew G. Knepley     PetscFormKey key;
2406528b96dSMatthew G. Knepley 
24107218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
2426528b96dSMatthew G. Knepley     key.value = 0;
2436528b96dSMatthew G. Knepley     key.field = 0;
24406ad1575SMatthew G. Knepley     key.part  = 0;
2456528b96dSMatthew G. Knepley     if (!key.label) {
2469566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
2476528b96dSMatthew G. Knepley       cellIS = allcellIS;
2486528b96dSMatthew G. Knepley     } else {
2496528b96dSMatthew G. Knepley       IS pointIS;
2506528b96dSMatthew G. Knepley 
2516528b96dSMatthew G. Knepley       key.value = 1;
2529566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
2539566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
2549566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
2556528b96dSMatthew G. Knepley     }
2566528b96dSMatthew G. Knepley     if (!s) {
2579566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
2589566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
2599566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
2609566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
2616528b96dSMatthew G. Knepley     }
2629566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user));
2639566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
2646528b96dSMatthew G. Knepley   }
2659566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2669566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
2673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
268756a1f44SMatthew G. Knepley }
269756a1f44SMatthew G. Knepley 
270cb36c0f9SMatthew G. Knepley /*@
27120f4b53cSBarry Smith   DMPlexTSComputeRHSFunctionFEM - Form the local residual `locG` from the local input `locX` using pointwise functions specified by the user
272cb36c0f9SMatthew G. Knepley 
273cb36c0f9SMatthew G. Knepley   Input Parameters:
274cb36c0f9SMatthew G. Knepley + dm   - The mesh
275b43aa488SJacob Faibussowitsch . time - The time
276cb36c0f9SMatthew G. Knepley . locX - Local solution
277cb36c0f9SMatthew G. Knepley - user - The user context
278cb36c0f9SMatthew G. Knepley 
279cb36c0f9SMatthew G. Knepley   Output Parameter:
280cb36c0f9SMatthew G. Knepley . locG - Local output vector
281cb36c0f9SMatthew G. Knepley 
282cb36c0f9SMatthew G. Knepley   Level: developer
283cb36c0f9SMatthew G. Knepley 
2841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeIJacobianFEM()`
285cb36c0f9SMatthew G. Knepley @*/
286d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locG, void *user)
287d71ae5a4SJacob Faibussowitsch {
288cb36c0f9SMatthew G. Knepley   DM       plex;
289cb36c0f9SMatthew G. Knepley   IS       allcellIS;
290cb36c0f9SMatthew G. Knepley   PetscInt Nds, s;
291cb36c0f9SMatthew G. Knepley 
292cb36c0f9SMatthew G. Knepley   PetscFunctionBegin;
2939566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
2949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2959566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
296cb36c0f9SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
297cb36c0f9SMatthew G. Knepley     PetscDS      ds;
298cb36c0f9SMatthew G. Knepley     IS           cellIS;
299cb36c0f9SMatthew G. Knepley     PetscFormKey key;
300cb36c0f9SMatthew G. Knepley 
30107218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
302cb36c0f9SMatthew G. Knepley     key.value = 0;
303cb36c0f9SMatthew G. Knepley     key.field = 0;
304cb36c0f9SMatthew G. Knepley     key.part  = 100;
305cb36c0f9SMatthew G. Knepley     if (!key.label) {
3069566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
307cb36c0f9SMatthew G. Knepley       cellIS = allcellIS;
308cb36c0f9SMatthew G. Knepley     } else {
309cb36c0f9SMatthew G. Knepley       IS pointIS;
310cb36c0f9SMatthew G. Knepley 
311cb36c0f9SMatthew G. Knepley       key.value = 1;
3129566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
3139566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
3149566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
315cb36c0f9SMatthew G. Knepley     }
3169566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locG, user));
3179566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
318cb36c0f9SMatthew G. Knepley   }
3199566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
3209566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
3213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
322cb36c0f9SMatthew G. Knepley }
323cb36c0f9SMatthew G. Knepley 
324bee9a294SMatthew G. Knepley /*@C
325f2cacb80SMatthew G. Knepley   DMTSCheckResidual - Check the residual of the exact solution
326f2cacb80SMatthew G. Knepley 
327f2cacb80SMatthew G. Knepley   Input Parameters:
328bcf0153eSBarry Smith + ts  - the `TS` object
329bcf0153eSBarry Smith . dm  - the `DM`
330f2cacb80SMatthew G. Knepley . t   - the time
331bcf0153eSBarry Smith . u   - a `DM` vector
332bcf0153eSBarry Smith . u_t - a `DM` vector
333f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
334f2cacb80SMatthew G. Knepley 
3352fe279fdSBarry Smith   Output Parameter:
33620f4b53cSBarry Smith . residual - The residual norm of the exact solution, or `NULL`
337f2cacb80SMatthew G. Knepley 
338f2cacb80SMatthew G. Knepley   Level: developer
339f2cacb80SMatthew G. Knepley 
3401cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCheckFromOptions()`, `DMTSCheckJacobian()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
341f2cacb80SMatthew G. Knepley @*/
342d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual)
343d71ae5a4SJacob Faibussowitsch {
344f2cacb80SMatthew G. Knepley   MPI_Comm  comm;
345f2cacb80SMatthew G. Knepley   Vec       r;
346f2cacb80SMatthew G. Knepley   PetscReal res;
347f2cacb80SMatthew G. Knepley 
348f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
349f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
350f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
351064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
3524f572ea9SToby Isaac   if (residual) PetscAssertPointer(residual, 7);
3539566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
3549566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
3559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
3569566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
3579566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
358f2cacb80SMatthew G. Knepley   if (tol >= 0.0) {
3593c633725SBarry Smith     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
360f2cacb80SMatthew G. Knepley   } else if (residual) {
361f2cacb80SMatthew G. Knepley     *residual = res;
362f2cacb80SMatthew G. Knepley   } else {
3639566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
36493ec0da9SPierre Jolivet     PetscCall(VecFilter(r, 1.0e-10));
3659566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)dm));
3669566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
3679566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
3689566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
3699566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL));
370f2cacb80SMatthew G. Knepley   }
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
373f2cacb80SMatthew G. Knepley }
374f2cacb80SMatthew G. Knepley 
375f2cacb80SMatthew G. Knepley /*@C
376f2cacb80SMatthew G. Knepley   DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
377f2cacb80SMatthew G. Knepley 
378f2cacb80SMatthew G. Knepley   Input Parameters:
37920f4b53cSBarry Smith + ts  - the `TS` object
38020f4b53cSBarry Smith . dm  - the `DM`
381f2cacb80SMatthew G. Knepley . t   - the time
38220f4b53cSBarry Smith . u   - a `DM` vector
38320f4b53cSBarry Smith . u_t - a `DM` vector
384f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
385f2cacb80SMatthew G. Knepley 
386f2cacb80SMatthew G. Knepley   Output Parameters:
38720f4b53cSBarry Smith + isLinear - Flag indicaing that the function looks linear, or `NULL`
38820f4b53cSBarry Smith - convRate - The rate of convergence of the linear model, or `NULL`
389f2cacb80SMatthew G. Knepley 
390f2cacb80SMatthew G. Knepley   Level: developer
391f2cacb80SMatthew G. Knepley 
3921cc06b55SBarry Smith .seealso: [](ch_ts), `DNTSCheckFromOptions()`, `DMTSCheckResidual()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
393f2cacb80SMatthew G. Knepley @*/
394d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
395d71ae5a4SJacob Faibussowitsch {
396f2cacb80SMatthew G. Knepley   MPI_Comm     comm;
397f2cacb80SMatthew G. Knepley   PetscDS      ds;
398f2cacb80SMatthew G. Knepley   Mat          J, M;
399f2cacb80SMatthew G. Knepley   MatNullSpace nullspace;
400f2cacb80SMatthew G. Knepley   PetscReal    dt, shift, slope, intercept;
401f2cacb80SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
402f2cacb80SMatthew G. Knepley 
403f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
404f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
405f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
406064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
4074f572ea9SToby Isaac   if (isLinear) PetscAssertPointer(isLinear, 7);
4084f572ea9SToby Isaac   if (convRate) PetscAssertPointer(convRate, 8);
4099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
4109566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
411f2cacb80SMatthew G. Knepley   /* Create and view matrices */
4129566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
413f2cacb80SMatthew G. Knepley   shift = 1.0 / dt;
4149566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
4159566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
4169566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
4179566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
418f2cacb80SMatthew G. Knepley   if (hasJac && hasPrec) {
4199566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
4209566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE));
4219566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix"));
4229566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
4239566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
4249566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
425f2cacb80SMatthew G. Knepley   } else {
4269566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE));
427f2cacb80SMatthew G. Knepley   }
4289566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
4299566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
4309566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
431f2cacb80SMatthew G. Knepley   /* Check nullspace */
4329566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
433f2cacb80SMatthew G. Knepley   if (nullspace) {
434f2cacb80SMatthew G. Knepley     PetscBool isNull;
4359566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
4363c633725SBarry Smith     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
437f2cacb80SMatthew G. Knepley   }
438f2cacb80SMatthew G. Knepley   /* Taylor test */
439f2cacb80SMatthew G. Knepley   {
440f2cacb80SMatthew G. Knepley     PetscRandom rand;
441f2cacb80SMatthew G. Knepley     Vec         du, uhat, uhat_t, r, rhat, df;
442f2cacb80SMatthew G. Knepley     PetscReal   h;
443f2cacb80SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
444f2cacb80SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
445f2cacb80SMatthew G. Knepley     PetscInt    Nv, v;
446f2cacb80SMatthew G. Knepley 
447f2cacb80SMatthew G. Knepley     /* Choose a perturbation direction */
4489566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
4499566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
4509566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
4519566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
4529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
4539566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
454f2cacb80SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
4559566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
4569566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
457f2cacb80SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
4589371c9d4SSatish Balay     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv)
4599371c9d4SSatish Balay       ;
4609566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
4619566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
4629566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat_t));
4639566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
464f2cacb80SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
4659566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
4669566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat_t, h * shift, du, u_t));
467f2cacb80SMatthew G. Knepley       /* F(\hat u, \hat u_t) \approx F(u, u_t) + J(u, u_t) (uhat - u) + J_t(u, u_t) (uhat_t - u_t) = F(u) + h * J(u) du + h * shift * J_t(u) du = F(u) + h F' du */
4689566063dSJacob Faibussowitsch       PetscCall(TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE));
4699566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
4709566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
471f2cacb80SMatthew G. Knepley 
472f2cacb80SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
473f2cacb80SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
474f2cacb80SMatthew G. Knepley     }
4759566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
4769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat_t));
4779566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
4789566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
4799566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
4809566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
481f2cacb80SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
482f2cacb80SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
483f2cacb80SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL) break;
484f2cacb80SMatthew G. Knepley     }
485f2cacb80SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
4869566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
4879566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
488f2cacb80SMatthew G. Knepley     /* Slope should be about 2 */
489f2cacb80SMatthew G. Knepley     if (tol >= 0) {
4903c633725SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
491f2cacb80SMatthew G. Knepley     } else if (isLinear || convRate) {
492f2cacb80SMatthew G. Knepley       if (isLinear) *isLinear = isLin;
493f2cacb80SMatthew G. Knepley       if (convRate) *convRate = slope;
494f2cacb80SMatthew G. Knepley     } else {
4959566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
4969566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
497f2cacb80SMatthew G. Knepley     }
498f2cacb80SMatthew G. Knepley   }
4999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
5003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
501f2cacb80SMatthew G. Knepley }
502f2cacb80SMatthew G. Knepley 
503f2cacb80SMatthew G. Knepley /*@C
50420f4b53cSBarry Smith   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information based on
50520f4b53cSBarry Smith   values in the options database
506bee9a294SMatthew G. Knepley 
507bee9a294SMatthew G. Knepley   Input Parameters:
508bcf0153eSBarry Smith + ts - the `TS` object
509bcf0153eSBarry Smith - u  - representative `TS` vector
5107f96f943SMatthew G. Knepley 
51120f4b53cSBarry Smith   Level: developer
51220f4b53cSBarry Smith 
513bcf0153eSBarry Smith   Note:
514bcf0153eSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
515bee9a294SMatthew G. Knepley 
516b43aa488SJacob Faibussowitsch   Developer Notes:
51720f4b53cSBarry Smith   What is the purpose of `u`, does it need to already have a solution or some other value in it?
51820f4b53cSBarry Smith 
51920f4b53cSBarry Smith .seealso: `DMTS`
520bee9a294SMatthew G. Knepley @*/
521d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u)
522d71ae5a4SJacob Faibussowitsch {
5237cdb2a12SMatthew G. Knepley   DM        dm;
5247cdb2a12SMatthew G. Knepley   SNES      snes;
525f2cacb80SMatthew G. Knepley   Vec       sol, u_t;
526f2cacb80SMatthew G. Knepley   PetscReal t;
5277cdb2a12SMatthew G. Knepley   PetscBool check;
5287cdb2a12SMatthew G. Knepley 
5297cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
5309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-dmts_check", &check));
5313ba16761SJacob Faibussowitsch   if (!check) PetscFunctionReturn(PETSC_SUCCESS);
5329566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
5339566063dSJacob Faibussowitsch   PetscCall(VecCopy(u, sol));
5349566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, u));
5359566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
5369566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
5379566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
5389566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, u));
539f2cacb80SMatthew G. Knepley 
5409566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
5419566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL));
5429566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u_t));
5439566063dSJacob Faibussowitsch   PetscCall(DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL));
5449566063dSJacob Faibussowitsch   PetscCall(DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL));
5459566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u_t));
546f2cacb80SMatthew G. Knepley 
5479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5497cdb2a12SMatthew G. Knepley }
550