xref: /petsc/src/ts/utils/dmplexts.c (revision d70f29a362ae60d541be4e9a72e9494be00f9e3d)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h>     /*I "petscts.h" I*/
337c2070cSMatthew G. Knepley #include <petsc/private/snesimpl.h>
4924a1b8fSMatthew G. Knepley #include <petscds.h>
56dbbd306SMatthew G. Knepley #include <petscfv.h>
66dbbd306SMatthew G. Knepley 
76da023fcSToby Isaac static PetscErrorCode DMTSConvertPlex(DM dm, DM *plex, PetscBool copy)
86da023fcSToby Isaac {
96da023fcSToby Isaac   PetscBool      isPlex;
106da023fcSToby Isaac   PetscErrorCode ierr;
116da023fcSToby Isaac 
126da023fcSToby Isaac   PetscFunctionBegin;
136da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
146da023fcSToby Isaac   if (isPlex) {
156da023fcSToby Isaac     *plex = dm;
166da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
17f7148743SMatthew G. Knepley   } else {
18f7148743SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
19f7148743SMatthew G. Knepley     if (!*plex) {
206da023fcSToby Isaac       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
21f7148743SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
226da023fcSToby Isaac       if (copy) {
236da023fcSToby Isaac         ierr = DMCopyDMTS(dm, *plex);CHKERRQ(ierr);
246da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
259a2a23afSMatthew G. Knepley         ierr = DMCopyAuxiliaryVec(dm, *plex);CHKERRQ(ierr);
266da023fcSToby Isaac       }
27f7148743SMatthew G. Knepley     } else {
28f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
29f7148743SMatthew G. Knepley     }
306da023fcSToby Isaac   }
316da023fcSToby Isaac   PetscFunctionReturn(0);
326da023fcSToby Isaac }
336da023fcSToby Isaac 
3408449791SMatthew G. Knepley /*@
3508449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
3608449791SMatthew G. Knepley 
3708449791SMatthew G. Knepley   Input Parameters:
3808449791SMatthew G. Knepley + dm - The mesh
3908449791SMatthew G. Knepley . t - The time
4008449791SMatthew G. Knepley . locX  - Local solution
4108449791SMatthew G. Knepley - user - The user context
4208449791SMatthew G. Knepley 
4308449791SMatthew G. Knepley   Output Parameter:
443b16df42SMatthew G. Knepley . F  - Global output vector
4508449791SMatthew G. Knepley 
4608449791SMatthew G. Knepley   Level: developer
4708449791SMatthew G. Knepley 
4808449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
4908449791SMatthew G. Knepley @*/
503b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
51254c1ad2SMatthew G. Knepley {
523b16df42SMatthew G. Knepley   Vec            locF;
534a3e9fdbSToby Isaac   IS             cellIS;
546da023fcSToby Isaac   DM             plex;
554a3e9fdbSToby Isaac   PetscInt       depth;
56*d70f29a3SPierre Jolivet   PetscFormKey key = {NULL, 0, 0, 0};
57254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
58254c1ad2SMatthew G. Knepley 
59254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
606da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
614a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
62aeadca18SToby Isaac   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
634a3e9fdbSToby Isaac   if (!cellIS) {
644a3e9fdbSToby Isaac     ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);
654a3e9fdbSToby Isaac   }
666da023fcSToby Isaac   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
673b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
686528b96dSMatthew G. Knepley   ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locF, user);CHKERRQ(ierr);
699bda831aSToby Isaac   ierr = DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
709bda831aSToby Isaac   ierr = DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
716da023fcSToby Isaac   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
724a3e9fdbSToby Isaac   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
73b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
74254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
75254c1ad2SMatthew G. Knepley }
76254c1ad2SMatthew G. Knepley 
77c5d70e09SMatthew G. Knepley /*@
78ef68eab9SMatthew 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
79c5d70e09SMatthew G. Knepley 
80c5d70e09SMatthew G. Knepley   Input Parameters:
81c5d70e09SMatthew G. Knepley + dm - The mesh
82a40652d4SToby Isaac . t - The time
83a40652d4SToby Isaac . locX  - Local solution
84a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL
85c5d70e09SMatthew G. Knepley - user - The user context
86c5d70e09SMatthew G. Knepley 
87c5d70e09SMatthew G. Knepley   Level: developer
88c5d70e09SMatthew G. Knepley 
89a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM()
90c5d70e09SMatthew G. Knepley @*/
91ef68eab9SMatthew G. Knepley PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
92c5d70e09SMatthew G. Knepley {
93c5d70e09SMatthew G. Knepley   DM             plex;
94ef68eab9SMatthew G. Knepley   Vec            faceGeometryFVM = NULL;
95ef68eab9SMatthew G. Knepley   PetscInt       Nf, f;
96c5d70e09SMatthew G. Knepley   PetscErrorCode ierr;
97c5d70e09SMatthew G. Knepley 
98c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
99c5d70e09SMatthew G. Knepley   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
100ef68eab9SMatthew G. Knepley   ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr);
1019586001cSMatthew G. Knepley   if (!locX_t) {
1029586001cSMatthew G. Knepley     /* This is the RHS part */
103ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
104ef68eab9SMatthew G. Knepley       PetscObject  obj;
105ef68eab9SMatthew G. Knepley       PetscClassId id;
106ef68eab9SMatthew G. Knepley 
10744a7f3ddSMatthew G. Knepley       ierr = DMGetField(plex, f, NULL, &obj);CHKERRQ(ierr);
108ef68eab9SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
109ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
1103e9753d6SMatthew G. Knepley         ierr = DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
111ef68eab9SMatthew G. Knepley         break;
112ef68eab9SMatthew G. Knepley       }
113ef68eab9SMatthew G. Knepley     }
1149586001cSMatthew G. Knepley   }
1159586001cSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
11656cf3b9cSMatthew G. Knepley   ierr = DMPlexInsertTimeDerivativeBoundaryValues(plex, PETSC_TRUE, locX_t, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
117c5d70e09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
118c5d70e09SMatthew G. Knepley   PetscFunctionReturn(0);
119c5d70e09SMatthew G. Knepley }
120c5d70e09SMatthew G. Knepley 
12124cdb843SMatthew G. Knepley /*@
12224cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
12324cdb843SMatthew G. Knepley 
12424cdb843SMatthew G. Knepley   Input Parameters:
12524cdb843SMatthew G. Knepley + dm - The mesh
12624cdb843SMatthew G. Knepley . t - The time
12708449791SMatthew G. Knepley . locX  - Local solution
12808449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
12924cdb843SMatthew G. Knepley - user - The user context
13024cdb843SMatthew G. Knepley 
13124cdb843SMatthew G. Knepley   Output Parameter:
13208449791SMatthew G. Knepley . locF  - Local output vector
13324cdb843SMatthew G. Knepley 
13424cdb843SMatthew G. Knepley   Level: developer
13524cdb843SMatthew G. Knepley 
136d64986aaSMatthew G. Knepley .seealso: DMPlexTSComputeIFunctionFEM(), DMPlexTSComputeRHSFunctionFEM()
13724cdb843SMatthew G. Knepley @*/
13808449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
13924cdb843SMatthew G. Knepley {
1406da023fcSToby Isaac   DM             plex;
1416528b96dSMatthew G. Knepley   IS             allcellIS;
1426528b96dSMatthew G. Knepley   PetscInt       Nds, s;
14324cdb843SMatthew G. Knepley   PetscErrorCode ierr;
14424cdb843SMatthew G. Knepley 
14524cdb843SMatthew G. Knepley   PetscFunctionBegin;
1466da023fcSToby Isaac   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
1476528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
1486528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1496528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1506528b96dSMatthew G. Knepley     PetscDS          ds;
1516528b96dSMatthew G. Knepley     IS               cellIS;
15206ad1575SMatthew G. Knepley     PetscFormKey key;
1536528b96dSMatthew G. Knepley 
1546528b96dSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
1556528b96dSMatthew G. Knepley     key.value = 0;
1566528b96dSMatthew G. Knepley     key.field = 0;
15706ad1575SMatthew G. Knepley     key.part  = 0;
1586528b96dSMatthew G. Knepley     if (!key.label) {
1596528b96dSMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1606528b96dSMatthew G. Knepley       cellIS = allcellIS;
1616528b96dSMatthew G. Knepley     } else {
1626528b96dSMatthew G. Knepley       IS pointIS;
1636528b96dSMatthew G. Knepley 
1646528b96dSMatthew G. Knepley       key.value = 1;
1656528b96dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
1666528b96dSMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1676528b96dSMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1684a3e9fdbSToby Isaac     }
1696528b96dSMatthew G. Knepley     ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, locX_t, time, locF, user);CHKERRQ(ierr);
1704a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1716528b96dSMatthew G. Knepley   }
1726528b96dSMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
173b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
17424cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
17524cdb843SMatthew G. Knepley }
1767cdb2a12SMatthew G. Knepley 
177756a1f44SMatthew G. Knepley /*@
178756a1f44SMatthew G. Knepley   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
179756a1f44SMatthew G. Knepley 
180756a1f44SMatthew G. Knepley   Input Parameters:
181756a1f44SMatthew G. Knepley + dm - The mesh
182756a1f44SMatthew G. Knepley . t - The time
183756a1f44SMatthew G. Knepley . locX  - Local solution
184756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
185756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t
186756a1f44SMatthew G. Knepley - user - The user context
187756a1f44SMatthew G. Knepley 
188756a1f44SMatthew G. Knepley   Output Parameter:
189756a1f44SMatthew G. Knepley . locF  - Local output vector
190756a1f44SMatthew G. Knepley 
191756a1f44SMatthew G. Knepley   Level: developer
192756a1f44SMatthew G. Knepley 
193d64986aaSMatthew G. Knepley .seealso: DMPlexTSComputeIFunctionFEM(), DMPlexTSComputeRHSFunctionFEM()
194756a1f44SMatthew G. Knepley @*/
195756a1f44SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
196756a1f44SMatthew G. Knepley {
197756a1f44SMatthew G. Knepley   DM             plex;
1986528b96dSMatthew G. Knepley   IS             allcellIS;
199f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
2006528b96dSMatthew G. Knepley   PetscInt       Nds, s;
201756a1f44SMatthew G. Knepley   PetscErrorCode ierr;
202756a1f44SMatthew G. Knepley 
203756a1f44SMatthew G. Knepley   PetscFunctionBegin;
204756a1f44SMatthew G. Knepley   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
2056528b96dSMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
2066528b96dSMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
2076528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
2086528b96dSMatthew G. Knepley     PetscDS          ds;
2096528b96dSMatthew G. Knepley     IS               cellIS;
21006ad1575SMatthew G. Knepley     PetscFormKey key;
2116528b96dSMatthew G. Knepley 
2126528b96dSMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
2136528b96dSMatthew G. Knepley     key.value = 0;
2146528b96dSMatthew G. Knepley     key.field = 0;
21506ad1575SMatthew G. Knepley     key.part  = 0;
2166528b96dSMatthew G. Knepley     if (!key.label) {
2176528b96dSMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
2186528b96dSMatthew G. Knepley       cellIS = allcellIS;
2196528b96dSMatthew G. Knepley     } else {
2206528b96dSMatthew G. Knepley       IS pointIS;
2216528b96dSMatthew G. Knepley 
2226528b96dSMatthew G. Knepley       key.value = 1;
2236528b96dSMatthew G. Knepley       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
2246528b96dSMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
2256528b96dSMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
2266528b96dSMatthew G. Knepley     }
2276528b96dSMatthew G. Knepley     if (!s) {
2286528b96dSMatthew G. Knepley       ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
2296528b96dSMatthew G. Knepley       ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
230f04eb4edSMatthew G. Knepley       if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
231f04eb4edSMatthew G. Knepley       ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
2326528b96dSMatthew G. Knepley     }
2336528b96dSMatthew G. Knepley     ierr = DMPlexComputeJacobian_Internal(plex, key, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
2344a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
2356528b96dSMatthew G. Knepley   }
2366528b96dSMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
237756a1f44SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
238756a1f44SMatthew G. Knepley   PetscFunctionReturn(0);
239756a1f44SMatthew G. Knepley }
240756a1f44SMatthew G. Knepley 
241cb36c0f9SMatthew G. Knepley /*@
242d64986aaSMatthew G. Knepley   DMPlexTSComputeRHSFunctionFEM - Form the local residual G from the local input X using pointwise functions specified by the user
243cb36c0f9SMatthew G. Knepley 
244cb36c0f9SMatthew G. Knepley   Input Parameters:
245cb36c0f9SMatthew G. Knepley + dm - The mesh
246cb36c0f9SMatthew G. Knepley . t - The time
247cb36c0f9SMatthew G. Knepley . locX  - Local solution
248cb36c0f9SMatthew G. Knepley - user - The user context
249cb36c0f9SMatthew G. Knepley 
250cb36c0f9SMatthew G. Knepley   Output Parameter:
251cb36c0f9SMatthew G. Knepley . locG  - Local output vector
252cb36c0f9SMatthew G. Knepley 
253cb36c0f9SMatthew G. Knepley   Level: developer
254cb36c0f9SMatthew G. Knepley 
255cb36c0f9SMatthew G. Knepley .seealso: DMPlexTSComputeIFunctionFEM(), DMPlexTSComputeIJacobianFEM()
256cb36c0f9SMatthew G. Knepley @*/
257d64986aaSMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locG, void *user)
258cb36c0f9SMatthew G. Knepley {
259cb36c0f9SMatthew G. Knepley   DM             plex;
260cb36c0f9SMatthew G. Knepley   IS             allcellIS;
261cb36c0f9SMatthew G. Knepley   PetscInt       Nds, s;
262cb36c0f9SMatthew G. Knepley   PetscErrorCode ierr;
263cb36c0f9SMatthew G. Knepley 
264cb36c0f9SMatthew G. Knepley   PetscFunctionBegin;
265cb36c0f9SMatthew G. Knepley   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
266cb36c0f9SMatthew G. Knepley   ierr = DMPlexGetAllCells_Internal(plex, &allcellIS);CHKERRQ(ierr);
267cb36c0f9SMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
268cb36c0f9SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
269cb36c0f9SMatthew G. Knepley     PetscDS          ds;
270cb36c0f9SMatthew G. Knepley     IS               cellIS;
271cb36c0f9SMatthew G. Knepley     PetscFormKey key;
272cb36c0f9SMatthew G. Knepley 
273cb36c0f9SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);CHKERRQ(ierr);
274cb36c0f9SMatthew G. Knepley     key.value = 0;
275cb36c0f9SMatthew G. Knepley     key.field = 0;
276cb36c0f9SMatthew G. Knepley     key.part  = 100;
277cb36c0f9SMatthew G. Knepley     if (!key.label) {
278cb36c0f9SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
279cb36c0f9SMatthew G. Knepley       cellIS = allcellIS;
280cb36c0f9SMatthew G. Knepley     } else {
281cb36c0f9SMatthew G. Knepley       IS pointIS;
282cb36c0f9SMatthew G. Knepley 
283cb36c0f9SMatthew G. Knepley       key.value = 1;
284cb36c0f9SMatthew G. Knepley       ierr = DMLabelGetStratumIS(key.label, key.value, &pointIS);CHKERRQ(ierr);
285cb36c0f9SMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
286cb36c0f9SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
287cb36c0f9SMatthew G. Knepley     }
288cb36c0f9SMatthew G. Knepley     ierr = DMPlexComputeResidual_Internal(plex, key, cellIS, time, locX, NULL, time, locG, user);CHKERRQ(ierr);
289cb36c0f9SMatthew G. Knepley     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
290cb36c0f9SMatthew G. Knepley   }
291cb36c0f9SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
292cb36c0f9SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
293cb36c0f9SMatthew G. Knepley   PetscFunctionReturn(0);
294cb36c0f9SMatthew G. Knepley }
295cb36c0f9SMatthew G. Knepley 
296bee9a294SMatthew G. Knepley /*@C
297f2cacb80SMatthew G. Knepley   DMTSCheckResidual - Check the residual of the exact solution
298f2cacb80SMatthew G. Knepley 
299f2cacb80SMatthew G. Knepley   Input Parameters:
300f2cacb80SMatthew G. Knepley + ts  - the TS object
301f2cacb80SMatthew G. Knepley . dm  - the DM
302f2cacb80SMatthew G. Knepley . t   - the time
303f2cacb80SMatthew G. Knepley . u   - a DM vector
304f2cacb80SMatthew G. Knepley . u_t - a DM vector
305f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
306f2cacb80SMatthew G. Knepley 
307f2cacb80SMatthew G. Knepley   Output Parameters:
308f2cacb80SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
309f2cacb80SMatthew G. Knepley 
310f2cacb80SMatthew G. Knepley   Level: developer
311f2cacb80SMatthew G. Knepley 
312f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckJacobian(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
313f2cacb80SMatthew G. Knepley @*/
314f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual)
315f2cacb80SMatthew G. Knepley {
316f2cacb80SMatthew G. Knepley   MPI_Comm       comm;
317f2cacb80SMatthew G. Knepley   Vec            r;
318f2cacb80SMatthew G. Knepley   PetscReal      res;
319f2cacb80SMatthew G. Knepley   PetscErrorCode ierr;
320f2cacb80SMatthew G. Knepley 
321f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
322f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
323f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
324064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
325064a246eSJacob Faibussowitsch   if (residual) PetscValidRealPointer(residual, 7);
326f2cacb80SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
327f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr);
328f2cacb80SMatthew G. Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
329f2cacb80SMatthew G. Knepley   ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr);
330f2cacb80SMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
331f2cacb80SMatthew G. Knepley   if (tol >= 0.0) {
3322c71b3e2SJacob Faibussowitsch     PetscCheckFalse(res > tol,comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
333f2cacb80SMatthew G. Knepley   } else if (residual) {
334f2cacb80SMatthew G. Knepley     *residual = res;
335f2cacb80SMatthew G. Knepley   } else {
336f2cacb80SMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
337f2cacb80SMatthew G. Knepley     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
3380c364540SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", (PetscObject) dm);CHKERRQ(ierr);
339f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
340f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
341f2cacb80SMatthew G. Knepley     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
3420c364540SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", NULL);CHKERRQ(ierr);
343f2cacb80SMatthew G. Knepley   }
344f2cacb80SMatthew G. Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
345f2cacb80SMatthew G. Knepley   PetscFunctionReturn(0);
346f2cacb80SMatthew G. Knepley }
347f2cacb80SMatthew G. Knepley 
348f2cacb80SMatthew G. Knepley /*@C
349f2cacb80SMatthew G. Knepley   DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
350f2cacb80SMatthew G. Knepley 
351f2cacb80SMatthew G. Knepley   Input Parameters:
352f2cacb80SMatthew G. Knepley + ts  - the TS object
353f2cacb80SMatthew G. Knepley . dm  - the DM
354f2cacb80SMatthew G. Knepley . t   - the time
355f2cacb80SMatthew G. Knepley . u   - a DM vector
356f2cacb80SMatthew G. Knepley . u_t - a DM vector
357f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
358f2cacb80SMatthew G. Knepley 
359f2cacb80SMatthew G. Knepley   Output Parameters:
360f2cacb80SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
361f2cacb80SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
362f2cacb80SMatthew G. Knepley 
363f2cacb80SMatthew G. Knepley   Level: developer
364f2cacb80SMatthew G. Knepley 
365f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckResidual(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
366f2cacb80SMatthew G. Knepley @*/
367f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
368f2cacb80SMatthew G. Knepley {
369f2cacb80SMatthew G. Knepley   MPI_Comm       comm;
370f2cacb80SMatthew G. Knepley   PetscDS        ds;
371f2cacb80SMatthew G. Knepley   Mat            J, M;
372f2cacb80SMatthew G. Knepley   MatNullSpace   nullspace;
373f2cacb80SMatthew G. Knepley   PetscReal      dt, shift, slope, intercept;
374f2cacb80SMatthew G. Knepley   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
375f2cacb80SMatthew G. Knepley   PetscErrorCode ierr;
376f2cacb80SMatthew G. Knepley 
377f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
378f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
379f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
380064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
381064a246eSJacob Faibussowitsch   if (isLinear) PetscValidBoolPointer(isLinear, 7);
382064a246eSJacob Faibussowitsch   if (convRate) PetscValidRealPointer(convRate, 8);
383f2cacb80SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
384f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr);
385f2cacb80SMatthew G. Knepley   /* Create and view matrices */
386f2cacb80SMatthew G. Knepley   ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr);
387f2cacb80SMatthew G. Knepley   shift = 1.0/dt;
388f2cacb80SMatthew G. Knepley   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
389f2cacb80SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
390f2cacb80SMatthew G. Knepley   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
391f2cacb80SMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
392f2cacb80SMatthew G. Knepley   if (hasJac && hasPrec) {
393f2cacb80SMatthew G. Knepley     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
394f2cacb80SMatthew G. Knepley     ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE);CHKERRQ(ierr);
395f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
396f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
397f2cacb80SMatthew G. Knepley     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
398f2cacb80SMatthew G. Knepley     ierr = MatDestroy(&M);CHKERRQ(ierr);
399f2cacb80SMatthew G. Knepley   } else {
400f2cacb80SMatthew G. Knepley     ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE);CHKERRQ(ierr);
401f2cacb80SMatthew G. Knepley   }
402f2cacb80SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
403f2cacb80SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
404f2cacb80SMatthew G. Knepley   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
405f2cacb80SMatthew G. Knepley   /* Check nullspace */
406f2cacb80SMatthew G. Knepley   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
407f2cacb80SMatthew G. Knepley   if (nullspace) {
408f2cacb80SMatthew G. Knepley     PetscBool isNull;
409f2cacb80SMatthew G. Knepley     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
4102c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!isNull,comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
411f2cacb80SMatthew G. Knepley   }
412f2cacb80SMatthew G. Knepley   /* Taylor test */
413f2cacb80SMatthew G. Knepley   {
414f2cacb80SMatthew G. Knepley     PetscRandom rand;
415f2cacb80SMatthew G. Knepley     Vec         du, uhat, uhat_t, r, rhat, df;
416f2cacb80SMatthew G. Knepley     PetscReal   h;
417f2cacb80SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
418f2cacb80SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
419f2cacb80SMatthew G. Knepley     PetscInt    Nv, v;
420f2cacb80SMatthew G. Knepley 
421f2cacb80SMatthew G. Knepley     /* Choose a perturbation direction */
422f2cacb80SMatthew G. Knepley     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
423f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
424f2cacb80SMatthew G. Knepley     ierr = VecSetRandom(du, rand);CHKERRQ(ierr);
425f2cacb80SMatthew G. Knepley     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
426f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
427f2cacb80SMatthew G. Knepley     ierr = MatMult(J, du, df);CHKERRQ(ierr);
428f2cacb80SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
429f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
430f2cacb80SMatthew G. Knepley     ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr);
431f2cacb80SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
432f2cacb80SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
433f2cacb80SMatthew G. Knepley     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
434f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
435f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat_t);CHKERRQ(ierr);
436f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
437f2cacb80SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
438f2cacb80SMatthew G. Knepley       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
439f2cacb80SMatthew G. Knepley       ierr = VecWAXPY(uhat_t, h*shift, du, u_t);CHKERRQ(ierr);
440f2cacb80SMatthew 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 */
441f2cacb80SMatthew G. Knepley       ierr = TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE);CHKERRQ(ierr);
442f2cacb80SMatthew G. Knepley       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
443f2cacb80SMatthew G. Knepley       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
444f2cacb80SMatthew G. Knepley 
445f2cacb80SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
446f2cacb80SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
447f2cacb80SMatthew G. Knepley     }
448f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
449f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&uhat_t);CHKERRQ(ierr);
450f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
451f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&df);CHKERRQ(ierr);
452f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&r);CHKERRQ(ierr);
453f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&du);CHKERRQ(ierr);
454f2cacb80SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
455f2cacb80SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
456f2cacb80SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL)    break;
457f2cacb80SMatthew G. Knepley     }
458f2cacb80SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
459f2cacb80SMatthew G. Knepley     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
460f2cacb80SMatthew G. Knepley     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
461f2cacb80SMatthew G. Knepley     /* Slope should be about 2 */
462f2cacb80SMatthew G. Knepley     if (tol >= 0) {
4632c71b3e2SJacob Faibussowitsch       PetscCheckFalse(!isLin && PetscAbsReal(2 - slope) > tol,comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
464f2cacb80SMatthew G. Knepley     } else if (isLinear || convRate) {
465f2cacb80SMatthew G. Knepley       if (isLinear) *isLinear = isLin;
466f2cacb80SMatthew G. Knepley       if (convRate) *convRate = slope;
467f2cacb80SMatthew G. Knepley     } else {
468f2cacb80SMatthew G. Knepley       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
469f2cacb80SMatthew G. Knepley       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
470f2cacb80SMatthew G. Knepley     }
471f2cacb80SMatthew G. Knepley   }
472f2cacb80SMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
473f2cacb80SMatthew G. Knepley   PetscFunctionReturn(0);
474f2cacb80SMatthew G. Knepley }
475f2cacb80SMatthew G. Knepley 
476f2cacb80SMatthew G. Knepley /*@C
477bee9a294SMatthew G. Knepley   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
478bee9a294SMatthew G. Knepley 
479bee9a294SMatthew G. Knepley   Input Parameters:
480bee9a294SMatthew G. Knepley + ts - the TS object
4817f96f943SMatthew G. Knepley - u  - representative TS vector
4827f96f943SMatthew G. Knepley 
4837f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
484bee9a294SMatthew G. Knepley 
485bee9a294SMatthew G. Knepley   Level: developer
486bee9a294SMatthew G. Knepley @*/
4877f96f943SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u)
4887cdb2a12SMatthew G. Knepley {
4897cdb2a12SMatthew G. Knepley   DM             dm;
4907cdb2a12SMatthew G. Knepley   SNES           snes;
491f2cacb80SMatthew G. Knepley   Vec            sol, u_t;
492f2cacb80SMatthew G. Knepley   PetscReal      t;
4937cdb2a12SMatthew G. Knepley   PetscBool      check;
4947cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
4957cdb2a12SMatthew G. Knepley 
4967cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
497c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
4987cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
49972fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
5006fb91e28SMatthew G. Knepley   ierr = VecCopy(u, sol);CHKERRQ(ierr);
5016fb91e28SMatthew G. Knepley   ierr = TSSetSolution(ts, u);CHKERRQ(ierr);
5027cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
5037cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
5047cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
5056fb91e28SMatthew G. Knepley   ierr = SNESSetSolution(snes, u);CHKERRQ(ierr);
506f2cacb80SMatthew G. Knepley 
507f2cacb80SMatthew G. Knepley   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
508f2cacb80SMatthew G. Knepley   ierr = DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL);CHKERRQ(ierr);
509f2cacb80SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &u_t);CHKERRQ(ierr);
510f2cacb80SMatthew G. Knepley   ierr = DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL);CHKERRQ(ierr);
511f2cacb80SMatthew G. Knepley   ierr = DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL);CHKERRQ(ierr);
512f2cacb80SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &u_t);CHKERRQ(ierr);
513f2cacb80SMatthew G. Knepley 
51472fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
5157cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
5167cdb2a12SMatthew G. Knepley }
517