xref: /petsc/src/ts/utils/dmplexts.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
7*9371c9d4SSatish Balay static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy) {
86da023fcSToby Isaac   PetscBool isPlex;
96da023fcSToby Isaac 
106da023fcSToby Isaac   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex));
126da023fcSToby Isaac   if (isPlex) {
136da023fcSToby Isaac     *plex = dm;
149566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
15f7148743SMatthew G. Knepley   } else {
169566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex));
17f7148743SMatthew G. Knepley     if (!*plex) {
189566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm, DMPLEX, plex));
199566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex));
206da023fcSToby Isaac       if (copy) {
219566063dSJacob Faibussowitsch         PetscCall(DMCopyDMTS(dm, *plex));
229566063dSJacob Faibussowitsch         PetscCall(DMCopyDMSNES(dm, *plex));
239566063dSJacob Faibussowitsch         PetscCall(DMCopyAuxiliaryVec(dm, *plex));
246da023fcSToby Isaac       }
25f7148743SMatthew G. Knepley     } else {
269566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)*plex));
27f7148743SMatthew G. Knepley     }
286da023fcSToby Isaac   }
296da023fcSToby Isaac   PetscFunctionReturn(0);
306da023fcSToby Isaac }
316da023fcSToby Isaac 
3208449791SMatthew G. Knepley /*@
3308449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
3408449791SMatthew G. Knepley 
3508449791SMatthew G. Knepley   Input Parameters:
3608449791SMatthew G. Knepley + dm - The mesh
3708449791SMatthew G. Knepley . t - The time
3808449791SMatthew G. Knepley . locX  - Local solution
3908449791SMatthew G. Knepley - user - The user context
4008449791SMatthew G. Knepley 
4108449791SMatthew G. Knepley   Output Parameter:
423b16df42SMatthew G. Knepley . F  - Global output vector
4308449791SMatthew G. Knepley 
4408449791SMatthew G. Knepley   Level: developer
4508449791SMatthew G. Knepley 
46db781477SPatrick Sanan .seealso: `DMPlexComputeJacobianActionFEM()`
4708449791SMatthew G. Knepley @*/
48*9371c9d4SSatish Balay PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user) {
493b16df42SMatthew G. Knepley   Vec          locF;
504a3e9fdbSToby Isaac   IS           cellIS;
516da023fcSToby Isaac   DM           plex;
524a3e9fdbSToby Isaac   PetscInt     depth;
53d70f29a3SPierre Jolivet   PetscFormKey key = {NULL, 0, 0, 0};
54254c1ad2SMatthew G. Knepley 
55254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
569566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
579566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
589566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(plex, "dim", depth, &cellIS));
599566063dSJacob Faibussowitsch   if (!cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, &cellIS));
609566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(plex, &locF));
619566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(locF));
629566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locF, user));
639566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F));
649566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F));
659566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(plex, &locF));
669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellIS));
679566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
68254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
69254c1ad2SMatthew G. Knepley }
70254c1ad2SMatthew G. Knepley 
71c5d70e09SMatthew G. Knepley /*@
72ef68eab9SMatthew G. Knepley   DMPlexTSComputeBoundary - Insert the essential boundary values for the local input X and/or its time derivative X_t using pointwise functions specified by the user
73c5d70e09SMatthew G. Knepley 
74c5d70e09SMatthew G. Knepley   Input Parameters:
75c5d70e09SMatthew G. Knepley + dm - The mesh
76a40652d4SToby Isaac . t - The time
77a40652d4SToby Isaac . locX  - Local solution
78a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL
79c5d70e09SMatthew G. Knepley - user - The user context
80c5d70e09SMatthew G. Knepley 
81c5d70e09SMatthew G. Knepley   Level: developer
82c5d70e09SMatthew G. Knepley 
83db781477SPatrick Sanan .seealso: `DMPlexComputeJacobianActionFEM()`
84c5d70e09SMatthew G. Knepley @*/
85*9371c9d4SSatish Balay PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user) {
86c5d70e09SMatthew G. Knepley   DM       plex;
87ef68eab9SMatthew G. Knepley   Vec      faceGeometryFVM = NULL;
88ef68eab9SMatthew G. Knepley   PetscInt Nf, f;
89c5d70e09SMatthew G. Knepley 
90c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
919566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
929566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(plex, &Nf));
939586001cSMatthew G. Knepley   if (!locX_t) {
949586001cSMatthew G. Knepley     /* This is the RHS part */
95ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
96ef68eab9SMatthew G. Knepley       PetscObject  obj;
97ef68eab9SMatthew G. Knepley       PetscClassId id;
98ef68eab9SMatthew G. Knepley 
999566063dSJacob Faibussowitsch       PetscCall(DMGetField(plex, f, NULL, &obj));
1009566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
101ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
1029566063dSJacob Faibussowitsch         PetscCall(DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL));
103ef68eab9SMatthew G. Knepley         break;
104ef68eab9SMatthew G. Knepley       }
105ef68eab9SMatthew G. Knepley     }
1069586001cSMatthew G. Knepley   }
1079566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL));
1089566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL));
1099566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
110c5d70e09SMatthew G. Knepley   PetscFunctionReturn(0);
111c5d70e09SMatthew G. Knepley }
112c5d70e09SMatthew G. Knepley 
11324cdb843SMatthew G. Knepley /*@
11424cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
11524cdb843SMatthew G. Knepley 
11624cdb843SMatthew G. Knepley   Input Parameters:
11724cdb843SMatthew G. Knepley + dm - The mesh
11824cdb843SMatthew G. Knepley . t - The time
11908449791SMatthew G. Knepley . locX  - Local solution
12008449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
12124cdb843SMatthew G. Knepley - user - The user context
12224cdb843SMatthew G. Knepley 
12324cdb843SMatthew G. Knepley   Output Parameter:
12408449791SMatthew G. Knepley . locF  - Local output vector
12524cdb843SMatthew G. Knepley 
12624cdb843SMatthew G. Knepley   Level: developer
12724cdb843SMatthew G. Knepley 
128db781477SPatrick Sanan .seealso: `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeRHSFunctionFEM()`
12924cdb843SMatthew G. Knepley @*/
130*9371c9d4SSatish Balay PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user) {
1316da023fcSToby Isaac   DM       plex;
1326528b96dSMatthew G. Knepley   IS       allcellIS;
1336528b96dSMatthew G. Knepley   PetscInt Nds, s;
13424cdb843SMatthew G. Knepley 
13524cdb843SMatthew G. Knepley   PetscFunctionBegin;
1369566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
1379566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
1389566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1396528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1406528b96dSMatthew G. Knepley     PetscDS      ds;
1416528b96dSMatthew G. Knepley     IS           cellIS;
14206ad1575SMatthew G. Knepley     PetscFormKey key;
1436528b96dSMatthew G. Knepley 
1449566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
1456528b96dSMatthew G. Knepley     key.value = 0;
1466528b96dSMatthew G. Knepley     key.field = 0;
14706ad1575SMatthew G. Knepley     key.part  = 0;
1486528b96dSMatthew G. Knepley     if (!key.label) {
1499566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
1506528b96dSMatthew G. Knepley       cellIS = allcellIS;
1516528b96dSMatthew G. Knepley     } else {
1526528b96dSMatthew G. Knepley       IS pointIS;
1536528b96dSMatthew G. Knepley 
1546528b96dSMatthew G. Knepley       key.value = 1;
1559566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
1569566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
1579566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1584a3e9fdbSToby Isaac     }
1599566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, locX_t, time, locF, user));
1609566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
1616528b96dSMatthew G. Knepley   }
1629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
1639566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
16424cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
16524cdb843SMatthew G. Knepley }
1667cdb2a12SMatthew G. Knepley 
167756a1f44SMatthew G. Knepley /*@
168756a1f44SMatthew G. Knepley   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
169756a1f44SMatthew G. Knepley 
170756a1f44SMatthew G. Knepley   Input Parameters:
171756a1f44SMatthew G. Knepley + dm - The mesh
172756a1f44SMatthew G. Knepley . t - The time
173756a1f44SMatthew G. Knepley . locX  - Local solution
174756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
175756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t
176756a1f44SMatthew G. Knepley - user - The user context
177756a1f44SMatthew G. Knepley 
178756a1f44SMatthew G. Knepley   Output Parameter:
179756a1f44SMatthew G. Knepley . locF  - Local output vector
180756a1f44SMatthew G. Knepley 
181756a1f44SMatthew G. Knepley   Level: developer
182756a1f44SMatthew G. Knepley 
183db781477SPatrick Sanan .seealso: `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeRHSFunctionFEM()`
184756a1f44SMatthew G. Knepley @*/
185*9371c9d4SSatish Balay PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user) {
186756a1f44SMatthew G. Knepley   DM        plex;
1876528b96dSMatthew G. Knepley   IS        allcellIS;
188f04eb4edSMatthew G. Knepley   PetscBool hasJac, hasPrec;
1896528b96dSMatthew G. Knepley   PetscInt  Nds, s;
190756a1f44SMatthew G. Knepley 
191756a1f44SMatthew G. Knepley   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
1939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
1949566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1956528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1966528b96dSMatthew G. Knepley     PetscDS      ds;
1976528b96dSMatthew G. Knepley     IS           cellIS;
19806ad1575SMatthew G. Knepley     PetscFormKey key;
1996528b96dSMatthew G. Knepley 
2009566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
2016528b96dSMatthew G. Knepley     key.value = 0;
2026528b96dSMatthew G. Knepley     key.field = 0;
20306ad1575SMatthew G. Knepley     key.part  = 0;
2046528b96dSMatthew G. Knepley     if (!key.label) {
2059566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
2066528b96dSMatthew G. Knepley       cellIS = allcellIS;
2076528b96dSMatthew G. Knepley     } else {
2086528b96dSMatthew G. Knepley       IS pointIS;
2096528b96dSMatthew G. Knepley 
2106528b96dSMatthew G. Knepley       key.value = 1;
2119566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
2129566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
2139566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
2146528b96dSMatthew G. Knepley     }
2156528b96dSMatthew G. Knepley     if (!s) {
2169566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
2179566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
2189566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
2199566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
2206528b96dSMatthew G. Knepley     }
2219566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user));
2229566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
2236528b96dSMatthew G. Knepley   }
2249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2259566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
226756a1f44SMatthew G. Knepley   PetscFunctionReturn(0);
227756a1f44SMatthew G. Knepley }
228756a1f44SMatthew G. Knepley 
229cb36c0f9SMatthew G. Knepley /*@
230d64986aaSMatthew G. Knepley   DMPlexTSComputeRHSFunctionFEM - Form the local residual G from the local input X using pointwise functions specified by the user
231cb36c0f9SMatthew G. Knepley 
232cb36c0f9SMatthew G. Knepley   Input Parameters:
233cb36c0f9SMatthew G. Knepley + dm - The mesh
234cb36c0f9SMatthew G. Knepley . t - The time
235cb36c0f9SMatthew G. Knepley . locX  - Local solution
236cb36c0f9SMatthew G. Knepley - user - The user context
237cb36c0f9SMatthew G. Knepley 
238cb36c0f9SMatthew G. Knepley   Output Parameter:
239cb36c0f9SMatthew G. Knepley . locG  - Local output vector
240cb36c0f9SMatthew G. Knepley 
241cb36c0f9SMatthew G. Knepley   Level: developer
242cb36c0f9SMatthew G. Knepley 
243db781477SPatrick Sanan .seealso: `DMPlexTSComputeIFunctionFEM()`, `DMPlexTSComputeIJacobianFEM()`
244cb36c0f9SMatthew G. Knepley @*/
245*9371c9d4SSatish Balay PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locG, void *user) {
246cb36c0f9SMatthew G. Knepley   DM       plex;
247cb36c0f9SMatthew G. Knepley   IS       allcellIS;
248cb36c0f9SMatthew G. Knepley   PetscInt Nds, s;
249cb36c0f9SMatthew G. Knepley 
250cb36c0f9SMatthew G. Knepley   PetscFunctionBegin;
2519566063dSJacob Faibussowitsch   PetscCall(DMTSConvertPlex(dm, &plex, PETSC_TRUE));
2529566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
2539566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
254cb36c0f9SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
255cb36c0f9SMatthew G. Knepley     PetscDS      ds;
256cb36c0f9SMatthew G. Knepley     IS           cellIS;
257cb36c0f9SMatthew G. Knepley     PetscFormKey key;
258cb36c0f9SMatthew G. Knepley 
2599566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
260cb36c0f9SMatthew G. Knepley     key.value = 0;
261cb36c0f9SMatthew G. Knepley     key.field = 0;
262cb36c0f9SMatthew G. Knepley     key.part  = 100;
263cb36c0f9SMatthew G. Knepley     if (!key.label) {
2649566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)allcellIS));
265cb36c0f9SMatthew G. Knepley       cellIS = allcellIS;
266cb36c0f9SMatthew G. Knepley     } else {
267cb36c0f9SMatthew G. Knepley       IS pointIS;
268cb36c0f9SMatthew G. Knepley 
269cb36c0f9SMatthew G. Knepley       key.value = 1;
2709566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
2719566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
2729566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
273cb36c0f9SMatthew G. Knepley     }
2749566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locG, user));
2759566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
276cb36c0f9SMatthew G. Knepley   }
2779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
2789566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
279cb36c0f9SMatthew G. Knepley   PetscFunctionReturn(0);
280cb36c0f9SMatthew G. Knepley }
281cb36c0f9SMatthew G. Knepley 
282bee9a294SMatthew G. Knepley /*@C
283f2cacb80SMatthew G. Knepley   DMTSCheckResidual - Check the residual of the exact solution
284f2cacb80SMatthew G. Knepley 
285f2cacb80SMatthew G. Knepley   Input Parameters:
286f2cacb80SMatthew G. Knepley + ts  - the TS object
287f2cacb80SMatthew G. Knepley . dm  - the DM
288f2cacb80SMatthew G. Knepley . t   - the time
289f2cacb80SMatthew G. Knepley . u   - a DM vector
290f2cacb80SMatthew G. Knepley . u_t - a DM vector
291f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
292f2cacb80SMatthew G. Knepley 
293f2cacb80SMatthew G. Knepley   Output Parameters:
294f2cacb80SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
295f2cacb80SMatthew G. Knepley 
296f2cacb80SMatthew G. Knepley   Level: developer
297f2cacb80SMatthew G. Knepley 
298db781477SPatrick Sanan .seealso: `DNTSCheckFromOptions()`, `DMTSCheckJacobian()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
299f2cacb80SMatthew G. Knepley @*/
300*9371c9d4SSatish Balay PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual) {
301f2cacb80SMatthew G. Knepley   MPI_Comm  comm;
302f2cacb80SMatthew G. Knepley   Vec       r;
303f2cacb80SMatthew G. Knepley   PetscReal res;
304f2cacb80SMatthew G. Knepley 
305f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
306f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
307f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
308064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
309064a246eSJacob Faibussowitsch   if (residual) PetscValidRealPointer(residual, 7);
3109566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
3119566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
3129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
3139566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
3149566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
315f2cacb80SMatthew G. Knepley   if (tol >= 0.0) {
3163c633725SBarry Smith     PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol);
317f2cacb80SMatthew G. Knepley   } else if (residual) {
318f2cacb80SMatthew G. Knepley     *residual = res;
319f2cacb80SMatthew G. Knepley   } else {
3209566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
3219566063dSJacob Faibussowitsch     PetscCall(VecChop(r, 1.0e-10));
3229566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)dm));
3239566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual"));
3249566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_"));
3259566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
3269566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL));
327f2cacb80SMatthew G. Knepley   }
3289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
329f2cacb80SMatthew G. Knepley   PetscFunctionReturn(0);
330f2cacb80SMatthew G. Knepley }
331f2cacb80SMatthew G. Knepley 
332f2cacb80SMatthew G. Knepley /*@C
333f2cacb80SMatthew G. Knepley   DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
334f2cacb80SMatthew G. Knepley 
335f2cacb80SMatthew G. Knepley   Input Parameters:
336f2cacb80SMatthew G. Knepley + ts  - the TS object
337f2cacb80SMatthew G. Knepley . dm  - the DM
338f2cacb80SMatthew G. Knepley . t   - the time
339f2cacb80SMatthew G. Knepley . u   - a DM vector
340f2cacb80SMatthew G. Knepley . u_t - a DM vector
341f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
342f2cacb80SMatthew G. Knepley 
343f2cacb80SMatthew G. Knepley   Output Parameters:
344f2cacb80SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
345f2cacb80SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
346f2cacb80SMatthew G. Knepley 
347f2cacb80SMatthew G. Knepley   Level: developer
348f2cacb80SMatthew G. Knepley 
349db781477SPatrick Sanan .seealso: `DNTSCheckFromOptions()`, `DMTSCheckResidual()`, `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
350f2cacb80SMatthew G. Knepley @*/
351*9371c9d4SSatish Balay PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) {
352f2cacb80SMatthew G. Knepley   MPI_Comm     comm;
353f2cacb80SMatthew G. Knepley   PetscDS      ds;
354f2cacb80SMatthew G. Knepley   Mat          J, M;
355f2cacb80SMatthew G. Knepley   MatNullSpace nullspace;
356f2cacb80SMatthew G. Knepley   PetscReal    dt, shift, slope, intercept;
357f2cacb80SMatthew G. Knepley   PetscBool    hasJac, hasPrec, isLin = PETSC_FALSE;
358f2cacb80SMatthew G. Knepley 
359f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
360f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
361f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
362064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
363064a246eSJacob Faibussowitsch   if (isLinear) PetscValidBoolPointer(isLinear, 7);
364064a246eSJacob Faibussowitsch   if (convRate) PetscValidRealPointer(convRate, 8);
3659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
3669566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, u_t));
367f2cacb80SMatthew G. Knepley   /* Create and view matrices */
3689566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
369f2cacb80SMatthew G. Knepley   shift = 1.0 / dt;
3709566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
3719566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
3729566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
3739566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
374f2cacb80SMatthew G. Knepley   if (hasJac && hasPrec) {
3759566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
3769566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE));
3779566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix"));
3789566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_"));
3799566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
3809566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
381f2cacb80SMatthew G. Knepley   } else {
3829566063dSJacob Faibussowitsch     PetscCall(TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE));
383f2cacb80SMatthew G. Knepley   }
3849566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
3859566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_"));
3869566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
387f2cacb80SMatthew G. Knepley   /* Check nullspace */
3889566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
389f2cacb80SMatthew G. Knepley   if (nullspace) {
390f2cacb80SMatthew G. Knepley     PetscBool isNull;
3919566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
3923c633725SBarry Smith     PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
393f2cacb80SMatthew G. Knepley   }
394f2cacb80SMatthew G. Knepley   /* Taylor test */
395f2cacb80SMatthew G. Knepley   {
396f2cacb80SMatthew G. Knepley     PetscRandom rand;
397f2cacb80SMatthew G. Knepley     Vec         du, uhat, uhat_t, r, rhat, df;
398f2cacb80SMatthew G. Knepley     PetscReal   h;
399f2cacb80SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
400f2cacb80SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
401f2cacb80SMatthew G. Knepley     PetscInt    Nv, v;
402f2cacb80SMatthew G. Knepley 
403f2cacb80SMatthew G. Knepley     /* Choose a perturbation direction */
4049566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
4059566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
4069566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
4079566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
4089566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
4099566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
410f2cacb80SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
4119566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
4129566063dSJacob Faibussowitsch     PetscCall(TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE));
413f2cacb80SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
414*9371c9d4SSatish Balay     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv)
415*9371c9d4SSatish Balay       ;
4169566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
4179566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
4189566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat_t));
4199566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
420f2cacb80SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
4219566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
4229566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat_t, h * shift, du, u_t));
423f2cacb80SMatthew 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 */
4249566063dSJacob Faibussowitsch       PetscCall(TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE));
4259566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
4269566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
427f2cacb80SMatthew G. Knepley 
428f2cacb80SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
429f2cacb80SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
430f2cacb80SMatthew G. Knepley     }
4319566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
4329566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat_t));
4339566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
4349566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
4359566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
4369566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
437f2cacb80SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
438f2cacb80SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
439f2cacb80SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL) break;
440f2cacb80SMatthew G. Knepley     }
441f2cacb80SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
4429566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
4439566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
444f2cacb80SMatthew G. Knepley     /* Slope should be about 2 */
445f2cacb80SMatthew G. Knepley     if (tol >= 0) {
4463c633725SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope);
447f2cacb80SMatthew G. Knepley     } else if (isLinear || convRate) {
448f2cacb80SMatthew G. Knepley       if (isLinear) *isLinear = isLin;
449f2cacb80SMatthew G. Knepley       if (convRate) *convRate = slope;
450f2cacb80SMatthew G. Knepley     } else {
4519566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope));
4529566063dSJacob Faibussowitsch       else PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
453f2cacb80SMatthew G. Knepley     }
454f2cacb80SMatthew G. Knepley   }
4559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
456f2cacb80SMatthew G. Knepley   PetscFunctionReturn(0);
457f2cacb80SMatthew G. Knepley }
458f2cacb80SMatthew G. Knepley 
459f2cacb80SMatthew G. Knepley /*@C
460bee9a294SMatthew G. Knepley   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
461bee9a294SMatthew G. Knepley 
462bee9a294SMatthew G. Knepley   Input Parameters:
463bee9a294SMatthew G. Knepley + ts - the TS object
4647f96f943SMatthew G. Knepley - u  - representative TS vector
4657f96f943SMatthew G. Knepley 
4667f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
467bee9a294SMatthew G. Knepley 
468bee9a294SMatthew G. Knepley   Level: developer
469bee9a294SMatthew G. Knepley @*/
470*9371c9d4SSatish Balay PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u) {
4717cdb2a12SMatthew G. Knepley   DM        dm;
4727cdb2a12SMatthew G. Knepley   SNES      snes;
473f2cacb80SMatthew G. Knepley   Vec       sol, u_t;
474f2cacb80SMatthew G. Knepley   PetscReal t;
4757cdb2a12SMatthew G. Knepley   PetscBool check;
4767cdb2a12SMatthew G. Knepley 
4777cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
4789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-dmts_check", &check));
4797cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
4809566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
4819566063dSJacob Faibussowitsch   PetscCall(VecCopy(u, sol));
4829566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, u));
4839566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
4849566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
4859566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
4869566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, u));
487f2cacb80SMatthew G. Knepley 
4889566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
4899566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL));
4909566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u_t));
4919566063dSJacob Faibussowitsch   PetscCall(DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL));
4929566063dSJacob Faibussowitsch   PetscCall(DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL));
4939566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u_t));
494f2cacb80SMatthew G. Knepley 
4959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
4967cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
4977cdb2a12SMatthew G. Knepley }
498