xref: /petsc/src/snes/utils/dmplexsnes.c (revision f2cacb80554376ba60c1c183c00cd40ea8a0a80d)
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>
4a925c78cSMatthew G. Knepley #include <petscblaslapack.h>
5af0996ceSBarry Smith #include <petsc/private/petscimpl.h>
6af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h>
7552f7358SJed Brown 
824cdb843SMatthew G. Knepley /************************** Interpolation *******************************/
924cdb843SMatthew G. Knepley 
106da023fcSToby Isaac static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
116da023fcSToby Isaac {
126da023fcSToby Isaac   PetscBool      isPlex;
136da023fcSToby Isaac   PetscErrorCode ierr;
146da023fcSToby Isaac 
156da023fcSToby Isaac   PetscFunctionBegin;
166da023fcSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);CHKERRQ(ierr);
176da023fcSToby Isaac   if (isPlex) {
186da023fcSToby Isaac     *plex = dm;
196da023fcSToby Isaac     ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
20f7148743SMatthew G. Knepley   } else {
21f7148743SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);CHKERRQ(ierr);
22f7148743SMatthew G. Knepley     if (!*plex) {
236da023fcSToby Isaac       ierr = DMConvert(dm,DMPLEX,plex);CHKERRQ(ierr);
24f7148743SMatthew G. Knepley       ierr = PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);CHKERRQ(ierr);
256da023fcSToby Isaac       if (copy) {
266da023fcSToby Isaac         PetscInt    i;
276da023fcSToby Isaac         PetscObject obj;
286da023fcSToby Isaac         const char *comps[3] = {"A","dmAux","dmCh"};
296da023fcSToby Isaac 
306da023fcSToby Isaac         ierr = DMCopyDMSNES(dm, *plex);CHKERRQ(ierr);
316da023fcSToby Isaac         for (i = 0; i < 3; i++) {
326da023fcSToby Isaac           ierr = PetscObjectQuery((PetscObject) dm, comps[i], &obj);CHKERRQ(ierr);
336da023fcSToby Isaac           ierr = PetscObjectCompose((PetscObject) *plex, comps[i], obj);CHKERRQ(ierr);
346da023fcSToby Isaac         }
356da023fcSToby Isaac       }
36f7148743SMatthew G. Knepley     } else {
37f7148743SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) *plex);CHKERRQ(ierr);
38f7148743SMatthew G. Knepley     }
396da023fcSToby Isaac   }
406da023fcSToby Isaac   PetscFunctionReturn(0);
416da023fcSToby Isaac }
426da023fcSToby Isaac 
434267b1a3SMatthew G. Knepley /*@C
444267b1a3SMatthew G. Knepley   DMInterpolationCreate - Creates a DMInterpolationInfo context
454267b1a3SMatthew G. Knepley 
46d083f849SBarry Smith   Collective
474267b1a3SMatthew G. Knepley 
484267b1a3SMatthew G. Knepley   Input Parameter:
494267b1a3SMatthew G. Knepley . comm - the communicator
504267b1a3SMatthew G. Knepley 
514267b1a3SMatthew G. Knepley   Output Parameter:
524267b1a3SMatthew G. Knepley . ctx - the context
534267b1a3SMatthew G. Knepley 
544267b1a3SMatthew G. Knepley   Level: beginner
554267b1a3SMatthew G. Knepley 
564267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
574267b1a3SMatthew G. Knepley @*/
580adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
590adebc6cSBarry Smith {
60552f7358SJed Brown   PetscErrorCode ierr;
61552f7358SJed Brown 
62552f7358SJed Brown   PetscFunctionBegin;
63552f7358SJed Brown   PetscValidPointer(ctx, 2);
6495dccacaSBarry Smith   ierr = PetscNew(ctx);CHKERRQ(ierr);
651aa26658SKarl Rupp 
66552f7358SJed Brown   (*ctx)->comm   = comm;
67552f7358SJed Brown   (*ctx)->dim    = -1;
68552f7358SJed Brown   (*ctx)->nInput = 0;
690298fd71SBarry Smith   (*ctx)->points = NULL;
700298fd71SBarry Smith   (*ctx)->cells  = NULL;
71552f7358SJed Brown   (*ctx)->n      = -1;
720298fd71SBarry Smith   (*ctx)->coords = NULL;
73552f7358SJed Brown   PetscFunctionReturn(0);
74552f7358SJed Brown }
75552f7358SJed Brown 
764267b1a3SMatthew G. Knepley /*@C
774267b1a3SMatthew G. Knepley   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
784267b1a3SMatthew G. Knepley 
794267b1a3SMatthew G. Knepley   Not collective
804267b1a3SMatthew G. Knepley 
814267b1a3SMatthew G. Knepley   Input Parameters:
824267b1a3SMatthew G. Knepley + ctx - the context
834267b1a3SMatthew G. Knepley - dim - the spatial dimension
844267b1a3SMatthew G. Knepley 
854267b1a3SMatthew G. Knepley   Level: intermediate
864267b1a3SMatthew G. Knepley 
874267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
884267b1a3SMatthew G. Knepley @*/
890adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
900adebc6cSBarry Smith {
91552f7358SJed Brown   PetscFunctionBegin;
92ff1e0c32SBarry Smith   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim);
93552f7358SJed Brown   ctx->dim = dim;
94552f7358SJed Brown   PetscFunctionReturn(0);
95552f7358SJed Brown }
96552f7358SJed Brown 
974267b1a3SMatthew G. Knepley /*@C
984267b1a3SMatthew G. Knepley   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
994267b1a3SMatthew G. Knepley 
1004267b1a3SMatthew G. Knepley   Not collective
1014267b1a3SMatthew G. Knepley 
1024267b1a3SMatthew G. Knepley   Input Parameter:
1034267b1a3SMatthew G. Knepley . ctx - the context
1044267b1a3SMatthew G. Knepley 
1054267b1a3SMatthew G. Knepley   Output Parameter:
1064267b1a3SMatthew G. Knepley . dim - the spatial dimension
1074267b1a3SMatthew G. Knepley 
1084267b1a3SMatthew G. Knepley   Level: intermediate
1094267b1a3SMatthew G. Knepley 
1104267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
1114267b1a3SMatthew G. Knepley @*/
1120adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
1130adebc6cSBarry Smith {
114552f7358SJed Brown   PetscFunctionBegin;
115552f7358SJed Brown   PetscValidIntPointer(dim, 2);
116552f7358SJed Brown   *dim = ctx->dim;
117552f7358SJed Brown   PetscFunctionReturn(0);
118552f7358SJed Brown }
119552f7358SJed Brown 
1204267b1a3SMatthew G. Knepley /*@C
1214267b1a3SMatthew G. Knepley   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
1224267b1a3SMatthew G. Knepley 
1234267b1a3SMatthew G. Knepley   Not collective
1244267b1a3SMatthew G. Knepley 
1254267b1a3SMatthew G. Knepley   Input Parameters:
1264267b1a3SMatthew G. Knepley + ctx - the context
1274267b1a3SMatthew G. Knepley - dof - the number of fields
1284267b1a3SMatthew G. Knepley 
1294267b1a3SMatthew G. Knepley   Level: intermediate
1304267b1a3SMatthew G. Knepley 
1314267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
1324267b1a3SMatthew G. Knepley @*/
1330adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
1340adebc6cSBarry Smith {
135552f7358SJed Brown   PetscFunctionBegin;
136ff1e0c32SBarry Smith   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof);
137552f7358SJed Brown   ctx->dof = dof;
138552f7358SJed Brown   PetscFunctionReturn(0);
139552f7358SJed Brown }
140552f7358SJed Brown 
1414267b1a3SMatthew G. Knepley /*@C
1424267b1a3SMatthew G. Knepley   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
1434267b1a3SMatthew G. Knepley 
1444267b1a3SMatthew G. Knepley   Not collective
1454267b1a3SMatthew G. Knepley 
1464267b1a3SMatthew G. Knepley   Input Parameter:
1474267b1a3SMatthew G. Knepley . ctx - the context
1484267b1a3SMatthew G. Knepley 
1494267b1a3SMatthew G. Knepley   Output Parameter:
1504267b1a3SMatthew G. Knepley . dof - the number of fields
1514267b1a3SMatthew G. Knepley 
1524267b1a3SMatthew G. Knepley   Level: intermediate
1534267b1a3SMatthew G. Knepley 
1544267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
1554267b1a3SMatthew G. Knepley @*/
1560adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
1570adebc6cSBarry Smith {
158552f7358SJed Brown   PetscFunctionBegin;
159552f7358SJed Brown   PetscValidIntPointer(dof, 2);
160552f7358SJed Brown   *dof = ctx->dof;
161552f7358SJed Brown   PetscFunctionReturn(0);
162552f7358SJed Brown }
163552f7358SJed Brown 
1644267b1a3SMatthew G. Knepley /*@C
1654267b1a3SMatthew G. Knepley   DMInterpolationAddPoints - Add points at which we will interpolate the fields
1664267b1a3SMatthew G. Knepley 
1674267b1a3SMatthew G. Knepley   Not collective
1684267b1a3SMatthew G. Knepley 
1694267b1a3SMatthew G. Knepley   Input Parameters:
1704267b1a3SMatthew G. Knepley + ctx    - the context
1714267b1a3SMatthew G. Knepley . n      - the number of points
1724267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim
1734267b1a3SMatthew G. Knepley 
1744267b1a3SMatthew G. Knepley   Note: The coordinate information is copied.
1754267b1a3SMatthew G. Knepley 
1764267b1a3SMatthew G. Knepley   Level: intermediate
1774267b1a3SMatthew G. Knepley 
1784267b1a3SMatthew G. Knepley .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
1794267b1a3SMatthew G. Knepley @*/
1800adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
1810adebc6cSBarry Smith {
182552f7358SJed Brown   PetscErrorCode ierr;
183552f7358SJed Brown 
184552f7358SJed Brown   PetscFunctionBegin;
1850adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
1860adebc6cSBarry Smith   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
187552f7358SJed Brown   ctx->nInput = n;
1881aa26658SKarl Rupp 
189785e854fSJed Brown   ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr);
190580bdb30SBarry Smith   ierr = PetscArraycpy(ctx->points, points, n*ctx->dim);CHKERRQ(ierr);
191552f7358SJed Brown   PetscFunctionReturn(0);
192552f7358SJed Brown }
193552f7358SJed Brown 
1944267b1a3SMatthew G. Knepley /*@C
1954267b1a3SMatthew G. Knepley   DMInterpolationSetUp - Computea spatial indices that add in point location during interpolation
1964267b1a3SMatthew G. Knepley 
1974267b1a3SMatthew G. Knepley   Collective on ctx
1984267b1a3SMatthew G. Knepley 
1994267b1a3SMatthew G. Knepley   Input Parameters:
2004267b1a3SMatthew G. Knepley + ctx - the context
2014267b1a3SMatthew G. Knepley . dm  - the DM for the function space used for interpolation
2024267b1a3SMatthew G. Knepley - redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
2034267b1a3SMatthew G. Knepley 
2044267b1a3SMatthew G. Knepley   Level: intermediate
2054267b1a3SMatthew G. Knepley 
2064267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
2074267b1a3SMatthew G. Knepley @*/
2080adebc6cSBarry Smith PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
2090adebc6cSBarry Smith {
210552f7358SJed Brown   MPI_Comm          comm = ctx->comm;
211552f7358SJed Brown   PetscScalar       *a;
212552f7358SJed Brown   PetscInt          p, q, i;
213552f7358SJed Brown   PetscMPIInt       rank, size;
214552f7358SJed Brown   PetscErrorCode    ierr;
215552f7358SJed Brown   Vec               pointVec;
2163a93e3b7SToby Isaac   PetscSF           cellSF;
217552f7358SJed Brown   PetscLayout       layout;
218552f7358SJed Brown   PetscReal         *globalPoints;
219cb313848SJed Brown   PetscScalar       *globalPointsScalar;
220552f7358SJed Brown   const PetscInt    *ranges;
221552f7358SJed Brown   PetscMPIInt       *counts, *displs;
2223a93e3b7SToby Isaac   const PetscSFNode *foundCells;
2233a93e3b7SToby Isaac   const PetscInt    *foundPoints;
224552f7358SJed Brown   PetscMPIInt       *foundProcs, *globalProcs;
2253a93e3b7SToby Isaac   PetscInt          n, N, numFound;
226552f7358SJed Brown 
22719436ca2SJed Brown   PetscFunctionBegin;
22819436ca2SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
22919436ca2SJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
23019436ca2SJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2310adebc6cSBarry Smith   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
23219436ca2SJed Brown   /* Locate points */
23319436ca2SJed Brown   n = ctx->nInput;
234552f7358SJed Brown   if (!redundantPoints) {
235552f7358SJed Brown     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
236552f7358SJed Brown     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
237552f7358SJed Brown     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
238552f7358SJed Brown     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
239552f7358SJed Brown     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
240552f7358SJed Brown     /* Communicate all points to all processes */
241dcca6d9dSJed Brown     ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr);
242552f7358SJed Brown     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
243552f7358SJed Brown     for (p = 0; p < size; ++p) {
244552f7358SJed Brown       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
245552f7358SJed Brown       displs[p] = ranges[p]*ctx->dim;
246552f7358SJed Brown     }
247552f7358SJed Brown     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr);
248552f7358SJed Brown   } else {
249552f7358SJed Brown     N = n;
250552f7358SJed Brown     globalPoints = ctx->points;
25138ea73c8SJed Brown     counts = displs = NULL;
25238ea73c8SJed Brown     layout = NULL;
253552f7358SJed Brown   }
254552f7358SJed Brown #if 0
255dcca6d9dSJed Brown   ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
25619436ca2SJed Brown   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
257552f7358SJed Brown #else
258cb313848SJed Brown #if defined(PETSC_USE_COMPLEX)
25946b3086cSToby Isaac   ierr = PetscMalloc1(N*ctx->dim,&globalPointsScalar);CHKERRQ(ierr);
26046b3086cSToby Isaac   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
261cb313848SJed Brown #else
262cb313848SJed Brown   globalPointsScalar = globalPoints;
263cb313848SJed Brown #endif
26404706141SMatthew G Knepley   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
265dcca6d9dSJed Brown   ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
2667b5f2079SMatthew G. Knepley   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
2673a93e3b7SToby Isaac   cellSF = NULL;
2682d1fa6caSMatthew G. Knepley   ierr = DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);CHKERRQ(ierr);
2693a93e3b7SToby Isaac   ierr = PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);CHKERRQ(ierr);
270552f7358SJed Brown #endif
2713a93e3b7SToby Isaac   for (p = 0; p < numFound; ++p) {
2723a93e3b7SToby Isaac     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
273552f7358SJed Brown   }
274552f7358SJed Brown   /* Let the lowest rank process own each point */
275b2566f29SBarry Smith   ierr   = MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
276552f7358SJed Brown   ctx->n = 0;
277552f7358SJed Brown   for (p = 0; p < N; ++p) {
278e5b268a4SToby Isaac     if (globalProcs[p] == size) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %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));
2791aa26658SKarl Rupp     else if (globalProcs[p] == rank) ctx->n++;
280552f7358SJed Brown   }
281552f7358SJed Brown   /* Create coordinates vector and array of owned cells */
282785e854fSJed Brown   ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr);
283552f7358SJed Brown   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
284552f7358SJed Brown   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
285552f7358SJed Brown   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
286c0dedaeaSBarry Smith   ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr);
287552f7358SJed Brown   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
288552f7358SJed Brown   for (p = 0, q = 0, i = 0; p < N; ++p) {
289552f7358SJed Brown     if (globalProcs[p] == rank) {
290552f7358SJed Brown       PetscInt d;
291552f7358SJed Brown 
2921aa26658SKarl Rupp       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
293f317b747SMatthew G. Knepley       ctx->cells[q] = foundCells[q].index;
294f317b747SMatthew G. Knepley       ++q;
295552f7358SJed Brown     }
296552f7358SJed Brown   }
297552f7358SJed Brown   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
298552f7358SJed Brown #if 0
299552f7358SJed Brown   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
300552f7358SJed Brown #else
301552f7358SJed Brown   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
3023a93e3b7SToby Isaac   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
303552f7358SJed Brown   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
304552f7358SJed Brown #endif
305cb313848SJed Brown   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
306d343d804SMatthew G. Knepley   if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);}
307552f7358SJed Brown   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
308552f7358SJed Brown   PetscFunctionReturn(0);
309552f7358SJed Brown }
310552f7358SJed Brown 
3114267b1a3SMatthew G. Knepley /*@C
3124267b1a3SMatthew G. Knepley   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
3134267b1a3SMatthew G. Knepley 
3144267b1a3SMatthew G. Knepley   Collective on ctx
3154267b1a3SMatthew G. Knepley 
3164267b1a3SMatthew G. Knepley   Input Parameter:
3174267b1a3SMatthew G. Knepley . ctx - the context
3184267b1a3SMatthew G. Knepley 
3194267b1a3SMatthew G. Knepley   Output Parameter:
3204267b1a3SMatthew G. Knepley . coordinates  - the coordinates of interpolation points
3214267b1a3SMatthew G. Knepley 
3224267b1a3SMatthew 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.
3234267b1a3SMatthew G. Knepley 
3244267b1a3SMatthew G. Knepley   Level: intermediate
3254267b1a3SMatthew G. Knepley 
3264267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
3274267b1a3SMatthew G. Knepley @*/
3280adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
3290adebc6cSBarry Smith {
330552f7358SJed Brown   PetscFunctionBegin;
331552f7358SJed Brown   PetscValidPointer(coordinates, 2);
3320adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
333552f7358SJed Brown   *coordinates = ctx->coords;
334552f7358SJed Brown   PetscFunctionReturn(0);
335552f7358SJed Brown }
336552f7358SJed Brown 
3374267b1a3SMatthew G. Knepley /*@C
3384267b1a3SMatthew G. Knepley   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
3394267b1a3SMatthew G. Knepley 
3404267b1a3SMatthew G. Knepley   Collective on ctx
3414267b1a3SMatthew G. Knepley 
3424267b1a3SMatthew G. Knepley   Input Parameter:
3434267b1a3SMatthew G. Knepley . ctx - the context
3444267b1a3SMatthew G. Knepley 
3454267b1a3SMatthew G. Knepley   Output Parameter:
3464267b1a3SMatthew G. Knepley . v  - a vector capable of holding the interpolated field values
3474267b1a3SMatthew G. Knepley 
3484267b1a3SMatthew G. Knepley   Note: This vector should be returned using DMInterpolationRestoreVector().
3494267b1a3SMatthew G. Knepley 
3504267b1a3SMatthew G. Knepley   Level: intermediate
3514267b1a3SMatthew G. Knepley 
3524267b1a3SMatthew G. Knepley .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
3534267b1a3SMatthew G. Knepley @*/
3540adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
3550adebc6cSBarry Smith {
356552f7358SJed Brown   PetscErrorCode ierr;
357552f7358SJed Brown 
358552f7358SJed Brown   PetscFunctionBegin;
359552f7358SJed Brown   PetscValidPointer(v, 2);
3600adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
361552f7358SJed Brown   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
362552f7358SJed Brown   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
363552f7358SJed Brown   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
364c0dedaeaSBarry Smith   ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr);
365552f7358SJed Brown   PetscFunctionReturn(0);
366552f7358SJed Brown }
367552f7358SJed Brown 
3684267b1a3SMatthew G. Knepley /*@C
3694267b1a3SMatthew G. Knepley   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
3704267b1a3SMatthew G. Knepley 
3714267b1a3SMatthew G. Knepley   Collective on ctx
3724267b1a3SMatthew G. Knepley 
3734267b1a3SMatthew G. Knepley   Input Parameters:
3744267b1a3SMatthew G. Knepley + ctx - the context
3754267b1a3SMatthew G. Knepley - v  - a vector capable of holding the interpolated field values
3764267b1a3SMatthew G. Knepley 
3774267b1a3SMatthew G. Knepley   Level: intermediate
3784267b1a3SMatthew G. Knepley 
3794267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
3804267b1a3SMatthew G. Knepley @*/
3810adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
3820adebc6cSBarry Smith {
383552f7358SJed Brown   PetscErrorCode ierr;
384552f7358SJed Brown 
385552f7358SJed Brown   PetscFunctionBegin;
386552f7358SJed Brown   PetscValidPointer(v, 2);
3870adebc6cSBarry Smith   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
388552f7358SJed Brown   ierr = VecDestroy(v);CHKERRQ(ierr);
389552f7358SJed Brown   PetscFunctionReturn(0);
390552f7358SJed Brown }
391552f7358SJed Brown 
3927a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
393a6dfd86eSKarl Rupp {
394552f7358SJed Brown   PetscReal      *v0, *J, *invJ, detJ;
39556044e6dSMatthew G. Knepley   const PetscScalar *coords;
39656044e6dSMatthew G. Knepley   PetscScalar    *a;
397552f7358SJed Brown   PetscInt       p;
398552f7358SJed Brown   PetscErrorCode ierr;
399552f7358SJed Brown 
400552f7358SJed Brown   PetscFunctionBegin;
401dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
40256044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
403552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
404552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
405552f7358SJed Brown     PetscInt     c = ctx->cells[p];
406a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
407552f7358SJed Brown     PetscReal    xi[4];
408552f7358SJed Brown     PetscInt     d, f, comp;
409552f7358SJed Brown 
4108e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
411ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
4120298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
4131aa26658SKarl Rupp     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
4141aa26658SKarl Rupp 
415552f7358SJed Brown     for (d = 0; d < ctx->dim; ++d) {
416552f7358SJed Brown       xi[d] = 0.0;
4171aa26658SKarl 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]);
4181aa26658SKarl 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];
419552f7358SJed Brown     }
4200298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
421552f7358SJed Brown   }
422552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
42356044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
424552f7358SJed Brown   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
425552f7358SJed Brown   PetscFunctionReturn(0);
426552f7358SJed Brown }
427552f7358SJed Brown 
4287a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
4297a1931ceSMatthew G. Knepley {
4307a1931ceSMatthew G. Knepley   PetscReal      *v0, *J, *invJ, detJ;
43156044e6dSMatthew G. Knepley   const PetscScalar *coords;
43256044e6dSMatthew G. Knepley   PetscScalar    *a;
4337a1931ceSMatthew G. Knepley   PetscInt       p;
4347a1931ceSMatthew G. Knepley   PetscErrorCode ierr;
4357a1931ceSMatthew G. Knepley 
4367a1931ceSMatthew G. Knepley   PetscFunctionBegin;
437dcca6d9dSJed Brown   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
43856044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
4397a1931ceSMatthew G. Knepley   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
4407a1931ceSMatthew G. Knepley   for (p = 0; p < ctx->n; ++p) {
4417a1931ceSMatthew G. Knepley     PetscInt       c = ctx->cells[p];
4427a1931ceSMatthew G. Knepley     const PetscInt order[3] = {2, 1, 3};
4432584bbe8SMatthew G. Knepley     PetscScalar   *x = NULL;
4447a1931ceSMatthew G. Knepley     PetscReal      xi[4];
4457a1931ceSMatthew G. Knepley     PetscInt       d, f, comp;
4467a1931ceSMatthew G. Knepley 
4478e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
448ff1e0c32SBarry Smith     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
4497a1931ceSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
4507a1931ceSMatthew G. Knepley     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
4517a1931ceSMatthew G. Knepley 
4527a1931ceSMatthew G. Knepley     for (d = 0; d < ctx->dim; ++d) {
4537a1931ceSMatthew G. Knepley       xi[d] = 0.0;
4547a1931ceSMatthew 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]);
4557a1931ceSMatthew 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];
4567a1931ceSMatthew G. Knepley     }
4577a1931ceSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
4587a1931ceSMatthew G. Knepley   }
4597a1931ceSMatthew G. Knepley   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
46056044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
4617a1931ceSMatthew G. Knepley   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
4627a1931ceSMatthew G. Knepley   PetscFunctionReturn(0);
4637a1931ceSMatthew G. Knepley }
4647a1931ceSMatthew G. Knepley 
4655820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
466552f7358SJed Brown {
467552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
468552f7358SJed Brown   const PetscScalar x0        = vertices[0];
469552f7358SJed Brown   const PetscScalar y0        = vertices[1];
470552f7358SJed Brown   const PetscScalar x1        = vertices[2];
471552f7358SJed Brown   const PetscScalar y1        = vertices[3];
472552f7358SJed Brown   const PetscScalar x2        = vertices[4];
473552f7358SJed Brown   const PetscScalar y2        = vertices[5];
474552f7358SJed Brown   const PetscScalar x3        = vertices[6];
475552f7358SJed Brown   const PetscScalar y3        = vertices[7];
476552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
477552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
478552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
479552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
480552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
481552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
48256044e6dSMatthew G. Knepley   const PetscScalar *ref;
48356044e6dSMatthew G. Knepley   PetscScalar       *real;
484552f7358SJed Brown   PetscErrorCode    ierr;
485552f7358SJed Brown 
486552f7358SJed Brown   PetscFunctionBegin;
48756044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
488552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
489552f7358SJed Brown   {
490552f7358SJed Brown     const PetscScalar p0 = ref[0];
491552f7358SJed Brown     const PetscScalar p1 = ref[1];
492552f7358SJed Brown 
493552f7358SJed Brown     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
494552f7358SJed Brown     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
495552f7358SJed Brown   }
496552f7358SJed Brown   ierr = PetscLogFlops(28);CHKERRQ(ierr);
49756044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
498552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
499552f7358SJed Brown   PetscFunctionReturn(0);
500552f7358SJed Brown }
501552f7358SJed Brown 
502af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
503d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
504552f7358SJed Brown {
505552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
506552f7358SJed Brown   const PetscScalar x0        = vertices[0];
507552f7358SJed Brown   const PetscScalar y0        = vertices[1];
508552f7358SJed Brown   const PetscScalar x1        = vertices[2];
509552f7358SJed Brown   const PetscScalar y1        = vertices[3];
510552f7358SJed Brown   const PetscScalar x2        = vertices[4];
511552f7358SJed Brown   const PetscScalar y2        = vertices[5];
512552f7358SJed Brown   const PetscScalar x3        = vertices[6];
513552f7358SJed Brown   const PetscScalar y3        = vertices[7];
514552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
515552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
51656044e6dSMatthew G. Knepley   const PetscScalar *ref;
517552f7358SJed Brown   PetscErrorCode    ierr;
518552f7358SJed Brown 
519552f7358SJed Brown   PetscFunctionBegin;
52056044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
521552f7358SJed Brown   {
522552f7358SJed Brown     const PetscScalar x       = ref[0];
523552f7358SJed Brown     const PetscScalar y       = ref[1];
524552f7358SJed Brown     const PetscInt    rows[2] = {0, 1};
525da80777bSKarl Rupp     PetscScalar       values[4];
526da80777bSKarl Rupp 
527da80777bSKarl Rupp     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
528da80777bSKarl Rupp     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
52994ab13aaSBarry Smith     ierr      = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
530552f7358SJed Brown   }
531552f7358SJed Brown   ierr = PetscLogFlops(30);CHKERRQ(ierr);
53256044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
53394ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
53494ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
535552f7358SJed Brown   PetscFunctionReturn(0);
536552f7358SJed Brown }
537552f7358SJed Brown 
538a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
539a6dfd86eSKarl Rupp {
540fafc0619SMatthew G Knepley   DM             dmCoord;
541de73a395SMatthew G. Knepley   PetscFE        fem = NULL;
542552f7358SJed Brown   SNES           snes;
543552f7358SJed Brown   KSP            ksp;
544552f7358SJed Brown   PC             pc;
545552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
546552f7358SJed Brown   Mat            J;
547ef0bb6c7SMatthew G. Knepley   PetscTabulation    T;
54856044e6dSMatthew G. Knepley   const PetscScalar *coords;
54956044e6dSMatthew G. Knepley   PetscScalar    *a;
550ef0bb6c7SMatthew G. Knepley   PetscReal      xir[2];
551de73a395SMatthew G. Knepley   PetscInt       Nf, p;
5525509d985SMatthew G. Knepley   const PetscInt dof = ctx->dof;
553552f7358SJed Brown   PetscErrorCode ierr;
554552f7358SJed Brown 
555552f7358SJed Brown   PetscFunctionBegin;
556de73a395SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
55744a7f3ddSMatthew G. Knepley   if (Nf) {ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fem);CHKERRQ(ierr);}
558552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
559fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
560552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
561552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
562552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
563552f7358SJed Brown   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
564c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
565552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
566552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
567552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
568552f7358SJed Brown   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
569552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
570552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
5710298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
5720298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
573552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
574552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
575552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
576552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
577552f7358SJed Brown 
57856044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
579552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
580ef0bb6c7SMatthew G. Knepley   ierr = PetscFECreateTabulation(fem, 1, 1, xir, 0, &T);CHKERRQ(ierr);
581552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
582a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
583552f7358SJed Brown     PetscScalar *xi;
584552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
585552f7358SJed Brown 
586552f7358SJed Brown     /* Can make this do all points at once */
5870298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
588ff1e0c32SBarry Smith     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2);
5890298fd71SBarry Smith     ierr   = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
5900298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
5910298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
592552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
593552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
594552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
595552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
596552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
597552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
598cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
599cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
6005509d985SMatthew G. Knepley     if (4*dof != xSize) {
6015509d985SMatthew G. Knepley       PetscInt d;
6021aa26658SKarl Rupp 
6035509d985SMatthew G. Knepley       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
604ef0bb6c7SMatthew G. Knepley       ierr = PetscFEComputeTabulation(fem, 1, xir, 0, T);CHKERRQ(ierr);
6055509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp) {
6065509d985SMatthew G. Knepley         a[p*dof+comp] = 0.0;
6075509d985SMatthew G. Knepley         for (d = 0; d < xSize/dof; ++d) {
608ef0bb6c7SMatthew G. Knepley           a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
6095509d985SMatthew G. Knepley         }
6105509d985SMatthew G. Knepley       }
6115509d985SMatthew G. Knepley     } else {
6125509d985SMatthew G. Knepley       for (comp = 0; comp < dof; ++comp)
6135509d985SMatthew 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];
6145509d985SMatthew G. Knepley     }
615552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
6160298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
6170298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
618552f7358SJed Brown   }
619ef0bb6c7SMatthew G. Knepley   ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
620552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
62156044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
622552f7358SJed Brown 
623552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
624552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
625552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
626552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
627552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
628552f7358SJed Brown   PetscFunctionReturn(0);
629552f7358SJed Brown }
630552f7358SJed Brown 
6315820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
632552f7358SJed Brown {
633552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
634552f7358SJed Brown   const PetscScalar x0        = vertices[0];
635552f7358SJed Brown   const PetscScalar y0        = vertices[1];
636552f7358SJed Brown   const PetscScalar z0        = vertices[2];
6377a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
6387a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
6397a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
640552f7358SJed Brown   const PetscScalar x2        = vertices[6];
641552f7358SJed Brown   const PetscScalar y2        = vertices[7];
642552f7358SJed Brown   const PetscScalar z2        = vertices[8];
6437a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
6447a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
6457a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
646552f7358SJed Brown   const PetscScalar x4        = vertices[12];
647552f7358SJed Brown   const PetscScalar y4        = vertices[13];
648552f7358SJed Brown   const PetscScalar z4        = vertices[14];
649552f7358SJed Brown   const PetscScalar x5        = vertices[15];
650552f7358SJed Brown   const PetscScalar y5        = vertices[16];
651552f7358SJed Brown   const PetscScalar z5        = vertices[17];
652552f7358SJed Brown   const PetscScalar x6        = vertices[18];
653552f7358SJed Brown   const PetscScalar y6        = vertices[19];
654552f7358SJed Brown   const PetscScalar z6        = vertices[20];
655552f7358SJed Brown   const PetscScalar x7        = vertices[21];
656552f7358SJed Brown   const PetscScalar y7        = vertices[22];
657552f7358SJed Brown   const PetscScalar z7        = vertices[23];
658552f7358SJed Brown   const PetscScalar f_1       = x1 - x0;
659552f7358SJed Brown   const PetscScalar g_1       = y1 - y0;
660552f7358SJed Brown   const PetscScalar h_1       = z1 - z0;
661552f7358SJed Brown   const PetscScalar f_3       = x3 - x0;
662552f7358SJed Brown   const PetscScalar g_3       = y3 - y0;
663552f7358SJed Brown   const PetscScalar h_3       = z3 - z0;
664552f7358SJed Brown   const PetscScalar f_4       = x4 - x0;
665552f7358SJed Brown   const PetscScalar g_4       = y4 - y0;
666552f7358SJed Brown   const PetscScalar h_4       = z4 - z0;
667552f7358SJed Brown   const PetscScalar f_01      = x2 - x1 - x3 + x0;
668552f7358SJed Brown   const PetscScalar g_01      = y2 - y1 - y3 + y0;
669552f7358SJed Brown   const PetscScalar h_01      = z2 - z1 - z3 + z0;
670552f7358SJed Brown   const PetscScalar f_12      = x7 - x3 - x4 + x0;
671552f7358SJed Brown   const PetscScalar g_12      = y7 - y3 - y4 + y0;
672552f7358SJed Brown   const PetscScalar h_12      = z7 - z3 - z4 + z0;
673552f7358SJed Brown   const PetscScalar f_02      = x5 - x1 - x4 + x0;
674552f7358SJed Brown   const PetscScalar g_02      = y5 - y1 - y4 + y0;
675552f7358SJed Brown   const PetscScalar h_02      = z5 - z1 - z4 + z0;
676552f7358SJed Brown   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
677552f7358SJed Brown   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
678552f7358SJed Brown   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
67956044e6dSMatthew G. Knepley   const PetscScalar *ref;
68056044e6dSMatthew G. Knepley   PetscScalar       *real;
681552f7358SJed Brown   PetscErrorCode    ierr;
682552f7358SJed Brown 
683552f7358SJed Brown   PetscFunctionBegin;
68456044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
685552f7358SJed Brown   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
686552f7358SJed Brown   {
687552f7358SJed Brown     const PetscScalar p0 = ref[0];
688552f7358SJed Brown     const PetscScalar p1 = ref[1];
689552f7358SJed Brown     const PetscScalar p2 = ref[2];
690552f7358SJed Brown 
691552f7358SJed 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;
692552f7358SJed 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;
693552f7358SJed 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;
694552f7358SJed Brown   }
695552f7358SJed Brown   ierr = PetscLogFlops(114);CHKERRQ(ierr);
69656044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
697552f7358SJed Brown   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
698552f7358SJed Brown   PetscFunctionReturn(0);
699552f7358SJed Brown }
700552f7358SJed Brown 
701d1e9a80fSBarry Smith PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
702552f7358SJed Brown {
703552f7358SJed Brown   const PetscScalar *vertices = (const PetscScalar*) ctx;
704552f7358SJed Brown   const PetscScalar x0        = vertices[0];
705552f7358SJed Brown   const PetscScalar y0        = vertices[1];
706552f7358SJed Brown   const PetscScalar z0        = vertices[2];
7077a1931ceSMatthew G. Knepley   const PetscScalar x1        = vertices[9];
7087a1931ceSMatthew G. Knepley   const PetscScalar y1        = vertices[10];
7097a1931ceSMatthew G. Knepley   const PetscScalar z1        = vertices[11];
710552f7358SJed Brown   const PetscScalar x2        = vertices[6];
711552f7358SJed Brown   const PetscScalar y2        = vertices[7];
712552f7358SJed Brown   const PetscScalar z2        = vertices[8];
7137a1931ceSMatthew G. Knepley   const PetscScalar x3        = vertices[3];
7147a1931ceSMatthew G. Knepley   const PetscScalar y3        = vertices[4];
7157a1931ceSMatthew G. Knepley   const PetscScalar z3        = vertices[5];
716552f7358SJed Brown   const PetscScalar x4        = vertices[12];
717552f7358SJed Brown   const PetscScalar y4        = vertices[13];
718552f7358SJed Brown   const PetscScalar z4        = vertices[14];
719552f7358SJed Brown   const PetscScalar x5        = vertices[15];
720552f7358SJed Brown   const PetscScalar y5        = vertices[16];
721552f7358SJed Brown   const PetscScalar z5        = vertices[17];
722552f7358SJed Brown   const PetscScalar x6        = vertices[18];
723552f7358SJed Brown   const PetscScalar y6        = vertices[19];
724552f7358SJed Brown   const PetscScalar z6        = vertices[20];
725552f7358SJed Brown   const PetscScalar x7        = vertices[21];
726552f7358SJed Brown   const PetscScalar y7        = vertices[22];
727552f7358SJed Brown   const PetscScalar z7        = vertices[23];
728552f7358SJed Brown   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
729552f7358SJed Brown   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
730552f7358SJed Brown   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
731552f7358SJed Brown   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
732552f7358SJed Brown   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
733552f7358SJed Brown   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
734552f7358SJed Brown   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
735552f7358SJed Brown   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
736552f7358SJed Brown   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
737552f7358SJed Brown   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
738552f7358SJed Brown   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
739552f7358SJed Brown   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
74056044e6dSMatthew G. Knepley   const PetscScalar *ref;
741552f7358SJed Brown   PetscErrorCode    ierr;
742552f7358SJed Brown 
743552f7358SJed Brown   PetscFunctionBegin;
74456044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(Xref,  &ref);CHKERRQ(ierr);
745552f7358SJed Brown   {
746552f7358SJed Brown     const PetscScalar x       = ref[0];
747552f7358SJed Brown     const PetscScalar y       = ref[1];
748552f7358SJed Brown     const PetscScalar z       = ref[2];
749552f7358SJed Brown     const PetscInt    rows[3] = {0, 1, 2};
750da80777bSKarl Rupp     PetscScalar       values[9];
751da80777bSKarl Rupp 
752da80777bSKarl Rupp     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
753da80777bSKarl Rupp     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
754da80777bSKarl Rupp     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
755da80777bSKarl Rupp     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
756da80777bSKarl Rupp     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
757da80777bSKarl Rupp     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
758da80777bSKarl Rupp     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
759da80777bSKarl Rupp     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
760da80777bSKarl Rupp     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
7611aa26658SKarl Rupp 
76294ab13aaSBarry Smith     ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
763552f7358SJed Brown   }
764552f7358SJed Brown   ierr = PetscLogFlops(152);CHKERRQ(ierr);
76556044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(Xref,  &ref);CHKERRQ(ierr);
76694ab13aaSBarry Smith   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76794ab13aaSBarry Smith   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
768552f7358SJed Brown   PetscFunctionReturn(0);
769552f7358SJed Brown }
770552f7358SJed Brown 
771a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
772a6dfd86eSKarl Rupp {
773fafc0619SMatthew G Knepley   DM             dmCoord;
774552f7358SJed Brown   SNES           snes;
775552f7358SJed Brown   KSP            ksp;
776552f7358SJed Brown   PC             pc;
777552f7358SJed Brown   Vec            coordsLocal, r, ref, real;
778552f7358SJed Brown   Mat            J;
77956044e6dSMatthew G. Knepley   const PetscScalar *coords;
78056044e6dSMatthew G. Knepley   PetscScalar    *a;
781552f7358SJed Brown   PetscInt       p;
782552f7358SJed Brown   PetscErrorCode ierr;
783552f7358SJed Brown 
784552f7358SJed Brown   PetscFunctionBegin;
785552f7358SJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
786fafc0619SMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
787552f7358SJed Brown   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
788552f7358SJed Brown   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
789552f7358SJed Brown   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
790552f7358SJed Brown   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
791c0dedaeaSBarry Smith   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
792552f7358SJed Brown   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
793552f7358SJed Brown   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
794552f7358SJed Brown   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
795552f7358SJed Brown   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
796552f7358SJed Brown   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
797552f7358SJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
7980298fd71SBarry Smith   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
7990298fd71SBarry Smith   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
800552f7358SJed Brown   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
801552f7358SJed Brown   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
802552f7358SJed Brown   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
803552f7358SJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
804552f7358SJed Brown 
80556044e6dSMatthew G. Knepley   ierr = VecGetArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
806552f7358SJed Brown   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
807552f7358SJed Brown   for (p = 0; p < ctx->n; ++p) {
808a1e44745SMatthew G. Knepley     PetscScalar *x = NULL, *vertices = NULL;
809552f7358SJed Brown     PetscScalar *xi;
810cb313848SJed Brown     PetscReal    xir[3];
811552f7358SJed Brown     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
812552f7358SJed Brown 
813552f7358SJed Brown     /* Can make this do all points at once */
8140298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
815ff1e0c32SBarry Smith     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3);
8160298fd71SBarry Smith     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
817ff1e0c32SBarry Smith     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof);
8180298fd71SBarry Smith     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
8190298fd71SBarry Smith     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
820552f7358SJed Brown     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
821552f7358SJed Brown     xi[0]  = coords[p*ctx->dim+0];
822552f7358SJed Brown     xi[1]  = coords[p*ctx->dim+1];
823552f7358SJed Brown     xi[2]  = coords[p*ctx->dim+2];
824552f7358SJed Brown     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
825552f7358SJed Brown     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
826552f7358SJed Brown     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
827cb313848SJed Brown     xir[0] = PetscRealPart(xi[0]);
828cb313848SJed Brown     xir[1] = PetscRealPart(xi[1]);
829cb313848SJed Brown     xir[2] = PetscRealPart(xi[2]);
830552f7358SJed Brown     for (comp = 0; comp < ctx->dof; ++comp) {
831552f7358SJed Brown       a[p*ctx->dof+comp] =
832cb313848SJed Brown         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
8337a1931ceSMatthew G. Knepley         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
834cb313848SJed Brown         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
8357a1931ceSMatthew G. Knepley         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
836cb313848SJed Brown         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
837cb313848SJed Brown         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
838cb313848SJed Brown         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
839cb313848SJed Brown         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
840552f7358SJed Brown     }
841552f7358SJed Brown     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
8420298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
8430298fd71SBarry Smith     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
844552f7358SJed Brown   }
845552f7358SJed Brown   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
84656044e6dSMatthew G. Knepley   ierr = VecRestoreArrayRead(ctx->coords, &coords);CHKERRQ(ierr);
847552f7358SJed Brown 
848552f7358SJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
849552f7358SJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
850552f7358SJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
851552f7358SJed Brown   ierr = VecDestroy(&real);CHKERRQ(ierr);
852552f7358SJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
853552f7358SJed Brown   PetscFunctionReturn(0);
854552f7358SJed Brown }
855552f7358SJed Brown 
8564267b1a3SMatthew G. Knepley /*@C
8574267b1a3SMatthew G. Knepley   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
8584267b1a3SMatthew G. Knepley 
859552f7358SJed Brown   Input Parameters:
860552f7358SJed Brown + ctx - The DMInterpolationInfo context
861552f7358SJed Brown . dm  - The DM
862552f7358SJed Brown - x   - The local vector containing the field to be interpolated
863552f7358SJed Brown 
864552f7358SJed Brown   Output Parameters:
865552f7358SJed Brown . v   - The vector containing the interpolated values
8664267b1a3SMatthew G. Knepley 
8674267b1a3SMatthew G. Knepley   Note: A suitable v can be obtained using DMInterpolationGetVector().
8684267b1a3SMatthew G. Knepley 
8694267b1a3SMatthew G. Knepley   Level: beginner
8704267b1a3SMatthew G. Knepley 
8714267b1a3SMatthew G. Knepley .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
8724267b1a3SMatthew G. Knepley @*/
8730adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
8740adebc6cSBarry Smith {
875552f7358SJed Brown   PetscInt       dim, coneSize, n;
876552f7358SJed Brown   PetscErrorCode ierr;
877552f7358SJed Brown 
878552f7358SJed Brown   PetscFunctionBegin;
879552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
880552f7358SJed Brown   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
881552f7358SJed Brown   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
882552f7358SJed Brown   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
883ff1e0c32SBarry Smith   if (n != ctx->n*ctx->dof) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %D should be %D", n, ctx->n*ctx->dof);
884552f7358SJed Brown   if (n) {
885c73cfb54SMatthew G. Knepley     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
886552f7358SJed Brown     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
887552f7358SJed Brown     if (dim == 2) {
888552f7358SJed Brown       if (coneSize == 3) {
8897a1931ceSMatthew G. Knepley         ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr);
890552f7358SJed Brown       } else if (coneSize == 4) {
891552f7358SJed Brown         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
892ff1e0c32SBarry Smith       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %D for point interpolation", dim);
893552f7358SJed Brown     } else if (dim == 3) {
894552f7358SJed Brown       if (coneSize == 4) {
8957a1931ceSMatthew G. Knepley         ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr);
896552f7358SJed Brown       } else {
897552f7358SJed Brown         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
898552f7358SJed Brown       }
899ff1e0c32SBarry Smith     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %D for point interpolation", dim);
900552f7358SJed Brown   }
901552f7358SJed Brown   PetscFunctionReturn(0);
902552f7358SJed Brown }
903552f7358SJed Brown 
9044267b1a3SMatthew G. Knepley /*@C
9054267b1a3SMatthew G. Knepley   DMInterpolationDestroy - Destroys a DMInterpolationInfo context
9064267b1a3SMatthew G. Knepley 
9074267b1a3SMatthew G. Knepley   Collective on ctx
9084267b1a3SMatthew G. Knepley 
9094267b1a3SMatthew G. Knepley   Input Parameter:
9104267b1a3SMatthew G. Knepley . ctx - the context
9114267b1a3SMatthew G. Knepley 
9124267b1a3SMatthew G. Knepley   Level: beginner
9134267b1a3SMatthew G. Knepley 
9144267b1a3SMatthew G. Knepley .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
9154267b1a3SMatthew G. Knepley @*/
9160adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
9170adebc6cSBarry Smith {
918552f7358SJed Brown   PetscErrorCode ierr;
919552f7358SJed Brown 
920552f7358SJed Brown   PetscFunctionBegin;
921552f7358SJed Brown   PetscValidPointer(ctx, 2);
922552f7358SJed Brown   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
923552f7358SJed Brown   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
924552f7358SJed Brown   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
925552f7358SJed Brown   ierr = PetscFree(*ctx);CHKERRQ(ierr);
9260298fd71SBarry Smith   *ctx = NULL;
927552f7358SJed Brown   PetscFunctionReturn(0);
928552f7358SJed Brown }
929cc0c4584SMatthew G. Knepley 
930cc0c4584SMatthew G. Knepley /*@C
931cc0c4584SMatthew G. Knepley   SNESMonitorFields - Monitors the residual for each field separately
932cc0c4584SMatthew G. Knepley 
933cc0c4584SMatthew G. Knepley   Collective on SNES
934cc0c4584SMatthew G. Knepley 
935cc0c4584SMatthew G. Knepley   Input Parameters:
936cc0c4584SMatthew G. Knepley + snes   - the SNES context
937cc0c4584SMatthew G. Knepley . its    - iteration number
938cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual
939d43b4f6eSBarry Smith - vf  - PetscViewerAndFormat of type ASCII
940cc0c4584SMatthew G. Knepley 
941cc0c4584SMatthew G. Knepley   Notes:
942cc0c4584SMatthew G. Knepley   This routine prints the residual norm at each iteration.
943cc0c4584SMatthew G. Knepley 
944cc0c4584SMatthew G. Knepley   Level: intermediate
945cc0c4584SMatthew G. Knepley 
946cc0c4584SMatthew G. Knepley .seealso: SNESMonitorSet(), SNESMonitorDefault()
947cc0c4584SMatthew G. Knepley @*/
948d43b4f6eSBarry Smith PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
949cc0c4584SMatthew G. Knepley {
950d43b4f6eSBarry Smith   PetscViewer        viewer = vf->viewer;
951cc0c4584SMatthew G. Knepley   Vec                res;
952cc0c4584SMatthew G. Knepley   DM                 dm;
953cc0c4584SMatthew G. Knepley   PetscSection       s;
954cc0c4584SMatthew G. Knepley   const PetscScalar *r;
955cc0c4584SMatthew G. Knepley   PetscReal         *lnorms, *norms;
956cc0c4584SMatthew G. Knepley   PetscInt           numFields, f, pStart, pEnd, p;
957cc0c4584SMatthew G. Knepley   PetscErrorCode     ierr;
958cc0c4584SMatthew G. Knepley 
959cc0c4584SMatthew G. Knepley   PetscFunctionBegin;
9604d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
9619e5d0892SLisandro Dalcin   ierr = SNESGetFunction(snes, &res, NULL, NULL);CHKERRQ(ierr);
962cc0c4584SMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
96392fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr);
964cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
965cc0c4584SMatthew G. Knepley   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
966cc0c4584SMatthew G. Knepley   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
967cc0c4584SMatthew G. Knepley   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
968cc0c4584SMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
969cc0c4584SMatthew G. Knepley     for (f = 0; f < numFields; ++f) {
970cc0c4584SMatthew G. Knepley       PetscInt fdof, foff, d;
971cc0c4584SMatthew G. Knepley 
972cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
973cc0c4584SMatthew G. Knepley       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
974cc0c4584SMatthew G. Knepley       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
975cc0c4584SMatthew G. Knepley     }
976cc0c4584SMatthew G. Knepley   }
977cc0c4584SMatthew G. Knepley   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
978b2566f29SBarry Smith   ierr = MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
979d43b4f6eSBarry Smith   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
980cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
981cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
982cc0c4584SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
983cc0c4584SMatthew G. Knepley     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
984cc0c4584SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
985cc0c4584SMatthew G. Knepley   }
986cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
987cc0c4584SMatthew G. Knepley   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
988d43b4f6eSBarry Smith   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
989cc0c4584SMatthew G. Knepley   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
990cc0c4584SMatthew G. Knepley   PetscFunctionReturn(0);
991cc0c4584SMatthew G. Knepley }
99224cdb843SMatthew G. Knepley 
99324cdb843SMatthew G. Knepley /********************* Residual Computation **************************/
99424cdb843SMatthew G. Knepley 
99524cdb843SMatthew G. Knepley /*@
99624cdb843SMatthew G. Knepley   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
99724cdb843SMatthew G. Knepley 
99824cdb843SMatthew G. Knepley   Input Parameters:
99924cdb843SMatthew G. Knepley + dm - The mesh
100024cdb843SMatthew G. Knepley . X  - Local solution
100124cdb843SMatthew G. Knepley - user - The user context
100224cdb843SMatthew G. Knepley 
100324cdb843SMatthew G. Knepley   Output Parameter:
100424cdb843SMatthew G. Knepley . F  - Local output vector
100524cdb843SMatthew G. Knepley 
100624cdb843SMatthew G. Knepley   Level: developer
100724cdb843SMatthew G. Knepley 
10087a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
100924cdb843SMatthew G. Knepley @*/
101024cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
101124cdb843SMatthew G. Knepley {
10126da023fcSToby Isaac   DM             plex;
1013083401c6SMatthew G. Knepley   IS             allcellIS;
1014083401c6SMatthew G. Knepley   PetscInt       Nds, s, depth;
101524cdb843SMatthew G. Knepley   PetscErrorCode ierr;
101624cdb843SMatthew G. Knepley 
101724cdb843SMatthew G. Knepley   PetscFunctionBegin;
1018083401c6SMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
10196da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
10204a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
1021083401c6SMatthew G. Knepley   ierr = DMGetStratumIS(plex, "dim", depth, &allcellIS);CHKERRQ(ierr);
1022083401c6SMatthew G. Knepley   if (!allcellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &allcellIS);CHKERRQ(ierr);}
1023083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1024083401c6SMatthew G. Knepley     PetscDS ds;
1025083401c6SMatthew G. Knepley     DMLabel label;
1026083401c6SMatthew G. Knepley     IS      cellIS;
1027083401c6SMatthew G. Knepley 
1028083401c6SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1029083401c6SMatthew G. Knepley     if (!label) {
1030083401c6SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1031083401c6SMatthew G. Knepley       cellIS = allcellIS;
1032083401c6SMatthew G. Knepley     } else {
1033083401c6SMatthew G. Knepley       IS pointIS;
1034083401c6SMatthew G. Knepley 
1035083401c6SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1036083401c6SMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1037083401c6SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
10384a3e9fdbSToby Isaac     }
10394bee2e38SMatthew G. Knepley     ierr = DMPlexComputeResidual_Internal(plex, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);CHKERRQ(ierr);
10404a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1041083401c6SMatthew G. Knepley   }
1042083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
10439a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
104424cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
104524cdb843SMatthew G. Knepley }
104624cdb843SMatthew G. Knepley 
1047bdd6f66aSToby Isaac /*@
1048bdd6f66aSToby Isaac   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1049bdd6f66aSToby Isaac 
1050bdd6f66aSToby Isaac   Input Parameters:
1051bdd6f66aSToby Isaac + dm - The mesh
1052bdd6f66aSToby Isaac - user - The user context
1053bdd6f66aSToby Isaac 
1054bdd6f66aSToby Isaac   Output Parameter:
1055bdd6f66aSToby Isaac . X  - Local solution
1056bdd6f66aSToby Isaac 
1057bdd6f66aSToby Isaac   Level: developer
1058bdd6f66aSToby Isaac 
10597a73cf09SMatthew G. Knepley .seealso: DMPlexComputeJacobianAction()
1060bdd6f66aSToby Isaac @*/
1061bdd6f66aSToby Isaac PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1062bdd6f66aSToby Isaac {
1063bdd6f66aSToby Isaac   DM             plex;
1064bdd6f66aSToby Isaac   PetscErrorCode ierr;
1065bdd6f66aSToby Isaac 
1066bdd6f66aSToby Isaac   PetscFunctionBegin;
1067bdd6f66aSToby Isaac   ierr = DMSNESConvertPlex(dm,&plex,PETSC_TRUE);CHKERRQ(ierr);
1068bdd6f66aSToby Isaac   ierr = DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);CHKERRQ(ierr);
1069bdd6f66aSToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1070bdd6f66aSToby Isaac   PetscFunctionReturn(0);
1071bdd6f66aSToby Isaac }
1072bdd6f66aSToby Isaac 
10737a73cf09SMatthew G. Knepley /*@
10747a73cf09SMatthew G. Knepley   DMPlexComputeJacobianAction - Form the local portion of the Jacobian action Z = J(X) Y at the local solution X using pointwise functions specified by the user.
10757a73cf09SMatthew G. Knepley 
10767a73cf09SMatthew G. Knepley   Input Parameters:
10777a73cf09SMatthew G. Knepley + dm - The mesh
10787a73cf09SMatthew G. Knepley . cellIS -
10797a73cf09SMatthew G. Knepley . t  - The time
10807a73cf09SMatthew G. Knepley . X_tShift - The multiplier for the Jacobian with repsect to X_t
10817a73cf09SMatthew G. Knepley . X  - Local solution vector
10827a73cf09SMatthew G. Knepley . X_t  - Time-derivative of the local solution vector
10837a73cf09SMatthew G. Knepley . Y  - Local input vector
10847a73cf09SMatthew G. Knepley - user - The user context
10857a73cf09SMatthew G. Knepley 
10867a73cf09SMatthew G. Knepley   Output Parameter:
10877a73cf09SMatthew G. Knepley . Z - Local output vector
10887a73cf09SMatthew G. Knepley 
10897a73cf09SMatthew G. Knepley   Note:
10907a73cf09SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
10917a73cf09SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
10927a73cf09SMatthew G. Knepley 
10937a73cf09SMatthew G. Knepley   Level: developer
10947a73cf09SMatthew G. Knepley 
10957a73cf09SMatthew G. Knepley .seealso: FormFunctionLocal()
10967a73cf09SMatthew G. Knepley @*/
10977a73cf09SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
1098a925c78cSMatthew G. Knepley {
1099a925c78cSMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
1100a925c78cSMatthew G. Knepley   const char       *name  = "Jacobian";
11017a73cf09SMatthew G. Knepley   DM                dmAux, plex, plexAux = NULL;
1102a6e0b375SMatthew G. Knepley   DMEnclosureType   encAux;
1103c330f8ffSToby Isaac   Vec               A;
1104a925c78cSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
1105a925c78cSMatthew G. Knepley   PetscQuadrature   quad;
1106a925c78cSMatthew G. Knepley   PetscSection      section, globalSection, sectionAux;
1107a925c78cSMatthew G. Knepley   PetscScalar      *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
11089044fa66SMatthew G. Knepley   PetscInt          Nf, fieldI, fieldJ;
11094d0b9603SSander Arens   PetscInt          totDim, totDimAux = 0;
11109044fa66SMatthew G. Knepley   const PetscInt   *cells;
11119044fa66SMatthew G. Knepley   PetscInt          cStart, cEnd, numCells, c;
111275b37a90SMatthew G. Knepley   PetscBool         hasDyn;
1113c330f8ffSToby Isaac   DMField           coordField;
1114a925c78cSMatthew G. Knepley   PetscErrorCode    ierr;
1115a925c78cSMatthew G. Knepley 
1116a925c78cSMatthew G. Knepley   PetscFunctionBegin;
1117a925c78cSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
11187a73cf09SMatthew G. Knepley   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
11197a73cf09SMatthew G. Knepley   if (!cellIS) {
11207a73cf09SMatthew G. Knepley     PetscInt depth;
11217a73cf09SMatthew G. Knepley 
11227a73cf09SMatthew G. Knepley     ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
11237a73cf09SMatthew G. Knepley     ierr = DMGetStratumIS(plex, "dim", depth, &cellIS);CHKERRQ(ierr);
11247a73cf09SMatthew G. Knepley     if (!cellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &cellIS);CHKERRQ(ierr);}
11257a73cf09SMatthew G. Knepley   } else {
11267a73cf09SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) cellIS);CHKERRQ(ierr);
11277a73cf09SMatthew G. Knepley   }
112892fd8e1eSJed Brown   ierr = DMGetLocalSection(dm, &section);CHKERRQ(ierr);
1129e87a4003SBarry Smith   ierr = DMGetGlobalSection(dm, &globalSection);CHKERRQ(ierr);
1130d28747ceSMatthew G. Knepley   ierr = ISGetLocalSize(cellIS, &numCells);CHKERRQ(ierr);
1131d28747ceSMatthew G. Knepley   ierr = ISGetPointRange(cellIS, &cStart, &cEnd, &cells);CHKERRQ(ierr);
1132d28747ceSMatthew G. Knepley   ierr = DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob);CHKERRQ(ierr);
1133d28747ceSMatthew G. Knepley   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
1134a925c78cSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1135a925c78cSMatthew G. Knepley   ierr = PetscDSHasDynamicJacobian(prob, &hasDyn);CHKERRQ(ierr);
1136a925c78cSMatthew G. Knepley   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
1137a925c78cSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
1138a925c78cSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
1139a925c78cSMatthew G. Knepley   if (dmAux) {
1140a6e0b375SMatthew G. Knepley     ierr = DMGetEnclosureRelation(dmAux, dm, &encAux);CHKERRQ(ierr);
11417a73cf09SMatthew G. Knepley     ierr = DMConvert(dmAux, DMPLEX, &plexAux);CHKERRQ(ierr);
114292fd8e1eSJed Brown     ierr = DMGetLocalSection(plexAux, &sectionAux);CHKERRQ(ierr);
1143a925c78cSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
1144a925c78cSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
1145a925c78cSMatthew G. Knepley   }
1146a925c78cSMatthew G. Knepley   ierr = VecSet(Z, 0.0);CHKERRQ(ierr);
1147a925c78cSMatthew G. Knepley   ierr = PetscMalloc6(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim*totDim,&elemMat,hasDyn ? numCells*totDim*totDim : 0, &elemMatD,numCells*totDim,&y,totDim,&z);CHKERRQ(ierr);
1148a925c78cSMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
11494a3e9fdbSToby Isaac   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
1150a925c78cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
11519044fa66SMatthew G. Knepley     const PetscInt cell = cells ? cells[c] : c;
11529044fa66SMatthew G. Knepley     const PetscInt cind = c - cStart;
1153a925c78cSMatthew G. Knepley     PetscScalar   *x = NULL,  *x_t = NULL;
1154a925c78cSMatthew G. Knepley     PetscInt       i;
1155a925c78cSMatthew G. Knepley 
11569044fa66SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
11579044fa66SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
11589044fa66SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);CHKERRQ(ierr);
1159a925c78cSMatthew G. Knepley     if (X_t) {
11609044fa66SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
11619044fa66SMatthew G. Knepley       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
11629044fa66SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);CHKERRQ(ierr);
1163a925c78cSMatthew G. Knepley     }
1164a925c78cSMatthew G. Knepley     if (dmAux) {
116544171101SMatthew G. Knepley       PetscInt subcell;
1166a6e0b375SMatthew G. Knepley       ierr = DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);CHKERRQ(ierr);
116744171101SMatthew G. Knepley       ierr = DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
11689044fa66SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
116944171101SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);CHKERRQ(ierr);
1170a925c78cSMatthew G. Knepley     }
11719044fa66SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
11729044fa66SMatthew G. Knepley     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
11739044fa66SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);CHKERRQ(ierr);
1174a925c78cSMatthew G. Knepley   }
1175580bdb30SBarry Smith   ierr = PetscArrayzero(elemMat, numCells*totDim*totDim);CHKERRQ(ierr);
1176580bdb30SBarry Smith   if (hasDyn)  {ierr = PetscArrayzero(elemMatD, numCells*totDim*totDim);CHKERRQ(ierr);}
1177a925c78cSMatthew G. Knepley   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1178a925c78cSMatthew G. Knepley     PetscFE  fe;
1179c330f8ffSToby Isaac     PetscInt Nb;
1180a925c78cSMatthew G. Knepley     /* Conforming batches */
1181a925c78cSMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1182a925c78cSMatthew G. Knepley     /* Remainder */
1183c330f8ffSToby Isaac     PetscInt Nr, offset, Nq;
1184c330f8ffSToby Isaac     PetscQuadrature qGeom = NULL;
1185b7260050SToby Isaac     PetscInt    maxDegree;
1186c330f8ffSToby Isaac     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
1187a925c78cSMatthew G. Knepley 
1188a925c78cSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
1189a925c78cSMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
1190a925c78cSMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1191a925c78cSMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
1192b7260050SToby Isaac     ierr = DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);CHKERRQ(ierr);
1193b7260050SToby Isaac     if (maxDegree <= 1) {ierr = DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);CHKERRQ(ierr);}
1194c330f8ffSToby Isaac     if (!qGeom) {
1195c330f8ffSToby Isaac       ierr = PetscFEGetQuadrature(fe,&qGeom);CHKERRQ(ierr);
1196c330f8ffSToby Isaac       ierr = PetscObjectReference((PetscObject)qGeom);CHKERRQ(ierr);
1197c330f8ffSToby Isaac     }
1198c330f8ffSToby Isaac     ierr = PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
11994a3e9fdbSToby Isaac     ierr = DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1200c330f8ffSToby Isaac     blockSize = Nb;
1201a925c78cSMatthew G. Knepley     batchSize = numBlocks * blockSize;
1202a925c78cSMatthew G. Knepley     ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
1203a925c78cSMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
1204a925c78cSMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
1205a925c78cSMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
1206a925c78cSMatthew G. Knepley     offset    = numCells - Nr;
1207c330f8ffSToby Isaac     ierr = PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
1208c330f8ffSToby Isaac     ierr = PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
1209a925c78cSMatthew G. Knepley     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
12104bee2e38SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);CHKERRQ(ierr);
12114bee2e38SMatthew G. Knepley       ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);CHKERRQ(ierr);
1212a925c78cSMatthew G. Knepley       if (hasDyn) {
12134bee2e38SMatthew G. Knepley         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);CHKERRQ(ierr);
12144bee2e38SMatthew G. Knepley         ierr = PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);CHKERRQ(ierr);
1215a925c78cSMatthew G. Knepley       }
1216a925c78cSMatthew G. Knepley     }
1217c330f8ffSToby Isaac     ierr = PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);CHKERRQ(ierr);
1218c330f8ffSToby Isaac     ierr = PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);CHKERRQ(ierr);
12194a3e9fdbSToby Isaac     ierr = DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);CHKERRQ(ierr);
1220c330f8ffSToby Isaac     ierr = PetscQuadratureDestroy(&qGeom);CHKERRQ(ierr);
1221a925c78cSMatthew G. Knepley   }
1222a925c78cSMatthew G. Knepley   if (hasDyn) {
12239044fa66SMatthew G. Knepley     for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
1224a925c78cSMatthew G. Knepley   }
1225a925c78cSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12269044fa66SMatthew G. Knepley     const PetscInt     cell = cells ? cells[c] : c;
12279044fa66SMatthew G. Knepley     const PetscInt     cind = c - cStart;
1228a925c78cSMatthew G. Knepley     const PetscBLASInt M = totDim, one = 1;
1229a925c78cSMatthew G. Knepley     const PetscScalar  a = 1.0, b = 0.0;
1230a925c78cSMatthew G. Knepley 
12319044fa66SMatthew G. Knepley     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
1232a925c78cSMatthew G. Knepley     if (mesh->printFEM > 1) {
12339044fa66SMatthew G. Knepley       ierr = DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);CHKERRQ(ierr);
12349044fa66SMatthew G. Knepley       ierr = DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);CHKERRQ(ierr);
1235a925c78cSMatthew G. Knepley       ierr = DMPrintCellVector(c, "Z",  totDim, z);CHKERRQ(ierr);
1236a925c78cSMatthew G. Knepley     }
12379044fa66SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);CHKERRQ(ierr);
1238a925c78cSMatthew G. Knepley   }
1239a925c78cSMatthew G. Knepley   ierr = PetscFree6(u,u_t,elemMat,elemMatD,y,z);CHKERRQ(ierr);
1240a925c78cSMatthew G. Knepley   if (mesh->printFEM) {
1241ea13f565SStefano Zampini     ierr = PetscPrintf(PetscObjectComm((PetscObject)Z), "Z:\n");CHKERRQ(ierr);
1242a8c5bd36SStefano Zampini     ierr = VecView(Z, NULL);CHKERRQ(ierr);
1243a925c78cSMatthew G. Knepley   }
12447a73cf09SMatthew G. Knepley   ierr = PetscFree(a);CHKERRQ(ierr);
12457a73cf09SMatthew G. Knepley   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
12467a73cf09SMatthew G. Knepley   ierr = DMDestroy(&plexAux);CHKERRQ(ierr);
12477a73cf09SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1248a925c78cSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr);
1249a925c78cSMatthew G. Knepley   PetscFunctionReturn(0);
1250a925c78cSMatthew G. Knepley }
1251a925c78cSMatthew G. Knepley 
125224cdb843SMatthew G. Knepley /*@
125324cdb843SMatthew G. Knepley   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
125424cdb843SMatthew G. Knepley 
125524cdb843SMatthew G. Knepley   Input Parameters:
125624cdb843SMatthew G. Knepley + dm - The mesh
125724cdb843SMatthew G. Knepley . X  - Local input vector
125824cdb843SMatthew G. Knepley - user - The user context
125924cdb843SMatthew G. Knepley 
126024cdb843SMatthew G. Knepley   Output Parameter:
126124cdb843SMatthew G. Knepley . Jac  - Jacobian matrix
126224cdb843SMatthew G. Knepley 
126324cdb843SMatthew G. Knepley   Note:
126424cdb843SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
126524cdb843SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
126624cdb843SMatthew G. Knepley 
126724cdb843SMatthew G. Knepley   Level: developer
126824cdb843SMatthew G. Knepley 
126924cdb843SMatthew G. Knepley .seealso: FormFunctionLocal()
127024cdb843SMatthew G. Knepley @*/
127124cdb843SMatthew G. Knepley PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
127224cdb843SMatthew G. Knepley {
12736da023fcSToby Isaac   DM             plex;
1274083401c6SMatthew G. Knepley   IS             allcellIS;
1275f04eb4edSMatthew G. Knepley   PetscBool      hasJac, hasPrec;
1276083401c6SMatthew G. Knepley   PetscInt       Nds, s, depth;
127724cdb843SMatthew G. Knepley   PetscErrorCode ierr;
127824cdb843SMatthew G. Knepley 
127924cdb843SMatthew G. Knepley   PetscFunctionBegin;
1280083401c6SMatthew G. Knepley   ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
12816da023fcSToby Isaac   ierr = DMSNESConvertPlex(dm, &plex, PETSC_TRUE);CHKERRQ(ierr);
12824a3e9fdbSToby Isaac   ierr = DMPlexGetDepth(plex, &depth);CHKERRQ(ierr);
1283083401c6SMatthew G. Knepley   ierr = DMGetStratumIS(plex, "dim", depth, &allcellIS);CHKERRQ(ierr);
1284083401c6SMatthew G. Knepley   if (!allcellIS) {ierr = DMGetStratumIS(plex, "depth", depth, &allcellIS);CHKERRQ(ierr);}
1285083401c6SMatthew G. Knepley   for (s = 0; s < Nds; ++s) {
1286083401c6SMatthew G. Knepley     PetscDS ds;
1287083401c6SMatthew G. Knepley     DMLabel label;
1288083401c6SMatthew G. Knepley     IS      cellIS;
1289083401c6SMatthew G. Knepley 
1290083401c6SMatthew G. Knepley     ierr = DMGetRegionNumDS(dm, s, &label, NULL, &ds);CHKERRQ(ierr);
1291083401c6SMatthew G. Knepley     if (!label) {
1292083401c6SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) allcellIS);CHKERRQ(ierr);
1293083401c6SMatthew G. Knepley       cellIS = allcellIS;
1294083401c6SMatthew G. Knepley     } else {
1295083401c6SMatthew G. Knepley       IS pointIS;
1296083401c6SMatthew G. Knepley 
1297083401c6SMatthew G. Knepley       ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr);
1298083401c6SMatthew G. Knepley       ierr = ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);CHKERRQ(ierr);
1299083401c6SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1300083401c6SMatthew G. Knepley     }
1301083401c6SMatthew G. Knepley     if (!s) {
1302083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1303083401c6SMatthew G. Knepley       ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1304f04eb4edSMatthew G. Knepley       if (hasJac && hasPrec) {ierr = MatZeroEntries(Jac);CHKERRQ(ierr);}
1305f04eb4edSMatthew G. Knepley       ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
1306083401c6SMatthew G. Knepley     }
13074a3e9fdbSToby Isaac     ierr = DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);CHKERRQ(ierr);
13084a3e9fdbSToby Isaac     ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1309083401c6SMatthew G. Knepley   }
1310083401c6SMatthew G. Knepley   ierr = ISDestroy(&allcellIS);CHKERRQ(ierr);
13119a81d013SToby Isaac   ierr = DMDestroy(&plex);CHKERRQ(ierr);
131224cdb843SMatthew G. Knepley   PetscFunctionReturn(0);
131324cdb843SMatthew G. Knepley }
13141878804aSMatthew G. Knepley 
131538cfc46eSPierre Jolivet /*
131638cfc46eSPierre Jolivet      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
131738cfc46eSPierre Jolivet 
131838cfc46eSPierre Jolivet    Input Parameters:
131938cfc46eSPierre Jolivet +     X - SNES linearization point
132038cfc46eSPierre Jolivet .     ovl - index set of overlapping subdomains
132138cfc46eSPierre Jolivet 
132238cfc46eSPierre Jolivet    Output Parameter:
132338cfc46eSPierre Jolivet .     J - unassembled (Neumann) local matrix
132438cfc46eSPierre Jolivet 
132538cfc46eSPierre Jolivet    Level: intermediate
132638cfc46eSPierre Jolivet 
132738cfc46eSPierre Jolivet .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
132838cfc46eSPierre Jolivet */
132938cfc46eSPierre Jolivet static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
133038cfc46eSPierre Jolivet {
133138cfc46eSPierre Jolivet   SNES           snes;
133238cfc46eSPierre Jolivet   Mat            pJ;
133338cfc46eSPierre Jolivet   DM             ovldm,origdm;
133438cfc46eSPierre Jolivet   DMSNES         sdm;
133538cfc46eSPierre Jolivet   PetscErrorCode (*bfun)(DM,Vec,void*);
133638cfc46eSPierre Jolivet   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
133738cfc46eSPierre Jolivet   void           *bctx,*jctx;
133838cfc46eSPierre Jolivet   PetscErrorCode ierr;
133938cfc46eSPierre Jolivet 
134038cfc46eSPierre Jolivet   PetscFunctionBegin;
134138cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);CHKERRQ(ierr);
134238cfc46eSPierre Jolivet   if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");CHKERRQ(ierr);
134338cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);CHKERRQ(ierr);
134438cfc46eSPierre Jolivet   if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");CHKERRQ(ierr);
134538cfc46eSPierre Jolivet   ierr = MatGetDM(pJ,&ovldm);CHKERRQ(ierr);
134638cfc46eSPierre Jolivet   ierr = DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);CHKERRQ(ierr);
134738cfc46eSPierre Jolivet   ierr = DMSNESSetBoundaryLocal(ovldm,bfun,bctx);CHKERRQ(ierr);
134838cfc46eSPierre Jolivet   ierr = DMSNESGetJacobianLocal(origdm,&jfun,&jctx);CHKERRQ(ierr);
134938cfc46eSPierre Jolivet   ierr = DMSNESSetJacobianLocal(ovldm,jfun,jctx);CHKERRQ(ierr);
135038cfc46eSPierre Jolivet   ierr = PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);CHKERRQ(ierr);
135138cfc46eSPierre Jolivet   if (!snes) {
135238cfc46eSPierre Jolivet     ierr = SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);CHKERRQ(ierr);
135338cfc46eSPierre Jolivet     ierr = SNESSetDM(snes,ovldm);CHKERRQ(ierr);
135438cfc46eSPierre Jolivet     ierr = PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);CHKERRQ(ierr);
135538cfc46eSPierre Jolivet     ierr = PetscObjectDereference((PetscObject)snes);CHKERRQ(ierr);
135638cfc46eSPierre Jolivet   }
135738cfc46eSPierre Jolivet   ierr = DMGetDMSNES(ovldm,&sdm);CHKERRQ(ierr);
135838cfc46eSPierre Jolivet   ierr = VecLockReadPush(X);CHKERRQ(ierr);
135938cfc46eSPierre Jolivet   PetscStackPush("SNES user Jacobian function");
136038cfc46eSPierre Jolivet   ierr = (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);CHKERRQ(ierr);
136138cfc46eSPierre Jolivet   PetscStackPop;
136238cfc46eSPierre Jolivet   ierr = VecLockReadPop(X);CHKERRQ(ierr);
136338cfc46eSPierre Jolivet   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
136438cfc46eSPierre Jolivet   {
136538cfc46eSPierre Jolivet     Mat locpJ;
136638cfc46eSPierre Jolivet 
136738cfc46eSPierre Jolivet     ierr = MatISGetLocalMat(pJ,&locpJ);CHKERRQ(ierr);
136838cfc46eSPierre Jolivet     ierr = MatCopy(locpJ,J,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
136938cfc46eSPierre Jolivet   }
137038cfc46eSPierre Jolivet   PetscFunctionReturn(0);
137138cfc46eSPierre Jolivet }
137238cfc46eSPierre Jolivet 
1373a925c78cSMatthew G. Knepley /*@
13749f520fc2SToby Isaac   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
13759f520fc2SToby Isaac 
13769f520fc2SToby Isaac   Input Parameters:
13779f520fc2SToby Isaac + dm - The DM object
1378dff059c6SToby Isaac . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
13799f520fc2SToby Isaac . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
13809f520fc2SToby Isaac - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
13811a244344SSatish Balay 
13821a244344SSatish Balay   Level: developer
13839f520fc2SToby Isaac @*/
13849f520fc2SToby Isaac PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
13859f520fc2SToby Isaac {
13869f520fc2SToby Isaac   PetscErrorCode ierr;
13879f520fc2SToby Isaac 
13889f520fc2SToby Isaac   PetscFunctionBegin;
13899f520fc2SToby Isaac   ierr = DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);CHKERRQ(ierr);
13909f520fc2SToby Isaac   ierr = DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);CHKERRQ(ierr);
13919f520fc2SToby Isaac   ierr = DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);CHKERRQ(ierr);
139238cfc46eSPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);CHKERRQ(ierr);
13939f520fc2SToby Isaac   PetscFunctionReturn(0);
13949f520fc2SToby Isaac }
13959f520fc2SToby Isaac 
1396f017bcb6SMatthew G. Knepley /*@C
1397f017bcb6SMatthew G. Knepley   DMSNESCheckDiscretization - Check the discretization error of the exact solution
1398f017bcb6SMatthew G. Knepley 
1399f017bcb6SMatthew G. Knepley   Input Parameters:
1400f017bcb6SMatthew G. Knepley + snes - the SNES object
1401f017bcb6SMatthew G. Knepley . dm   - the DM
1402*f2cacb80SMatthew G. Knepley . t    - the time
1403f017bcb6SMatthew G. Knepley . u    - a DM vector
1404f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1405f017bcb6SMatthew G. Knepley 
1406f3c94560SMatthew G. Knepley   Output Parameters:
1407f3c94560SMatthew G. Knepley . error - An array which holds the discretization error in each field, or NULL
1408f3c94560SMatthew G. Knepley 
14097f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
14107f96f943SMatthew G. Knepley 
1411f017bcb6SMatthew G. Knepley   Level: developer
1412f017bcb6SMatthew G. Knepley 
14137f96f943SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian(), PetscDSSetExactSolution()
1414f017bcb6SMatthew G. Knepley @*/
1415*f2cacb80SMatthew G. Knepley PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
14161878804aSMatthew G. Knepley {
1417f017bcb6SMatthew G. Knepley   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1418f017bcb6SMatthew G. Knepley   void            **ectxs;
1419f3c94560SMatthew G. Knepley   PetscReal        *err;
14207f96f943SMatthew G. Knepley   MPI_Comm          comm;
14217f96f943SMatthew G. Knepley   PetscInt          Nf, f;
14221878804aSMatthew G. Knepley   PetscErrorCode    ierr;
14231878804aSMatthew G. Knepley 
14241878804aSMatthew G. Knepley   PetscFunctionBegin;
1425f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1426f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1427f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1428534a8f05SLisandro Dalcin   if (error) PetscValidRealPointer(error, 6);
14297f96f943SMatthew G. Knepley 
1430*f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
14317f96f943SMatthew G. Knepley   ierr = VecViewFromOptions(u, NULL, "-vec_view");CHKERRQ(ierr);
14327f96f943SMatthew G. Knepley 
1433f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1434f017bcb6SMatthew G. Knepley   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
1435083401c6SMatthew G. Knepley   ierr = PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);CHKERRQ(ierr);
14367f96f943SMatthew G. Knepley   {
14377f96f943SMatthew G. Knepley     PetscInt Nds, s;
14387f96f943SMatthew G. Knepley 
1439083401c6SMatthew G. Knepley     ierr = DMGetNumDS(dm, &Nds);CHKERRQ(ierr);
1440083401c6SMatthew G. Knepley     for (s = 0; s < Nds; ++s) {
14417f96f943SMatthew G. Knepley       PetscDS         ds;
1442083401c6SMatthew G. Knepley       DMLabel         label;
1443083401c6SMatthew G. Knepley       IS              fieldIS;
14447f96f943SMatthew G. Knepley       const PetscInt *fields;
14457f96f943SMatthew G. Knepley       PetscInt        dsNf, f;
1446083401c6SMatthew G. Knepley 
1447083401c6SMatthew G. Knepley       ierr = DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);CHKERRQ(ierr);
1448083401c6SMatthew G. Knepley       ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr);
1449083401c6SMatthew G. Knepley       ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);
1450083401c6SMatthew G. Knepley       for (f = 0; f < dsNf; ++f) {
1451083401c6SMatthew G. Knepley         const PetscInt field = fields[f];
1452083401c6SMatthew G. Knepley         ierr = PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]);CHKERRQ(ierr);
1453083401c6SMatthew G. Knepley       }
1454083401c6SMatthew G. Knepley       ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);
1455083401c6SMatthew G. Knepley     }
1456083401c6SMatthew G. Knepley   }
1457f017bcb6SMatthew G. Knepley   if (Nf > 1) {
1458*f2cacb80SMatthew G. Knepley     ierr = DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err);CHKERRQ(ierr);
1459f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1460f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1461f3c94560SMatthew G. Knepley         if (err[f] > tol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %D exceeds tolerance %g", (double) err[f], f, (double) tol);
14621878804aSMatthew G. Knepley       }
1463b878532fSMatthew G. Knepley     } else if (error) {
1464b878532fSMatthew G. Knepley       for (f = 0; f < Nf; ++f) error[f] = err[f];
14651878804aSMatthew G. Knepley     } else {
1466f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr);
1467f017bcb6SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
1468f017bcb6SMatthew G. Knepley         if (f) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
1469f3c94560SMatthew G. Knepley         ierr = PetscPrintf(comm, "%g", (double)err[f]);CHKERRQ(ierr);
14701878804aSMatthew G. Knepley       }
1471f017bcb6SMatthew G. Knepley       ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr);
1472f017bcb6SMatthew G. Knepley     }
1473f017bcb6SMatthew G. Knepley   } else {
1474*f2cacb80SMatthew G. Knepley     ierr = DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]);CHKERRQ(ierr);
1475f017bcb6SMatthew G. Knepley     if (tol >= 0.0) {
1476f3c94560SMatthew G. Knepley       if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1477b878532fSMatthew G. Knepley     } else if (error) {
1478b878532fSMatthew G. Knepley       error[0] = err[0];
1479f017bcb6SMatthew G. Knepley     } else {
1480f3c94560SMatthew G. Knepley       ierr = PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]);CHKERRQ(ierr);
1481f017bcb6SMatthew G. Knepley     }
1482f017bcb6SMatthew G. Knepley   }
1483f3c94560SMatthew G. Knepley   ierr = PetscFree3(exacts, ectxs, err);CHKERRQ(ierr);
1484f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1485f017bcb6SMatthew G. Knepley }
1486f017bcb6SMatthew G. Knepley 
1487f017bcb6SMatthew G. Knepley /*@C
1488f017bcb6SMatthew G. Knepley   DMSNESCheckResidual - Check the residual of the exact solution
1489f017bcb6SMatthew G. Knepley 
1490f017bcb6SMatthew G. Knepley   Input Parameters:
1491f017bcb6SMatthew G. Knepley + snes - the SNES object
1492f017bcb6SMatthew G. Knepley . dm   - the DM
1493f017bcb6SMatthew G. Knepley . u    - a DM vector
1494f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1495f017bcb6SMatthew G. Knepley 
1496f3c94560SMatthew G. Knepley   Output Parameters:
1497f3c94560SMatthew G. Knepley . residual - The residual norm of the exact solution, or NULL
1498f3c94560SMatthew G. Knepley 
1499f017bcb6SMatthew G. Knepley   Level: developer
1500f017bcb6SMatthew G. Knepley 
1501f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
1502f017bcb6SMatthew G. Knepley @*/
1503f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1504f017bcb6SMatthew G. Knepley {
1505f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1506f017bcb6SMatthew G. Knepley   Vec            r;
1507f017bcb6SMatthew G. Knepley   PetscReal      res;
1508f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1509f017bcb6SMatthew G. Knepley 
1510f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1511f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1512f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1513f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1514534a8f05SLisandro Dalcin   if (residual) PetscValidRealPointer(residual, 5);
1515f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1516*f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1517f017bcb6SMatthew G. Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
15181878804aSMatthew G. Knepley   ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
15191878804aSMatthew G. Knepley   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
1520f017bcb6SMatthew G. Knepley   if (tol >= 0.0) {
1521f017bcb6SMatthew G. Knepley     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1522b878532fSMatthew G. Knepley   } else if (residual) {
1523b878532fSMatthew G. Knepley     *residual = res;
1524f017bcb6SMatthew G. Knepley   } else {
1525f017bcb6SMatthew G. Knepley     ierr = PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
15261878804aSMatthew G. Knepley     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
15271878804aSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) r, "Initial Residual");CHKERRQ(ierr);
1528685405a1SBarry Smith     ierr = PetscObjectSetOptionsPrefix((PetscObject)r,"res_");CHKERRQ(ierr);
1529685405a1SBarry Smith     ierr = VecViewFromOptions(r, NULL, "-vec_view");CHKERRQ(ierr);
1530f017bcb6SMatthew G. Knepley   }
1531f017bcb6SMatthew G. Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
1532f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1533f017bcb6SMatthew G. Knepley }
1534f017bcb6SMatthew G. Knepley 
1535f017bcb6SMatthew G. Knepley /*@C
1536f3c94560SMatthew G. Knepley   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1537f017bcb6SMatthew G. Knepley 
1538f017bcb6SMatthew G. Knepley   Input Parameters:
1539f017bcb6SMatthew G. Knepley + snes - the SNES object
1540f017bcb6SMatthew G. Knepley . dm   - the DM
1541f017bcb6SMatthew G. Knepley . u    - a DM vector
1542f0fc11ceSJed Brown - tol  - A tolerance for the check, or -1 to print the results instead
1543f017bcb6SMatthew G. Knepley 
1544f3c94560SMatthew G. Knepley   Output Parameters:
1545f3c94560SMatthew G. Knepley + isLinear - Flag indicaing that the function looks linear, or NULL
1546f3c94560SMatthew G. Knepley - convRate - The rate of convergence of the linear model, or NULL
1547f3c94560SMatthew G. Knepley 
1548f017bcb6SMatthew G. Knepley   Level: developer
1549f017bcb6SMatthew G. Knepley 
1550f017bcb6SMatthew G. Knepley .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
1551f017bcb6SMatthew G. Knepley @*/
1552f3c94560SMatthew G. Knepley PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1553f017bcb6SMatthew G. Knepley {
1554f017bcb6SMatthew G. Knepley   MPI_Comm       comm;
1555f017bcb6SMatthew G. Knepley   PetscDS        ds;
1556f017bcb6SMatthew G. Knepley   Mat            J, M;
1557f017bcb6SMatthew G. Knepley   MatNullSpace   nullspace;
1558f3c94560SMatthew G. Knepley   PetscReal      slope, intercept;
15597f96f943SMatthew G. Knepley   PetscBool      hasJac, hasPrec, isLin = PETSC_FALSE;
1560f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1561f017bcb6SMatthew G. Knepley 
1562f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1563f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1564f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1565f3c94560SMatthew G. Knepley   PetscValidHeaderSpecific(u, VEC_CLASSID, 3);
1566534a8f05SLisandro Dalcin   if (isLinear) PetscValidBoolPointer(isLinear, 5);
1567534a8f05SLisandro Dalcin   if (convRate) PetscValidRealPointer(convRate, 5);
1568f017bcb6SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) snes, &comm);CHKERRQ(ierr);
1569*f2cacb80SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, 0.0, u, NULL);CHKERRQ(ierr);
1570f017bcb6SMatthew G. Knepley   /* Create and view matrices */
1571f017bcb6SMatthew G. Knepley   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
15727f96f943SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1573f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobian(ds, &hasJac);CHKERRQ(ierr);
1574f017bcb6SMatthew G. Knepley   ierr = PetscDSHasJacobianPreconditioner(ds, &hasPrec);CHKERRQ(ierr);
1575282e7bb4SMatthew G. Knepley   if (hasJac && hasPrec) {
1576282e7bb4SMatthew G. Knepley     ierr = DMCreateMatrix(dm, &M);CHKERRQ(ierr);
1577282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, M);CHKERRQ(ierr);
1578f017bcb6SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");CHKERRQ(ierr);
1579282e7bb4SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");CHKERRQ(ierr);
1580282e7bb4SMatthew G. Knepley     ierr = MatViewFromOptions(M, NULL, "-mat_view");CHKERRQ(ierr);
1581282e7bb4SMatthew G. Knepley     ierr = MatDestroy(&M);CHKERRQ(ierr);
1582282e7bb4SMatthew G. Knepley   } else {
1583282e7bb4SMatthew G. Knepley     ierr = SNESComputeJacobian(snes, u, J, J);CHKERRQ(ierr);
1584282e7bb4SMatthew G. Knepley   }
1585f017bcb6SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
1586282e7bb4SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");CHKERRQ(ierr);
1587282e7bb4SMatthew G. Knepley   ierr = MatViewFromOptions(J, NULL, "-mat_view");CHKERRQ(ierr);
1588f017bcb6SMatthew G. Knepley   /* Check nullspace */
1589f017bcb6SMatthew G. Knepley   ierr = MatGetNullSpace(J, &nullspace);CHKERRQ(ierr);
1590f017bcb6SMatthew G. Knepley   if (nullspace) {
15911878804aSMatthew G. Knepley     PetscBool isNull;
1592f017bcb6SMatthew G. Knepley     ierr = MatNullSpaceTest(nullspace, J, &isNull);CHKERRQ(ierr);
1593f017bcb6SMatthew G. Knepley     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
15941878804aSMatthew G. Knepley   }
1595f017bcb6SMatthew G. Knepley   ierr = MatNullSpaceDestroy(&nullspace);CHKERRQ(ierr);
1596f017bcb6SMatthew G. Knepley   /* Taylor test */
1597f017bcb6SMatthew G. Knepley   {
1598f017bcb6SMatthew G. Knepley     PetscRandom rand;
1599f017bcb6SMatthew G. Knepley     Vec         du, uhat, r, rhat, df;
1600f3c94560SMatthew G. Knepley     PetscReal   h;
1601f017bcb6SMatthew G. Knepley     PetscReal  *es, *hs, *errors;
1602f017bcb6SMatthew G. Knepley     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1603f017bcb6SMatthew G. Knepley     PetscInt    Nv, v;
1604f017bcb6SMatthew G. Knepley 
1605f017bcb6SMatthew G. Knepley     /* Choose a perturbation direction */
1606f017bcb6SMatthew G. Knepley     ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr);
1607f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &du);CHKERRQ(ierr);
1608f017bcb6SMatthew G. Knepley     ierr = VecSetRandom(du, rand); CHKERRQ(ierr);
1609f017bcb6SMatthew G. Knepley     ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
1610f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &df);CHKERRQ(ierr);
1611f017bcb6SMatthew G. Knepley     ierr = MatMult(J, du, df);CHKERRQ(ierr);
1612f017bcb6SMatthew G. Knepley     /* Evaluate residual at u, F(u), save in vector r */
1613f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
1614f017bcb6SMatthew G. Knepley     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
1615f017bcb6SMatthew G. Knepley     /* Look at the convergence of our Taylor approximation as we approach u */
1616f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
1617f017bcb6SMatthew G. Knepley     ierr = PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);CHKERRQ(ierr);
1618f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &uhat);CHKERRQ(ierr);
1619f017bcb6SMatthew G. Knepley     ierr = VecDuplicate(u, &rhat);CHKERRQ(ierr);
1620f017bcb6SMatthew G. Knepley     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
1621f017bcb6SMatthew G. Knepley       ierr = VecWAXPY(uhat, h, du, u);CHKERRQ(ierr);
1622f017bcb6SMatthew G. Knepley       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
1623f017bcb6SMatthew G. Knepley       ierr = SNESComputeFunction(snes, uhat, rhat);CHKERRQ(ierr);
1624f017bcb6SMatthew G. Knepley       ierr = VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);CHKERRQ(ierr);
1625f017bcb6SMatthew G. Knepley       ierr = VecNorm(rhat, NORM_2, &errors[Nv]);CHKERRQ(ierr);
1626f017bcb6SMatthew G. Knepley 
1627f017bcb6SMatthew G. Knepley       es[Nv] = PetscLog10Real(errors[Nv]);
1628f017bcb6SMatthew G. Knepley       hs[Nv] = PetscLog10Real(h);
1629f017bcb6SMatthew G. Knepley     }
1630f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&uhat);CHKERRQ(ierr);
1631f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&rhat);CHKERRQ(ierr);
1632f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&df);CHKERRQ(ierr);
16331878804aSMatthew G. Knepley     ierr = VecDestroy(&r);CHKERRQ(ierr);
1634f017bcb6SMatthew G. Knepley     ierr = VecDestroy(&du);CHKERRQ(ierr);
1635f017bcb6SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
1636f017bcb6SMatthew G. Knepley       if ((tol >= 0) && (errors[v] > tol)) break;
1637f017bcb6SMatthew G. Knepley       else if (errors[v] > PETSC_SMALL)    break;
1638f017bcb6SMatthew G. Knepley     }
1639f3c94560SMatthew G. Knepley     if (v == Nv) isLin = PETSC_TRUE;
1640f017bcb6SMatthew G. Knepley     ierr = PetscLinearRegression(Nv, hs, es, &slope, &intercept);CHKERRQ(ierr);
1641f017bcb6SMatthew G. Knepley     ierr = PetscFree3(es, hs, errors);CHKERRQ(ierr);
1642f017bcb6SMatthew G. Knepley     /* Slope should be about 2 */
1643f017bcb6SMatthew G. Knepley     if (tol >= 0) {
1644b878532fSMatthew G. Knepley       if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
1645b878532fSMatthew G. Knepley     } else if (isLinear || convRate) {
1646b878532fSMatthew G. Knepley       if (isLinear) *isLinear = isLin;
1647b878532fSMatthew G. Knepley       if (convRate) *convRate = slope;
1648f017bcb6SMatthew G. Knepley     } else {
1649f3c94560SMatthew G. Knepley       if (!isLin) {ierr = PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);CHKERRQ(ierr);}
1650f017bcb6SMatthew G. Knepley       else        {ierr = PetscPrintf(comm, "Function appears to be linear\n");CHKERRQ(ierr);}
1651f017bcb6SMatthew G. Knepley     }
1652f017bcb6SMatthew G. Knepley   }
16531878804aSMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
16541878804aSMatthew G. Knepley   PetscFunctionReturn(0);
16551878804aSMatthew G. Knepley }
16561878804aSMatthew G. Knepley 
16577f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1658f017bcb6SMatthew G. Knepley {
1659f017bcb6SMatthew G. Knepley   PetscErrorCode ierr;
1660f017bcb6SMatthew G. Knepley 
1661f017bcb6SMatthew G. Knepley   PetscFunctionBegin;
1662*f2cacb80SMatthew G. Knepley   ierr = DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL);CHKERRQ(ierr);
1663f3c94560SMatthew G. Knepley   ierr = DMSNESCheckResidual(snes, dm, u, -1.0, NULL);CHKERRQ(ierr);
1664f3c94560SMatthew G. Knepley   ierr = DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);CHKERRQ(ierr);
1665f017bcb6SMatthew G. Knepley   PetscFunctionReturn(0);
1666f017bcb6SMatthew G. Knepley }
1667f017bcb6SMatthew G. Knepley 
1668bee9a294SMatthew G. Knepley /*@C
1669bee9a294SMatthew G. Knepley   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1670bee9a294SMatthew G. Knepley 
1671bee9a294SMatthew G. Knepley   Input Parameters:
1672bee9a294SMatthew G. Knepley + snes - the SNES object
16737f96f943SMatthew G. Knepley - u    - representative SNES vector
16747f96f943SMatthew G. Knepley 
16757f96f943SMatthew G. Knepley   Note: The user must call PetscDSSetExactSolution() beforehand
1676bee9a294SMatthew G. Knepley 
1677bee9a294SMatthew G. Knepley   Level: developer
1678bee9a294SMatthew G. Knepley @*/
16797f96f943SMatthew G. Knepley PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
16801878804aSMatthew G. Knepley {
16811878804aSMatthew G. Knepley   DM             dm;
16821878804aSMatthew G. Knepley   Vec            sol;
16831878804aSMatthew G. Knepley   PetscBool      check;
16841878804aSMatthew G. Knepley   PetscErrorCode ierr;
16851878804aSMatthew G. Knepley 
16861878804aSMatthew G. Knepley   PetscFunctionBegin;
1687c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);CHKERRQ(ierr);
16881878804aSMatthew G. Knepley   if (!check) PetscFunctionReturn(0);
16891878804aSMatthew G. Knepley   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
16901878804aSMatthew G. Knepley   ierr = VecDuplicate(u, &sol);CHKERRQ(ierr);
16911878804aSMatthew G. Knepley   ierr = SNESSetSolution(snes, sol);CHKERRQ(ierr);
16927f96f943SMatthew G. Knepley   ierr = DMSNESCheck_Internal(snes, dm, sol);CHKERRQ(ierr);
16931878804aSMatthew G. Knepley   ierr = VecDestroy(&sol);CHKERRQ(ierr);
1694552f7358SJed Brown   PetscFunctionReturn(0);
1695552f7358SJed Brown }
1696