xref: /petsc/src/ts/utils/dmplexts.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 
DMTSConvertPlex(DM dm,DM * plex,PetscBool copy)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));
21cbf8eb3cSStefano Zampini     } else {
22cbf8eb3cSStefano Zampini       PetscCall(PetscObjectReference((PetscObject)*plex));
23cbf8eb3cSStefano Zampini     }
246da023fcSToby Isaac     if (copy) {
259566063dSJacob Faibussowitsch       PetscCall(DMCopyDMTS(dm, *plex));
269566063dSJacob Faibussowitsch       PetscCall(DMCopyDMSNES(dm, *plex));
279566063dSJacob Faibussowitsch       PetscCall(DMCopyAuxiliaryVec(dm, *plex));
286da023fcSToby Isaac     }
296da023fcSToby Isaac   }
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
316da023fcSToby Isaac }
326da023fcSToby Isaac 
3308449791SMatthew G. Knepley /*@
3420f4b53cSBarry Smith   DMPlexTSComputeRHSFunctionFVM - Form the forcing `F` from the local input `locX` using pointwise functions specified by the user
3508449791SMatthew G. Knepley 
3608449791SMatthew G. Knepley   Input Parameters:
3708449791SMatthew G. Knepley + dm   - The mesh
38b43aa488SJacob Faibussowitsch . time - The time
3908449791SMatthew G. Knepley . locX - Local solution
40*2a8381b2SBarry Smith - ctx  - The application context
4108449791SMatthew G. Knepley 
4208449791SMatthew G. Knepley   Output Parameter:
433b16df42SMatthew G. Knepley . F - Global output vector
4408449791SMatthew G. Knepley 
4508449791SMatthew G. Knepley   Level: developer
4608449791SMatthew G. Knepley 
471cc06b55SBarry Smith .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexComputeJacobianActionFEM()`
4808449791SMatthew G. Knepley @*/
DMPlexTSComputeRHSFunctionFVM(DM dm,PetscReal time,Vec locX,Vec F,PetscCtx ctx)49*2a8381b2SBarry Smith PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, PetscCtx ctx)
50d71ae5a4SJacob Faibussowitsch {
513b16df42SMatthew G. Knepley   Vec          locF;
524a3e9fdbSToby Isaac   IS           cellIS;
536da023fcSToby Isaac   DM           plex;
544a3e9fdbSToby Isaac   PetscInt     depth;
55d70f29a3SPierre Jolivet   PetscFormKey key = {NULL, 0, 0, 0};
56254c1ad2SMatthew G. Knepley 
57254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
609566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(plex, "dim", depth, &cellIS));
619566063dSJacob Faibussowitsch   if (!cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &cellIS));
629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(plex, &locF));
639566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
64*2a8381b2SBarry Smith   PetscCall(DMPlexComputeResidualByKey(plex, key, cellIS, time, locX, NULL, time, locF, ctx));
659566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F));
669566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F));
679566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(plex, &locF));
689566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
699566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
705962854dSMatthew G. Knepley   PetscCall(VecViewFromOptions(F, NULL, "-fv_rhs_view"));
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72254c1ad2SMatthew G. Knepley }
73254c1ad2SMatthew G. Knepley 
74c5d70e09SMatthew G. Knepley /*@
7520f4b53cSBarry 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
76c5d70e09SMatthew G. Knepley 
77c5d70e09SMatthew G. Knepley   Input Parameters:
78c5d70e09SMatthew G. Knepley + dm     - The mesh
79b43aa488SJacob Faibussowitsch . time   - The time
80a40652d4SToby Isaac . locX   - Local solution
8120f4b53cSBarry Smith . locX_t - Local solution time derivative, or `NULL`
82*2a8381b2SBarry Smith - ctx    - The application context
83c5d70e09SMatthew G. Knepley 
84c5d70e09SMatthew G. Knepley   Level: developer
85c5d70e09SMatthew G. Knepley 
861cc06b55SBarry Smith .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexComputeJacobianActionFEM()`
87c5d70e09SMatthew G. Knepley @*/
DMPlexTSComputeBoundary(DM dm,PetscReal time,Vec locX,Vec locX_t,PetscCtx ctx)88*2a8381b2SBarry Smith PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscCtx ctx)
89d71ae5a4SJacob Faibussowitsch {
90c5d70e09SMatthew G. Knepley   DM       plex;
91ef68eab9SMatthew G. Knepley   Vec      faceGeometryFVM = NULL;
92ef68eab9SMatthew G. Knepley   PetscInt Nf, f;
93c5d70e09SMatthew G. Knepley 
94c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
969566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(plex, &Nf));
979586001cSMatthew G. Knepley   if (!locX_t) {
989586001cSMatthew G. Knepley     /* This is the RHS part */
99ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
100ef68eab9SMatthew G. Knepley       PetscObject  obj;
101ef68eab9SMatthew G. Knepley       PetscClassId id;
102ef68eab9SMatthew G. Knepley 
1039566063dSJacob Faibussowitsch       PetscCall(DMGetField(plex, f, NULL, &obj));
1049566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
105ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
1069566063dSJacob Faibussowitsch         PetscCall(DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL));
107ef68eab9SMatthew G. Knepley         break;
108ef68eab9SMatthew G. Knepley       }
109ef68eab9SMatthew G. Knepley     }
1109586001cSMatthew G. Knepley   }
1119566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL));
1129566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL));
1139566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115c5d70e09SMatthew G. Knepley }
116c5d70e09SMatthew G. Knepley 
11724cdb843SMatthew G. Knepley /*@
11820f4b53cSBarry Smith   DMPlexTSComputeIFunctionFEM - Form the local residual `locF` from the local input `locX` using pointwise functions specified by the user
11924cdb843SMatthew G. Knepley 
12024cdb843SMatthew G. Knepley   Input Parameters:
12124cdb843SMatthew G. Knepley + dm     - The mesh
122b43aa488SJacob Faibussowitsch . time   - The time
12308449791SMatthew G. Knepley . locX   - Local solution
12420f4b53cSBarry Smith . locX_t - Local solution time derivative, or `NULL`
125*2a8381b2SBarry Smith - ctx    - The application context
12624cdb843SMatthew G. Knepley 
12724cdb843SMatthew G. Knepley   Output Parameter:
12808449791SMatthew G. Knepley . locF - Local output vector
12924cdb843SMatthew G. Knepley 
13024cdb843SMatthew G. Knepley   Level: developer
13124cdb843SMatthew G. Knepley 
13242747ad1SJacob Faibussowitsch .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexTSComputeRHSFunctionFEM()`
13324cdb843SMatthew G. Knepley @*/
DMPlexTSComputeIFunctionFEM(DM dm,PetscReal time,Vec locX,Vec locX_t,Vec locF,PetscCtx ctx)134*2a8381b2SBarry Smith PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, PetscCtx ctx)
135d71ae5a4SJacob Faibussowitsch {
1366da023fcSToby Isaac   DM       plex;
1376528b96dSMatthew G. Knepley   IS       allcellIS;
1386528b96dSMatthew G. Knepley   PetscInt Nds, s;
13924cdb843SMatthew G. Knepley 
14024cdb843SMatthew G. Knepley   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
1429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
1439566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1446528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1456528b96dSMatthew G. Knepley     PetscDS      ds;
1466528b96dSMatthew G. Knepley     IS           cellIS;
14706ad1575SMatthew G. Knepley     PetscFormKey key;
1486528b96dSMatthew G. Knepley 
14907218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
1506528b96dSMatthew G. Knepley     key.value = 0;
1516528b96dSMatthew G. Knepley     key.field = 0;
15206ad1575SMatthew G. Knepley     key.part  = 0;
1536528b96dSMatthew G. Knepley     if (!key.label) {
1549566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
1556528b96dSMatthew G. Knepley       cellIS = allcellIS;
1566528b96dSMatthew G. Knepley     } else {
1576528b96dSMatthew G. Knepley       IS pointIS;
1586528b96dSMatthew G. Knepley 
1596528b96dSMatthew G. Knepley       key.value = 1;
1609566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
1619566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
1629566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1634a3e9fdbSToby Isaac     }
164*2a8381b2SBarry Smith     PetscCall(DMPlexComputeResidualByKey(plex, key, cellIS, time, locX, locX_t, time, locF, ctx));
1659566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
1666528b96dSMatthew G. Knepley   }
1679566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
1689566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17024cdb843SMatthew G. Knepley }
1717cdb2a12SMatthew G. Knepley 
172756a1f44SMatthew G. Knepley /*@
17320f4b53cSBarry Smith   DMPlexTSComputeIJacobianFEM - Form the Jacobian `Jac` from the local input `locX` using pointwise functions specified by the user
174756a1f44SMatthew G. Knepley 
175756a1f44SMatthew G. Knepley   Input Parameters:
176756a1f44SMatthew G. Knepley + dm       - The mesh
177b43aa488SJacob Faibussowitsch . time     - The time
178756a1f44SMatthew G. Knepley . locX     - Local solution
17920f4b53cSBarry Smith . locX_t   - Local solution time derivative, or `NULL`
18020f4b53cSBarry Smith . X_tShift - The multiplicative parameter for dF/du_t
181*2a8381b2SBarry Smith - ctx      - The application context
182756a1f44SMatthew G. Knepley 
18320f4b53cSBarry Smith   Output Parameters:
18420f4b53cSBarry Smith + Jac  - the Jacobian
18520f4b53cSBarry Smith - JacP - an additional approximation for the Jacobian to be used to compute the preconditioner (often is `Jac`)
186756a1f44SMatthew G. Knepley 
187756a1f44SMatthew G. Knepley   Level: developer
188756a1f44SMatthew G. Knepley 
1891cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeRHSFunctionFEM()`
190756a1f44SMatthew G. Knepley @*/
DMPlexTSComputeIJacobianFEM(DM dm,PetscReal time,Vec locX,Vec locX_t,PetscReal X_tShift,Mat Jac,Mat JacP,PetscCtx ctx)191*2a8381b2SBarry Smith PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, PetscCtx ctx)
192d71ae5a4SJacob Faibussowitsch {
193756a1f44SMatthew G. Knepley   DM        plex;
1946528b96dSMatthew G. Knepley   IS        allcellIS;
195f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
1966528b96dSMatthew G. Knepley   PetscInt  Nds, s;
197756a1f44SMatthew G. Knepley 
198756a1f44SMatthew G. Knepley   PetscFunctionBegin;
1999566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
2009566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2019566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
2026528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
2036528b96dSMatthew G. Knepley     PetscDS      ds;
2046528b96dSMatthew G. Knepley     IS           cellIS;
20506ad1575SMatthew G. Knepley     PetscFormKey key;
2066528b96dSMatthew G. Knepley 
20707218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
2086528b96dSMatthew G. Knepley     key.value = 0;
2096528b96dSMatthew G. Knepley     key.field = 0;
21006ad1575SMatthew G. Knepley     key.part  = 0;
2116528b96dSMatthew G. Knepley     if (!key.label) {
2129566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
2136528b96dSMatthew G. Knepley       cellIS = allcellIS;
2146528b96dSMatthew G. Knepley     } else {
2156528b96dSMatthew G. Knepley       IS pointIS;
2166528b96dSMatthew G. Knepley 
2176528b96dSMatthew G. Knepley       key.value = 1;
2189566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
2199566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
2209566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
2216528b96dSMatthew G. Knepley     }
2226528b96dSMatthew G. Knepley     if (!s) {
2239566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
2249566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
2259566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
2269566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
2276528b96dSMatthew G. Knepley     }
228*2a8381b2SBarry Smith     PetscCall(DMPlexComputeJacobianByKey(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, ctx));
2299566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
2306528b96dSMatthew G. Knepley   }
2319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2329566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234756a1f44SMatthew G. Knepley }
235756a1f44SMatthew G. Knepley 
236cb36c0f9SMatthew G. Knepley /*@
23720f4b53cSBarry Smith   DMPlexTSComputeRHSFunctionFEM - Form the local residual `locG` from the local input `locX` using pointwise functions specified by the user
238cb36c0f9SMatthew G. Knepley 
239cb36c0f9SMatthew G. Knepley   Input Parameters:
240cb36c0f9SMatthew G. Knepley + dm   - The mesh
241b43aa488SJacob Faibussowitsch . time - The time
242cb36c0f9SMatthew G. Knepley . locX - Local solution
243*2a8381b2SBarry Smith - ctx  - The application context
244cb36c0f9SMatthew G. Knepley 
245cb36c0f9SMatthew G. Knepley   Output Parameter:
246cb36c0f9SMatthew G. Knepley . locG - Local output vector
247cb36c0f9SMatthew G. Knepley 
248cb36c0f9SMatthew G. Knepley   Level: developer
249cb36c0f9SMatthew G. Knepley 
2501cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeIJacobianFEM()`
251cb36c0f9SMatthew G. Knepley @*/
DMPlexTSComputeRHSFunctionFEM(DM dm,PetscReal time,Vec locX,Vec locG,PetscCtx ctx)252*2a8381b2SBarry Smith PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locG, PetscCtx ctx)
253d71ae5a4SJacob Faibussowitsch {
254cb36c0f9SMatthew G. Knepley   DM       plex;
255cb36c0f9SMatthew G. Knepley   IS       allcellIS;
256cb36c0f9SMatthew G. Knepley   PetscInt Nds, s;
257cb36c0f9SMatthew G. Knepley 
258cb36c0f9SMatthew G. Knepley   PetscFunctionBegin;
2599566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
2609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2619566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
262cb36c0f9SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
263cb36c0f9SMatthew G. Knepley     PetscDS      ds;
264cb36c0f9SMatthew G. Knepley     IS           cellIS;
265cb36c0f9SMatthew G. Knepley     PetscFormKey key;
266cb36c0f9SMatthew G. Knepley 
26707218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
268cb36c0f9SMatthew G. Knepley     key.value = 0;
269cb36c0f9SMatthew G. Knepley     key.field = 0;
270cb36c0f9SMatthew G. Knepley     key.part  = 100;
271cb36c0f9SMatthew G. Knepley     if (!key.label) {
2729566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
273cb36c0f9SMatthew G. Knepley       cellIS = allcellIS;
274cb36c0f9SMatthew G. Knepley     } else {
275cb36c0f9SMatthew G. Knepley       IS pointIS;
276cb36c0f9SMatthew G. Knepley 
277cb36c0f9SMatthew G. Knepley       key.value = 1;
2789566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
2799566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
2809566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
281cb36c0f9SMatthew G. Knepley     }
282*2a8381b2SBarry Smith     PetscCall(DMPlexComputeResidualByKey(plex, key, cellIS, time, locX, NULL, time, locG, ctx));
2839566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
284cb36c0f9SMatthew G. Knepley   }
2859566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
288cb36c0f9SMatthew G. Knepley }
289cb36c0f9SMatthew G. Knepley 
290cc4c1da9SBarry Smith /*@
291f2cacb80SMatthew G. Knepley   DMTSCheckResidual - Check the residual of the exact solution
292f2cacb80SMatthew G. Knepley 
293f2cacb80SMatthew G. Knepley   Input Parameters:
294bcf0153eSBarry Smith + ts  - the `TS` object
295bcf0153eSBarry Smith . dm  - the `DM`
296f2cacb80SMatthew G. Knepley . t   - the time
297bcf0153eSBarry Smith . u   - a `DM` vector
298bcf0153eSBarry Smith . u_t - a `DM` vector
299f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
300f2cacb80SMatthew G. Knepley 
3012fe279fdSBarry Smith   Output Parameter:
30220f4b53cSBarry Smith . residual - The residual norm of the exact solution, or `NULL`
303f2cacb80SMatthew G. Knepley 
304f2cacb80SMatthew G. Knepley   Level: developer
305f2cacb80SMatthew G. Knepley 
3061cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCheckFromOptions()`, `DMTSCheckJacobian()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
307f2cacb80SMatthew G. Knepley @*/
DMTSCheckResidual(TS ts,DM dm,PetscReal t,Vec u,Vec u_t,PetscReal tol,PetscReal * residual)308d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual)
309d71ae5a4SJacob Faibussowitsch {
310f2cacb80SMatthew G. Knepley   MPI_Comm  comm;
311f2cacb80SMatthew G. Knepley   Vec       r;
312f2cacb80SMatthew G. Knepley   PetscReal res;
313f2cacb80SMatthew G. Knepley 
314f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
315f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
316f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
317064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
3184f572ea9SToby Isaac   if (residual) PetscAssertPointer(residual, 7);
3199566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
3209566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
3219566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
3229566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
3239566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
324f2cacb80SMatthew G. Knepley   if (tol >= 0.0) {
3253c633725SBarry Smith     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
326f2cacb80SMatthew G. Knepley   } else if (residual) {
327f2cacb80SMatthew G. Knepley     *residual = res;
328f2cacb80SMatthew G. Knepley   } else {
3299566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
33093ec0da9SPierre Jolivet     PetscCall(VecFilter(r, 1.0e-10));
3319566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)dm));
3329566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
3339566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
3349566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
3359566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL));
336f2cacb80SMatthew G. Knepley   }
3379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
339f2cacb80SMatthew G. Knepley }
340f2cacb80SMatthew G. Knepley 
341cc4c1da9SBarry Smith /*@
342f2cacb80SMatthew G. Knepley   DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
343f2cacb80SMatthew G. Knepley 
344f2cacb80SMatthew G. Knepley   Input Parameters:
34520f4b53cSBarry Smith + ts  - the `TS` object
34620f4b53cSBarry Smith . dm  - the `DM`
347f2cacb80SMatthew G. Knepley . t   - the time
34820f4b53cSBarry Smith . u   - a `DM` vector
34920f4b53cSBarry Smith . u_t - a `DM` vector
350f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
351f2cacb80SMatthew G. Knepley 
352f2cacb80SMatthew G. Knepley   Output Parameters:
35343e72464SBrad Aagaard + isLinear - Flag indicating that the function looks linear, or `NULL`
35420f4b53cSBarry Smith - convRate - The rate of convergence of the linear model, or `NULL`
355f2cacb80SMatthew G. Knepley 
356f2cacb80SMatthew G. Knepley   Level: developer
357f2cacb80SMatthew G. Knepley 
3581cc06b55SBarry Smith .seealso: [](ch_ts), `DNTSCheckFromOptions()`, `DMTSCheckResidual()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
359f2cacb80SMatthew G. Knepley @*/
DMTSCheckJacobian(TS ts,DM dm,PetscReal t,Vec u,Vec u_t,PetscReal tol,PetscBool * isLinear,PetscReal * convRate)360d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
361d71ae5a4SJacob Faibussowitsch {
362f2cacb80SMatthew G. Knepley   MPI_Comm     comm;
363f2cacb80SMatthew G. Knepley   PetscDS      ds;
364f2cacb80SMatthew G. Knepley   Mat          J, M;
365f2cacb80SMatthew G. Knepley   MatNullSpace nullspace;
366f2cacb80SMatthew G. Knepley   PetscReal    dt, shift, slope, intercept;
367f2cacb80SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
368f2cacb80SMatthew G. Knepley 
369f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
370f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
371f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
372064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
3734f572ea9SToby Isaac   if (isLinear) PetscAssertPointer(isLinear, 7);
3744f572ea9SToby Isaac   if (convRate) PetscAssertPointer(convRate, 8);
3759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
3769566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
377f2cacb80SMatthew G. Knepley   /* Create and view matrices */
3789566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
379f2cacb80SMatthew G. Knepley   shift = 1.0 / dt;
3809566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
3819566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
3829566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
3839566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
384f2cacb80SMatthew G. Knepley   if (hasJac && hasPrec) {
3859566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
3869566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE));
3877addb90fSBarry Smith     PetscCall(PetscObjectSetName((PetscObject)M, "Matrix used to construct the preconditioner"));
3889566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
3899566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
3909566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
391f2cacb80SMatthew G. Knepley   } else {
3929566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE));
393f2cacb80SMatthew G. Knepley   }
3949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
3959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
3969566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
397f2cacb80SMatthew G. Knepley   /* Check nullspace */
3989566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
399f2cacb80SMatthew G. Knepley   if (nullspace) {
400f2cacb80SMatthew G. Knepley     PetscBool isNull;
4019566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
4023c633725SBarry Smith     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
403f2cacb80SMatthew G. Knepley   }
404f2cacb80SMatthew G. Knepley   /* Taylor test */
405f2cacb80SMatthew G. Knepley   {
406f2cacb80SMatthew G. Knepley     PetscRandom rand;
407f2cacb80SMatthew G. Knepley     Vec         du, uhat, uhat_t, r, rhat, df;
408f2cacb80SMatthew G. Knepley     PetscReal   h;
409f2cacb80SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
410f2cacb80SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
411f2cacb80SMatthew G. Knepley     PetscInt    Nv, v;
412f2cacb80SMatthew G. Knepley 
413f2cacb80SMatthew G. Knepley     /* Choose a perturbation direction */
4149566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
4159566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
4169566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
4179566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
4189566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
4199566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
420f2cacb80SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
4219566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
4229566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
423f2cacb80SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
424fbccb6d4SPierre Jolivet     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
4259566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
4269566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
4279566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat_t));
4289566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
429f2cacb80SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
4309566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
4319566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat_t, h * shift, du, u_t));
432f2cacb80SMatthew 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 */
4339566063dSJacob Faibussowitsch       PetscCall(TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE));
4349566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
4359566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
436f2cacb80SMatthew G. Knepley 
437f2cacb80SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
438f2cacb80SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
439f2cacb80SMatthew G. Knepley     }
4409566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
4419566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat_t));
4429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
4439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
4449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
4459566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
446f2cacb80SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
44743e72464SBrad Aagaard       if (tol >= 0) {
44843e72464SBrad Aagaard         if (errors[v] > tol) break;
44943e72464SBrad Aagaard       } else if (errors[v] > PETSC_SMALL) {
45043e72464SBrad Aagaard         break;
45143e72464SBrad Aagaard       }
452f2cacb80SMatthew G. Knepley     }
453f2cacb80SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
4549566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
4559566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
45643e72464SBrad Aagaard     if (isLinear || convRate) {
457f2cacb80SMatthew G. Knepley       if (isLinear) *isLinear = isLin;
458f2cacb80SMatthew G. Knepley       if (convRate) *convRate = slope;
45943e72464SBrad Aagaard     } else if (tol >= 0) {
46043e72464SBrad Aagaard       /* Slope should be about 2 */
46143e72464SBrad Aagaard       PetscCheck(PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
462f2cacb80SMatthew G. Knepley     } else {
4639566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
4649566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
465f2cacb80SMatthew G. Knepley     }
466f2cacb80SMatthew G. Knepley   }
4679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
4683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
469f2cacb80SMatthew G. Knepley }
470f2cacb80SMatthew G. Knepley 
471f2cacb80SMatthew G. Knepley /*@C
47220f4b53cSBarry Smith   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information based on
47320f4b53cSBarry Smith   values in the options database
474bee9a294SMatthew G. Knepley 
475bee9a294SMatthew G. Knepley   Input Parameters:
476bcf0153eSBarry Smith + ts - the `TS` object
477bcf0153eSBarry Smith - u  - representative `TS` vector
4787f96f943SMatthew G. Knepley 
47920f4b53cSBarry Smith   Level: developer
48020f4b53cSBarry Smith 
481bcf0153eSBarry Smith   Note:
482bcf0153eSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
483bee9a294SMatthew G. Knepley 
484b43aa488SJacob Faibussowitsch   Developer Notes:
48520f4b53cSBarry Smith   What is the purpose of `u`, does it need to already have a solution or some other value in it?
48620f4b53cSBarry Smith 
48720f4b53cSBarry Smith .seealso: `DMTS`
488bee9a294SMatthew G. Knepley @*/
DMTSCheckFromOptions(TS ts,Vec u)489d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u)
490d71ae5a4SJacob Faibussowitsch {
4917cdb2a12SMatthew G. Knepley   DM        dm;
4927cdb2a12SMatthew G. Knepley   SNES      snes;
493f2cacb80SMatthew G. Knepley   Vec       sol, u_t;
494f2cacb80SMatthew G. Knepley   PetscReal t;
4957cdb2a12SMatthew G. Knepley   PetscBool check;
4967cdb2a12SMatthew G. Knepley 
4977cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
4989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-dmts_check", &check));
4993ba16761SJacob Faibussowitsch   if (!check) PetscFunctionReturn(PETSC_SUCCESS);
5009566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
5019566063dSJacob Faibussowitsch   PetscCall(VecCopy(u, sol));
5029566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, u));
5039566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
5049566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
5059566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
5069566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, u));
507f2cacb80SMatthew G. Knepley 
5089566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
5099566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL));
5109566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u_t));
5119566063dSJacob Faibussowitsch   PetscCall(DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL));
5129566063dSJacob Faibussowitsch   PetscCall(DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL));
5139566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u_t));
514f2cacb80SMatthew G. Knepley 
5159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5177cdb2a12SMatthew G. Knepley }
518