xref: /petsc/src/snes/utils/dmplexsnes.c (revision 800f99ff9e85495c69e9e5819c0be0dbd8cbc57c)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I "petscdmplex.h" I*/
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>     /*I "petscsnes.h"   I*/
324cdb843SMatthew G. Knepley #include <petscds.h>
4af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
5af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
6552f7358SJed Brown 
7649ef022SMatthew Knepley static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux,
8649ef022SMatthew Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
9649ef022SMatthew Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10649ef022SMatthew Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
11649ef022SMatthew Knepley {
12649ef022SMatthew Knepley   p[0] = u[uOff[1]];
13649ef022SMatthew Knepley }
14649ef022SMatthew Knepley 
15649ef022SMatthew Knepley /*
16649ef022SMatthew Knepley   SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately.
17649ef022SMatthew Knepley 
18649ef022SMatthew Knepley   Collective on SNES
19649ef022SMatthew Knepley 
20649ef022SMatthew Knepley   Input Parameters:
21649ef022SMatthew Knepley + snes      - The SNES
22649ef022SMatthew Knepley . pfield    - The field number for pressure
23649ef022SMatthew Knepley . nullspace - The pressure nullspace
24649ef022SMatthew Knepley . u         - The solution vector
25649ef022SMatthew Knepley - ctx       - An optional user context
26649ef022SMatthew Knepley 
27649ef022SMatthew Knepley   Output Parameter:
28649ef022SMatthew Knepley . u         - The solution with a continuum pressure integral of zero
29649ef022SMatthew Knepley 
30649ef022SMatthew Knepley   Notes:
31649ef022SMatthew Knepley   If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0. We assume that the nullspace is a single vector given explicitly.
32649ef022SMatthew Knepley 
33649ef022SMatthew Knepley   Level: developer
34649ef022SMatthew Knepley 
35db781477SPatrick Sanan .seealso: `SNESConvergedCorrectPressure()`
36649ef022SMatthew Knepley */
37649ef022SMatthew Knepley static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx)
38649ef022SMatthew Knepley {
39649ef022SMatthew Knepley   DM             dm;
40649ef022SMatthew Knepley   PetscDS        ds;
41649ef022SMatthew Knepley   const Vec     *nullvecs;
42649ef022SMatthew Knepley   PetscScalar    pintd, *intc, *intn;
43649ef022SMatthew Knepley   MPI_Comm       comm;
44649ef022SMatthew Knepley   PetscInt       Nf, Nv;
45649ef022SMatthew Knepley 
46649ef022SMatthew Knepley   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) snes, &comm));
489566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
4928b400f6SJacob Faibussowitsch   PetscCheck(dm,comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
5028b400f6SJacob Faibussowitsch   PetscCheck(nullspace,comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
519566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
529566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(ds, pfield, pressure_Private));
539566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs));
5463a3b9bcSJacob Faibussowitsch   PetscCheck(Nv == 1,comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %" PetscInt_FMT, Nv);
559566063dSJacob Faibussowitsch   PetscCall(VecDot(nullvecs[0], u, &pintd));
5608401ef6SPierre Jolivet   PetscCheck(PetscAbsScalar(pintd) <= PETSC_SMALL,comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g", (double) PetscRealPart(pintd));
579566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &intc, Nf, &intn));
599566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx));
609566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx));
619566063dSJacob Faibussowitsch   PetscCall(VecAXPY(u, -intc[pfield]/intn[pfield], nullvecs[0]));
62649ef022SMatthew Knepley #if defined (PETSC_USE_DEBUG)
639566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx));
6408401ef6SPierre Jolivet   PetscCheck(PetscAbsScalar(intc[pfield]) <= PETSC_SMALL,comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g", (double) PetscRealPart(intc[pfield]));
65649ef022SMatthew Knepley #endif
669566063dSJacob Faibussowitsch   PetscCall(PetscFree2(intc, intn));
67649ef022SMatthew Knepley   PetscFunctionReturn(0);
68649ef022SMatthew Knepley }
69649ef022SMatthew Knepley 
70649ef022SMatthew Knepley /*@C
71649ef022SMatthew Knepley    SNESConvergedCorrectPressure - Convergence test that adds a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately. The convergence test itself just mimics SNESConvergedDefault().
72649ef022SMatthew Knepley 
73649ef022SMatthew Knepley    Logically Collective on SNES
74649ef022SMatthew Knepley 
75649ef022SMatthew Knepley    Input Parameters:
76649ef022SMatthew Knepley +  snes - the SNES context
77649ef022SMatthew Knepley .  it - the iteration (0 indicates before any Newton steps)
78649ef022SMatthew Knepley .  xnorm - 2-norm of current iterate
79649ef022SMatthew Knepley .  snorm - 2-norm of current step
80649ef022SMatthew Knepley .  fnorm - 2-norm of function at current iterate
81649ef022SMatthew Knepley -  ctx   - Optional user context
82649ef022SMatthew Knepley 
83649ef022SMatthew Knepley    Output Parameter:
84649ef022SMatthew Knepley .  reason  - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
85649ef022SMatthew Knepley 
86649ef022SMatthew Knepley    Notes:
87f362779dSJed Brown    In order to use this convergence test, you must set up several PETSc structures. First fields must be added to the DM, and a PetscDS must be created with discretizations of those fields. We currently assume that the pressure field has index 1. The pressure field must have a nullspace, likely created using the DMSetNullSpaceConstructor() interface. Last we must be able to integrate the pressure over the domain, so the DM attached to the SNES must be a Plex at this time.
88f362779dSJed Brown 
89f362779dSJed Brown    Options Database Keys:
90f362779dSJed Brown .  -snes_convergence_test correct_pressure - see SNESSetFromOptions()
91649ef022SMatthew Knepley 
92649ef022SMatthew Knepley    Level: advanced
93649ef022SMatthew Knepley 
94db781477SPatrick Sanan .seealso: `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `DMSetNullSpaceConstructor()`
95649ef022SMatthew Knepley @*/
96649ef022SMatthew Knepley PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx)
97649ef022SMatthew Knepley {
98649ef022SMatthew Knepley   PetscBool      monitorIntegral = PETSC_FALSE;
99649ef022SMatthew Knepley 
100649ef022SMatthew Knepley   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx));
102649ef022SMatthew Knepley   if (monitorIntegral) {
103649ef022SMatthew Knepley     Mat          J;
104649ef022SMatthew Knepley     Vec          u;
105649ef022SMatthew Knepley     MatNullSpace nullspace;
106649ef022SMatthew Knepley     const Vec   *nullvecs;
107649ef022SMatthew Knepley     PetscScalar  pintd;
108649ef022SMatthew Knepley 
1099566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
1109566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
1119566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(J, &nullspace));
1129566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs));
1139566063dSJacob Faibussowitsch     PetscCall(VecDot(nullvecs[0], u, &pintd));
1149566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "SNES: Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd)));
115649ef022SMatthew Knepley   }
116649ef022SMatthew Knepley   if (*reason > 0) {
117649ef022SMatthew Knepley     Mat          J;
118649ef022SMatthew Knepley     Vec          u;
119649ef022SMatthew Knepley     MatNullSpace nullspace;
120649ef022SMatthew Knepley     PetscInt     pfield = 1;
121649ef022SMatthew Knepley 
1229566063dSJacob Faibussowitsch     PetscCall(SNESGetSolution(snes, &u));
1239566063dSJacob Faibussowitsch     PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL));
1249566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(J, &nullspace));
1259566063dSJacob Faibussowitsch     PetscCall(SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx));
126649ef022SMatthew Knepley   }
127649ef022SMatthew Knepley   PetscFunctionReturn(0);
128649ef022SMatthew Knepley }
129649ef022SMatthew Knepley 
13024cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
13124cdb843SMatthew G. Knepley 
1326da023fcSToby Isaac static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
1336da023fcSToby Isaac {
1346da023fcSToby Isaac   PetscBool      isPlex;
1356da023fcSToby Isaac 
1366da023fcSToby Isaac   PetscFunctionBegin;
1379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex));
1386da023fcSToby Isaac   if (isPlex) {
1396da023fcSToby Isaac     *plex = dm;
1409566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject) dm));
141f7148743SMatthew G. Knepley   } else {
1429566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex));
143f7148743SMatthew G. Knepley     if (!*plex) {
1449566063dSJacob Faibussowitsch       PetscCall(DMConvert(dm,DMPLEX,plex));
1459566063dSJacob Faibussowitsch       PetscCall(PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex));
1466da023fcSToby Isaac       if (copy) {
1479566063dSJacob Faibussowitsch         PetscCall(DMCopyDMSNES(dm, *plex));
1489566063dSJacob Faibussowitsch         PetscCall(DMCopyAuxiliaryVec(dm, *plex));
1496da023fcSToby Isaac       }
150f7148743SMatthew G. Knepley     } else {
1519566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject) *plex));
152f7148743SMatthew G. Knepley     }
1536da023fcSToby Isaac   }
1546da023fcSToby Isaac   PetscFunctionReturn(0);
1556da023fcSToby Isaac }
1566da023fcSToby Isaac 
1574267b1a3SMatthew G. Knepley /*@C
1584267b1a3SMatthew G. Knepley   DMInterpolationCreate - Creates a DMInterpolationInfo context
1594267b1a3SMatthew G. Knepley 
160d083f849SBarry Smith   Collective
1614267b1a3SMatthew G. Knepley 
1624267b1a3SMatthew G. Knepley   Input Parameter:
1634267b1a3SMatthew G. Knepley . comm - the communicator
1644267b1a3SMatthew G. Knepley 
1654267b1a3SMatthew G. Knepley   Output Parameter:
1664267b1a3SMatthew G. Knepley . ctx - the context
1674267b1a3SMatthew G. Knepley 
1684267b1a3SMatthew G. Knepley   Level: beginner
1694267b1a3SMatthew G. Knepley 
170db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationDestroy()`
1714267b1a3SMatthew G. Knepley @*/
1720adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
1730adebc6cSBarry Smith {
174552f7358SJed Brown   PetscFunctionBegin;
175552f7358SJed Brown   PetscValidPointer(ctx, 2);
1769566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
1771aa26658SKarl Rupp 
178552f7358SJed Brown   (*ctx)->comm   = comm;
179552f7358SJed Brown   (*ctx)->dim    = -1;
180552f7358SJed Brown   (*ctx)->nInput = 0;
1810298fd71SBarry Smith   (*ctx)->points = NULL;
1820298fd71SBarry Smith   (*ctx)->cells  = NULL;
183552f7358SJed Brown   (*ctx)->n      = -1;
1840298fd71SBarry Smith   (*ctx)->coords = NULL;
185552f7358SJed Brown   PetscFunctionReturn(0);
186552f7358SJed Brown }
187552f7358SJed Brown 
1884267b1a3SMatthew G. Knepley /*@C
1894267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
1904267b1a3SMatthew G. Knepley 
1914267b1a3SMatthew G. Knepley   Not collective
1924267b1a3SMatthew G. Knepley 
1934267b1a3SMatthew G. Knepley   Input Parameters:
1944267b1a3SMatthew G. Knepley + ctx - the context
1954267b1a3SMatthew G. Knepley - dim - the spatial dimension
1964267b1a3SMatthew G. Knepley 
1974267b1a3SMatthew G. Knepley   Level: intermediate
1984267b1a3SMatthew G. Knepley 
199db781477SPatrick Sanan .seealso: `DMInterpolationGetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2004267b1a3SMatthew G. Knepley @*/
2010adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
2020adebc6cSBarry Smith {
203552f7358SJed Brown   PetscFunctionBegin;
20463a3b9bcSJacob Faibussowitsch   PetscCheck(!(dim < 1) && !(dim > 3),ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %" PetscInt_FMT, dim);
205552f7358SJed Brown   ctx->dim = dim;
206552f7358SJed Brown   PetscFunctionReturn(0);
207552f7358SJed Brown }
208552f7358SJed Brown 
2094267b1a3SMatthew G. Knepley /*@C
2104267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
2114267b1a3SMatthew G. Knepley 
2124267b1a3SMatthew G. Knepley   Not collective
2134267b1a3SMatthew G. Knepley 
2144267b1a3SMatthew G. Knepley   Input Parameter:
2154267b1a3SMatthew G. Knepley . ctx - the context
2164267b1a3SMatthew G. Knepley 
2174267b1a3SMatthew G. Knepley   Output Parameter:
2184267b1a3SMatthew G. Knepley . dim - the spatial dimension
2194267b1a3SMatthew G. Knepley 
2204267b1a3SMatthew G. Knepley   Level: intermediate
2214267b1a3SMatthew G. Knepley 
222db781477SPatrick Sanan .seealso: `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2234267b1a3SMatthew G. Knepley @*/
2240adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
2250adebc6cSBarry Smith {
226552f7358SJed Brown   PetscFunctionBegin;
227552f7358SJed Brown   PetscValidIntPointer(dim, 2);
228552f7358SJed Brown   *dim = ctx->dim;
229552f7358SJed Brown   PetscFunctionReturn(0);
230552f7358SJed Brown }
231552f7358SJed Brown 
2324267b1a3SMatthew G. Knepley /*@C
2334267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
2344267b1a3SMatthew G. Knepley 
2354267b1a3SMatthew G. Knepley   Not collective
2364267b1a3SMatthew G. Knepley 
2374267b1a3SMatthew G. Knepley   Input Parameters:
2384267b1a3SMatthew G. Knepley + ctx - the context
2394267b1a3SMatthew G. Knepley - dof - the number of fields
2404267b1a3SMatthew G. Knepley 
2414267b1a3SMatthew G. Knepley   Level: intermediate
2424267b1a3SMatthew G. Knepley 
243db781477SPatrick Sanan .seealso: `DMInterpolationGetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2444267b1a3SMatthew G. Knepley @*/
2450adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
2460adebc6cSBarry Smith {
247552f7358SJed Brown   PetscFunctionBegin;
24863a3b9bcSJacob Faibussowitsch   PetscCheck(dof >= 1,ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %" PetscInt_FMT, dof);
249552f7358SJed Brown   ctx->dof = dof;
250552f7358SJed Brown   PetscFunctionReturn(0);
251552f7358SJed Brown }
252552f7358SJed Brown 
2534267b1a3SMatthew G. Knepley /*@C
2544267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
2554267b1a3SMatthew G. Knepley 
2564267b1a3SMatthew G. Knepley   Not collective
2574267b1a3SMatthew G. Knepley 
2584267b1a3SMatthew G. Knepley   Input Parameter:
2594267b1a3SMatthew G. Knepley . ctx - the context
2604267b1a3SMatthew G. Knepley 
2614267b1a3SMatthew G. Knepley   Output Parameter:
2624267b1a3SMatthew G. Knepley . dof - the number of fields
2634267b1a3SMatthew G. Knepley 
2644267b1a3SMatthew G. Knepley   Level: intermediate
2654267b1a3SMatthew G. Knepley 
266db781477SPatrick Sanan .seealso: `DMInterpolationSetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
2674267b1a3SMatthew G. Knepley @*/
2680adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
2690adebc6cSBarry Smith {
270552f7358SJed Brown   PetscFunctionBegin;
271552f7358SJed Brown   PetscValidIntPointer(dof, 2);
272552f7358SJed Brown   *dof = ctx->dof;
273552f7358SJed Brown   PetscFunctionReturn(0);
274552f7358SJed Brown }
275552f7358SJed Brown 
2764267b1a3SMatthew G. Knepley /*@C
2774267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
2784267b1a3SMatthew G. Knepley 
2794267b1a3SMatthew G. Knepley   Not collective
2804267b1a3SMatthew G. Knepley 
2814267b1a3SMatthew G. Knepley   Input Parameters:
2824267b1a3SMatthew G. Knepley + ctx    - the context
2834267b1a3SMatthew G. Knepley . n      - the number of points
2844267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
2854267b1a3SMatthew G. Knepley 
2864267b1a3SMatthew G. Knepley   Note: The coordinate information is copied.
2874267b1a3SMatthew G. Knepley 
2884267b1a3SMatthew G. Knepley   Level: intermediate
2894267b1a3SMatthew G. Knepley 
290db781477SPatrick Sanan .seealso: `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationCreate()`
2914267b1a3SMatthew G. Knepley @*/
2920adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
2930adebc6cSBarry Smith {
294552f7358SJed Brown   PetscFunctionBegin;
29508401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0,ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
29628b400f6SJacob Faibussowitsch   PetscCheck(!ctx->points,ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
297552f7358SJed Brown   ctx->nInput = n;
2981aa26658SKarl Rupp 
2999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n*ctx->dim, &ctx->points));
3009566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ctx->points, points, n*ctx->dim));
301552f7358SJed Brown   PetscFunctionReturn(0);
302552f7358SJed Brown }
303552f7358SJed Brown 
3044267b1a3SMatthew G. Knepley /*@C
30552aa1562SMatthew G. Knepley   DMInterpolationSetUp - Compute spatial indices for point location during interpolation
3064267b1a3SMatthew G. Knepley 
3074267b1a3SMatthew G. Knepley   Collective on ctx
3084267b1a3SMatthew G. Knepley 
3094267b1a3SMatthew G. Knepley   Input Parameters:
3104267b1a3SMatthew G. Knepley + ctx - the context
3114267b1a3SMatthew G. Knepley . dm  - the DM for the function space used for interpolation
31252aa1562SMatthew G. Knepley . redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
3137deeda50SSatish Balay - ignoreOutsideDomain - If PETSC_TRUE, ignore points outside the domain, otherwise return an error
3144267b1a3SMatthew G. Knepley 
3154267b1a3SMatthew G. Knepley   Level: intermediate
3164267b1a3SMatthew G. Knepley 
317db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
3184267b1a3SMatthew G. Knepley @*/
31952aa1562SMatthew G. Knepley PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain)
3200adebc6cSBarry Smith {
321552f7358SJed Brown   MPI_Comm          comm = ctx->comm;
322552f7358SJed Brown   PetscScalar       *a;
323552f7358SJed Brown   PetscInt          p, q, i;
324552f7358SJed Brown   PetscMPIInt       rank, size;
325552f7358SJed Brown   Vec               pointVec;
3263a93e3b7SToby Isaac   PetscSF           cellSF;
327552f7358SJed Brown   PetscLayout       layout;
328552f7358SJed Brown   PetscReal         *globalPoints;
329cb313848SJed Brown   PetscScalar       *globalPointsScalar;
330552f7358SJed Brown   const PetscInt    *ranges;
331552f7358SJed Brown   PetscMPIInt       *counts, *displs;
3323a93e3b7SToby Isaac   const PetscSFNode *foundCells;
3333a93e3b7SToby Isaac   const PetscInt    *foundPoints;
334552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
3353a93e3b7SToby Isaac   PetscInt          n, N, numFound;
336552f7358SJed Brown 
33719436ca2SJed Brown   PetscFunctionBegin;
338064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
3399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
34108401ef6SPierre Jolivet   PetscCheck(ctx->dim >= 0,comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
34219436ca2SJed Brown   /* Locate points */
34319436ca2SJed Brown   n = ctx->nInput;
344552f7358SJed Brown   if (!redundantPoints) {
3459566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreate(comm, &layout));
3469566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetBlockSize(layout, 1));
3479566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetLocalSize(layout, n));
3489566063dSJacob Faibussowitsch     PetscCall(PetscLayoutSetUp(layout));
3499566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(layout, &N));
350552f7358SJed Brown     /* Communicate all points to all processes */
3519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs));
3529566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(layout, &ranges));
353552f7358SJed Brown     for (p = 0; p < size; ++p) {
354552f7358SJed Brown       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
355552f7358SJed Brown       displs[p] = ranges[p]*ctx->dim;
356552f7358SJed Brown     }
3579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm));
358552f7358SJed Brown   } else {
359552f7358SJed Brown     N = n;
360552f7358SJed Brown     globalPoints = ctx->points;
36138ea73c8SJed Brown     counts = displs = NULL;
36238ea73c8SJed Brown     layout = NULL;
363552f7358SJed Brown   }
364552f7358SJed Brown #if 0
3659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs));
36619436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
367552f7358SJed Brown #else
368cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
3699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N*ctx->dim,&globalPointsScalar));
37046b3086cSToby Isaac   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
371cb313848SJed Brown #else
372cb313848SJed Brown   globalPointsScalar = globalPoints;
373cb313848SJed Brown #endif
3749566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec));
3759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(N,&foundProcs,N,&globalProcs));
3767b5f2079SMatthew G. Knepley   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
3773a93e3b7SToby Isaac   cellSF = NULL;
3789566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF));
3799566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells));
380552f7358SJed Brown #endif
3813a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
3823a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
383552f7358SJed Brown   }
384552f7358SJed Brown   /* Let the lowest rank process own each point */
3851c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm));
386552f7358SJed Brown   ctx->n = 0;
387552f7358SJed Brown   for (p = 0; p < N; ++p) {
38852aa1562SMatthew G. Knepley     if (globalProcs[p] == size) {
38963a3b9bcSJacob Faibussowitsch       PetscCheck(ignoreOutsideDomain,comm, PETSC_ERR_PLIB, "Point %" PetscInt_FMT ": %g %g %g not located in mesh", p, (double)globalPoints[p*ctx->dim+0], (double)(ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0), (double)(ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0));
390f7d195e4SLawrence Mitchell       if (rank == 0) ++ctx->n;
39152aa1562SMatthew G. Knepley     } else if (globalProcs[p] == rank) ++ctx->n;
392552f7358SJed Brown   }
393552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
3949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ctx->n, &ctx->cells));
3959566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &ctx->coords));
3969566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE));
3979566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(ctx->coords, ctx->dim));
3989566063dSJacob Faibussowitsch   PetscCall(VecSetType(ctx->coords,VECSTANDARD));
3999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->coords, &a));
400552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
401552f7358SJed Brown     if (globalProcs[p] == rank) {
402552f7358SJed Brown       PetscInt d;
403552f7358SJed Brown 
4041aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
405f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
406f317b747SMatthew G. Knepley       ++q;
407552f7358SJed Brown     }
408dd400576SPatrick Sanan     if (globalProcs[p] == size && rank == 0) {
40952aa1562SMatthew G. Knepley       PetscInt d;
41052aa1562SMatthew G. Knepley 
41152aa1562SMatthew G. Knepley       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
41252aa1562SMatthew G. Knepley       ctx->cells[q] = -1;
41352aa1562SMatthew G. Knepley       ++q;
41452aa1562SMatthew G. Knepley     }
415552f7358SJed Brown   }
4169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->coords, &a));
417552f7358SJed Brown #if 0
4189566063dSJacob Faibussowitsch   PetscCall(PetscFree3(foundCells,foundProcs,globalProcs));
419552f7358SJed Brown #else
4209566063dSJacob Faibussowitsch   PetscCall(PetscFree2(foundProcs,globalProcs));
4219566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
4229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pointVec));
423552f7358SJed Brown #endif
4249566063dSJacob Faibussowitsch   if ((void*)globalPointsScalar != (void*)globalPoints) PetscCall(PetscFree(globalPointsScalar));
4259566063dSJacob Faibussowitsch   if (!redundantPoints) PetscCall(PetscFree3(globalPoints,counts,displs));
4269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&layout));
427552f7358SJed Brown   PetscFunctionReturn(0);
428552f7358SJed Brown }
429552f7358SJed Brown 
4304267b1a3SMatthew G. Knepley /*@C
4314267b1a3SMatthew G. Knepley   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
4324267b1a3SMatthew G. Knepley 
4334267b1a3SMatthew G. Knepley   Collective on ctx
4344267b1a3SMatthew G. Knepley 
4354267b1a3SMatthew G. Knepley   Input Parameter:
4364267b1a3SMatthew G. Knepley . ctx - the context
4374267b1a3SMatthew G. Knepley 
4384267b1a3SMatthew G. Knepley   Output Parameter:
4394267b1a3SMatthew G. Knepley . coordinates  - the coordinates of interpolation points
4404267b1a3SMatthew G. Knepley 
4414267b1a3SMatthew G. Knepley   Note: The local vector entries correspond to interpolation points lying on this process, according to the associated DM. This is a borrowed vector that the user should not destroy.
4424267b1a3SMatthew G. Knepley 
4434267b1a3SMatthew G. Knepley   Level: intermediate
4444267b1a3SMatthew G. Knepley 
445db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4464267b1a3SMatthew G. Knepley @*/
4470adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
4480adebc6cSBarry Smith {
449552f7358SJed Brown   PetscFunctionBegin;
450552f7358SJed Brown   PetscValidPointer(coordinates, 2);
45128b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords,ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
452552f7358SJed Brown   *coordinates = ctx->coords;
453552f7358SJed Brown   PetscFunctionReturn(0);
454552f7358SJed Brown }
455552f7358SJed Brown 
4564267b1a3SMatthew G. Knepley /*@C
4574267b1a3SMatthew G. Knepley   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
4584267b1a3SMatthew G. Knepley 
4594267b1a3SMatthew G. Knepley   Collective on ctx
4604267b1a3SMatthew G. Knepley 
4614267b1a3SMatthew G. Knepley   Input Parameter:
4624267b1a3SMatthew G. Knepley . ctx - the context
4634267b1a3SMatthew G. Knepley 
4644267b1a3SMatthew G. Knepley   Output Parameter:
4654267b1a3SMatthew G. Knepley . v  - a vector capable of holding the interpolated field values
4664267b1a3SMatthew G. Knepley 
4674267b1a3SMatthew G. Knepley   Note: This vector should be returned using DMInterpolationRestoreVector().
4684267b1a3SMatthew G. Knepley 
4694267b1a3SMatthew G. Knepley   Level: intermediate
4704267b1a3SMatthew G. Knepley 
471db781477SPatrick Sanan .seealso: `DMInterpolationRestoreVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4724267b1a3SMatthew G. Knepley @*/
4730adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
4740adebc6cSBarry Smith {
475552f7358SJed Brown   PetscFunctionBegin;
476552f7358SJed Brown   PetscValidPointer(v, 2);
47728b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords,ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
4789566063dSJacob Faibussowitsch   PetscCall(VecCreate(ctx->comm, v));
4799566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE));
4809566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*v, ctx->dof));
4819566063dSJacob Faibussowitsch   PetscCall(VecSetType(*v,VECSTANDARD));
482552f7358SJed Brown   PetscFunctionReturn(0);
483552f7358SJed Brown }
484552f7358SJed Brown 
4854267b1a3SMatthew G. Knepley /*@C
4864267b1a3SMatthew G. Knepley   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
4874267b1a3SMatthew G. Knepley 
4884267b1a3SMatthew G. Knepley   Collective on ctx
4894267b1a3SMatthew G. Knepley 
4904267b1a3SMatthew G. Knepley   Input Parameters:
4914267b1a3SMatthew G. Knepley + ctx - the context
4924267b1a3SMatthew G. Knepley - v  - a vector capable of holding the interpolated field values
4934267b1a3SMatthew G. Knepley 
4944267b1a3SMatthew G. Knepley   Level: intermediate
4954267b1a3SMatthew G. Knepley 
496db781477SPatrick Sanan .seealso: `DMInterpolationGetVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
4974267b1a3SMatthew G. Knepley @*/
4980adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
4990adebc6cSBarry Smith {
500552f7358SJed Brown   PetscFunctionBegin;
501552f7358SJed Brown   PetscValidPointer(v, 2);
50228b400f6SJacob Faibussowitsch   PetscCheck(ctx->coords,ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
5039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(v));
504552f7358SJed Brown   PetscFunctionReturn(0);
505552f7358SJed Brown }
506552f7358SJed Brown 
507cfe5744fSMatthew G. Knepley static inline PetscErrorCode DMInterpolate_Segment_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
508cfe5744fSMatthew G. Knepley {
509cfe5744fSMatthew G. Knepley   PetscReal          v0, J, invJ, detJ;
510cfe5744fSMatthew G. Knepley   const PetscInt     dof = ctx->dof;
511cfe5744fSMatthew G. Knepley   const PetscScalar *coords;
512cfe5744fSMatthew G. Knepley   PetscScalar       *a;
513cfe5744fSMatthew G. Knepley   PetscInt           p;
514cfe5744fSMatthew G. Knepley 
515cfe5744fSMatthew G. Knepley   PetscFunctionBegin;
5169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
518cfe5744fSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
519cfe5744fSMatthew G. Knepley     PetscInt     c = ctx->cells[p];
520cfe5744fSMatthew G. Knepley     PetscScalar *x = NULL;
521cfe5744fSMatthew G. Knepley     PetscReal    xir[1];
522cfe5744fSMatthew G. Knepley     PetscInt     xSize, comp;
523cfe5744fSMatthew G. Knepley 
5249566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ));
525cfe5744fSMatthew G. Knepley     PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double) detJ, c);
526cfe5744fSMatthew G. Knepley     xir[0] = invJ*PetscRealPart(coords[p] - v0);
5279566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
528cfe5744fSMatthew G. Knepley     if (2*dof == xSize) {
529cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) a[p*dof+comp] = x[0*dof+comp]*(1 - xir[0]) + x[1*dof+comp]*xir[0];
530cfe5744fSMatthew G. Knepley     } else if (dof == xSize) {
531cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) a[p*dof+comp] = x[0*dof+comp];
532cfe5744fSMatthew G. Knepley     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input closure size %" PetscInt_FMT " must be either %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 2*dof, dof);
5339566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
534cfe5744fSMatthew G. Knepley   }
5359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
537cfe5744fSMatthew G. Knepley   PetscFunctionReturn(0);
538cfe5744fSMatthew G. Knepley }
539cfe5744fSMatthew G. Knepley 
5409fbee547SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
541a6dfd86eSKarl Rupp {
542552f7358SJed Brown   PetscReal      *v0, *J, *invJ, detJ;
54356044e6dSMatthew G. Knepley   const PetscScalar *coords;
54456044e6dSMatthew G. Knepley   PetscScalar    *a;
545552f7358SJed Brown   PetscInt       p;
546552f7358SJed Brown 
547552f7358SJed Brown   PetscFunctionBegin;
5489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ));
5499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
551552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
552552f7358SJed Brown     PetscInt     c = ctx->cells[p];
553a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
554552f7358SJed Brown     PetscReal    xi[4];
555552f7358SJed Brown     PetscInt     d, f, comp;
556552f7358SJed Brown 
5579566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
55863a3b9bcSJacob Faibussowitsch     PetscCheck(detJ > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5599566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
5601aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5611aa26658SKarl Rupp 
562552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
563552f7358SJed Brown       xi[d] = 0.0;
5641aa26658SKarl Rupp       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
5651aa26658SKarl Rupp       for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
566552f7358SJed Brown     }
5679566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
568552f7358SJed Brown   }
5699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
5709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
5719566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
572552f7358SJed Brown   PetscFunctionReturn(0);
573552f7358SJed Brown }
574552f7358SJed Brown 
5759fbee547SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
5767a1931ceSMatthew G. Knepley {
5777a1931ceSMatthew G. Knepley   PetscReal      *v0, *J, *invJ, detJ;
57856044e6dSMatthew G. Knepley   const PetscScalar *coords;
57956044e6dSMatthew G. Knepley   PetscScalar    *a;
5807a1931ceSMatthew G. Knepley   PetscInt       p;
5817a1931ceSMatthew G. Knepley 
5827a1931ceSMatthew G. Knepley   PetscFunctionBegin;
5839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ));
5849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
5859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
5867a1931ceSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
5877a1931ceSMatthew G. Knepley     PetscInt       c = ctx->cells[p];
5887a1931ceSMatthew G. Knepley     const PetscInt order[3] = {2, 1, 3};
5892584bbe8SMatthew G. Knepley     PetscScalar   *x = NULL;
5907a1931ceSMatthew G. Knepley     PetscReal      xi[4];
5917a1931ceSMatthew G. Knepley     PetscInt       d, f, comp;
5927a1931ceSMatthew G. Knepley 
5939566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
59463a3b9bcSJacob Faibussowitsch     PetscCheck(detJ > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c);
5959566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x));
5967a1931ceSMatthew G. Knepley     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
5977a1931ceSMatthew G. Knepley 
5987a1931ceSMatthew G. Knepley     for (d = 0; d < ctx->dim; ++d) {
5997a1931ceSMatthew G. Knepley       xi[d] = 0.0;
6007a1931ceSMatthew G. Knepley       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
6017a1931ceSMatthew G. Knepley       for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[order[d]*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
6027a1931ceSMatthew G. Knepley     }
6039566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x));
6047a1931ceSMatthew G. Knepley   }
6059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
6069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
6079566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v0, J, invJ));
6087a1931ceSMatthew G. Knepley   PetscFunctionReturn(0);
6097a1931ceSMatthew G. Knepley }
6107a1931ceSMatthew G. Knepley 
6119fbee547SJacob Faibussowitsch static inline PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
612552f7358SJed Brown {
613552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
614552f7358SJed Brown   const PetscScalar x0        = vertices[0];
615552f7358SJed Brown   const PetscScalar y0        = vertices[1];
616552f7358SJed Brown   const PetscScalar x1        = vertices[2];
617552f7358SJed Brown   const PetscScalar y1        = vertices[3];
618552f7358SJed Brown   const PetscScalar x2        = vertices[4];
619552f7358SJed Brown   const PetscScalar y2        = vertices[5];
620552f7358SJed Brown   const PetscScalar x3        = vertices[6];
621552f7358SJed Brown   const PetscScalar y3        = vertices[7];
622552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
623552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
624552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
625552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
626552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
627552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
62856044e6dSMatthew G. Knepley   const PetscScalar *ref;
62956044e6dSMatthew G. Knepley   PetscScalar       *real;
630552f7358SJed Brown 
631552f7358SJed Brown   PetscFunctionBegin;
6329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref,  &ref));
6339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
634552f7358SJed Brown   {
635552f7358SJed Brown     const PetscScalar p0 = ref[0];
636552f7358SJed Brown     const PetscScalar p1 = ref[1];
637552f7358SJed Brown 
638552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
639552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
640552f7358SJed Brown   }
6419566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(28));
6429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref,  &ref));
6439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
644552f7358SJed Brown   PetscFunctionReturn(0);
645552f7358SJed Brown }
646552f7358SJed Brown 
647af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
6489fbee547SJacob Faibussowitsch static inline PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
649552f7358SJed Brown {
650552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
651552f7358SJed Brown   const PetscScalar x0        = vertices[0];
652552f7358SJed Brown   const PetscScalar y0        = vertices[1];
653552f7358SJed Brown   const PetscScalar x1        = vertices[2];
654552f7358SJed Brown   const PetscScalar y1        = vertices[3];
655552f7358SJed Brown   const PetscScalar x2        = vertices[4];
656552f7358SJed Brown   const PetscScalar y2        = vertices[5];
657552f7358SJed Brown   const PetscScalar x3        = vertices[6];
658552f7358SJed Brown   const PetscScalar y3        = vertices[7];
659552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
660552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
66156044e6dSMatthew G. Knepley   const PetscScalar *ref;
662552f7358SJed Brown 
663552f7358SJed Brown   PetscFunctionBegin;
6649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref,  &ref));
665552f7358SJed Brown   {
666552f7358SJed Brown     const PetscScalar x       = ref[0];
667552f7358SJed Brown     const PetscScalar y       = ref[1];
668552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
669da80777bSKarl Rupp     PetscScalar       values[4];
670da80777bSKarl Rupp 
671da80777bSKarl Rupp     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
672da80777bSKarl Rupp     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
6739566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES));
674552f7358SJed Brown   }
6759566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(30));
6769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref,  &ref));
6779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
6789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
679552f7358SJed Brown   PetscFunctionReturn(0);
680552f7358SJed Brown }
681552f7358SJed Brown 
6829fbee547SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
683a6dfd86eSKarl Rupp {
684fafc0619SMatthew G Knepley   DM                 dmCoord;
685de73a395SMatthew G. Knepley   PetscFE            fem = NULL;
686552f7358SJed Brown   SNES               snes;
687552f7358SJed Brown   KSP                ksp;
688552f7358SJed Brown   PC                 pc;
689552f7358SJed Brown   Vec                coordsLocal, r, ref, real;
690552f7358SJed Brown   Mat                J;
691716009bfSMatthew G. Knepley   PetscTabulation    T = NULL;
69256044e6dSMatthew G. Knepley   const PetscScalar *coords;
69356044e6dSMatthew G. Knepley   PetscScalar        *a;
694716009bfSMatthew G. Knepley   PetscReal          xir[2] = {0., 0.};
695de73a395SMatthew G. Knepley   PetscInt           Nf, p;
6965509d985SMatthew G. Knepley   const PetscInt     dof = ctx->dof;
697552f7358SJed Brown 
698552f7358SJed Brown   PetscFunctionBegin;
6999566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
700716009bfSMatthew G. Knepley   if (Nf) {
701cfe5744fSMatthew G. Knepley     PetscObject  obj;
702cfe5744fSMatthew G. Knepley     PetscClassId id;
703cfe5744fSMatthew G. Knepley 
7049566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, 0, NULL, &obj));
7059566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
7069566063dSJacob Faibussowitsch     if (id == PETSCFE_CLASSID) {fem = (PetscFE) obj; PetscCall(PetscFECreateTabulation(fem, 1, 1, xir, 0, &T));}
707716009bfSMatthew G. Knepley   }
7089566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
7099566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
7109566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
7119566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "quad_interp_"));
7129566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
7139566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 2, 2));
7149566063dSJacob Faibussowitsch   PetscCall(VecSetType(r,dm->vectype));
7159566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
7169566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
7179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
7189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 2, 2, 2, 2));
7199566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
7209566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
7219566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, QuadMap_Private, NULL));
7229566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL));
7239566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
7249566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
7259566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
7269566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
727552f7358SJed Brown 
7289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
7299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
730552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
731a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
732552f7358SJed Brown     PetscScalar *xi;
733552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
734552f7358SJed Brown 
735552f7358SJed Brown     /* Can make this do all points at once */
7369566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
73763a3b9bcSJacob Faibussowitsch     PetscCheck(4*2 == coordSize,ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %" PetscInt_FMT " should be %d", coordSize, 4*2);
7389566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
7399566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
7409566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
7419566063dSJacob Faibussowitsch     PetscCall(VecGetArray(real, &xi));
742552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
743552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
7449566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(real, &xi));
7459566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, real, ref));
7469566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ref, &xi));
747cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
748cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
749cfe5744fSMatthew G. Knepley     if (4*dof == xSize) {
750cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp)
751cfe5744fSMatthew G. Knepley         a[p*dof+comp] = x[0*dof+comp]*(1 - xir[0])*(1 - xir[1]) + x[1*dof+comp]*xir[0]*(1 - xir[1]) + x[2*dof+comp]*xir[0]*xir[1] + x[3*dof+comp]*(1 - xir[0])*xir[1];
752cfe5744fSMatthew G. Knepley     } else if (dof == xSize) {
753cfe5744fSMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) a[p*dof+comp] = x[0*dof+comp];
754cfe5744fSMatthew G. Knepley     } else {
7555509d985SMatthew G. Knepley       PetscInt d;
7561aa26658SKarl Rupp 
7572c71b3e2SJacob Faibussowitsch       PetscCheck(fem, ctx->comm, PETSC_ERR_ARG_WRONG, "Cannot have a higher order interpolant if the discretization is not PetscFE");
7585509d985SMatthew G. Knepley       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
7599566063dSJacob Faibussowitsch       PetscCall(PetscFEComputeTabulation(fem, 1, xir, 0, T));
7605509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) {
7615509d985SMatthew G. Knepley         a[p*dof+comp] = 0.0;
7625509d985SMatthew G. Knepley         for (d = 0; d < xSize/dof; ++d) {
763ef0bb6c7SMatthew G. Knepley           a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
7645509d985SMatthew G. Knepley         }
7655509d985SMatthew G. Knepley       }
7665509d985SMatthew G. Knepley     }
7679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ref, &xi));
7689566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
7699566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
770552f7358SJed Brown   }
7719566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&T));
7729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
7739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
774552f7358SJed Brown 
7759566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
7769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
7779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
7789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
7799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
780552f7358SJed Brown   PetscFunctionReturn(0);
781552f7358SJed Brown }
782552f7358SJed Brown 
7839fbee547SJacob Faibussowitsch static inline PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
784552f7358SJed Brown {
785552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
786552f7358SJed Brown   const PetscScalar x0        = vertices[0];
787552f7358SJed Brown   const PetscScalar y0        = vertices[1];
788552f7358SJed Brown   const PetscScalar z0        = vertices[2];
7897a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
7907a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
7917a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
792552f7358SJed Brown   const PetscScalar x2        = vertices[6];
793552f7358SJed Brown   const PetscScalar y2        = vertices[7];
794552f7358SJed Brown   const PetscScalar z2        = vertices[8];
7957a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
7967a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
7977a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
798552f7358SJed Brown   const PetscScalar x4        = vertices[12];
799552f7358SJed Brown   const PetscScalar y4        = vertices[13];
800552f7358SJed Brown   const PetscScalar z4        = vertices[14];
801552f7358SJed Brown   const PetscScalar x5        = vertices[15];
802552f7358SJed Brown   const PetscScalar y5        = vertices[16];
803552f7358SJed Brown   const PetscScalar z5        = vertices[17];
804552f7358SJed Brown   const PetscScalar x6        = vertices[18];
805552f7358SJed Brown   const PetscScalar y6        = vertices[19];
806552f7358SJed Brown   const PetscScalar z6        = vertices[20];
807552f7358SJed Brown   const PetscScalar x7        = vertices[21];
808552f7358SJed Brown   const PetscScalar y7        = vertices[22];
809552f7358SJed Brown   const PetscScalar z7        = vertices[23];
810552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
811552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
812552f7358SJed Brown   const PetscScalar h_1       = z1 - z0;
813552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
814552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
815552f7358SJed Brown   const PetscScalar h_3       = z3 - z0;
816552f7358SJed Brown   const PetscScalar f_4       = x4 - x0;
817552f7358SJed Brown   const PetscScalar g_4       = y4 - y0;
818552f7358SJed Brown   const PetscScalar h_4       = z4 - z0;
819552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
820552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
821552f7358SJed Brown   const PetscScalar h_01      = z2 - z1 - z3 + z0;
822552f7358SJed Brown   const PetscScalar f_12      = x7 - x3 - x4 + x0;
823552f7358SJed Brown   const PetscScalar g_12      = y7 - y3 - y4 + y0;
824552f7358SJed Brown   const PetscScalar h_12      = z7 - z3 - z4 + z0;
825552f7358SJed Brown   const PetscScalar f_02      = x5 - x1 - x4 + x0;
826552f7358SJed Brown   const PetscScalar g_02      = y5 - y1 - y4 + y0;
827552f7358SJed Brown   const PetscScalar h_02      = z5 - z1 - z4 + z0;
828552f7358SJed Brown   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
829552f7358SJed Brown   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
830552f7358SJed Brown   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
83156044e6dSMatthew G. Knepley   const PetscScalar *ref;
83256044e6dSMatthew G. Knepley   PetscScalar       *real;
833552f7358SJed Brown 
834552f7358SJed Brown   PetscFunctionBegin;
8359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref,  &ref));
8369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xreal, &real));
837552f7358SJed Brown   {
838552f7358SJed Brown     const PetscScalar p0 = ref[0];
839552f7358SJed Brown     const PetscScalar p1 = ref[1];
840552f7358SJed Brown     const PetscScalar p2 = ref[2];
841552f7358SJed Brown 
842552f7358SJed Brown     real[0] = x0 + f_1*p0 + f_3*p1 + f_4*p2 + f_01*p0*p1 + f_12*p1*p2 + f_02*p0*p2 + f_012*p0*p1*p2;
843552f7358SJed Brown     real[1] = y0 + g_1*p0 + g_3*p1 + g_4*p2 + g_01*p0*p1 + g_01*p0*p1 + g_12*p1*p2 + g_02*p0*p2 + g_012*p0*p1*p2;
844552f7358SJed Brown     real[2] = z0 + h_1*p0 + h_3*p1 + h_4*p2 + h_01*p0*p1 + h_01*p0*p1 + h_12*p1*p2 + h_02*p0*p2 + h_012*p0*p1*p2;
845552f7358SJed Brown   }
8469566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(114));
8479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref,  &ref));
8489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xreal, &real));
849552f7358SJed Brown   PetscFunctionReturn(0);
850552f7358SJed Brown }
851552f7358SJed Brown 
8529fbee547SJacob Faibussowitsch static inline PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
853552f7358SJed Brown {
854552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
855552f7358SJed Brown   const PetscScalar x0        = vertices[0];
856552f7358SJed Brown   const PetscScalar y0        = vertices[1];
857552f7358SJed Brown   const PetscScalar z0        = vertices[2];
8587a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
8597a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
8607a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
861552f7358SJed Brown   const PetscScalar x2        = vertices[6];
862552f7358SJed Brown   const PetscScalar y2        = vertices[7];
863552f7358SJed Brown   const PetscScalar z2        = vertices[8];
8647a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
8657a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
8667a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
867552f7358SJed Brown   const PetscScalar x4        = vertices[12];
868552f7358SJed Brown   const PetscScalar y4        = vertices[13];
869552f7358SJed Brown   const PetscScalar z4        = vertices[14];
870552f7358SJed Brown   const PetscScalar x5        = vertices[15];
871552f7358SJed Brown   const PetscScalar y5        = vertices[16];
872552f7358SJed Brown   const PetscScalar z5        = vertices[17];
873552f7358SJed Brown   const PetscScalar x6        = vertices[18];
874552f7358SJed Brown   const PetscScalar y6        = vertices[19];
875552f7358SJed Brown   const PetscScalar z6        = vertices[20];
876552f7358SJed Brown   const PetscScalar x7        = vertices[21];
877552f7358SJed Brown   const PetscScalar y7        = vertices[22];
878552f7358SJed Brown   const PetscScalar z7        = vertices[23];
879552f7358SJed Brown   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
880552f7358SJed Brown   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
881552f7358SJed Brown   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
882552f7358SJed Brown   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
883552f7358SJed Brown   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
884552f7358SJed Brown   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
885552f7358SJed Brown   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
886552f7358SJed Brown   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
887552f7358SJed Brown   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
888552f7358SJed Brown   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
889552f7358SJed Brown   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
890552f7358SJed Brown   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
89156044e6dSMatthew G. Knepley   const PetscScalar *ref;
892552f7358SJed Brown 
893552f7358SJed Brown   PetscFunctionBegin;
8949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xref,  &ref));
895552f7358SJed Brown   {
896552f7358SJed Brown     const PetscScalar x       = ref[0];
897552f7358SJed Brown     const PetscScalar y       = ref[1];
898552f7358SJed Brown     const PetscScalar z       = ref[2];
899552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
900da80777bSKarl Rupp     PetscScalar       values[9];
901da80777bSKarl Rupp 
902da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
903da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
904da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
905da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
906da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
907da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
908da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
909da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
910da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
9111aa26658SKarl Rupp 
9129566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES));
913552f7358SJed Brown   }
9149566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(152));
9159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xref,  &ref));
9169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
9179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
918552f7358SJed Brown   PetscFunctionReturn(0);
919552f7358SJed Brown }
920552f7358SJed Brown 
9219fbee547SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
922a6dfd86eSKarl Rupp {
923fafc0619SMatthew G Knepley   DM             dmCoord;
924552f7358SJed Brown   SNES           snes;
925552f7358SJed Brown   KSP            ksp;
926552f7358SJed Brown   PC             pc;
927552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
928552f7358SJed Brown   Mat            J;
92956044e6dSMatthew G. Knepley   const PetscScalar *coords;
93056044e6dSMatthew G. Knepley   PetscScalar    *a;
931552f7358SJed Brown   PetscInt       p;
932552f7358SJed Brown 
933552f7358SJed Brown   PetscFunctionBegin;
9349566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
9359566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &dmCoord));
9369566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
9379566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes, "hex_interp_"));
9389566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &r));
9399566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(r, 3, 3));
9409566063dSJacob Faibussowitsch   PetscCall(VecSetType(r,dm->vectype));
9419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &ref));
9429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(r, &real));
9439566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &J));
9449566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, 3, 3, 3, 3));
9459566063dSJacob Faibussowitsch   PetscCall(MatSetType(J, MATSEQDENSE));
9469566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
9479566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, r, HexMap_Private, NULL));
9489566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL));
9499566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
9509566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
9519566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCLU));
9529566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
953552f7358SJed Brown 
9549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->coords, &coords));
9559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(v, &a));
956552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
957a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
958552f7358SJed Brown     PetscScalar *xi;
959cb313848SJed Brown     PetscReal    xir[3];
960552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
961552f7358SJed Brown 
962552f7358SJed Brown     /* Can make this do all points at once */
9639566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
964cfe5744fSMatthew G. Knepley     PetscCheck(8*3 == coordSize,ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid coordinate closure size %" PetscInt_FMT " should be %d", coordSize, 8*3);
9659566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x));
966cfe5744fSMatthew G. Knepley     PetscCheck((8*ctx->dof == xSize) || (ctx->dof == xSize),ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input closure size %" PetscInt_FMT " should be %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 8*ctx->dof, ctx->dof);
9679566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, NULL, NULL, vertices));
9689566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices));
9699566063dSJacob Faibussowitsch     PetscCall(VecGetArray(real, &xi));
970552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
971552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
972552f7358SJed Brown     xi[2]  = coords[p*ctx->dim+2];
9739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(real, &xi));
9749566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, real, ref));
9759566063dSJacob Faibussowitsch     PetscCall(VecGetArray(ref, &xi));
976cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
977cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
978cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
979cfe5744fSMatthew G. Knepley     if (8*ctx->dof == xSize) {
980552f7358SJed Brown       for (comp = 0; comp < ctx->dof; ++comp) {
981552f7358SJed Brown         a[p*ctx->dof+comp] =
982cb313848SJed Brown           x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
9837a1931ceSMatthew G. Knepley           x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
984cb313848SJed Brown           x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
9857a1931ceSMatthew G. Knepley           x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
986cb313848SJed Brown           x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
987cb313848SJed Brown           x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
988cb313848SJed Brown           x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
989cb313848SJed Brown           x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
990552f7358SJed Brown       }
991cfe5744fSMatthew G. Knepley     } else {
992cfe5744fSMatthew G. Knepley       for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
993cfe5744fSMatthew G. Knepley     }
9949566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(ref, &xi));
9959566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices));
9969566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x));
997552f7358SJed Brown   }
9989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(v, &a));
9999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
1000552f7358SJed Brown 
10019566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
10029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
10039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
10049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&real));
10059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1006552f7358SJed Brown   PetscFunctionReturn(0);
1007552f7358SJed Brown }
1008552f7358SJed Brown 
10094267b1a3SMatthew G. Knepley /*@C
10104267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
10114267b1a3SMatthew G. Knepley 
1012552f7358SJed Brown   Input Parameters:
1013552f7358SJed Brown + ctx - The DMInterpolationInfo context
1014552f7358SJed Brown . dm  - The DM
1015552f7358SJed Brown - x   - The local vector containing the field to be interpolated
1016552f7358SJed Brown 
1017552f7358SJed Brown   Output Parameters:
1018552f7358SJed Brown . v   - The vector containing the interpolated values
10194267b1a3SMatthew G. Knepley 
10204267b1a3SMatthew G. Knepley   Note: A suitable v can be obtained using DMInterpolationGetVector().
10214267b1a3SMatthew G. Knepley 
10224267b1a3SMatthew G. Knepley   Level: beginner
10234267b1a3SMatthew G. Knepley 
1024db781477SPatrick Sanan .seealso: `DMInterpolationGetVector()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
10254267b1a3SMatthew G. Knepley @*/
10260adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
10270adebc6cSBarry Smith {
102866a0a883SMatthew G. Knepley   PetscDS        ds;
102966a0a883SMatthew G. Knepley   PetscInt       n, p, Nf, field;
103066a0a883SMatthew G. Knepley   PetscBool      useDS = PETSC_FALSE;
1031552f7358SJed Brown 
1032552f7358SJed Brown   PetscFunctionBegin;
1033552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1034552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
1035552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
10369566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
103763a3b9bcSJacob Faibussowitsch   PetscCheck(n == ctx->n*ctx->dof,ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %" PetscInt_FMT " should be %" PetscInt_FMT, n, ctx->n*ctx->dof);
103866a0a883SMatthew G. Knepley   if (!n) PetscFunctionReturn(0);
10399566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1040680d18d5SMatthew G. Knepley   if (ds) {
104166a0a883SMatthew G. Knepley     useDS = PETSC_TRUE;
10429566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(ds, &Nf));
1043680d18d5SMatthew G. Knepley     for (field = 0; field < Nf; ++field) {
104466a0a883SMatthew G. Knepley       PetscObject  obj;
1045680d18d5SMatthew G. Knepley       PetscClassId id;
1046680d18d5SMatthew G. Knepley 
10479566063dSJacob Faibussowitsch       PetscCall(PetscDSGetDiscretization(ds, field, &obj));
10489566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetClassId(obj, &id));
104966a0a883SMatthew G. Knepley       if (id != PETSCFE_CLASSID) {useDS = PETSC_FALSE; break;}
105066a0a883SMatthew G. Knepley     }
105166a0a883SMatthew G. Knepley   }
105266a0a883SMatthew G. Knepley   if (useDS) {
105366a0a883SMatthew G. Knepley     const PetscScalar *coords;
105466a0a883SMatthew G. Knepley     PetscScalar       *interpolant;
105566a0a883SMatthew G. Knepley     PetscInt           cdim, d;
105666a0a883SMatthew G. Knepley 
10579566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(dm, &cdim));
10589566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ctx->coords, &coords));
10599566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(v, &interpolant));
1060680d18d5SMatthew G. Knepley     for (p = 0; p < ctx->n; ++p) {
106166a0a883SMatthew G. Knepley       PetscReal    pcoords[3], xi[3];
1062680d18d5SMatthew G. Knepley       PetscScalar *xa   = NULL;
106366a0a883SMatthew G. Knepley       PetscInt     coff = 0, foff = 0, clSize;
1064680d18d5SMatthew G. Knepley 
106552aa1562SMatthew G. Knepley       if (ctx->cells[p] < 0) continue;
1066680d18d5SMatthew G. Knepley       for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p*cdim+d]);
10679566063dSJacob Faibussowitsch       PetscCall(DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi));
10689566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
106966a0a883SMatthew G. Knepley       for (field = 0; field < Nf; ++field) {
107066a0a883SMatthew G. Knepley         PetscTabulation T;
107166a0a883SMatthew G. Knepley         PetscFE         fe;
107266a0a883SMatthew G. Knepley 
10739566063dSJacob Faibussowitsch         PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *) &fe));
10749566063dSJacob Faibussowitsch         PetscCall(PetscFECreateTabulation(fe, 1, 1, xi, 0, &T));
1075680d18d5SMatthew G. Knepley         {
1076680d18d5SMatthew G. Knepley           const PetscReal *basis = T->T[0];
1077680d18d5SMatthew G. Knepley           const PetscInt   Nb    = T->Nb;
1078680d18d5SMatthew G. Knepley           const PetscInt   Nc    = T->Nc;
107966a0a883SMatthew G. Knepley           PetscInt         f, fc;
108066a0a883SMatthew G. Knepley 
1081680d18d5SMatthew G. Knepley           for (fc = 0; fc < Nc; ++fc) {
108266a0a883SMatthew G. Knepley             interpolant[p*ctx->dof+coff+fc] = 0.0;
1083680d18d5SMatthew G. Knepley             for (f = 0; f < Nb; ++f) {
108466a0a883SMatthew G. Knepley               interpolant[p*ctx->dof+coff+fc] += xa[foff+f]*basis[(0*Nb + f)*Nc + fc];
1085552f7358SJed Brown             }
1086680d18d5SMatthew G. Knepley           }
108766a0a883SMatthew G. Knepley           coff += Nc;
108866a0a883SMatthew G. Knepley           foff += Nb;
1089680d18d5SMatthew G. Knepley         }
10909566063dSJacob Faibussowitsch         PetscCall(PetscTabulationDestroy(&T));
1091680d18d5SMatthew G. Knepley       }
10929566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa));
109363a3b9bcSJacob Faibussowitsch       PetscCheck(coff == ctx->dof,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %" PetscInt_FMT " != %" PetscInt_FMT " dof specified for interpolation", coff, ctx->dof);
109463a3b9bcSJacob Faibussowitsch       PetscCheck(foff == clSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE space size %" PetscInt_FMT " != %" PetscInt_FMT " closure size", foff, clSize);
109566a0a883SMatthew G. Knepley     }
10969566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ctx->coords, &coords));
10979566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(v, &interpolant));
109866a0a883SMatthew G. Knepley   } else {
109966a0a883SMatthew G. Knepley     DMPolytopeType ct;
110066a0a883SMatthew G. Knepley 
1101680d18d5SMatthew G. Knepley     /* TODO Check each cell individually */
11029566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCellType(dm, ctx->cells[0], &ct));
1103680d18d5SMatthew G. Knepley     switch (ct) {
11049566063dSJacob Faibussowitsch       case DM_POLYTOPE_SEGMENT:       PetscCall(DMInterpolate_Segment_Private(ctx, dm, x, v));break;
11059566063dSJacob Faibussowitsch       case DM_POLYTOPE_TRIANGLE:      PetscCall(DMInterpolate_Triangle_Private(ctx, dm, x, v));break;
11069566063dSJacob Faibussowitsch       case DM_POLYTOPE_QUADRILATERAL: PetscCall(DMInterpolate_Quad_Private(ctx, dm, x, v));break;
11079566063dSJacob Faibussowitsch       case DM_POLYTOPE_TETRAHEDRON:   PetscCall(DMInterpolate_Tetrahedron_Private(ctx, dm, x, v));break;
11089566063dSJacob Faibussowitsch       case DM_POLYTOPE_HEXAHEDRON:    PetscCall(DMInterpolate_Hex_Private(ctx, dm, x, v));break;
1109cfe5744fSMatthew G. Knepley       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]);
1110680d18d5SMatthew G. Knepley     }
1111552f7358SJed Brown   }
1112552f7358SJed Brown   PetscFunctionReturn(0);
1113552f7358SJed Brown }
1114552f7358SJed Brown 
11154267b1a3SMatthew G. Knepley /*@C
11164267b1a3SMatthew G. Knepley   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
11174267b1a3SMatthew G. Knepley 
11184267b1a3SMatthew G. Knepley   Collective on ctx
11194267b1a3SMatthew G. Knepley 
11204267b1a3SMatthew G. Knepley   Input Parameter:
11214267b1a3SMatthew G. Knepley . ctx - the context
11224267b1a3SMatthew G. Knepley 
11234267b1a3SMatthew G. Knepley   Level: beginner
11244267b1a3SMatthew G. Knepley 
1125db781477SPatrick Sanan .seealso: `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`
11264267b1a3SMatthew G. Knepley @*/
11270adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
11280adebc6cSBarry Smith {
1129552f7358SJed Brown   PetscFunctionBegin;
1130064a246eSJacob Faibussowitsch   PetscValidPointer(ctx, 1);
11319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*ctx)->coords));
11329566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->points));
11339566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ctx)->cells));
11349566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ctx));
11350298fd71SBarry Smith   *ctx = NULL;
1136552f7358SJed Brown   PetscFunctionReturn(0);
1137552f7358SJed Brown }
1138cc0c4584SMatthew G. Knepley 
1139cc0c4584SMatthew G. Knepley /*@C
1140cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
1141cc0c4584SMatthew G. Knepley 
1142cc0c4584SMatthew G. Knepley   Collective on SNES
1143cc0c4584SMatthew G. Knepley 
1144cc0c4584SMatthew G. Knepley   Input Parameters:
1145cc0c4584SMatthew G. Knepley + snes   - the SNES context
1146cc0c4584SMatthew G. Knepley . its    - iteration number
1147cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
1148d43b4f6eSBarry Smith - vf  - PetscViewerAndFormat of type ASCII
1149cc0c4584SMatthew G. Knepley 
1150cc0c4584SMatthew G. Knepley   Notes:
1151cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
1152cc0c4584SMatthew G. Knepley 
1153cc0c4584SMatthew G. Knepley   Level: intermediate
1154cc0c4584SMatthew G. Knepley 
1155db781477SPatrick Sanan .seealso: `SNESMonitorSet()`, `SNESMonitorDefault()`
1156cc0c4584SMatthew G. Knepley @*/
1157d43b4f6eSBarry Smith PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
1158cc0c4584SMatthew G. Knepley {
1159d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
1160cc0c4584SMatthew G. Knepley   Vec                res;
1161cc0c4584SMatthew G. Knepley   DM                 dm;
1162cc0c4584SMatthew G. Knepley   PetscSection       s;
1163cc0c4584SMatthew G. Knepley   const PetscScalar *r;
1164cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
1165cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
1166cc0c4584SMatthew G. Knepley 
1167cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
11684d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
11699566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(snes, &res, NULL, NULL));
11709566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
11719566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &s));
11729566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(s, &numFields));
11739566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(s, &pStart, &pEnd));
11749566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms));
11759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(res, &r));
1176cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
1177cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
1178cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
1179cc0c4584SMatthew G. Knepley 
11809566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof));
11819566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff));
1182cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
1183cc0c4584SMatthew G. Knepley     }
1184cc0c4584SMatthew G. Knepley   }
11859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(res, &r));
11861c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm)));
11879566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(viewer,vf->format));
11889566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel));
118963a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double) fgnorm));
1190cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
11919566063dSJacob Faibussowitsch     if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
11929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f])));
1193cc0c4584SMatthew G. Knepley   }
11949566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "]\n"));
11959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel));
11969566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(viewer));
11979566063dSJacob Faibussowitsch   PetscCall(PetscFree2(lnorms, norms));
1198cc0c4584SMatthew G. Knepley   PetscFunctionReturn(0);
1199cc0c4584SMatthew G. Knepley }
120024cdb843SMatthew G. Knepley 
120124cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
120224cdb843SMatthew G. Knepley 
12036528b96dSMatthew G. Knepley PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS)
12046528b96dSMatthew G. Knepley {
12056528b96dSMatthew G. Knepley   PetscInt       depth;
12066528b96dSMatthew G. Knepley 
12076528b96dSMatthew G. Knepley   PetscFunctionBegin;
12089566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(plex, &depth));
12099566063dSJacob Faibussowitsch   PetscCall(DMGetStratumIS(plex, "dim", depth, cellIS));
12109566063dSJacob Faibussowitsch   if (!*cellIS) PetscCall(DMGetStratumIS(plex, "depth", depth, cellIS));
12116528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
12126528b96dSMatthew G. Knepley }
12136528b96dSMatthew G. Knepley 
121424cdb843SMatthew G. Knepley /*@
12158db23a53SJed Brown   DMPlexSNESComputeResidualFEM - Sums the local residual into vector F from the local input X using pointwise functions specified by the user
121624cdb843SMatthew G. Knepley 
121724cdb843SMatthew G. Knepley   Input Parameters:
121824cdb843SMatthew G. Knepley + dm - The mesh
121924cdb843SMatthew G. Knepley . X  - Local solution
122024cdb843SMatthew G. Knepley - user - The user context
122124cdb843SMatthew G. Knepley 
122224cdb843SMatthew G. Knepley   Output Parameter:
122324cdb843SMatthew G. Knepley . F  - Local output vector
122424cdb843SMatthew G. Knepley 
12258db23a53SJed Brown   Notes:
12268db23a53SJed Brown   The residual is summed into F; the caller is responsible for using VecZeroEntries() or otherwise ensuring that any data in F is intentional.
12278db23a53SJed Brown 
122824cdb843SMatthew G. Knepley   Level: developer
122924cdb843SMatthew G. Knepley 
1230db781477SPatrick Sanan .seealso: `DMPlexComputeJacobianAction()`
123124cdb843SMatthew G. Knepley @*/
123224cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
123324cdb843SMatthew G. Knepley {
12346da023fcSToby Isaac   DM             plex;
1235083401c6SMatthew G. Knepley   IS             allcellIS;
12366528b96dSMatthew G. Knepley   PetscInt       Nds, s;
123724cdb843SMatthew G. Knepley 
123824cdb843SMatthew G. Knepley   PetscFunctionBegin;
12399566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
12409566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
12419566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
12426528b96dSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
12436528b96dSMatthew G. Knepley     PetscDS          ds;
12446528b96dSMatthew G. Knepley     IS               cellIS;
124506ad1575SMatthew G. Knepley     PetscFormKey key;
12466528b96dSMatthew G. Knepley 
12479566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
12486528b96dSMatthew G. Knepley     key.value = 0;
12496528b96dSMatthew G. Knepley     key.field = 0;
125006ad1575SMatthew G. Knepley     key.part  = 0;
12516528b96dSMatthew G. Knepley     if (!key.label) {
12529566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject) allcellIS));
12536528b96dSMatthew G. Knepley       cellIS = allcellIS;
12546528b96dSMatthew G. Knepley     } else {
12556528b96dSMatthew G. Knepley       IS pointIS;
12566528b96dSMatthew G. Knepley 
12576528b96dSMatthew G. Knepley       key.value = 1;
12589566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
12599566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
12609566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
12616528b96dSMatthew G. Knepley     }
12629566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
12639566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
12646528b96dSMatthew G. Knepley   }
12659566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
12669566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
12676528b96dSMatthew G. Knepley   PetscFunctionReturn(0);
12686528b96dSMatthew G. Knepley }
12696528b96dSMatthew G. Knepley 
12706528b96dSMatthew G. Knepley PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user)
12716528b96dSMatthew G. Knepley {
12726528b96dSMatthew G. Knepley   DM             plex;
12736528b96dSMatthew G. Knepley   IS             allcellIS;
12746528b96dSMatthew G. Knepley   PetscInt       Nds, s;
12756528b96dSMatthew G. Knepley 
12766528b96dSMatthew G. Knepley   PetscFunctionBegin;
12779566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
12789566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
12799566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1280083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1281083401c6SMatthew G. Knepley     PetscDS ds;
1282083401c6SMatthew G. Knepley     DMLabel label;
1283083401c6SMatthew G. Knepley     IS      cellIS;
1284083401c6SMatthew G. Knepley 
12859566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds));
12866528b96dSMatthew G. Knepley     {
128706ad1575SMatthew G. Knepley       PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1};
12886528b96dSMatthew G. Knepley       PetscWeakForm     wf;
12896528b96dSMatthew G. Knepley       PetscInt          Nm = 2, m, Nk = 0, k, kp, off = 0;
129006ad1575SMatthew G. Knepley       PetscFormKey *reskeys;
12916528b96dSMatthew G. Knepley 
12926528b96dSMatthew G. Knepley       /* Get unique residual keys */
12936528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
12946528b96dSMatthew G. Knepley         PetscInt Nkm;
12959566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm));
12966528b96dSMatthew G. Knepley         Nk  += Nkm;
12976528b96dSMatthew G. Knepley       }
12989566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &reskeys));
12996528b96dSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
13009566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys));
13016528b96dSMatthew G. Knepley       }
130263a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
13039566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, reskeys));
13046528b96dSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
13056528b96dSMatthew G. Knepley         if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
13066528b96dSMatthew G. Knepley           ++k;
13076528b96dSMatthew G. Knepley           if (kp != k) reskeys[k] = reskeys[kp];
13086528b96dSMatthew G. Knepley         }
13096528b96dSMatthew G. Knepley       }
13106528b96dSMatthew G. Knepley       Nk = k;
13116528b96dSMatthew G. Knepley 
13129566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
13136528b96dSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
13146528b96dSMatthew G. Knepley         DMLabel  label = reskeys[k].label;
13156528b96dSMatthew G. Knepley         PetscInt val   = reskeys[k].value;
13166528b96dSMatthew G. Knepley 
1317083401c6SMatthew G. Knepley         if (!label) {
13189566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject) allcellIS));
1319083401c6SMatthew G. Knepley           cellIS = allcellIS;
1320083401c6SMatthew G. Knepley         } else {
1321083401c6SMatthew G. Knepley           IS pointIS;
1322083401c6SMatthew G. Knepley 
13239566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
13249566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
13259566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
13264a3e9fdbSToby Isaac         }
13279566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user));
13289566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
1329083401c6SMatthew G. Knepley       }
13309566063dSJacob Faibussowitsch       PetscCall(PetscFree(reskeys));
13316528b96dSMatthew G. Knepley     }
13326528b96dSMatthew G. Knepley   }
13339566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
13349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
133524cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
133624cdb843SMatthew G. Knepley }
133724cdb843SMatthew G. Knepley 
1338bdd6f66aSToby Isaac /*@
1339bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1340bdd6f66aSToby Isaac 
1341bdd6f66aSToby Isaac   Input Parameters:
1342bdd6f66aSToby Isaac + dm - The mesh
1343bdd6f66aSToby Isaac - user - The user context
1344bdd6f66aSToby Isaac 
1345bdd6f66aSToby Isaac   Output Parameter:
1346bdd6f66aSToby Isaac . X  - Local solution
1347bdd6f66aSToby Isaac 
1348bdd6f66aSToby Isaac   Level: developer
1349bdd6f66aSToby Isaac 
1350db781477SPatrick Sanan .seealso: `DMPlexComputeJacobianAction()`
1351bdd6f66aSToby Isaac @*/
1352bdd6f66aSToby Isaac PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1353bdd6f66aSToby Isaac {
1354bdd6f66aSToby Isaac   DM             plex;
1355bdd6f66aSToby Isaac 
1356bdd6f66aSToby Isaac   PetscFunctionBegin;
13579566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm,&plex,PETSC_TRUE));
13589566063dSJacob Faibussowitsch   PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL));
13599566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1360bdd6f66aSToby Isaac   PetscFunctionReturn(0);
1361bdd6f66aSToby Isaac }
1362bdd6f66aSToby Isaac 
13637a73cf09SMatthew G. Knepley /*@
13648e3a2eefSMatthew G. Knepley   DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y
13657a73cf09SMatthew G. Knepley 
13667a73cf09SMatthew G. Knepley   Input Parameters:
13678e3a2eefSMatthew G. Knepley + dm   - The DM
13687a73cf09SMatthew G. Knepley . X    - Local solution vector
13697a73cf09SMatthew G. Knepley . Y    - Local input vector
13707a73cf09SMatthew G. Knepley - user - The user context
13717a73cf09SMatthew G. Knepley 
13727a73cf09SMatthew G. Knepley   Output Parameter:
13738e3a2eefSMatthew G. Knepley . F    - lcoal output vector
13747a73cf09SMatthew G. Knepley 
13757a73cf09SMatthew G. Knepley   Level: developer
13767a73cf09SMatthew G. Knepley 
13778e3a2eefSMatthew G. Knepley   Notes:
13788e3a2eefSMatthew G. Knepley   Users will typically use DMSNESCreateJacobianMF() followed by MatMult() instead of calling this routine directly.
13798e3a2eefSMatthew G. Knepley 
1380db781477SPatrick Sanan .seealso: `DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()`
13817a73cf09SMatthew G. Knepley @*/
13828e3a2eefSMatthew G. Knepley PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user)
1383a925c78cSMatthew G. Knepley {
13848e3a2eefSMatthew G. Knepley   DM             plex;
13858e3a2eefSMatthew G. Knepley   IS             allcellIS;
13868e3a2eefSMatthew G. Knepley   PetscInt       Nds, s;
1387a925c78cSMatthew G. Knepley 
1388a925c78cSMatthew G. Knepley   PetscFunctionBegin;
13899566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
13909566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
13919566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
13928e3a2eefSMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
13938e3a2eefSMatthew G. Knepley     PetscDS ds;
13948e3a2eefSMatthew G. Knepley     DMLabel label;
13958e3a2eefSMatthew G. Knepley     IS      cellIS;
13967a73cf09SMatthew G. Knepley 
13979566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds));
13988e3a2eefSMatthew G. Knepley     {
139906ad1575SMatthew G. Knepley       PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3};
14008e3a2eefSMatthew G. Knepley       PetscWeakForm     wf;
14018e3a2eefSMatthew G. Knepley       PetscInt          Nm = 4, m, Nk = 0, k, kp, off = 0;
140206ad1575SMatthew G. Knepley       PetscFormKey *jackeys;
14038e3a2eefSMatthew G. Knepley 
14048e3a2eefSMatthew G. Knepley       /* Get unique Jacobian keys */
14058e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
14068e3a2eefSMatthew G. Knepley         PetscInt Nkm;
14079566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm));
14088e3a2eefSMatthew G. Knepley         Nk  += Nkm;
14098e3a2eefSMatthew G. Knepley       }
14109566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &jackeys));
14118e3a2eefSMatthew G. Knepley       for (m = 0; m < Nm; ++m) {
14129566063dSJacob Faibussowitsch         PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys));
14138e3a2eefSMatthew G. Knepley       }
141463a3b9bcSJacob Faibussowitsch       PetscCheck(off == Nk,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk);
14159566063dSJacob Faibussowitsch       PetscCall(PetscFormKeySort(Nk, jackeys));
14168e3a2eefSMatthew G. Knepley       for (k = 0, kp = 1; kp < Nk; ++kp) {
14178e3a2eefSMatthew G. Knepley         if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
14188e3a2eefSMatthew G. Knepley           ++k;
14198e3a2eefSMatthew G. Knepley           if (kp != k) jackeys[k] = jackeys[kp];
14208e3a2eefSMatthew G. Knepley         }
14218e3a2eefSMatthew G. Knepley       }
14228e3a2eefSMatthew G. Knepley       Nk = k;
14238e3a2eefSMatthew G. Knepley 
14249566063dSJacob Faibussowitsch       PetscCall(PetscDSGetWeakForm(ds, &wf));
14258e3a2eefSMatthew G. Knepley       for (k = 0; k < Nk; ++k) {
14268e3a2eefSMatthew G. Knepley         DMLabel  label = jackeys[k].label;
14278e3a2eefSMatthew G. Knepley         PetscInt val   = jackeys[k].value;
14288e3a2eefSMatthew G. Knepley 
14298e3a2eefSMatthew G. Knepley         if (!label) {
14309566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject) allcellIS));
14318e3a2eefSMatthew G. Knepley           cellIS = allcellIS;
14327a73cf09SMatthew G. Knepley         } else {
14338e3a2eefSMatthew G. Knepley           IS pointIS;
1434a925c78cSMatthew G. Knepley 
14359566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(label, val, &pointIS));
14369566063dSJacob Faibussowitsch           PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
14379566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&pointIS));
1438a925c78cSMatthew G. Knepley         }
14399566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user));
14409566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&cellIS));
14418e3a2eefSMatthew G. Knepley       }
14429566063dSJacob Faibussowitsch       PetscCall(PetscFree(jackeys));
14438e3a2eefSMatthew G. Knepley     }
14448e3a2eefSMatthew G. Knepley   }
14459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
14469566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
1447a925c78cSMatthew G. Knepley   PetscFunctionReturn(0);
1448a925c78cSMatthew G. Knepley }
1449a925c78cSMatthew G. Knepley 
145024cdb843SMatthew G. Knepley /*@
145124cdb843SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
145224cdb843SMatthew G. Knepley 
145324cdb843SMatthew G. Knepley   Input Parameters:
145424cdb843SMatthew G. Knepley + dm - The mesh
145524cdb843SMatthew G. Knepley . X  - Local input vector
145624cdb843SMatthew G. Knepley - user - The user context
145724cdb843SMatthew G. Knepley 
145824cdb843SMatthew G. Knepley   Output Parameter:
145924cdb843SMatthew G. Knepley . Jac  - Jacobian matrix
146024cdb843SMatthew G. Knepley 
146124cdb843SMatthew G. Knepley   Note:
146224cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
146324cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
146424cdb843SMatthew G. Knepley 
146524cdb843SMatthew G. Knepley   Level: developer
146624cdb843SMatthew G. Knepley 
1467db781477SPatrick Sanan .seealso: `FormFunctionLocal()`
146824cdb843SMatthew G. Knepley @*/
146924cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
147024cdb843SMatthew G. Knepley {
14716da023fcSToby Isaac   DM             plex;
1472083401c6SMatthew G. Knepley   IS             allcellIS;
1473f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
14746528b96dSMatthew G. Knepley   PetscInt       Nds, s;
147524cdb843SMatthew G. Knepley 
147624cdb843SMatthew G. Knepley   PetscFunctionBegin;
14779566063dSJacob Faibussowitsch   PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE));
14789566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS));
14799566063dSJacob Faibussowitsch   PetscCall(DMGetNumDS(dm, &Nds));
1480083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1481083401c6SMatthew G. Knepley     PetscDS          ds;
1482083401c6SMatthew G. Knepley     IS               cellIS;
148306ad1575SMatthew G. Knepley     PetscFormKey key;
1484083401c6SMatthew G. Knepley 
14859566063dSJacob Faibussowitsch     PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds));
14866528b96dSMatthew G. Knepley     key.value = 0;
14876528b96dSMatthew G. Knepley     key.field = 0;
148806ad1575SMatthew G. Knepley     key.part  = 0;
14896528b96dSMatthew G. Knepley     if (!key.label) {
14909566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject) allcellIS));
1491083401c6SMatthew G. Knepley       cellIS = allcellIS;
1492083401c6SMatthew G. Knepley     } else {
1493083401c6SMatthew G. Knepley       IS pointIS;
1494083401c6SMatthew G. Knepley 
14956528b96dSMatthew G. Knepley       key.value = 1;
14969566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS));
14979566063dSJacob Faibussowitsch       PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS));
14989566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&pointIS));
1499083401c6SMatthew G. Knepley     }
1500083401c6SMatthew G. Knepley     if (!s) {
15019566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobian(ds, &hasJac));
15029566063dSJacob Faibussowitsch       PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
15039566063dSJacob Faibussowitsch       if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac));
15049566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(JacP));
1505083401c6SMatthew G. Knepley     }
15069566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user));
15079566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&cellIS));
1508083401c6SMatthew G. Knepley   }
15099566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&allcellIS));
15109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
151124cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
151224cdb843SMatthew G. Knepley }
15131878804aSMatthew G. Knepley 
15148e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx
15158e3a2eefSMatthew G. Knepley {
15168e3a2eefSMatthew G. Knepley   DM    dm;
15178e3a2eefSMatthew G. Knepley   Vec   X;
15188e3a2eefSMatthew G. Knepley   void *ctx;
15198e3a2eefSMatthew G. Knepley };
15208e3a2eefSMatthew G. Knepley 
15218e3a2eefSMatthew G. Knepley static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A)
15228e3a2eefSMatthew G. Knepley {
15238e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15248e3a2eefSMatthew G. Knepley 
15258e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
15269566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
15279566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(A, NULL));
15289566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dm));
15299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->X));
15309566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
15318e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
15328e3a2eefSMatthew G. Knepley }
15338e3a2eefSMatthew G. Knepley 
15348e3a2eefSMatthew G. Knepley static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z)
15358e3a2eefSMatthew G. Knepley {
15368e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15378e3a2eefSMatthew G. Knepley 
15388e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
15399566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
15409566063dSJacob Faibussowitsch   PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx));
15418e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
15428e3a2eefSMatthew G. Knepley }
15438e3a2eefSMatthew G. Knepley 
15448e3a2eefSMatthew G. Knepley /*@
15458e3a2eefSMatthew G. Knepley   DMSNESCreateJacobianMF - Create a Mat which computes the action of the Jacobian matrix-free
15468e3a2eefSMatthew G. Knepley 
15478e3a2eefSMatthew G. Knepley   Collective on dm
15488e3a2eefSMatthew G. Knepley 
15498e3a2eefSMatthew G. Knepley   Input Parameters:
15508e3a2eefSMatthew G. Knepley + dm   - The DM
15518e3a2eefSMatthew G. Knepley . X    - The evaluation point for the Jacobian
15528e3a2eefSMatthew G. Knepley - user - A user context, or NULL
15538e3a2eefSMatthew G. Knepley 
15548e3a2eefSMatthew G. Knepley   Output Parameter:
15558e3a2eefSMatthew G. Knepley . J    - The Mat
15568e3a2eefSMatthew G. Knepley 
15578e3a2eefSMatthew G. Knepley   Level: advanced
15588e3a2eefSMatthew G. Knepley 
15598e3a2eefSMatthew G. Knepley   Notes:
15608e3a2eefSMatthew G. Knepley   Vec X is kept in Mat J, so updating X then updates the evaluation point.
15618e3a2eefSMatthew G. Knepley 
1562db781477SPatrick Sanan .seealso: `DMSNESComputeJacobianAction()`
15638e3a2eefSMatthew G. Knepley @*/
15648e3a2eefSMatthew G. Knepley PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J)
15658e3a2eefSMatthew G. Knepley {
15668e3a2eefSMatthew G. Knepley   struct _DMSNESJacobianMFCtx *ctx;
15678e3a2eefSMatthew G. Knepley   PetscInt                     n, N;
15688e3a2eefSMatthew G. Knepley 
15698e3a2eefSMatthew G. Knepley   PetscFunctionBegin;
15709566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject) dm), J));
15719566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATSHELL));
15729566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(X, &n));
15739566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
15749566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, n, n, N, N));
15759566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject) dm));
15769566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject) X));
15779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &ctx));
15788e3a2eefSMatthew G. Knepley   ctx->dm  = dm;
15798e3a2eefSMatthew G. Knepley   ctx->X   = X;
15808e3a2eefSMatthew G. Knepley   ctx->ctx = user;
15819566063dSJacob Faibussowitsch   PetscCall(MatShellSetContext(*J, ctx));
15829566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void)) DMSNESJacobianMF_Destroy_Private));
15839566063dSJacob Faibussowitsch   PetscCall(MatShellSetOperation(*J, MATOP_MULT,    (void (*)(void)) DMSNESJacobianMF_Mult_Private));
15848e3a2eefSMatthew G. Knepley   PetscFunctionReturn(0);
15858e3a2eefSMatthew G. Knepley }
15868e3a2eefSMatthew G. Knepley 
158738cfc46eSPierre Jolivet /*
158838cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
158938cfc46eSPierre Jolivet 
159038cfc46eSPierre Jolivet    Input Parameters:
159138cfc46eSPierre Jolivet +     X - SNES linearization point
159238cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
159338cfc46eSPierre Jolivet 
159438cfc46eSPierre Jolivet    Output Parameter:
159538cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
159638cfc46eSPierre Jolivet 
159738cfc46eSPierre Jolivet    Level: intermediate
159838cfc46eSPierre Jolivet 
1599db781477SPatrick Sanan .seealso: `DMCreateNeumannOverlap()`, `MATIS`, `PCHPDDMSetAuxiliaryMat()`
160038cfc46eSPierre Jolivet */
160138cfc46eSPierre Jolivet static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
160238cfc46eSPierre Jolivet {
160338cfc46eSPierre Jolivet   SNES           snes;
160438cfc46eSPierre Jolivet   Mat            pJ;
160538cfc46eSPierre Jolivet   DM             ovldm,origdm;
160638cfc46eSPierre Jolivet   DMSNES         sdm;
160738cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM,Vec,void*);
160838cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
160938cfc46eSPierre Jolivet   void           *bctx,*jctx;
161038cfc46eSPierre Jolivet 
161138cfc46eSPierre Jolivet   PetscFunctionBegin;
16129566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ));
161328b400f6SJacob Faibussowitsch   PetscCheck(pJ,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");
16149566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm));
161528b400f6SJacob Faibussowitsch   PetscCheck(origdm,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");
16169566063dSJacob Faibussowitsch   PetscCall(MatGetDM(pJ,&ovldm));
16179566063dSJacob Faibussowitsch   PetscCall(DMSNESGetBoundaryLocal(origdm,&bfun,&bctx));
16189566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(ovldm,bfun,bctx));
16199566063dSJacob Faibussowitsch   PetscCall(DMSNESGetJacobianLocal(origdm,&jfun,&jctx));
16209566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(ovldm,jfun,jctx));
16219566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes));
162238cfc46eSPierre Jolivet   if (!snes) {
16239566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl),&snes));
16249566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes,ovldm));
16259566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes));
16269566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)snes));
162738cfc46eSPierre Jolivet   }
16289566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(ovldm,&sdm));
16299566063dSJacob Faibussowitsch   PetscCall(VecLockReadPush(X));
1630*800f99ffSJeremy L Thompson   {
1631*800f99ffSJeremy L Thompson     void *ctx;
1632*800f99ffSJeremy L Thompson     PetscErrorCode (*J)(SNES,Vec,Mat,Mat,void*);
1633*800f99ffSJeremy L Thompson     PetscCall(DMSNESGetJacobian(ovldm,&J,&ctx));
1634*800f99ffSJeremy L Thompson     PetscCallBack("SNES callback Jacobian",(*J)(snes,X,pJ,pJ,ctx));
1635*800f99ffSJeremy L Thompson   }
16369566063dSJacob Faibussowitsch   PetscCall(VecLockReadPop(X));
163738cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
163838cfc46eSPierre Jolivet   {
163938cfc46eSPierre Jolivet     Mat locpJ;
164038cfc46eSPierre Jolivet 
16419566063dSJacob Faibussowitsch     PetscCall(MatISGetLocalMat(pJ,&locpJ));
16429566063dSJacob Faibussowitsch     PetscCall(MatCopy(locpJ,J,SAME_NONZERO_PATTERN));
164338cfc46eSPierre Jolivet   }
164438cfc46eSPierre Jolivet   PetscFunctionReturn(0);
164538cfc46eSPierre Jolivet }
164638cfc46eSPierre Jolivet 
1647a925c78cSMatthew G. Knepley /*@
16489f520fc2SToby Isaac   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
16499f520fc2SToby Isaac 
16509f520fc2SToby Isaac   Input Parameters:
16519f520fc2SToby Isaac + dm - The DM object
1652dff059c6SToby Isaac . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
16539f520fc2SToby Isaac . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
16549f520fc2SToby Isaac - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
16551a244344SSatish Balay 
16561a244344SSatish Balay   Level: developer
16579f520fc2SToby Isaac @*/
16589f520fc2SToby Isaac PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
16599f520fc2SToby Isaac {
16609f520fc2SToby Isaac   PetscFunctionBegin;
16619566063dSJacob Faibussowitsch   PetscCall(DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx));
16629566063dSJacob Faibussowitsch   PetscCall(DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx));
16639566063dSJacob Faibussowitsch   PetscCall(DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx));
16649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex));
16659f520fc2SToby Isaac   PetscFunctionReturn(0);
16669f520fc2SToby Isaac }
16679f520fc2SToby Isaac 
1668f017bcb6SMatthew G. Knepley /*@C
1669f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1670f017bcb6SMatthew G. Knepley 
1671f017bcb6SMatthew G. Knepley   Input Parameters:
1672f017bcb6SMatthew G. Knepley + snes - the SNES object
1673f017bcb6SMatthew G. Knepley . dm   - the DM
1674f2cacb80SMatthew G. Knepley . t    - the time
1675f017bcb6SMatthew G. Knepley . u    - a DM vector
1676f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1677f017bcb6SMatthew G. Knepley 
1678f3c94560SMatthew G. Knepley   Output Parameters:
1679f3c94560SMatthew G. Knepley . error - An array which holds the discretization error in each field, or NULL
1680f3c94560SMatthew G. Knepley 
16817f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
16827f96f943SMatthew G. Knepley 
1683f017bcb6SMatthew G. Knepley   Level: developer
1684f017bcb6SMatthew G. Knepley 
1685db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()`, `PetscDSSetExactSolution()`
1686f017bcb6SMatthew G. Knepley @*/
1687f2cacb80SMatthew G. Knepley PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
16881878804aSMatthew G. Knepley {
1689f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1690f017bcb6SMatthew G. Knepley   void            **ectxs;
1691f3c94560SMatthew G. Knepley   PetscReal        *err;
16927f96f943SMatthew G. Knepley   MPI_Comm          comm;
16937f96f943SMatthew G. Knepley   PetscInt          Nf, f;
16941878804aSMatthew G. Knepley 
16951878804aSMatthew G. Knepley   PetscFunctionBegin;
1696f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1697f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1698064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(u, VEC_CLASSID, 4);
1699534a8f05SLisandro Dalcin   if (error) PetscValidRealPointer(error, 6);
17007f96f943SMatthew G. Knepley 
17019566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
17029566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-vec_view"));
17037f96f943SMatthew G. Knepley 
17049566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) snes, &comm));
17059566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
17069566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err));
17077f96f943SMatthew G. Knepley   {
17087f96f943SMatthew G. Knepley     PetscInt Nds, s;
17097f96f943SMatthew G. Knepley 
17109566063dSJacob Faibussowitsch     PetscCall(DMGetNumDS(dm, &Nds));
1711083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
17127f96f943SMatthew G. Knepley       PetscDS         ds;
1713083401c6SMatthew G. Knepley       DMLabel         label;
1714083401c6SMatthew G. Knepley       IS              fieldIS;
17157f96f943SMatthew G. Knepley       const PetscInt *fields;
17167f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1717083401c6SMatthew G. Knepley 
17189566063dSJacob Faibussowitsch       PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds));
17199566063dSJacob Faibussowitsch       PetscCall(PetscDSGetNumFields(ds, &dsNf));
17209566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fieldIS, &fields));
1721083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1722083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
17239566063dSJacob Faibussowitsch         PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]));
1724083401c6SMatthew G. Knepley       }
17259566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fieldIS, &fields));
1726083401c6SMatthew G. Knepley     }
1727083401c6SMatthew G. Knepley   }
1728f017bcb6SMatthew G. Knepley   if (Nf > 1) {
17299566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err));
1730f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1731f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
173263a3b9bcSJacob Faibussowitsch         PetscCheck(err[f] <= tol,comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %" PetscInt_FMT " exceeds tolerance %g", (double) err[f], f, (double) tol);
17331878804aSMatthew G. Knepley       }
1734b878532fSMatthew G. Knepley     } else if (error) {
1735b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
17361878804aSMatthew G. Knepley     } else {
17379566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: ["));
1738f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
17399566063dSJacob Faibussowitsch         if (f) PetscCall(PetscPrintf(comm, ", "));
17409566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(comm, "%g", (double)err[f]));
17411878804aSMatthew G. Knepley       }
17429566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "]\n"));
1743f017bcb6SMatthew G. Knepley     }
1744f017bcb6SMatthew G. Knepley   } else {
17459566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]));
1746f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
174708401ef6SPierre Jolivet       PetscCheck(err[0] <= tol,comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1748b878532fSMatthew G. Knepley     } else if (error) {
1749b878532fSMatthew G. Knepley       error[0] = err[0];
1750f017bcb6SMatthew G. Knepley     } else {
17519566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]));
1752f017bcb6SMatthew G. Knepley     }
1753f017bcb6SMatthew G. Knepley   }
17549566063dSJacob Faibussowitsch   PetscCall(PetscFree3(exacts, ectxs, err));
1755f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1756f017bcb6SMatthew G. Knepley }
1757f017bcb6SMatthew G. Knepley 
1758f017bcb6SMatthew G. Knepley /*@C
1759f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1760f017bcb6SMatthew G. Knepley 
1761f017bcb6SMatthew G. Knepley   Input Parameters:
1762f017bcb6SMatthew G. Knepley + snes - the SNES object
1763f017bcb6SMatthew G. Knepley . dm   - the DM
1764f017bcb6SMatthew G. Knepley . u    - a DM vector
1765f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1766f017bcb6SMatthew G. Knepley 
1767f3c94560SMatthew G. Knepley   Output Parameters:
1768f3c94560SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
1769f3c94560SMatthew G. Knepley 
1770f017bcb6SMatthew G. Knepley   Level: developer
1771f017bcb6SMatthew G. Knepley 
1772db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()`
1773f017bcb6SMatthew G. Knepley @*/
1774f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1775f017bcb6SMatthew G. Knepley {
1776f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1777f017bcb6SMatthew G. Knepley   Vec            r;
1778f017bcb6SMatthew G. Knepley   PetscReal      res;
1779f017bcb6SMatthew G. Knepley 
1780f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1781f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1782f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1783f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1784534a8f05SLisandro Dalcin   if (residual) PetscValidRealPointer(residual, 5);
17859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) snes, &comm));
17869566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
17879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
17889566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, u, r));
17899566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
1790f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
179108401ef6SPierre Jolivet     PetscCheck(res <= tol,comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1792b878532fSMatthew G. Knepley   } else if (residual) {
1793b878532fSMatthew G. Knepley     *residual = res;
1794f017bcb6SMatthew G. Knepley   } else {
17959566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res));
17969566063dSJacob Faibussowitsch     PetscCall(VecChop(r, 1.0e-10));
17979566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) r, "Initial Residual"));
17989566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r,"res_"));
17999566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-vec_view"));
1800f017bcb6SMatthew G. Knepley   }
18019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
1802f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1803f017bcb6SMatthew G. Knepley }
1804f017bcb6SMatthew G. Knepley 
1805f017bcb6SMatthew G. Knepley /*@C
1806f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1807f017bcb6SMatthew G. Knepley 
1808f017bcb6SMatthew G. Knepley   Input Parameters:
1809f017bcb6SMatthew G. Knepley + snes - the SNES object
1810f017bcb6SMatthew G. Knepley . dm   - the DM
1811f017bcb6SMatthew G. Knepley . u    - a DM vector
1812f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1813f017bcb6SMatthew G. Knepley 
1814f3c94560SMatthew G. Knepley   Output Parameters:
1815f3c94560SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
1816f3c94560SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
1817f3c94560SMatthew G. Knepley 
1818f017bcb6SMatthew G. Knepley   Level: developer
1819f017bcb6SMatthew G. Knepley 
1820db781477SPatrick Sanan .seealso: `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()`
1821f017bcb6SMatthew G. Knepley @*/
1822f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1823f017bcb6SMatthew G. Knepley {
1824f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1825f017bcb6SMatthew G. Knepley   PetscDS        ds;
1826f017bcb6SMatthew G. Knepley   Mat            J, M;
1827f017bcb6SMatthew G. Knepley   MatNullSpace   nullspace;
1828f3c94560SMatthew G. Knepley   PetscReal      slope, intercept;
18297f96f943SMatthew G. Knepley   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
1830f017bcb6SMatthew G. Knepley 
1831f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1832f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1833f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1834f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1835534a8f05SLisandro Dalcin   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1836064a246eSJacob Faibussowitsch   if (convRate) PetscValidRealPointer(convRate, 6);
18379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) snes, &comm));
18389566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL));
1839f017bcb6SMatthew G. Knepley   /* Create and view matrices */
18409566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &J));
18419566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
18429566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobian(ds, &hasJac));
18439566063dSJacob Faibussowitsch   PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec));
1844282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
18459566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &M));
18469566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, M));
18479566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) M, "Preconditioning Matrix"));
18489566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_"));
18499566063dSJacob Faibussowitsch     PetscCall(MatViewFromOptions(M, NULL, "-mat_view"));
18509566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M));
1851282e7bb4SMatthew G. Knepley   } else {
18529566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, u, J, J));
1853282e7bb4SMatthew G. Knepley   }
18549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) J, "Jacobian"));
18559566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) J, "jac_"));
18569566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(J, NULL, "-mat_view"));
1857f017bcb6SMatthew G. Knepley   /* Check nullspace */
18589566063dSJacob Faibussowitsch   PetscCall(MatGetNullSpace(J, &nullspace));
1859f017bcb6SMatthew G. Knepley   if (nullspace) {
18601878804aSMatthew G. Knepley     PetscBool isNull;
18619566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceTest(nullspace, J, &isNull));
186228b400f6SJacob Faibussowitsch     PetscCheck(isNull,comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
18631878804aSMatthew G. Knepley   }
1864f017bcb6SMatthew G. Knepley   /* Taylor test */
1865f017bcb6SMatthew G. Knepley   {
1866f017bcb6SMatthew G. Knepley     PetscRandom rand;
1867f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1868f3c94560SMatthew G. Knepley     PetscReal   h;
1869f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1870f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1871f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1872f017bcb6SMatthew G. Knepley 
1873f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
18749566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(comm, &rand));
18759566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &du));
18769566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(du, rand));
18779566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&rand));
18789566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &df));
18799566063dSJacob Faibussowitsch     PetscCall(MatMult(J, du, df));
1880f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
18819566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &r));
18829566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, u, r));
1883f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
1884f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
18859566063dSJacob Faibussowitsch     PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors));
18869566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &uhat));
18879566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &rhat));
1888f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
18899566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(uhat, h, du, u));
1890f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
18919566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, uhat, rhat));
18929566063dSJacob Faibussowitsch       PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df));
18939566063dSJacob Faibussowitsch       PetscCall(VecNorm(rhat, NORM_2, &errors[Nv]));
1894f017bcb6SMatthew G. Knepley 
1895f017bcb6SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
1896f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1897f017bcb6SMatthew G. Knepley     }
18989566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&uhat));
18999566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhat));
19009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&df));
19019566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
19029566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&du));
1903f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1904f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1905f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL)    break;
1906f017bcb6SMatthew G. Knepley     }
1907f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
19089566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept));
19099566063dSJacob Faibussowitsch     PetscCall(PetscFree3(es, hs, errors));
1910f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1911f017bcb6SMatthew G. Knepley     if (tol >= 0) {
19120b121fc5SBarry Smith       PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol,comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
1913b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1914b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1915b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1916f017bcb6SMatthew G. Knepley     } else {
19179566063dSJacob Faibussowitsch       if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope));
19189566063dSJacob Faibussowitsch       else        PetscCall(PetscPrintf(comm, "Function appears to be linear\n"));
1919f017bcb6SMatthew G. Knepley     }
1920f017bcb6SMatthew G. Knepley   }
19219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
19221878804aSMatthew G. Knepley   PetscFunctionReturn(0);
19231878804aSMatthew G. Knepley }
19241878804aSMatthew G. Knepley 
19257f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1926f017bcb6SMatthew G. Knepley {
1927f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
19289566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL));
19299566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL));
19309566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL));
1931f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1932f017bcb6SMatthew G. Knepley }
1933f017bcb6SMatthew G. Knepley 
1934bee9a294SMatthew G. Knepley /*@C
1935bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1936bee9a294SMatthew G. Knepley 
1937bee9a294SMatthew G. Knepley   Input Parameters:
1938bee9a294SMatthew G. Knepley + snes - the SNES object
19397f96f943SMatthew G. Knepley - u    - representative SNES vector
19407f96f943SMatthew G. Knepley 
19417f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
1942bee9a294SMatthew G. Knepley 
1943bee9a294SMatthew G. Knepley   Level: developer
1944bee9a294SMatthew G. Knepley @*/
19457f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
19461878804aSMatthew G. Knepley {
19471878804aSMatthew G. Knepley   DM             dm;
19481878804aSMatthew G. Knepley   Vec            sol;
19491878804aSMatthew G. Knepley   PetscBool      check;
19501878804aSMatthew G. Knepley 
19511878804aSMatthew G. Knepley   PetscFunctionBegin;
19529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check));
19531878804aSMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
19549566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
19559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &sol));
19569566063dSJacob Faibussowitsch   PetscCall(SNESSetSolution(snes, sol));
19579566063dSJacob Faibussowitsch   PetscCall(DMSNESCheck_Internal(snes, dm, sol));
19589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sol));
1959552f7358SJed Brown   PetscFunctionReturn(0);
1960552f7358SJed Brown }
1961