xref: /petsc/src/ts/utils/dmplexts.c (revision 0c36454096925ed0dc33600f014db4c058649b2b)
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         PetscInt    i;
246da023fcSToby Isaac         PetscObject obj;
256da023fcSToby Isaac         const char *comps[3] = {"A","dmAux","dmCh"};
266da023fcSToby Isaac 
276da023fcSToby Isaac         ierr = DMCopyDMTS(dm, *plex);CHKERRQ(ierr);
286da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
296da023fcSToby Isaac         for (i = 0; i < 3; i++) {
306da023fcSToby Isaac           ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr);
316da023fcSToby Isaac           ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr);
326da023fcSToby Isaac         }
336da023fcSToby Isaac       }
34f7148743SMatthew G. Knepley     } else {
35f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
36f7148743SMatthew G. Knepley     }
376da023fcSToby Isaac   }
386da023fcSToby Isaac   PetscFunctionReturn(0);
396da023fcSToby Isaac }
406da023fcSToby Isaac 
4108449791SMatthew G. Knepley /*@
4208449791SMatthew G. Knepley   DMPlexTSComputeRHSFunctionFVM - Form the local forcing F from the local input X using pointwise functions specified by the user
4308449791SMatthew G. Knepley 
4408449791SMatthew G. Knepley   Input Parameters:
4508449791SMatthew G. Knepley + dm - The mesh
4608449791SMatthew G. Knepley . t - The time
4708449791SMatthew G. Knepley . locX  - Local solution
4808449791SMatthew G. Knepley - user - The user context
4908449791SMatthew G. Knepley 
5008449791SMatthew G. Knepley   Output Parameter:
513b16df42SMatthew G. Knepley . F  - Global output vector
5208449791SMatthew G. Knepley 
5308449791SMatthew G. Knepley   Level: developer
5408449791SMatthew G. Knepley 
5508449791SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
5608449791SMatthew G. Knepley @*/
573b16df42SMatthew G. Knepley PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM dm, PetscReal time, Vec locX, Vec F, void *user)
58254c1ad2SMatthew G. Knepley {
593b16df42SMatthew G. Knepley   Vec            locF;
604a3e9fdbSToby Isaac   IS             cellIS;
616da023fcSToby Isaac   DM             plex;
624a3e9fdbSToby Isaac   PetscInt       depth;
63254c1ad2SMatthew G. Knepley   PetscErrorCode ierr;
64254c1ad2SMatthew G. Knepley 
65254c1ad2SMatthew G. Knepley   PetscFunctionBegin;
666da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
674a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
68aeadca18SToby Isaac   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
694a3e9fdbSToby Isaac   if (!cellIS) {
704a3e9fdbSToby Isaac     ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);
714a3e9fdbSToby Isaac   }
726da023fcSToby Isaac   ierr = DMGetLocalVector(plex, &locF);CHKERRQ(ierr);
733b16df42SMatthew G. Knepley   ierr = VecZeroEntries(locF);CHKERRQ(ierr);
744a3e9fdbSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cellIS, time, locX, NULL, time, locF, user);CHKERRQ(ierr);
759bda831aSToby Isaac   ierr = DMLocalToGlobalBegin(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
769bda831aSToby Isaac   ierr = DMLocalToGlobalEnd(plex, locF, ADD_VALUES, F);CHKERRQ(ierr);
776da023fcSToby Isaac   ierr = DMRestoreLocalVector(plex, &locF);CHKERRQ(ierr);
784a3e9fdbSToby Isaac   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
79b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
80254c1ad2SMatthew G. Knepley   PetscFunctionReturn(0);
81254c1ad2SMatthew G. Knepley }
82254c1ad2SMatthew G. Knepley 
83c5d70e09SMatthew G. Knepley /*@
84ef68eab9SMatthew 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
85c5d70e09SMatthew G. Knepley 
86c5d70e09SMatthew G. Knepley   Input Parameters:
87c5d70e09SMatthew G. Knepley + dm - The mesh
88a40652d4SToby Isaac . t - The time
89a40652d4SToby Isaac . locX  - Local solution
90a40652d4SToby Isaac . locX_t - Local solution time derivative, or NULL
91c5d70e09SMatthew G. Knepley - user - The user context
92c5d70e09SMatthew G. Knepley 
93c5d70e09SMatthew G. Knepley   Level: developer
94c5d70e09SMatthew G. Knepley 
95a40652d4SToby Isaac .seealso: DMPlexComputeJacobianActionFEM()
96c5d70e09SMatthew G. Knepley @*/
97ef68eab9SMatthew G. Knepley PetscErrorCode DMPlexTSComputeBoundary(DM dm, PetscReal time, Vec locX, Vec locX_t, void *user)
98c5d70e09SMatthew G. Knepley {
99c5d70e09SMatthew G. Knepley   DM             plex;
100ef68eab9SMatthew G. Knepley   Vec            faceGeometryFVM = NULL;
101ef68eab9SMatthew G. Knepley   PetscInt       Nf, f;
102c5d70e09SMatthew G. Knepley   PetscErrorCode ierr;
103c5d70e09SMatthew G. Knepley 
104c5d70e09SMatthew G. Knepley   PetscFunctionBegin;
105c5d70e09SMatthew G. Knepley   ierr = DMTSConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
106ef68eab9SMatthew G. Knepley   ierr = DMGetNumFields(plex, &Nf);CHKERRQ(ierr);
1079586001cSMatthew G. Knepley   if (!locX_t) {
1089586001cSMatthew G. Knepley     /* This is the RHS part */
109ef68eab9SMatthew G. Knepley     for (f = 0; f < Nf; f++) {
110ef68eab9SMatthew G. Knepley       PetscObject  obj;
111ef68eab9SMatthew G. Knepley       PetscClassId id;
112ef68eab9SMatthew G. Knepley 
11344a7f3ddSMatthew G. Knepley       ierr = DMGetField(plex, f, NULL, &obj);CHKERRQ(ierr);
114ef68eab9SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
115ef68eab9SMatthew G. Knepley       if (id == PETSCFV_CLASSID) {
1163e9753d6SMatthew G. Knepley         ierr = DMPlexGetGeometryFVM(plex, &faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
117ef68eab9SMatthew G. Knepley         break;
118ef68eab9SMatthew G. Knepley       }
119ef68eab9SMatthew G. Knepley     }
1209586001cSMatthew G. Knepley   }
1219586001cSMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, time, faceGeometryFVM, NULL, NULL);CHKERRQ(ierr);
122a40652d4SToby Isaac   /* TODO: locX_t */
123c5d70e09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
124c5d70e09SMatthew G. Knepley   PetscFunctionReturn(0);
125c5d70e09SMatthew G. Knepley }
126c5d70e09SMatthew G. Knepley 
12724cdb843SMatthew G. Knepley /*@
12824cdb843SMatthew G. Knepley   DMPlexTSComputeIFunctionFEM - Form the local residual F from the local input X using pointwise functions specified by the user
12924cdb843SMatthew G. Knepley 
13024cdb843SMatthew G. Knepley   Input Parameters:
13124cdb843SMatthew G. Knepley + dm - The mesh
13224cdb843SMatthew G. Knepley . t - The time
13308449791SMatthew G. Knepley . locX  - Local solution
13408449791SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
13524cdb843SMatthew G. Knepley - user - The user context
13624cdb843SMatthew G. Knepley 
13724cdb843SMatthew G. Knepley   Output Parameter:
13808449791SMatthew G. Knepley . locF  - Local output vector
13924cdb843SMatthew G. Knepley 
14024cdb843SMatthew G. Knepley   Level: developer
14124cdb843SMatthew G. Knepley 
14224cdb843SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
14324cdb843SMatthew G. Knepley @*/
14408449791SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIFunctionFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
14524cdb843SMatthew G. Knepley {
1466da023fcSToby Isaac   DM             plex;
1474a3e9fdbSToby Isaac   IS             cellIS;
1484a3e9fdbSToby Isaac   PetscInt       depth;
14924cdb843SMatthew G. Knepley   PetscErrorCode ierr;
15024cdb843SMatthew G. Knepley 
15124cdb843SMatthew G. Knepley   PetscFunctionBegin;
1526da023fcSToby Isaac   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1534a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
154aeadca18SToby Isaac   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
1554a3e9fdbSToby Isaac   if (!cellIS) {
1564a3e9fdbSToby Isaac     ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);
1574a3e9fdbSToby Isaac   }
1584a3e9fdbSToby Isaac   ierr = DMPlexComputeResidual_Internal(plex, cellIS, time, locX, locX_t, time, locF, user);CHKERRQ(ierr);
1594a3e9fdbSToby Isaac   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
160b2338dbcSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
16124cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
16224cdb843SMatthew G. Knepley }
1637cdb2a12SMatthew G. Knepley 
164756a1f44SMatthew G. Knepley /*@
165756a1f44SMatthew G. Knepley   DMPlexTSComputeIJacobianFEM - Form the local Jacobian J from the local input X using pointwise functions specified by the user
166756a1f44SMatthew G. Knepley 
167756a1f44SMatthew G. Knepley   Input Parameters:
168756a1f44SMatthew G. Knepley + dm - The mesh
169756a1f44SMatthew G. Knepley . t - The time
170756a1f44SMatthew G. Knepley . locX  - Local solution
171756a1f44SMatthew G. Knepley . locX_t - Local solution time derivative, or NULL
172756a1f44SMatthew G. Knepley . X_tshift - The multiplicative parameter for dF/du_t
173756a1f44SMatthew G. Knepley - user - The user context
174756a1f44SMatthew G. Knepley 
175756a1f44SMatthew G. Knepley   Output Parameter:
176756a1f44SMatthew G. Knepley . locF  - Local output vector
177756a1f44SMatthew G. Knepley 
178756a1f44SMatthew G. Knepley   Level: developer
179756a1f44SMatthew G. Knepley 
180756a1f44SMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM()
181756a1f44SMatthew G. Knepley @*/
182756a1f44SMatthew G. Knepley PetscErrorCode DMPlexTSComputeIJacobianFEM(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
183756a1f44SMatthew G. Knepley {
184756a1f44SMatthew G. Knepley   DM             plex;
185f04eb4edSMatthew G. Knepley   PetscDS        prob;
186f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
1879044fa66SMatthew G. Knepley   IS             cellIS;
1884a3e9fdbSToby Isaac   PetscInt       depth;
189756a1f44SMatthew G. Knepley   PetscErrorCode ierr;
190756a1f44SMatthew G. Knepley 
191756a1f44SMatthew G. Knepley   PetscFunctionBegin;
192756a1f44SMatthew G. Knepley   ierr = DMTSConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
193f04eb4edSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
194f04eb4edSMatthew G. Knepley   ierr = PetscDSHasJacobian(prob, &hasJac);CHKERRQ(ierr);
195f04eb4edSMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(prob, &hasPrec);CHKERRQ(ierr);
196f04eb4edSMatthew G. Knepley   if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
197f04eb4edSMatthew G. Knepley   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1984a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex,&depth);CHKERRQ(ierr);
199aeadca18SToby Isaac   ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
2009044fa66SMatthew G. Knepley   if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);}
2014a3e9fdbSToby Isaac   ierr = DMPlexComputeJacobian_Internal(plex, cellIS, time, X_tShift, locX, locX_t, Jac, JacP, user);CHKERRQ(ierr);
2024a3e9fdbSToby Isaac   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
203756a1f44SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
204756a1f44SMatthew G. Knepley   PetscFunctionReturn(0);
205756a1f44SMatthew G. Knepley }
206756a1f44SMatthew G. Knepley 
207bee9a294SMatthew G. Knepley /*@C
208f2cacb80SMatthew G. Knepley   DMTSCheckResidual - Check the residual of the exact solution
209f2cacb80SMatthew G. Knepley 
210f2cacb80SMatthew G. Knepley   Input Parameters:
211f2cacb80SMatthew G. Knepley + ts  - the TS object
212f2cacb80SMatthew G. Knepley . dm  - the DM
213f2cacb80SMatthew G. Knepley . t   - the time
214f2cacb80SMatthew G. Knepley . u   - a DM vector
215f2cacb80SMatthew G. Knepley . u_t - a DM vector
216f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
217f2cacb80SMatthew G. Knepley 
218f2cacb80SMatthew G. Knepley   Output Parameters:
219f2cacb80SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
220f2cacb80SMatthew G. Knepley 
221f2cacb80SMatthew G. Knepley   Level: developer
222f2cacb80SMatthew G. Knepley 
223f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckJacobian(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
224f2cacb80SMatthew G. Knepley @*/
225f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckResidual(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscReal *residual)
226f2cacb80SMatthew G. Knepley {
227f2cacb80SMatthew G. Knepley   MPI_Comm       comm;
228f2cacb80SMatthew G. Knepley   Vec            r;
229f2cacb80SMatthew G. Knepley   PetscReal      res;
230f2cacb80SMatthew G. Knepley   PetscErrorCode ierr;
231f2cacb80SMatthew G. Knepley 
232f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
233f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
234f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
235f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
236f2cacb80SMatthew G. Knepley   if (residual) PetscValidRealPointer(residual, 5);
237f2cacb80SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
238f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr);
239f2cacb80SMatthew G. Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
240f2cacb80SMatthew G. Knepley   ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr);
241f2cacb80SMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
242f2cacb80SMatthew G. Knepley   if (tol >= 0.0) {
243f2cacb80SMatthew G. Knepley     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
244f2cacb80SMatthew G. Knepley   } else if (residual) {
245f2cacb80SMatthew G. Knepley     *residual = res;
246f2cacb80SMatthew G. Knepley   } else {
247f2cacb80SMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
248f2cacb80SMatthew G. Knepley     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
249*0c364540SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", (PetscObject) dm);CHKERRQ(ierr);
250f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
251f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
252f2cacb80SMatthew G. Knepley     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
253*0c364540SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) r, "__Vec_bc_zero__", NULL);CHKERRQ(ierr);
254f2cacb80SMatthew G. Knepley   }
255f2cacb80SMatthew G. Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
256f2cacb80SMatthew G. Knepley   PetscFunctionReturn(0);
257f2cacb80SMatthew G. Knepley }
258f2cacb80SMatthew G. Knepley 
259f2cacb80SMatthew G. Knepley /*@C
260f2cacb80SMatthew G. Knepley   DMTSCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
261f2cacb80SMatthew G. Knepley 
262f2cacb80SMatthew G. Knepley   Input Parameters:
263f2cacb80SMatthew G. Knepley + ts  - the TS object
264f2cacb80SMatthew G. Knepley . dm  - the DM
265f2cacb80SMatthew G. Knepley . t   - the time
266f2cacb80SMatthew G. Knepley . u   - a DM vector
267f2cacb80SMatthew G. Knepley . u_t - a DM vector
268f2cacb80SMatthew G. Knepley - tol - A tolerance for the check, or -1 to print the results instead
269f2cacb80SMatthew G. Knepley 
270f2cacb80SMatthew G. Knepley   Output Parameters:
271f2cacb80SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
272f2cacb80SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
273f2cacb80SMatthew G. Knepley 
274f2cacb80SMatthew G. Knepley   Level: developer
275f2cacb80SMatthew G. Knepley 
276f2cacb80SMatthew G. Knepley .seealso: DNTSCheckFromOptions(), DMTSCheckResidual(), DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
277f2cacb80SMatthew G. Knepley @*/
278f2cacb80SMatthew G. Knepley PetscErrorCode DMTSCheckJacobian(TS ts, DM dm, PetscReal t, Vec u, Vec u_t, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
279f2cacb80SMatthew G. Knepley {
280f2cacb80SMatthew G. Knepley   MPI_Comm       comm;
281f2cacb80SMatthew G. Knepley   PetscDS        ds;
282f2cacb80SMatthew G. Knepley   Mat            J, M;
283f2cacb80SMatthew G. Knepley   MatNullSpace   nullspace;
284f2cacb80SMatthew G. Knepley   PetscReal      dt, shift, slope, intercept;
285f2cacb80SMatthew G. Knepley   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
286f2cacb80SMatthew G. Knepley   PetscErrorCode ierr;
287f2cacb80SMatthew G. Knepley 
288f2cacb80SMatthew G. Knepley   PetscFunctionBegin;
289f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
290f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
291f2cacb80SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
292f2cacb80SMatthew G. Knepley   if (isLinear) PetscValidBoolPointer(isLinear, 5);
293f2cacb80SMatthew G. Knepley   if (convRate) PetscValidRealPointer(convRate, 5);
294f2cacb80SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
295f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, u_t);CHKERRQ(ierr);
296f2cacb80SMatthew G. Knepley   /* Create and view matrices */
297f2cacb80SMatthew G. Knepley   ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr);
298f2cacb80SMatthew G. Knepley   shift = 1.0/dt;
299f2cacb80SMatthew G. Knepley   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
300f2cacb80SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
301f2cacb80SMatthew G. Knepley   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
302f2cacb80SMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
303f2cacb80SMatthew G. Knepley   if (hasJac && hasPrec) {
304f2cacb80SMatthew G. Knepley     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
305f2cacb80SMatthew G. Knepley     ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, M, PETSC_FALSE);CHKERRQ(ierr);
306f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
307f2cacb80SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
308f2cacb80SMatthew G. Knepley     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
309f2cacb80SMatthew G. Knepley     ierr = MatDestroy(&M);CHKERRQ(ierr);
310f2cacb80SMatthew G. Knepley   } else {
311f2cacb80SMatthew G. Knepley     ierr = TSComputeIJacobian(ts, t, u, u_t, shift, J, J, PETSC_FALSE);CHKERRQ(ierr);
312f2cacb80SMatthew G. Knepley   }
313f2cacb80SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
314f2cacb80SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
315f2cacb80SMatthew G. Knepley   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
316f2cacb80SMatthew G. Knepley   /* Check nullspace */
317f2cacb80SMatthew G. Knepley   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
318f2cacb80SMatthew G. Knepley   if (nullspace) {
319f2cacb80SMatthew G. Knepley     PetscBool isNull;
320f2cacb80SMatthew G. Knepley     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
321f2cacb80SMatthew G. Knepley     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
322f2cacb80SMatthew G. Knepley   }
323f2cacb80SMatthew G. Knepley   ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
324f2cacb80SMatthew G. Knepley   /* Taylor test */
325f2cacb80SMatthew G. Knepley   {
326f2cacb80SMatthew G. Knepley     PetscRandom rand;
327f2cacb80SMatthew G. Knepley     Vec         du, uhat, uhat_t, r, rhat, df;
328f2cacb80SMatthew G. Knepley     PetscReal   h;
329f2cacb80SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
330f2cacb80SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
331f2cacb80SMatthew G. Knepley     PetscInt    Nv, v;
332f2cacb80SMatthew G. Knepley 
333f2cacb80SMatthew G. Knepley     /* Choose a perturbation direction */
334f2cacb80SMatthew G. Knepley     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
335f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
336f2cacb80SMatthew G. Knepley     ierr = VecSetRandom(du, rand); CHKERRQ(ierr);
337f2cacb80SMatthew G. Knepley     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
338f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
339f2cacb80SMatthew G. Knepley     ierr = MatMult(J, du, df);CHKERRQ(ierr);
340f2cacb80SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
341f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
342f2cacb80SMatthew G. Knepley     ierr = TSComputeIFunction(ts, t, u, u_t, r, PETSC_FALSE);CHKERRQ(ierr);
343f2cacb80SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
344f2cacb80SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
345f2cacb80SMatthew G. Knepley     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
346f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
347f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat_t);CHKERRQ(ierr);
348f2cacb80SMatthew G. Knepley     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
349f2cacb80SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
350f2cacb80SMatthew G. Knepley       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
351f2cacb80SMatthew G. Knepley       ierr = VecWAXPY(uhat_t, h*shift, du, u_t);CHKERRQ(ierr);
352f2cacb80SMatthew 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 */
353f2cacb80SMatthew G. Knepley       ierr = TSComputeIFunction(ts, t, uhat, uhat_t, rhat, PETSC_FALSE);CHKERRQ(ierr);
354f2cacb80SMatthew G. Knepley       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
355f2cacb80SMatthew G. Knepley       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
356f2cacb80SMatthew G. Knepley 
357f2cacb80SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
358f2cacb80SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
359f2cacb80SMatthew G. Knepley     }
360f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
361f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&uhat_t);CHKERRQ(ierr);
362f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
363f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&df);CHKERRQ(ierr);
364f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&r);CHKERRQ(ierr);
365f2cacb80SMatthew G. Knepley     ierr = VecDestroy(&du);CHKERRQ(ierr);
366f2cacb80SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
367f2cacb80SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
368f2cacb80SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL)    break;
369f2cacb80SMatthew G. Knepley     }
370f2cacb80SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
371f2cacb80SMatthew G. Knepley     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
372f2cacb80SMatthew G. Knepley     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
373f2cacb80SMatthew G. Knepley     /* Slope should be about 2 */
374f2cacb80SMatthew G. Knepley     if (tol >= 0) {
375f2cacb80SMatthew G. Knepley       if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
376f2cacb80SMatthew G. Knepley     } else if (isLinear || convRate) {
377f2cacb80SMatthew G. Knepley       if (isLinear) *isLinear = isLin;
378f2cacb80SMatthew G. Knepley       if (convRate) *convRate = slope;
379f2cacb80SMatthew G. Knepley     } else {
380f2cacb80SMatthew G. Knepley       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
381f2cacb80SMatthew G. Knepley       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
382f2cacb80SMatthew G. Knepley     }
383f2cacb80SMatthew G. Knepley   }
384f2cacb80SMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
385f2cacb80SMatthew G. Knepley   PetscFunctionReturn(0);
386f2cacb80SMatthew G. Knepley }
387f2cacb80SMatthew G. Knepley 
388f2cacb80SMatthew G. Knepley /*@C
389bee9a294SMatthew G. Knepley   DMTSCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
390bee9a294SMatthew G. Knepley 
391bee9a294SMatthew G. Knepley   Input Parameters:
392bee9a294SMatthew G. Knepley + ts - the TS object
3937f96f943SMatthew G. Knepley - u  - representative TS vector
3947f96f943SMatthew G. Knepley 
3957f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
396bee9a294SMatthew G. Knepley 
397bee9a294SMatthew G. Knepley   Level: developer
398bee9a294SMatthew G. Knepley @*/
3997f96f943SMatthew G. Knepley PetscErrorCode DMTSCheckFromOptions(TS ts, Vec u)
4007cdb2a12SMatthew G. Knepley {
4017cdb2a12SMatthew G. Knepley   DM             dm;
4027cdb2a12SMatthew G. Knepley   SNES           snes;
403f2cacb80SMatthew G. Knepley   Vec            sol, u_t;
404f2cacb80SMatthew G. Knepley   PetscReal      t;
4057cdb2a12SMatthew G. Knepley   PetscBool      check;
4067cdb2a12SMatthew G. Knepley   PetscErrorCode ierr;
4077cdb2a12SMatthew G. Knepley 
4087cdb2a12SMatthew G. Knepley   PetscFunctionBegin;
409c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)ts)->options,((PetscObject)ts)->prefix, "-dmts_check", &check);CHKERRQ(ierr);
4107cdb2a12SMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
41172fd7fc8SMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
4126fb91e28SMatthew G. Knepley   ierr = VecCopy(u, sol);CHKERRQ(ierr);
4136fb91e28SMatthew G. Knepley   ierr = TSSetSolution(ts, u);CHKERRQ(ierr);
4147cdb2a12SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
4157cdb2a12SMatthew G. Knepley   ierr = TSSetUp(ts);CHKERRQ(ierr);
4167cdb2a12SMatthew G. Knepley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
4176fb91e28SMatthew G. Knepley   ierr = SNESSetSolution(snes, u);CHKERRQ(ierr);
418f2cacb80SMatthew G. Knepley 
419f2cacb80SMatthew G. Knepley   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
420f2cacb80SMatthew G. Knepley   ierr = DMSNESCheckDiscretization(snes, dm, t, sol, -1.0, NULL);CHKERRQ(ierr);
421f2cacb80SMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &u_t);CHKERRQ(ierr);
422f2cacb80SMatthew G. Knepley   ierr = DMTSCheckResidual(ts, dm, t, sol, u_t, -1.0, NULL);CHKERRQ(ierr);
423f2cacb80SMatthew G. Knepley   ierr = DMTSCheckJacobian(ts, dm, t, sol, u_t, -1.0, NULL, NULL);CHKERRQ(ierr);
424f2cacb80SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &u_t);CHKERRQ(ierr);
425f2cacb80SMatthew G. Knepley 
42672fd7fc8SMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
4277cdb2a12SMatthew G. Knepley   PetscFunctionReturn(0);
4287cdb2a12SMatthew G. Knepley }
429