xref: /petsc/src/ts/utils/dmplexts.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
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
3808449791SMatthew G. Knepley . t - The time
3908449791SMatthew G. Knepley . locX  - Local solution
4008449791SMatthew G. Knepley - user - The user 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 
47*1cc06b55SBarry Smith .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexComputeJacobianActionFEM()`
4808449791SMatthew G. Knepley @*/
49d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
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));
649566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locF, user));
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));
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71254c1ad2SMatthew G. Knepley }
72254c1ad2SMatthew G. Knepley 
73c5d70e09SMatthew G. Knepley /*@
7420f4b53cSBarry 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
75c5d70e09SMatthew G. Knepley 
76c5d70e09SMatthew G. Knepley   Input Parameters:
77c5d70e09SMatthew G. Knepley + dm - The mesh
78a40652d4SToby Isaac . t - The time
79a40652d4SToby Isaac . locX  - Local solution
8020f4b53cSBarry Smith . locX_t - Local solution time derivative, or `NULL`
81c5d70e09SMatthew G. Knepley - user - The user context
82c5d70e09SMatthew G. Knepley 
83c5d70e09SMatthew G. Knepley   Level: developer
84c5d70e09SMatthew G. Knepley 
85*1cc06b55SBarry Smith .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexComputeJacobianActionFEM()`
86c5d70e09SMatthew G. Knepley @*/
87d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
88d71ae5a4SJacob Faibussowitsch {
89c5d70e09SMatthew G. Knepley   DM       plex;
90ef68eab9SMatthew G. Knepley   Vec      faceGeometryFVM = NULL;
91ef68eab9SMatthew G. Knepley   PetscInt Nf, f;
92c5d70e09SMatthew G. Knepley 
93c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
959566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(plex, &Nf));
969586001cSMatthew G. Knepley   if (!locX_t) {
979586001cSMatthew G. Knepley     /* This is the RHS part */
98ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
99ef68eab9SMatthew G. Knepley       PetscObject  obj;
100ef68eab9SMatthew G. Knepley       PetscClassId id;
101ef68eab9SMatthew G. Knepley 
1029566063dSJacob Faibussowitsch       PetscCall(DMGetField(plex, f, NULL, &obj));
1039566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
104ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
1059566063dSJacob Faibussowitsch         PetscCall(DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL));
106ef68eab9SMatthew G. Knepley         break;
107ef68eab9SMatthew G. Knepley       }
108ef68eab9SMatthew G. Knepley     }
1099586001cSMatthew G. Knepley   }
1109566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL));
1119566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL));
1129566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114c5d70e09SMatthew G. Knepley }
115c5d70e09SMatthew G. Knepley 
11624cdb843SMatthew G. Knepley /*@
11720f4b53cSBarry Smith   DMPlexTSComputeIFunctionFEM - Form the local residual `locF` from the local input `locX` using pointwise functions specified by the user
11824cdb843SMatthew G. Knepley 
11924cdb843SMatthew G. Knepley   Input Parameters:
12024cdb843SMatthew G. Knepley + dm - The mesh
12124cdb843SMatthew G. Knepley . t - The time
12208449791SMatthew G. Knepley . locX  - Local solution
12320f4b53cSBarry Smith . locX_t - Local solution time derivative, or `NULL`
12424cdb843SMatthew G. Knepley - user - The user context
12524cdb843SMatthew G. Knepley 
12624cdb843SMatthew G. Knepley   Output Parameter:
12708449791SMatthew G. Knepley . locF  - Local output vector
12824cdb843SMatthew G. Knepley 
12924cdb843SMatthew G. Knepley   Level: developer
13024cdb843SMatthew G. Knepley 
131*1cc06b55SBarry Smith .seealso: [](ch_ts), `DMPLEX`, `TS`, `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeRHSFunctionFEM()`
13224cdb843SMatthew G. Knepley @*/
133d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
134d71ae5a4SJacob Faibussowitsch {
1356da023fcSToby Isaac   DM       plex;
1366528b96dSMatthew G. Knepley   IS       allcellIS;
1376528b96dSMatthew G. Knepley   PetscInt Nds, s;
13824cdb843SMatthew G. Knepley 
13924cdb843SMatthew G. Knepley   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
1419566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
1429566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1436528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1446528b96dSMatthew G. Knepley     PetscDS      ds;
1456528b96dSMatthew G. Knepley     IS           cellIS;
14606ad1575SMatthew G. Knepley     PetscFormKey key;
1476528b96dSMatthew G. Knepley 
14807218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
1496528b96dSMatthew G. Knepley     key.value = 0;
1506528b96dSMatthew G. Knepley     key.field = 0;
15106ad1575SMatthew G. Knepley     key.part  = 0;
1526528b96dSMatthew G. Knepley     if (!key.label) {
1539566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
1546528b96dSMatthew G. Knepley       cellIS = allcellIS;
1556528b96dSMatthew G. Knepley     } else {
1566528b96dSMatthew G. Knepley       IS pointIS;
1576528b96dSMatthew G. Knepley 
1586528b96dSMatthew G. Knepley       key.value = 1;
1599566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
1609566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
1619566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1624a3e9fdbSToby Isaac     }
1639566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, locX_t, time, locF, user));
1649566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
1656528b96dSMatthew G. Knepley   }
1669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
1679566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16924cdb843SMatthew G. Knepley }
1707cdb2a12SMatthew G. Knepley 
171756a1f44SMatthew G. Knepley /*@
17220f4b53cSBarry Smith   DMPlexTSComputeIJacobianFEM - Form the Jacobian `Jac` from the local input `locX` using pointwise functions specified by the user
173756a1f44SMatthew G. Knepley 
174756a1f44SMatthew G. Knepley   Input Parameters:
175756a1f44SMatthew G. Knepley + dm - The mesh
176756a1f44SMatthew G. Knepley . t - The time
177756a1f44SMatthew G. Knepley . locX  - Local solution
17820f4b53cSBarry Smith . locX_t - Local solution time derivative, or `NULL`
17920f4b53cSBarry Smith . X_tShift - The multiplicative parameter for dF/du_t
180756a1f44SMatthew G. Knepley - user - The user context
181756a1f44SMatthew G. Knepley 
18220f4b53cSBarry Smith   Output Parameters:
18320f4b53cSBarry Smith + Jac - the Jacobian
18420f4b53cSBarry Smith - JacP - an additional approximation for the Jacobian to be used to compute the preconditioner (often is `Jac`)
185756a1f44SMatthew G. Knepley 
186756a1f44SMatthew G. Knepley   Level: developer
187756a1f44SMatthew G. Knepley 
188*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeRHSFunctionFEM()`
189756a1f44SMatthew G. Knepley @*/
190d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
191d71ae5a4SJacob Faibussowitsch {
192756a1f44SMatthew G. Knepley   DM        plex;
1936528b96dSMatthew G. Knepley   IS        allcellIS;
194f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
1956528b96dSMatthew G. Knepley   PetscInt  Nds, s;
196756a1f44SMatthew G. Knepley 
197756a1f44SMatthew G. Knepley   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
1999566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2009566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
2016528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
2026528b96dSMatthew G. Knepley     PetscDS      ds;
2036528b96dSMatthew G. Knepley     IS           cellIS;
20406ad1575SMatthew G. Knepley     PetscFormKey key;
2056528b96dSMatthew G. Knepley 
20607218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
2076528b96dSMatthew G. Knepley     key.value = 0;
2086528b96dSMatthew G. Knepley     key.field = 0;
20906ad1575SMatthew G. Knepley     key.part  = 0;
2106528b96dSMatthew G. Knepley     if (!key.label) {
2119566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
2126528b96dSMatthew G. Knepley       cellIS = allcellIS;
2136528b96dSMatthew G. Knepley     } else {
2146528b96dSMatthew G. Knepley       IS pointIS;
2156528b96dSMatthew G. Knepley 
2166528b96dSMatthew G. Knepley       key.value = 1;
2179566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
2189566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
2199566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
2206528b96dSMatthew G. Knepley     }
2216528b96dSMatthew G. Knepley     if (!s) {
2229566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
2239566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
2249566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
2259566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
2266528b96dSMatthew G. Knepley     }
2279566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user));
2289566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
2296528b96dSMatthew G. Knepley   }
2309566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2319566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
233756a1f44SMatthew G. Knepley }
234756a1f44SMatthew G. Knepley 
235cb36c0f9SMatthew G. Knepley /*@
23620f4b53cSBarry Smith   DMPlexTSComputeRHSFunctionFEM - Form the local residual `locG` from the local input `locX` using pointwise functions specified by the user
237cb36c0f9SMatthew G. Knepley 
238cb36c0f9SMatthew G. Knepley   Input Parameters:
239cb36c0f9SMatthew G. Knepley + dm - The mesh
240cb36c0f9SMatthew G. Knepley . t - The time
241cb36c0f9SMatthew G. Knepley . locX  - Local solution
242cb36c0f9SMatthew G. Knepley - user - The user context
243cb36c0f9SMatthew G. Knepley 
244cb36c0f9SMatthew G. Knepley   Output Parameter:
245cb36c0f9SMatthew G. Knepley . locG  - Local output vector
246cb36c0f9SMatthew G. Knepley 
247cb36c0f9SMatthew G. Knepley   Level: developer
248cb36c0f9SMatthew G. Knepley 
249*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `DM`, `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeIJacobianFEM()`
250cb36c0f9SMatthew G. Knepley @*/
251d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locG, void *user)
252d71ae5a4SJacob Faibussowitsch {
253cb36c0f9SMatthew G. Knepley   DM       plex;
254cb36c0f9SMatthew G. Knepley   IS       allcellIS;
255cb36c0f9SMatthew G. Knepley   PetscInt Nds, s;
256cb36c0f9SMatthew G. Knepley 
257cb36c0f9SMatthew G. Knepley   PetscFunctionBegin;
2589566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
2599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2609566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
261cb36c0f9SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
262cb36c0f9SMatthew G. Knepley     PetscDS      ds;
263cb36c0f9SMatthew G. Knepley     IS           cellIS;
264cb36c0f9SMatthew G. Knepley     PetscFormKey key;
265cb36c0f9SMatthew G. Knepley 
26607218a29SMatthew G. Knepley     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL));
267cb36c0f9SMatthew G. Knepley     key.value = 0;
268cb36c0f9SMatthew G. Knepley     key.field = 0;
269cb36c0f9SMatthew G. Knepley     key.part  = 100;
270cb36c0f9SMatthew G. Knepley     if (!key.label) {
2719566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
272cb36c0f9SMatthew G. Knepley       cellIS = allcellIS;
273cb36c0f9SMatthew G. Knepley     } else {
274cb36c0f9SMatthew G. Knepley       IS pointIS;
275cb36c0f9SMatthew G. Knepley 
276cb36c0f9SMatthew G. Knepley       key.value = 1;
2779566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
2789566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
2799566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
280cb36c0f9SMatthew G. Knepley     }
2819566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locG, user));
2829566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
283cb36c0f9SMatthew G. Knepley   }
2849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
287cb36c0f9SMatthew G. Knepley }
288cb36c0f9SMatthew G. Knepley 
289bee9a294SMatthew G. Knepley /*@C
290f2cacb80SMatthew G. Knepley   DMTSCheckResidual - Check the residual of the exact solution
291f2cacb80SMatthew G. Knepley 
292f2cacb80SMatthew G. Knepley   Input Parameters:
293bcf0153eSBarry Smith + ts  - the `TS` object
294bcf0153eSBarry Smith . dm  - the `DM`
295f2cacb80SMatthew G. Knepley . t   - the time
296bcf0153eSBarry Smith . u   - a `DM` vector
297bcf0153eSBarry Smith . u_t - a `DM` vector
298f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
299f2cacb80SMatthew G. Knepley 
3002fe279fdSBarry Smith   Output Parameter:
30120f4b53cSBarry Smith . residual - The residual norm of the exact solution, or `NULL`
302f2cacb80SMatthew G. Knepley 
303f2cacb80SMatthew G. Knepley   Level: developer
304f2cacb80SMatthew G. Knepley 
305*1cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCheckFromOptions()`, `DMTSCheckJacobian()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
306f2cacb80SMatthew G. Knepley @*/
307d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual)
308d71ae5a4SJacob Faibussowitsch {
309f2cacb80SMatthew G. Knepley   MPI_Comm  comm;
310f2cacb80SMatthew G. Knepley   Vec       r;
311f2cacb80SMatthew G. Knepley   PetscReal res;
312f2cacb80SMatthew G. Knepley 
313f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
314f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
315f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
316064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
317064a246eSJacob Faibussowitsch   if (residual) PetscValidRealPointer(residual, 7);
3189566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
3199566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
3209566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
3219566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
3229566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
323f2cacb80SMatthew G. Knepley   if (tol >= 0.0) {
3243c633725SBarry Smith     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
325f2cacb80SMatthew G. Knepley   } else if (residual) {
326f2cacb80SMatthew G. Knepley     *residual = res;
327f2cacb80SMatthew G. Knepley   } else {
3289566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
3299566063dSJacob Faibussowitsch     PetscCall(VecChop(r, 1.0e-10));
3309566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)dm));
3319566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
3329566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
3339566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
3349566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL));
335f2cacb80SMatthew G. Knepley   }
3369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
338f2cacb80SMatthew G. Knepley }
339f2cacb80SMatthew G. Knepley 
340f2cacb80SMatthew G. Knepley /*@C
341f2cacb80SMatthew G. Knepley   DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
342f2cacb80SMatthew G. Knepley 
343f2cacb80SMatthew G. Knepley   Input Parameters:
34420f4b53cSBarry Smith + ts  - the `TS` object
34520f4b53cSBarry Smith . dm  - the `DM`
346f2cacb80SMatthew G. Knepley . t   - the time
34720f4b53cSBarry Smith . u   - a `DM` vector
34820f4b53cSBarry Smith . u_t - a `DM` vector
349f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
350f2cacb80SMatthew G. Knepley 
351f2cacb80SMatthew G. Knepley   Output Parameters:
35220f4b53cSBarry Smith + isLinear - Flag indicaing that the function looks linear, or `NULL`
35320f4b53cSBarry Smith - convRate - The rate of convergence of the linear model, or `NULL`
354f2cacb80SMatthew G. Knepley 
355f2cacb80SMatthew G. Knepley   Level: developer
356f2cacb80SMatthew G. Knepley 
357*1cc06b55SBarry Smith .seealso: [](ch_ts), `DNTSCheckFromOptions()`, `DMTSCheckResidual()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
358f2cacb80SMatthew G. Knepley @*/
359d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
360d71ae5a4SJacob Faibussowitsch {
361f2cacb80SMatthew G. Knepley   MPI_Comm     comm;
362f2cacb80SMatthew G. Knepley   PetscDS      ds;
363f2cacb80SMatthew G. Knepley   Mat          J, M;
364f2cacb80SMatthew G. Knepley   MatNullSpace nullspace;
365f2cacb80SMatthew G. Knepley   PetscReal    dt, shift, slope, intercept;
366f2cacb80SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
367f2cacb80SMatthew G. Knepley 
368f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
369f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
370f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
371064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
372064a246eSJacob Faibussowitsch   if (isLinear) PetscValidBoolPointer(isLinear, 7);
373064a246eSJacob Faibussowitsch   if (convRate) PetscValidRealPointer(convRate, 8);
3749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
3759566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
376f2cacb80SMatthew G. Knepley   /* Create and view matrices */
3779566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
378f2cacb80SMatthew G. Knepley   shift = 1.0 / dt;
3799566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
3809566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
3819566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
3829566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
383f2cacb80SMatthew G. Knepley   if (hasJac && hasPrec) {
3849566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
3859566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE));
3869566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix"));
3879566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
3889566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
3899566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
390f2cacb80SMatthew G. Knepley   } else {
3919566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE));
392f2cacb80SMatthew G. Knepley   }
3939566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
3949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
3959566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
396f2cacb80SMatthew G. Knepley   /* Check nullspace */
3979566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
398f2cacb80SMatthew G. Knepley   if (nullspace) {
399f2cacb80SMatthew G. Knepley     PetscBool isNull;
4009566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
4013c633725SBarry Smith     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
402f2cacb80SMatthew G. Knepley   }
403f2cacb80SMatthew G. Knepley   /* Taylor test */
404f2cacb80SMatthew G. Knepley   {
405f2cacb80SMatthew G. Knepley     PetscRandom rand;
406f2cacb80SMatthew G. Knepley     Vec         du, uhat, uhat_t, r, rhat, df;
407f2cacb80SMatthew G. Knepley     PetscReal   h;
408f2cacb80SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
409f2cacb80SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
410f2cacb80SMatthew G. Knepley     PetscInt    Nv, v;
411f2cacb80SMatthew G. Knepley 
412f2cacb80SMatthew G. Knepley     /* Choose a perturbation direction */
4139566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
4149566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
4159566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
4169566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
4179566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
4189566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
419f2cacb80SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
4209566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
4219566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
422f2cacb80SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
4239371c9d4SSatish Balay     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv)
4249371c9d4SSatish Balay       ;
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) {
447f2cacb80SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
448f2cacb80SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL) break;
449f2cacb80SMatthew G. Knepley     }
450f2cacb80SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
4519566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
4529566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
453f2cacb80SMatthew G. Knepley     /* Slope should be about 2 */
454f2cacb80SMatthew G. Knepley     if (tol >= 0) {
4553c633725SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
456f2cacb80SMatthew G. Knepley     } else if (isLinear || convRate) {
457f2cacb80SMatthew G. Knepley       if (isLinear) *isLinear = isLin;
458f2cacb80SMatthew G. Knepley       if (convRate) *convRate = slope;
459f2cacb80SMatthew G. Knepley     } else {
4609566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
4619566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
462f2cacb80SMatthew G. Knepley     }
463f2cacb80SMatthew G. Knepley   }
4649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
466f2cacb80SMatthew G. Knepley }
467f2cacb80SMatthew G. Knepley 
468f2cacb80SMatthew G. Knepley /*@C
46920f4b53cSBarry Smith   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information based on
47020f4b53cSBarry Smith   values in the options database
471bee9a294SMatthew G. Knepley 
472bee9a294SMatthew G. Knepley   Input Parameters:
473bcf0153eSBarry Smith + ts - the `TS` object
474bcf0153eSBarry Smith - u  - representative `TS` vector
4757f96f943SMatthew G. Knepley 
47620f4b53cSBarry Smith   Level: developer
47720f4b53cSBarry Smith 
478bcf0153eSBarry Smith   Note:
479bcf0153eSBarry Smith   The user must call `PetscDSSetExactSolution()` beforehand
480bee9a294SMatthew G. Knepley 
48120f4b53cSBarry Smith   Developer Note:
48220f4b53cSBarry Smith   What is the purpose of `u`, does it need to already have a solution or some other value in it?
48320f4b53cSBarry Smith 
48420f4b53cSBarry Smith .seealso: `DMTS`
485bee9a294SMatthew G. Knepley @*/
486d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u)
487d71ae5a4SJacob Faibussowitsch {
4887cdb2a12SMatthew G. Knepley   DM        dm;
4897cdb2a12SMatthew G. Knepley   SNES      snes;
490f2cacb80SMatthew G. Knepley   Vec       sol, u_t;
491f2cacb80SMatthew G. Knepley   PetscReal t;
4927cdb2a12SMatthew G. Knepley   PetscBool check;
4937cdb2a12SMatthew G. Knepley 
4947cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
4959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-dmts_check", &check));
4963ba16761SJacob Faibussowitsch   if (!check) PetscFunctionReturn(PETSC_SUCCESS);
4979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
4989566063dSJacob Faibussowitsch   PetscCall(VecCopy(u, sol));
4999566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, u));
5009566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
5019566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
5029566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
5039566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, u));
504f2cacb80SMatthew G. Knepley 
5059566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
5069566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL));
5079566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u_t));
5089566063dSJacob Faibussowitsch   PetscCall(DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL));
5099566063dSJacob Faibussowitsch   PetscCall(DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL));
5109566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u_t));
511f2cacb80SMatthew G. Knepley 
5129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5147cdb2a12SMatthew G. Knepley }
515