Actual source code: dmplexsnes.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1:  #include <petsc/private/dmpleximpl.h>
  2:  #include <petsc/private/snesimpl.h>
  3:  #include <petscds.h>
  4:  #include <petscblaslapack.h>
  5:  #include <petsc/private/petscimpl.h>
  6:  #include <petsc/private/petscfeimpl.h>

  8: /************************** Interpolation *******************************/

 10: static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
 11: {
 12:   PetscBool      isPlex;

 16:   PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);
 17:   if (isPlex) {
 18:     *plex = dm;
 19:     PetscObjectReference((PetscObject) dm);
 20:   } else {
 21:     PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);
 22:     if (!*plex) {
 23:       DMConvert(dm,DMPLEX,plex);
 24:       PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);
 25:       if (copy) {
 26:         PetscInt    i;
 27:         PetscObject obj;
 28:         const char *comps[3] = {"A","dmAux","dmCh"};

 30:         DMCopyDMSNES(dm, *plex);
 31:         for (i = 0; i < 3; i++) {
 32:           PetscObjectQuery((PetscObject) dm, comps[i], &obj);
 33:           PetscObjectCompose((PetscObject) *plex, comps[i], obj);
 34:         }
 35:       }
 36:     } else {
 37:       PetscObjectReference((PetscObject) *plex);
 38:     }
 39:   }
 40:   return(0);
 41: }

 43: PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
 44: {

 49:   PetscNew(ctx);

 51:   (*ctx)->comm   = comm;
 52:   (*ctx)->dim    = -1;
 53:   (*ctx)->nInput = 0;
 54:   (*ctx)->points = NULL;
 55:   (*ctx)->cells  = NULL;
 56:   (*ctx)->n      = -1;
 57:   (*ctx)->coords = NULL;
 58:   return(0);
 59: }

 61: PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
 62: {
 64:   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);
 65:   ctx->dim = dim;
 66:   return(0);
 67: }

 69: PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
 70: {
 73:   *dim = ctx->dim;
 74:   return(0);
 75: }

 77: PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
 78: {
 80:   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);
 81:   ctx->dof = dof;
 82:   return(0);
 83: }

 85: PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
 86: {
 89:   *dof = ctx->dof;
 90:   return(0);
 91: }

 93: PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
 94: {

 98:   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
 99:   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
100:   ctx->nInput = n;

102:   PetscMalloc1(n*ctx->dim, &ctx->points);
103:   PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));
104:   return(0);
105: }

107: PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
108: {
109:   MPI_Comm          comm = ctx->comm;
110:   PetscScalar       *a;
111:   PetscInt          p, q, i;
112:   PetscMPIInt       rank, size;
113:   PetscErrorCode    ierr;
114:   Vec               pointVec;
115:   PetscSF           cellSF;
116:   PetscLayout       layout;
117:   PetscReal         *globalPoints;
118:   PetscScalar       *globalPointsScalar;
119:   const PetscInt    *ranges;
120:   PetscMPIInt       *counts, *displs;
121:   const PetscSFNode *foundCells;
122:   const PetscInt    *foundPoints;
123:   PetscMPIInt       *foundProcs, *globalProcs;
124:   PetscInt          n, N, numFound;

128:   MPI_Comm_size(comm, &size);
129:   MPI_Comm_rank(comm, &rank);
130:   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
131:   /* Locate points */
132:   n = ctx->nInput;
133:   if (!redundantPoints) {
134:     PetscLayoutCreate(comm, &layout);
135:     PetscLayoutSetBlockSize(layout, 1);
136:     PetscLayoutSetLocalSize(layout, n);
137:     PetscLayoutSetUp(layout);
138:     PetscLayoutGetSize(layout, &N);
139:     /* Communicate all points to all processes */
140:     PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);
141:     PetscLayoutGetRanges(layout, &ranges);
142:     for (p = 0; p < size; ++p) {
143:       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
144:       displs[p] = ranges[p]*ctx->dim;
145:     }
146:     MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);
147:   } else {
148:     N = n;
149:     globalPoints = ctx->points;
150:     counts = displs = NULL;
151:     layout = NULL;
152:   }
153: #if 0
154:   PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);
155:   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
156: #else
157: #if defined(PETSC_USE_COMPLEX)
158:   PetscMalloc1(N*ctx->dim,&globalPointsScalar);
159:   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
160: #else
161:   globalPointsScalar = globalPoints;
162: #endif
163:   VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);
164:   PetscMalloc2(N,&foundProcs,N,&globalProcs);
165:   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
166:   cellSF = NULL;
167:   DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);
168:   PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);
169: #endif
170:   for (p = 0; p < numFound; ++p) {
171:     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
172:   }
173:   /* Let the lowest rank process own each point */
174:   MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);
175:   ctx->n = 0;
176:   for (p = 0; p < N; ++p) {
177:     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));
178:     else if (globalProcs[p] == rank) ctx->n++;
179:   }
180:   /* Create coordinates vector and array of owned cells */
181:   PetscMalloc1(ctx->n, &ctx->cells);
182:   VecCreate(comm, &ctx->coords);
183:   VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);
184:   VecSetBlockSize(ctx->coords, ctx->dim);
185:   VecSetType(ctx->coords,VECSTANDARD);
186:   VecGetArray(ctx->coords, &a);
187:   for (p = 0, q = 0, i = 0; p < N; ++p) {
188:     if (globalProcs[p] == rank) {
189:       PetscInt d;

191:       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
192:       ctx->cells[q] = foundCells[q].index;
193:       ++q;
194:     }
195:   }
196:   VecRestoreArray(ctx->coords, &a);
197: #if 0
198:   PetscFree3(foundCells,foundProcs,globalProcs);
199: #else
200:   PetscFree2(foundProcs,globalProcs);
201:   PetscSFDestroy(&cellSF);
202:   VecDestroy(&pointVec);
203: #endif
204:   if ((void*)globalPointsScalar != (void*)globalPoints) {PetscFree(globalPointsScalar);}
205:   if (!redundantPoints) {PetscFree3(globalPoints,counts,displs);}
206:   PetscLayoutDestroy(&layout);
207:   return(0);
208: }

210: PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
211: {
214:   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
215:   *coordinates = ctx->coords;
216:   return(0);
217: }

219: PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
220: {

225:   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
226:   VecCreate(ctx->comm, v);
227:   VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);
228:   VecSetBlockSize(*v, ctx->dof);
229:   VecSetType(*v,VECSTANDARD);
230:   return(0);
231: }

233: PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
234: {

239:   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
240:   VecDestroy(v);
241:   return(0);
242: }

244: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
245: {
246:   PetscReal      *v0, *J, *invJ, detJ;
247:   const PetscScalar *coords;
248:   PetscScalar    *a;
249:   PetscInt       p;

253:   PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);
254:   VecGetArrayRead(ctx->coords, &coords);
255:   VecGetArray(v, &a);
256:   for (p = 0; p < ctx->n; ++p) {
257:     PetscInt     c = ctx->cells[p];
258:     PetscScalar *x = NULL;
259:     PetscReal    xi[4];
260:     PetscInt     d, f, comp;

262:     DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);
263:     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", (double)detJ, c);
264:     DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);
265:     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];

267:     for (d = 0; d < ctx->dim; ++d) {
268:       xi[d] = 0.0;
269:       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
270:       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];
271:     }
272:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);
273:   }
274:   VecRestoreArray(v, &a);
275:   VecRestoreArrayRead(ctx->coords, &coords);
276:   PetscFree3(v0, J, invJ);
277:   return(0);
278: }

280: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
281: {
282:   PetscReal      *v0, *J, *invJ, detJ;
283:   const PetscScalar *coords;
284:   PetscScalar    *a;
285:   PetscInt       p;

289:   PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);
290:   VecGetArrayRead(ctx->coords, &coords);
291:   VecGetArray(v, &a);
292:   for (p = 0; p < ctx->n; ++p) {
293:     PetscInt       c = ctx->cells[p];
294:     const PetscInt order[3] = {2, 1, 3};
295:     PetscScalar   *x = NULL;
296:     PetscReal      xi[4];
297:     PetscInt       d, f, comp;

299:     DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);
300:     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", (double)detJ, c);
301:     DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);
302:     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];

304:     for (d = 0; d < ctx->dim; ++d) {
305:       xi[d] = 0.0;
306:       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
307:       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];
308:     }
309:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);
310:   }
311:   VecRestoreArray(v, &a);
312:   VecRestoreArrayRead(ctx->coords, &coords);
313:   PetscFree3(v0, J, invJ);
314:   return(0);
315: }

317: PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
318: {
319:   const PetscScalar *vertices = (const PetscScalar*) ctx;
320:   const PetscScalar x0        = vertices[0];
321:   const PetscScalar y0        = vertices[1];
322:   const PetscScalar x1        = vertices[2];
323:   const PetscScalar y1        = vertices[3];
324:   const PetscScalar x2        = vertices[4];
325:   const PetscScalar y2        = vertices[5];
326:   const PetscScalar x3        = vertices[6];
327:   const PetscScalar y3        = vertices[7];
328:   const PetscScalar f_1       = x1 - x0;
329:   const PetscScalar g_1       = y1 - y0;
330:   const PetscScalar f_3       = x3 - x0;
331:   const PetscScalar g_3       = y3 - y0;
332:   const PetscScalar f_01      = x2 - x1 - x3 + x0;
333:   const PetscScalar g_01      = y2 - y1 - y3 + y0;
334:   const PetscScalar *ref;
335:   PetscScalar       *real;
336:   PetscErrorCode    ierr;

339:   VecGetArrayRead(Xref,  &ref);
340:   VecGetArray(Xreal, &real);
341:   {
342:     const PetscScalar p0 = ref[0];
343:     const PetscScalar p1 = ref[1];

345:     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
346:     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
347:   }
348:   PetscLogFlops(28);
349:   VecRestoreArrayRead(Xref,  &ref);
350:   VecRestoreArray(Xreal, &real);
351:   return(0);
352: }

354:  #include <petsc/private/dmimpl.h>
355: PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
356: {
357:   const PetscScalar *vertices = (const PetscScalar*) ctx;
358:   const PetscScalar x0        = vertices[0];
359:   const PetscScalar y0        = vertices[1];
360:   const PetscScalar x1        = vertices[2];
361:   const PetscScalar y1        = vertices[3];
362:   const PetscScalar x2        = vertices[4];
363:   const PetscScalar y2        = vertices[5];
364:   const PetscScalar x3        = vertices[6];
365:   const PetscScalar y3        = vertices[7];
366:   const PetscScalar f_01      = x2 - x1 - x3 + x0;
367:   const PetscScalar g_01      = y2 - y1 - y3 + y0;
368:   const PetscScalar *ref;
369:   PetscErrorCode    ierr;

372:   VecGetArrayRead(Xref,  &ref);
373:   {
374:     const PetscScalar x       = ref[0];
375:     const PetscScalar y       = ref[1];
376:     const PetscInt    rows[2] = {0, 1};
377:     PetscScalar       values[4];

379:     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
380:     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
381:     MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);
382:   }
383:   PetscLogFlops(30);
384:   VecRestoreArrayRead(Xref,  &ref);
385:   MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
386:   MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
387:   return(0);
388: }

390: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
391: {
392:   DM             dmCoord;
393:   PetscFE        fem = NULL;
394:   SNES           snes;
395:   KSP            ksp;
396:   PC             pc;
397:   Vec            coordsLocal, r, ref, real;
398:   Mat            J;
399:   const PetscScalar *coords;
400:   PetscScalar    *a;
401:   PetscInt       Nf, p;
402:   const PetscInt dof = ctx->dof;

406:   DMGetNumFields(dm, &Nf);
407:   if (Nf) {DMGetField(dm, 0, (PetscObject *) &fem);}
408:   DMGetCoordinatesLocal(dm, &coordsLocal);
409:   DMGetCoordinateDM(dm, &dmCoord);
410:   SNESCreate(PETSC_COMM_SELF, &snes);
411:   SNESSetOptionsPrefix(snes, "quad_interp_");
412:   VecCreate(PETSC_COMM_SELF, &r);
413:   VecSetSizes(r, 2, 2);
414:   VecSetType(r,dm->vectype);
415:   VecDuplicate(r, &ref);
416:   VecDuplicate(r, &real);
417:   MatCreate(PETSC_COMM_SELF, &J);
418:   MatSetSizes(J, 2, 2, 2, 2);
419:   MatSetType(J, MATSEQDENSE);
420:   MatSetUp(J);
421:   SNESSetFunction(snes, r, QuadMap_Private, NULL);
422:   SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);
423:   SNESGetKSP(snes, &ksp);
424:   KSPGetPC(ksp, &pc);
425:   PCSetType(pc, PCLU);
426:   SNESSetFromOptions(snes);

428:   VecGetArrayRead(ctx->coords, &coords);
429:   VecGetArray(v, &a);
430:   for (p = 0; p < ctx->n; ++p) {
431:     PetscScalar *x = NULL, *vertices = NULL;
432:     PetscScalar *xi;
433:     PetscReal    xir[2];
434:     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;

436:     /* Can make this do all points at once */
437:     DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
438:     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);
439:     DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);
440:     SNESSetFunction(snes, NULL, NULL, (void*) vertices);
441:     SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);
442:     VecGetArray(real, &xi);
443:     xi[0]  = coords[p*ctx->dim+0];
444:     xi[1]  = coords[p*ctx->dim+1];
445:     VecRestoreArray(real, &xi);
446:     SNESSolve(snes, real, ref);
447:     VecGetArray(ref, &xi);
448:     xir[0] = PetscRealPart(xi[0]);
449:     xir[1] = PetscRealPart(xi[1]);
450:     if (4*dof != xSize) {
451:       PetscReal *B;
452:       PetscInt   d;

454:       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
455:       PetscFEGetTabulation(fem, 1, xir, &B, NULL, NULL);
456:       for (comp = 0; comp < dof; ++comp) {
457:         a[p*dof+comp] = 0.0;
458:         for (d = 0; d < xSize/dof; ++d) {
459:           a[p*dof+comp] += x[d*dof+comp]*B[d*dof+comp];
460:         }
461:       }
462:       PetscFERestoreTabulation(fem, 1, xir, &B, NULL, NULL);
463:     } else {
464:       for (comp = 0; comp < dof; ++comp)
465:         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];
466:     }
467:     VecRestoreArray(ref, &xi);
468:     DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
469:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);
470:   }
471:   VecRestoreArray(v, &a);
472:   VecRestoreArrayRead(ctx->coords, &coords);

474:   SNESDestroy(&snes);
475:   VecDestroy(&r);
476:   VecDestroy(&ref);
477:   VecDestroy(&real);
478:   MatDestroy(&J);
479:   return(0);
480: }

482: PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
483: {
484:   const PetscScalar *vertices = (const PetscScalar*) ctx;
485:   const PetscScalar x0        = vertices[0];
486:   const PetscScalar y0        = vertices[1];
487:   const PetscScalar z0        = vertices[2];
488:   const PetscScalar x1        = vertices[9];
489:   const PetscScalar y1        = vertices[10];
490:   const PetscScalar z1        = vertices[11];
491:   const PetscScalar x2        = vertices[6];
492:   const PetscScalar y2        = vertices[7];
493:   const PetscScalar z2        = vertices[8];
494:   const PetscScalar x3        = vertices[3];
495:   const PetscScalar y3        = vertices[4];
496:   const PetscScalar z3        = vertices[5];
497:   const PetscScalar x4        = vertices[12];
498:   const PetscScalar y4        = vertices[13];
499:   const PetscScalar z4        = vertices[14];
500:   const PetscScalar x5        = vertices[15];
501:   const PetscScalar y5        = vertices[16];
502:   const PetscScalar z5        = vertices[17];
503:   const PetscScalar x6        = vertices[18];
504:   const PetscScalar y6        = vertices[19];
505:   const PetscScalar z6        = vertices[20];
506:   const PetscScalar x7        = vertices[21];
507:   const PetscScalar y7        = vertices[22];
508:   const PetscScalar z7        = vertices[23];
509:   const PetscScalar f_1       = x1 - x0;
510:   const PetscScalar g_1       = y1 - y0;
511:   const PetscScalar h_1       = z1 - z0;
512:   const PetscScalar f_3       = x3 - x0;
513:   const PetscScalar g_3       = y3 - y0;
514:   const PetscScalar h_3       = z3 - z0;
515:   const PetscScalar f_4       = x4 - x0;
516:   const PetscScalar g_4       = y4 - y0;
517:   const PetscScalar h_4       = z4 - z0;
518:   const PetscScalar f_01      = x2 - x1 - x3 + x0;
519:   const PetscScalar g_01      = y2 - y1 - y3 + y0;
520:   const PetscScalar h_01      = z2 - z1 - z3 + z0;
521:   const PetscScalar f_12      = x7 - x3 - x4 + x0;
522:   const PetscScalar g_12      = y7 - y3 - y4 + y0;
523:   const PetscScalar h_12      = z7 - z3 - z4 + z0;
524:   const PetscScalar f_02      = x5 - x1 - x4 + x0;
525:   const PetscScalar g_02      = y5 - y1 - y4 + y0;
526:   const PetscScalar h_02      = z5 - z1 - z4 + z0;
527:   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
528:   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
529:   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
530:   const PetscScalar *ref;
531:   PetscScalar       *real;
532:   PetscErrorCode    ierr;

535:   VecGetArrayRead(Xref,  &ref);
536:   VecGetArray(Xreal, &real);
537:   {
538:     const PetscScalar p0 = ref[0];
539:     const PetscScalar p1 = ref[1];
540:     const PetscScalar p2 = ref[2];

542:     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;
543:     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;
544:     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;
545:   }
546:   PetscLogFlops(114);
547:   VecRestoreArrayRead(Xref,  &ref);
548:   VecRestoreArray(Xreal, &real);
549:   return(0);
550: }

552: PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
553: {
554:   const PetscScalar *vertices = (const PetscScalar*) ctx;
555:   const PetscScalar x0        = vertices[0];
556:   const PetscScalar y0        = vertices[1];
557:   const PetscScalar z0        = vertices[2];
558:   const PetscScalar x1        = vertices[9];
559:   const PetscScalar y1        = vertices[10];
560:   const PetscScalar z1        = vertices[11];
561:   const PetscScalar x2        = vertices[6];
562:   const PetscScalar y2        = vertices[7];
563:   const PetscScalar z2        = vertices[8];
564:   const PetscScalar x3        = vertices[3];
565:   const PetscScalar y3        = vertices[4];
566:   const PetscScalar z3        = vertices[5];
567:   const PetscScalar x4        = vertices[12];
568:   const PetscScalar y4        = vertices[13];
569:   const PetscScalar z4        = vertices[14];
570:   const PetscScalar x5        = vertices[15];
571:   const PetscScalar y5        = vertices[16];
572:   const PetscScalar z5        = vertices[17];
573:   const PetscScalar x6        = vertices[18];
574:   const PetscScalar y6        = vertices[19];
575:   const PetscScalar z6        = vertices[20];
576:   const PetscScalar x7        = vertices[21];
577:   const PetscScalar y7        = vertices[22];
578:   const PetscScalar z7        = vertices[23];
579:   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
580:   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
581:   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
582:   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
583:   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
584:   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
585:   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
586:   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
587:   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
588:   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
589:   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
590:   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
591:   const PetscScalar *ref;
592:   PetscErrorCode    ierr;

595:   VecGetArrayRead(Xref,  &ref);
596:   {
597:     const PetscScalar x       = ref[0];
598:     const PetscScalar y       = ref[1];
599:     const PetscScalar z       = ref[2];
600:     const PetscInt    rows[3] = {0, 1, 2};
601:     PetscScalar       values[9];

603:     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
604:     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
605:     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
606:     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
607:     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
608:     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
609:     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
610:     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
611:     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;

613:     MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);
614:   }
615:   PetscLogFlops(152);
616:   VecRestoreArrayRead(Xref,  &ref);
617:   MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
618:   MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
619:   return(0);
620: }

622: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
623: {
624:   DM             dmCoord;
625:   SNES           snes;
626:   KSP            ksp;
627:   PC             pc;
628:   Vec            coordsLocal, r, ref, real;
629:   Mat            J;
630:   const PetscScalar *coords;
631:   PetscScalar    *a;
632:   PetscInt       p;

636:   DMGetCoordinatesLocal(dm, &coordsLocal);
637:   DMGetCoordinateDM(dm, &dmCoord);
638:   SNESCreate(PETSC_COMM_SELF, &snes);
639:   SNESSetOptionsPrefix(snes, "hex_interp_");
640:   VecCreate(PETSC_COMM_SELF, &r);
641:   VecSetSizes(r, 3, 3);
642:   VecSetType(r,dm->vectype);
643:   VecDuplicate(r, &ref);
644:   VecDuplicate(r, &real);
645:   MatCreate(PETSC_COMM_SELF, &J);
646:   MatSetSizes(J, 3, 3, 3, 3);
647:   MatSetType(J, MATSEQDENSE);
648:   MatSetUp(J);
649:   SNESSetFunction(snes, r, HexMap_Private, NULL);
650:   SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);
651:   SNESGetKSP(snes, &ksp);
652:   KSPGetPC(ksp, &pc);
653:   PCSetType(pc, PCLU);
654:   SNESSetFromOptions(snes);

656:   VecGetArrayRead(ctx->coords, &coords);
657:   VecGetArray(v, &a);
658:   for (p = 0; p < ctx->n; ++p) {
659:     PetscScalar *x = NULL, *vertices = NULL;
660:     PetscScalar *xi;
661:     PetscReal    xir[3];
662:     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;

664:     /* Can make this do all points at once */
665:     DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
666:     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);
667:     DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);
668:     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);
669:     SNESSetFunction(snes, NULL, NULL, (void*) vertices);
670:     SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);
671:     VecGetArray(real, &xi);
672:     xi[0]  = coords[p*ctx->dim+0];
673:     xi[1]  = coords[p*ctx->dim+1];
674:     xi[2]  = coords[p*ctx->dim+2];
675:     VecRestoreArray(real, &xi);
676:     SNESSolve(snes, real, ref);
677:     VecGetArray(ref, &xi);
678:     xir[0] = PetscRealPart(xi[0]);
679:     xir[1] = PetscRealPart(xi[1]);
680:     xir[2] = PetscRealPart(xi[2]);
681:     for (comp = 0; comp < ctx->dof; ++comp) {
682:       a[p*ctx->dof+comp] =
683:         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
684:         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
685:         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
686:         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
687:         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
688:         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
689:         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
690:         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
691:     }
692:     VecRestoreArray(ref, &xi);
693:     DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
694:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);
695:   }
696:   VecRestoreArray(v, &a);
697:   VecRestoreArrayRead(ctx->coords, &coords);

699:   SNESDestroy(&snes);
700:   VecDestroy(&r);
701:   VecDestroy(&ref);
702:   VecDestroy(&real);
703:   MatDestroy(&J);
704:   return(0);
705: }

707: /*
708:   Input Parameters:
709: + ctx - The DMInterpolationInfo context
710: . dm  - The DM
711: - x   - The local vector containing the field to be interpolated

713:   Output Parameters:
714: . v   - The vector containing the interpolated values
715: */
716: PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
717: {
718:   PetscInt       dim, coneSize, n;

725:   VecGetLocalSize(v, &n);
726:   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);
727:   if (n) {
728:     DMGetDimension(dm, &dim);
729:     DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);
730:     if (dim == 2) {
731:       if (coneSize == 3) {
732:         DMInterpolate_Triangle_Private(ctx, dm, x, v);
733:       } else if (coneSize == 4) {
734:         DMInterpolate_Quad_Private(ctx, dm, x, v);
735:       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
736:     } else if (dim == 3) {
737:       if (coneSize == 4) {
738:         DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);
739:       } else {
740:         DMInterpolate_Hex_Private(ctx, dm, x, v);
741:       }
742:     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
743:   }
744:   return(0);
745: }

747: PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
748: {

753:   VecDestroy(&(*ctx)->coords);
754:   PetscFree((*ctx)->points);
755:   PetscFree((*ctx)->cells);
756:   PetscFree(*ctx);
757:   *ctx = NULL;
758:   return(0);
759: }

761: /*@C
762:   SNESMonitorFields - Monitors the residual for each field separately

764:   Collective on SNES

766:   Input Parameters:
767: + snes   - the SNES context
768: . its    - iteration number
769: . fgnorm - 2-norm of residual
770: - vf  - PetscViewerAndFormat of type ASCII

772:   Notes:
773:   This routine prints the residual norm at each iteration.

775:   Level: intermediate

777: .keywords: SNES, nonlinear, default, monitor, norm
778: .seealso: SNESMonitorSet(), SNESMonitorDefault()
779: @*/
780: PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
781: {
782:   PetscViewer        viewer = vf->viewer;
783:   Vec                res;
784:   DM                 dm;
785:   PetscSection       s;
786:   const PetscScalar *r;
787:   PetscReal         *lnorms, *norms;
788:   PetscInt           numFields, f, pStart, pEnd, p;
789:   PetscErrorCode     ierr;

793:   SNESGetFunction(snes, &res, 0, 0);
794:   SNESGetDM(snes, &dm);
795:   DMGetSection(dm, &s);
796:   PetscSectionGetNumFields(s, &numFields);
797:   PetscSectionGetChart(s, &pStart, &pEnd);
798:   PetscCalloc2(numFields, &lnorms, numFields, &norms);
799:   VecGetArrayRead(res, &r);
800:   for (p = pStart; p < pEnd; ++p) {
801:     for (f = 0; f < numFields; ++f) {
802:       PetscInt fdof, foff, d;

804:       PetscSectionGetFieldDof(s, p, f, &fdof);
805:       PetscSectionGetFieldOffset(s, p, f, &foff);
806:       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
807:     }
808:   }
809:   VecRestoreArrayRead(res, &r);
810:   MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));
811:   PetscViewerPushFormat(viewer,vf->format);
812:   PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);
813:   PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);
814:   for (f = 0; f < numFields; ++f) {
815:     if (f > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
816:     PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));
817:   }
818:   PetscViewerASCIIPrintf(viewer, "]\n");
819:   PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);
820:   PetscViewerPopFormat(viewer);
821:   PetscFree2(lnorms, norms);
822:   return(0);
823: }

825: /********************* Residual Computation **************************/


828: /*@
829:   DMPlexSNESGetGeometryFVM - Return precomputed geometric data

831:   Input Parameter:
832: . dm - The DM

834:   Output Parameters:
835: + facegeom - The values precomputed from face geometry
836: . cellgeom - The values precomputed from cell geometry
837: - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell

839:   Level: developer

841: .seealso: DMPlexTSSetRHSFunctionLocal()
842: @*/
843: PetscErrorCode DMPlexSNESGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
844: {
845:   DM             plex;

850:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
851:   DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL);
852:   if (minRadius) {DMPlexGetMinRadius(plex, minRadius);}
853:   DMDestroy(&plex);
854:   return(0);
855: }

857: /*@
858:   DMPlexSNESGetGradientDM - Return gradient data layout

860:   Input Parameters:
861: + dm - The DM
862: - fv - The PetscFV

864:   Output Parameter:
865: . dmGrad - The layout for gradient values

867:   Level: developer

869: .seealso: DMPlexSNESGetGeometryFVM()
870: @*/
871: PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
872: {
873:   DM             plex;
874:   PetscBool      computeGradients;

881:   PetscFVGetComputeGradients(fv, &computeGradients);
882:   if (!computeGradients) {*dmGrad = NULL; return(0);}
883:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
884:   DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad);
885:   DMDestroy(&plex);
886:   return(0);
887: }

889: /*@C
890:   DMPlexGetCellFields - Retrieve the field values values for a chunk of cells

892:   Input Parameters:
893: + dm     - The DM
894: . cellIS - The cells to include
895: . locX   - A local vector with the solution fields
896: . locX_t - A local vector with solution field time derivatives, or NULL
897: - locA   - A local vector with auxiliary fields, or NULL

899:   Output Parameters:
900: + u   - The field coefficients
901: . u_t - The fields derivative coefficients
902: - a   - The auxiliary field coefficients

904:   Level: developer

906: .seealso: DMPlexGetFaceFields()
907: @*/
908: PetscErrorCode DMPlexGetCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
909: {
910:   DM              plex, plexA = NULL;
911:   PetscSection    section, sectionAux;
912:   PetscDS         prob;
913:   const PetscInt *cells;
914:   PetscInt        cStart, cEnd, numCells, totDim, totDimAux, c;
915:   PetscErrorCode  ierr;

925:   DMSNESConvertPlex(dm, &plex, PETSC_FALSE);
926:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
927:   DMGetSection(dm, &section);
928:   DMGetDS(dm, &prob);
929:   PetscDSGetTotalDimension(prob, &totDim);
930:   if (locA) {
931:     DM      dmAux;
932:     PetscDS probAux;

934:     VecGetDM(locA, &dmAux);
935:     DMSNESConvertPlex(dmAux, &plexA, PETSC_FALSE);
936:     DMGetSection(dmAux, &sectionAux);
937:     DMGetDS(dmAux, &probAux);
938:     PetscDSGetTotalDimension(probAux, &totDimAux);
939:   }
940:   numCells = cEnd - cStart;
941:   DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, u);
942:   if (locX_t) {DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, u_t);} else {*u_t = NULL;}
943:   if (locA)   {DMGetWorkArray(dm, numCells*totDimAux, MPIU_SCALAR, a);} else {*a = NULL;}
944:   for (c = cStart; c < cEnd; ++c) {
945:     const PetscInt cell = cells ? cells[c] : c;
946:     const PetscInt cind = c - cStart;
947:     PetscScalar   *x = NULL, *x_t = NULL, *ul = *u, *ul_t = *u_t, *al = *a;
948:     PetscInt       i;

950:     DMPlexVecGetClosure(plex, section, locX, cell, NULL, &x);
951:     for (i = 0; i < totDim; ++i) ul[cind*totDim+i] = x[i];
952:     DMPlexVecRestoreClosure(plex, section, locX, cell, NULL, &x);
953:     if (locX_t) {
954:       DMPlexVecGetClosure(plex, section, locX_t, cell, NULL, &x_t);
955:       for (i = 0; i < totDim; ++i) ul_t[cind*totDim+i] = x_t[i];
956:       DMPlexVecRestoreClosure(plex, section, locX_t, cell, NULL, &x_t);
957:     }
958:     if (locA) {
959:       DMPlexVecGetClosure(plexA, sectionAux, locA, cell, NULL, &x);
960:       for (i = 0; i < totDimAux; ++i) al[cind*totDimAux+i] = x[i];
961:       DMPlexVecRestoreClosure(plexA, sectionAux, locA, cell, NULL, &x);
962:     }
963:   }
964:   DMDestroy(&plex);
965:   if (locA) {DMDestroy(&plexA);}
966:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
967:   return(0);
968: }

970: /*@C
971:   DMPlexRestoreCellFields - Restore the field values values for a chunk of cells

973:   Input Parameters:
974: + dm     - The DM
975: . cellIS - The cells to include
976: . locX   - A local vector with the solution fields
977: . locX_t - A local vector with solution field time derivatives, or NULL
978: - locA   - A local vector with auxiliary fields, or NULL

980:   Output Parameters:
981: + u   - The field coefficients
982: . u_t - The fields derivative coefficients
983: - a   - The auxiliary field coefficients

985:   Level: developer

987: .seealso: DMPlexGetFaceFields()
988: @*/
989: PetscErrorCode DMPlexRestoreCellFields(DM dm, IS cellIS, Vec locX, Vec locX_t, Vec locA, PetscScalar **u, PetscScalar **u_t, PetscScalar **a)
990: {

994:   DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u);
995:   if (locX_t) {DMRestoreWorkArray(dm, 0, MPIU_SCALAR, u_t);}
996:   if (locA)   {DMRestoreWorkArray(dm, 0, MPIU_SCALAR, a);}
997:   return(0);
998: }

1000: /*@C
1001:   DMPlexGetFaceFields - Retrieve the field values values for a chunk of faces

1003:   Input Parameters:
1004: + dm     - The DM
1005: . fStart - The first face to include
1006: . fEnd   - The first face to exclude
1007: . locX   - A local vector with the solution fields
1008: . locX_t - A local vector with solution field time derivatives, or NULL
1009: . faceGeometry - A local vector with face geometry
1010: . cellGeometry - A local vector with cell geometry
1011: - locaGrad - A local vector with field gradients, or NULL

1013:   Output Parameters:
1014: + Nface - The number of faces with field values
1015: . uL - The field values at the left side of the face
1016: - uR - The field values at the right side of the face

1018:   Level: developer

1020: .seealso: DMPlexGetCellFields()
1021: @*/
1022: PetscErrorCode DMPlexGetFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscInt *Nface, PetscScalar **uL, PetscScalar **uR)
1023: {
1024:   DM                 dmFace, dmCell, dmGrad = NULL;
1025:   PetscSection       section;
1026:   PetscDS            prob;
1027:   DMLabel            ghostLabel;
1028:   const PetscScalar *facegeom, *cellgeom, *x, *lgrad;
1029:   PetscBool         *isFE;
1030:   PetscInt           dim, Nf, f, Nc, numFaces = fEnd - fStart, iface, face;
1031:   PetscErrorCode     ierr;

1042:   DMGetDimension(dm, &dim);
1043:   DMGetDS(dm, &prob);
1044:   DMGetSection(dm, &section);
1045:   PetscDSGetNumFields(prob, &Nf);
1046:   PetscDSGetTotalComponents(prob, &Nc);
1047:   PetscMalloc1(Nf, &isFE);
1048:   for (f = 0; f < Nf; ++f) {
1049:     PetscObject  obj;
1050:     PetscClassId id;

1052:     PetscDSGetDiscretization(prob, f, &obj);
1053:     PetscObjectGetClassId(obj, &id);
1054:     if (id == PETSCFE_CLASSID)      {isFE[f] = PETSC_TRUE;}
1055:     else if (id == PETSCFV_CLASSID) {isFE[f] = PETSC_FALSE;}
1056:     else                            {isFE[f] = PETSC_FALSE;}
1057:   }
1058:   DMGetLabel(dm, "ghost", &ghostLabel);
1059:   VecGetArrayRead(locX, &x);
1060:   VecGetDM(faceGeometry, &dmFace);
1061:   VecGetArrayRead(faceGeometry, &facegeom);
1062:   VecGetDM(cellGeometry, &dmCell);
1063:   VecGetArrayRead(cellGeometry, &cellgeom);
1064:   if (locGrad) {
1065:     VecGetDM(locGrad, &dmGrad);
1066:     VecGetArrayRead(locGrad, &lgrad);
1067:   }
1068:   DMGetWorkArray(dm, numFaces*Nc, MPIU_SCALAR, uL);
1069:   DMGetWorkArray(dm, numFaces*Nc, MPIU_SCALAR, uR);
1070:   /* Right now just eat the extra work for FE (could make a cell loop) */
1071:   for (face = fStart, iface = 0; face < fEnd; ++face) {
1072:     const PetscInt        *cells;
1073:     PetscFVFaceGeom       *fg;
1074:     PetscFVCellGeom       *cgL, *cgR;
1075:     PetscScalar           *xL, *xR, *gL, *gR;
1076:     PetscScalar           *uLl = *uL, *uRl = *uR;
1077:     PetscInt               ghost, nsupp, nchild;

1079:     DMLabelGetValue(ghostLabel, face, &ghost);
1080:     DMPlexGetSupportSize(dm, face, &nsupp);
1081:     DMPlexGetTreeChildren(dm, face, &nchild, NULL);
1082:     if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
1083:     DMPlexPointLocalRead(dmFace, face, facegeom, &fg);
1084:     DMPlexGetSupport(dm, face, &cells);
1085:     DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);
1086:     DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);
1087:     for (f = 0; f < Nf; ++f) {
1088:       PetscInt off;

1090:       PetscDSGetComponentOffset(prob, f, &off);
1091:       if (isFE[f]) {
1092:         const PetscInt *cone;
1093:         PetscInt        comp, coneSizeL, coneSizeR, faceLocL, faceLocR, ldof, rdof, d;

1095:         xL = xR = NULL;
1096:         PetscSectionGetFieldComponents(section, f, &comp);
1097:         DMPlexVecGetClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);
1098:         DMPlexVecGetClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);
1099:         DMPlexGetCone(dm, cells[0], &cone);
1100:         DMPlexGetConeSize(dm, cells[0], &coneSizeL);
1101:         for (faceLocL = 0; faceLocL < coneSizeL; ++faceLocL) if (cone[faceLocL] == face) break;
1102:         DMPlexGetCone(dm, cells[1], &cone);
1103:         DMPlexGetConeSize(dm, cells[1], &coneSizeR);
1104:         for (faceLocR = 0; faceLocR < coneSizeR; ++faceLocR) if (cone[faceLocR] == face) break;
1105:         if (faceLocL == coneSizeL && faceLocR == coneSizeR) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %d in cone of cell %d or cell %d", face, cells[0], cells[1]);
1106:         /* Check that FEM field has values in the right cell (sometimes its an FV ghost cell) */
1107:         /* TODO: this is a hack that might not be right for nonconforming */
1108:         if (faceLocL < coneSizeL) {
1109:           EvaluateFaceFields(prob, f, faceLocL, xL, &uLl[iface*Nc+off]);
1110:           if (rdof == ldof && faceLocR < coneSizeR) {EvaluateFaceFields(prob, f, faceLocR, xR, &uRl[iface*Nc+off]);}
1111:           else              {for(d = 0; d < comp; ++d) uRl[iface*Nc+off+d] = uLl[iface*Nc+off+d];}
1112:         }
1113:         else {
1114:           EvaluateFaceFields(prob, f, faceLocR, xR, &uRl[iface*Nc+off]);
1115:           PetscSectionGetFieldComponents(section, f, &comp);
1116:           for(d = 0; d < comp; ++d) uLl[iface*Nc+off+d] = uRl[iface*Nc+off+d];
1117:         }
1118:         DMPlexVecRestoreClosure(dm, section, locX, cells[0], &ldof, (PetscScalar **) &xL);
1119:         DMPlexVecRestoreClosure(dm, section, locX, cells[1], &rdof, (PetscScalar **) &xR);
1120:       } else {
1121:         PetscFV  fv;
1122:         PetscInt numComp, c;

1124:         PetscDSGetDiscretization(prob, f, (PetscObject *) &fv);
1125:         PetscFVGetNumComponents(fv, &numComp);
1126:         DMPlexPointLocalFieldRead(dm, cells[0], f, x, &xL);
1127:         DMPlexPointLocalFieldRead(dm, cells[1], f, x, &xR);
1128:         if (dmGrad) {
1129:           PetscReal dxL[3], dxR[3];

1131:           DMPlexPointLocalRead(dmGrad, cells[0], lgrad, &gL);
1132:           DMPlexPointLocalRead(dmGrad, cells[1], lgrad, &gR);
1133:           DMPlex_WaxpyD_Internal(dim, -1, cgL->centroid, fg->centroid, dxL);
1134:           DMPlex_WaxpyD_Internal(dim, -1, cgR->centroid, fg->centroid, dxR);
1135:           for (c = 0; c < numComp; ++c) {
1136:             uLl[iface*Nc+off+c] = xL[c] + DMPlex_DotD_Internal(dim, &gL[c*dim], dxL);
1137:             uRl[iface*Nc+off+c] = xR[c] + DMPlex_DotD_Internal(dim, &gR[c*dim], dxR);
1138:           }
1139:         } else {
1140:           for (c = 0; c < numComp; ++c) {
1141:             uLl[iface*Nc+off+c] = xL[c];
1142:             uRl[iface*Nc+off+c] = xR[c];
1143:           }
1144:         }
1145:       }
1146:     }
1147:     ++iface;
1148:   }
1149:   *Nface = iface;
1150:   VecRestoreArrayRead(locX, &x);
1151:   VecRestoreArrayRead(faceGeometry, &facegeom);
1152:   VecRestoreArrayRead(cellGeometry, &cellgeom);
1153:   if (locGrad) {
1154:     VecRestoreArrayRead(locGrad, &lgrad);
1155:   }
1156:   PetscFree(isFE);
1157:   return(0);
1158: }

1160: /*@C
1161:   DMPlexRestoreFaceFields - Restore the field values values for a chunk of faces

1163:   Input Parameters:
1164: + dm     - The DM
1165: . fStart - The first face to include
1166: . fEnd   - The first face to exclude
1167: . locX   - A local vector with the solution fields
1168: . locX_t - A local vector with solution field time derivatives, or NULL
1169: . faceGeometry - A local vector with face geometry
1170: . cellGeometry - A local vector with cell geometry
1171: - locaGrad - A local vector with field gradients, or NULL

1173:   Output Parameters:
1174: + Nface - The number of faces with field values
1175: . uL - The field values at the left side of the face
1176: - uR - The field values at the right side of the face

1178:   Level: developer

1180: .seealso: DMPlexGetFaceFields()
1181: @*/
1182: PetscErrorCode DMPlexRestoreFaceFields(DM dm, PetscInt fStart, PetscInt fEnd, Vec locX, Vec locX_t, Vec faceGeometry, Vec cellGeometry, Vec locGrad, PetscInt *Nface, PetscScalar **uL, PetscScalar **uR)
1183: {

1187:   DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uL);
1188:   DMRestoreWorkArray(dm, 0, MPIU_SCALAR, uR);
1189:   return(0);
1190: }

1192: /*@C
1193:   DMPlexGetFaceGeometry - Retrieve the geometric values for a chunk of faces

1195:   Input Parameters:
1196: + dm     - The DM
1197: . fStart - The first face to include
1198: . fEnd   - The first face to exclude
1199: . faceGeometry - A local vector with face geometry
1200: - cellGeometry - A local vector with cell geometry

1202:   Output Parameters:
1203: + Nface - The number of faces with field values
1204: . fgeom - The extract the face centroid and normal
1205: - vol   - The cell volume

1207:   Level: developer

1209: .seealso: DMPlexGetCellFields()
1210: @*/
1211: PetscErrorCode DMPlexGetFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom **fgeom, PetscReal **vol)
1212: {
1213:   DM                 dmFace, dmCell;
1214:   DMLabel            ghostLabel;
1215:   const PetscScalar *facegeom, *cellgeom;
1216:   PetscInt           dim, numFaces = fEnd - fStart, iface, face;
1217:   PetscErrorCode     ierr;

1225:   DMGetDimension(dm, &dim);
1226:   DMGetLabel(dm, "ghost", &ghostLabel);
1227:   VecGetDM(faceGeometry, &dmFace);
1228:   VecGetArrayRead(faceGeometry, &facegeom);
1229:   VecGetDM(cellGeometry, &dmCell);
1230:   VecGetArrayRead(cellGeometry, &cellgeom);
1231:   PetscMalloc1(numFaces, fgeom);
1232:   DMGetWorkArray(dm, numFaces*2, MPIU_SCALAR, vol);
1233:   for (face = fStart, iface = 0; face < fEnd; ++face) {
1234:     const PetscInt        *cells;
1235:     PetscFVFaceGeom       *fg;
1236:     PetscFVCellGeom       *cgL, *cgR;
1237:     PetscFVFaceGeom       *fgeoml = *fgeom;
1238:     PetscReal             *voll   = *vol;
1239:     PetscInt               ghost, d, nchild, nsupp;

1241:     DMLabelGetValue(ghostLabel, face, &ghost);
1242:     DMPlexGetSupportSize(dm, face, &nsupp);
1243:     DMPlexGetTreeChildren(dm, face, &nchild, NULL);
1244:     if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
1245:     DMPlexPointLocalRead(dmFace, face, facegeom, &fg);
1246:     DMPlexGetSupport(dm, face, &cells);
1247:     DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cgL);
1248:     DMPlexPointLocalRead(dmCell, cells[1], cellgeom, &cgR);
1249:     for (d = 0; d < dim; ++d) {
1250:       fgeoml[iface].centroid[d] = fg->centroid[d];
1251:       fgeoml[iface].normal[d]   = fg->normal[d];
1252:     }
1253:     voll[iface*2+0] = cgL->volume;
1254:     voll[iface*2+1] = cgR->volume;
1255:     ++iface;
1256:   }
1257:   *Nface = iface;
1258:   VecRestoreArrayRead(faceGeometry, &facegeom);
1259:   VecRestoreArrayRead(cellGeometry, &cellgeom);
1260:   return(0);
1261: }

1263: /*@C
1264:   DMPlexRestoreFaceGeometry - Restore the field values values for a chunk of faces

1266:   Input Parameters:
1267: + dm     - The DM
1268: . fStart - The first face to include
1269: . fEnd   - The first face to exclude
1270: . faceGeometry - A local vector with face geometry
1271: - cellGeometry - A local vector with cell geometry

1273:   Output Parameters:
1274: + Nface - The number of faces with field values
1275: . fgeom - The extract the face centroid and normal
1276: - vol   - The cell volume

1278:   Level: developer

1280: .seealso: DMPlexGetFaceFields()
1281: @*/
1282: PetscErrorCode DMPlexRestoreFaceGeometry(DM dm, PetscInt fStart, PetscInt fEnd, Vec faceGeometry, Vec cellGeometry, PetscInt *Nface, PetscFVFaceGeom **fgeom, PetscReal **vol)
1283: {

1287:   PetscFree(*fgeom);
1288:   DMRestoreWorkArray(dm, 0, MPIU_REAL, vol);
1289:   return(0);
1290: }

1292: static PetscErrorCode DMPlexComputeBdResidual_Single_Internal(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, Vec locF, DMField coordField, IS facetIS)
1293: {
1294:   DM_Plex         *mesh = (DM_Plex *) dm->data;
1295:   DM               plex = NULL, plexA = NULL;
1296:   PetscDS          prob, probAux = NULL;
1297:   PetscSection     section, sectionAux = NULL;
1298:   Vec              locA = NULL;
1299:   PetscScalar     *u = NULL, *u_t = NULL, *a = NULL, *elemVec = NULL;
1300:   PetscInt         v;
1301:   PetscInt         totDim, totDimAux = 0;
1302:   PetscErrorCode   ierr;

1305:   DMConvert(dm, DMPLEX, &plex);
1306:   DMGetSection(dm, &section);
1307:   DMGetDS(dm, &prob);
1308:   PetscDSGetTotalDimension(prob, &totDim);
1309:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);
1310:   if (locA) {
1311:     DM dmAux;

1313:     VecGetDM(locA, &dmAux);
1314:     DMConvert(dmAux, DMPLEX, &plexA);
1315:     DMGetDS(plexA, &probAux);
1316:     PetscDSGetTotalDimension(probAux, &totDimAux);
1317:     DMGetSection(plexA, &sectionAux);
1318:   }
1319:   for (v = 0; v < numValues; ++v) {
1320:     PetscFEGeom    *fgeom;
1321:     PetscInt        maxDegree;
1322:     PetscQuadrature qGeom = NULL;
1323:     IS              pointIS;
1324:     const PetscInt *points;
1325:     PetscInt        numFaces, face, Nq;

1327:     DMLabelGetStratumIS(label, values[v], &pointIS);
1328:     if (!pointIS) continue; /* No points with that id on this process */
1329:     {
1330:       IS isectIS;

1332:       /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
1333:       ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);
1334:       ISDestroy(&pointIS);
1335:       pointIS = isectIS;
1336:     }
1337:     ISGetLocalSize(pointIS,&numFaces);
1338:     ISGetIndices(pointIS,&points);
1339:     PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim, &elemVec, locA ? numFaces*totDimAux : 0, &a);
1340:     DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);
1341:     if (maxDegree <= 1) {
1342:       DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);
1343:     }
1344:     if (!qGeom) {
1345:       PetscFE fe;

1347:       PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1348:       PetscFEGetFaceQuadrature(fe, &qGeom);
1349:       PetscObjectReference((PetscObject)qGeom);
1350:     }
1351:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
1352:     DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
1353:     for (face = 0; face < numFaces; ++face) {
1354:       const PetscInt point = points[face], *support, *cone;
1355:       PetscScalar   *x     = NULL;
1356:       PetscInt       i, coneSize, faceLoc;

1358:       DMPlexGetSupport(dm, point, &support);
1359:       DMPlexGetConeSize(dm, support[0], &coneSize);
1360:       DMPlexGetCone(dm, support[0], &cone);
1361:       for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break;
1362:       if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %D in cone of support[0] %D", point, support[0]);
1363:       fgeom->face[face][0] = faceLoc;
1364:       DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);
1365:       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
1366:       DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);
1367:       if (locX_t) {
1368:         DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);
1369:         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
1370:         DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);
1371:       }
1372:       if (locA) {
1373:         DMLabel  spmap;
1374:         PetscInt subp = support[0];

1376:         /* If dm is a submesh, do not get subpoint */
1377:         DMPlexGetSubpointMap(dm, &spmap);
1378:         if (!spmap) {DMPlexGetSubpoint(plexA, support[0], &subp);}
1379:         DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);
1380:         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
1381:         DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);
1382:       }
1383:     }
1384:     PetscMemzero(elemVec, numFaces*totDim * sizeof(PetscScalar));
1385:     {
1386:       PetscFE         fe;
1387:       PetscInt        Nb;
1388:       PetscFEGeom     *chunkGeom = NULL;
1389:       /* Conforming batches */
1390:       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1391:       /* Remainder */
1392:       PetscInt        Nr, offset;

1394:       PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1395:       PetscFEGetDimension(fe, &Nb);
1396:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1397:       /* TODO: documentation is unclear about what is going on with these numbers: how should Nb / Nq factor in ? */
1398:       blockSize = Nb;
1399:       batchSize = numBlocks * blockSize;
1400:        PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1401:       numChunks = numFaces / (numBatches*batchSize);
1402:       Ne        = numChunks*numBatches*batchSize;
1403:       Nr        = numFaces % (numBatches*batchSize);
1404:       offset    = numFaces - Nr;
1405:       PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);
1406:       PetscFEIntegrateBdResidual(fe, prob, field, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);
1407:       PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom);
1408:       PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);
1409:       PetscFEIntegrateBdResidual(fe, prob, field, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, &elemVec[offset*totDim]);
1410:       PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);
1411:     }
1412:     for (face = 0; face < numFaces; ++face) {
1413:       const PetscInt point = points[face], *support;

1415:       if (mesh->printFEM > 1) {DMPrintCellVector(point, "BdResidual", totDim, &elemVec[face*totDim]);}
1416:       DMPlexGetSupport(plex, point, &support);
1417:       DMPlexVecSetClosure(plex, NULL, locF, support[0], &elemVec[face*totDim], ADD_ALL_VALUES);
1418:     }
1419:     DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
1420:     PetscQuadratureDestroy(&qGeom);
1421:     ISRestoreIndices(pointIS, &points);
1422:     ISDestroy(&pointIS);
1423:     PetscFree4(u, u_t, elemVec, a);
1424:   }
1425:   if (plex)  {DMDestroy(&plex);}
1426:   if (plexA) {DMDestroy(&plexA);}
1427:   return(0);
1428: }

1430: PetscErrorCode DMPlexComputeBdResidualSingle(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, Vec locF)
1431: {
1432:   DMField        coordField;
1433:   DMLabel        depthLabel;
1434:   IS             facetIS;
1435:   PetscInt       dim;

1439:   DMGetDimension(dm, &dim);
1440:   DMPlexGetDepthLabel(dm, &depthLabel);
1441:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
1442:   DMGetCoordinateField(dm, &coordField);
1443:   DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);
1444:   return(0);
1445: }

1447: PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1448: {
1449:   PetscDS        prob;
1450:   PetscInt       dim, numBd, bd;
1451:   DMLabel        depthLabel;
1452:   DMField        coordField = NULL;
1453:   IS             facetIS;

1457:   DMGetDS(dm, &prob);
1458:   DMPlexGetDepthLabel(dm, &depthLabel);
1459:   DMGetDimension(dm, &dim);
1460:   DMLabelGetStratumIS(depthLabel,dim - 1,&facetIS);
1461:   PetscDSGetNumBoundary(prob, &numBd);
1462:   DMGetCoordinateField(dm, &coordField);
1463:   for (bd = 0; bd < numBd; ++bd) {
1464:     DMBoundaryConditionType type;
1465:     const char             *bdLabel;
1466:     DMLabel                 label;
1467:     const PetscInt         *values;
1468:     PetscInt                field, numValues;
1469:     PetscObject             obj;
1470:     PetscClassId            id;

1472:     PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);
1473:     PetscDSGetDiscretization(prob, field, &obj);
1474:     PetscObjectGetClassId(obj, &id);
1475:     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
1476:     DMGetLabel(dm, bdLabel, &label);
1477:     DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);
1478:   }
1479:   ISDestroy(&facetIS);
1480:   return(0);
1481: }

1483: PetscErrorCode DMPlexComputeResidual_Internal(DM dm, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1484: {
1485:   DM_Plex         *mesh       = (DM_Plex *) dm->data;
1486:   const char      *name       = "Residual";
1487:   DM               dmAux      = NULL;
1488:   DM               dmGrad     = NULL;
1489:   DMLabel          ghostLabel = NULL;
1490:   PetscDS          prob       = NULL;
1491:   PetscDS          probAux    = NULL;
1492:   PetscSection     section    = NULL;
1493:   PetscBool        useFEM     = PETSC_FALSE;
1494:   PetscBool        useFVM     = PETSC_FALSE;
1495:   PetscBool        isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
1496:   PetscFV          fvm        = NULL;
1497:   PetscFVCellGeom *cgeomFVM   = NULL;
1498:   PetscFVFaceGeom *fgeomFVM   = NULL;
1499:   DMField          coordField = NULL;
1500:   Vec              locA, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL;
1501:   PetscScalar     *u = NULL, *u_t, *a, *uL, *uR;
1502:   IS               chunkIS;
1503:   const PetscInt  *cells;
1504:   PetscInt         cStart, cEnd, numCells;
1505:   PetscInt         Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
1506:   PetscInt         maxDegree = PETSC_MAX_INT;
1507:   PetscQuadrature  affineQuad = NULL, *quads = NULL;
1508:   PetscFEGeom     *affineGeom = NULL, **geoms = NULL;
1509:   PetscErrorCode   ierr;

1512:   PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);
1513:   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
1514:   /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
1515:   /* FEM+FVM */
1516:   /* 1: Get sizes from dm and dmAux */
1517:   DMGetSection(dm, &section);
1518:   DMGetLabel(dm, "ghost", &ghostLabel);
1519:   DMGetDS(dm, &prob);
1520:   PetscDSGetNumFields(prob, &Nf);
1521:   PetscDSGetTotalDimension(prob, &totDim);
1522:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);
1523:   if (locA) {
1524:     VecGetDM(locA, &dmAux);
1525:     DMGetDS(dmAux, &probAux);
1526:     PetscDSGetTotalDimension(probAux, &totDimAux);
1527:   }
1528:   /* 2: Get geometric data */
1529:   for (f = 0; f < Nf; ++f) {
1530:     PetscObject  obj;
1531:     PetscClassId id;
1532:     PetscBool    fimp;

1534:     PetscDSGetImplicit(prob, f, &fimp);
1535:     if (isImplicit != fimp) continue;
1536:     PetscDSGetDiscretization(prob, f, &obj);
1537:     PetscObjectGetClassId(obj, &id);
1538:     if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;}
1539:     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
1540:   }
1541:   if (useFEM) {
1542:     DMGetCoordinateField(dm, &coordField);
1543:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
1544:     if (maxDegree <= 1) {
1545:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);
1546:       if (affineQuad) {
1547:         DMSNESGetFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
1548:       }
1549:     } else {
1550:       PetscCalloc2(Nf,&quads,Nf,&geoms);
1551:       for (f = 0; f < Nf; ++f) {
1552:         PetscObject  obj;
1553:         PetscClassId id;
1554:         PetscBool    fimp;

1556:         PetscDSGetImplicit(prob, f, &fimp);
1557:         if (isImplicit != fimp) continue;
1558:         PetscDSGetDiscretization(prob, f, &obj);
1559:         PetscObjectGetClassId(obj, &id);
1560:         if (id == PETSCFE_CLASSID) {
1561:           PetscFE fe = (PetscFE) obj;

1563:           PetscFEGetQuadrature(fe, &quads[f]);
1564:           PetscObjectReference((PetscObject)quads[f]);
1565:           DMSNESGetFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
1566:         }
1567:       }
1568:     }
1569:   }
1570:   if (useFVM) {
1571:     DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);
1572:     VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);
1573:     VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);
1574:     /* Reconstruct and limit cell gradients */
1575:     DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);
1576:     if (dmGrad) {
1577:       DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
1578:       DMGetGlobalVector(dmGrad, &grad);
1579:       DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);
1580:       /* Communicate gradient values */
1581:       DMGetLocalVector(dmGrad, &locGrad);
1582:       DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);
1583:       DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);
1584:       DMRestoreGlobalVector(dmGrad, &grad);
1585:     }
1586:     /* Handle non-essential (e.g. outflow) boundary values */
1587:     DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);
1588:   }
1589:   /* Loop over chunks */
1590:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
1591:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
1592:   if (useFEM) {ISCreate(PETSC_COMM_SELF, &chunkIS);}
1593:   numCells      = cEnd - cStart;
1594:   numChunks     = 1;
1595:   cellChunkSize = numCells/numChunks;
1596:   faceChunkSize = (fEnd - fStart)/numChunks;
1597:   numChunks     = PetscMin(1,numCells);
1598:   for (chunk = 0; chunk < numChunks; ++chunk) {
1599:     PetscScalar     *elemVec, *fluxL, *fluxR;
1600:     PetscReal       *vol;
1601:     PetscFVFaceGeom *fgeom;
1602:     PetscInt         cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c;
1603:     PetscInt         fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = 0, face;

1605:     /* Extract field coefficients */
1606:     if (useFEM) {
1607:       ISGetPointSubrange(chunkIS, cS, cE, cells);
1608:       DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
1609:       DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
1610:       PetscMemzero(elemVec, numCells*totDim * sizeof(PetscScalar));
1611:     }
1612:     if (useFVM) {
1613:       DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);
1614:       DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);
1615:       DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);
1616:       DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);
1617:       PetscMemzero(fluxL, numFaces*totDim * sizeof(PetscScalar));
1618:       PetscMemzero(fluxR, numFaces*totDim * sizeof(PetscScalar));
1619:     }
1620:     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
1621:     /* Loop over fields */
1622:     for (f = 0; f < Nf; ++f) {
1623:       PetscObject  obj;
1624:       PetscClassId id;
1625:       PetscBool    fimp;
1626:       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

1628:       PetscDSGetImplicit(prob, f, &fimp);
1629:       if (isImplicit != fimp) continue;
1630:       PetscDSGetDiscretization(prob, f, &obj);
1631:       PetscObjectGetClassId(obj, &id);
1632:       if (id == PETSCFE_CLASSID) {
1633:         PetscFE         fe = (PetscFE) obj;
1634:         PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[f];
1635:         PetscFEGeom    *chunkGeom = NULL;
1636:         PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
1637:         PetscInt        Nq, Nb;

1639:         PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1640:         PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
1641:         PetscFEGetDimension(fe, &Nb);
1642:         blockSize = Nb;
1643:         batchSize = numBlocks * blockSize;
1644:         PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1645:         numChunks = numCells / (numBatches*batchSize);
1646:         Ne        = numChunks*numBatches*batchSize;
1647:         Nr        = numCells % (numBatches*batchSize);
1648:         offset    = numCells - Nr;
1649:         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
1650:         /*   For FV, I think we use a P0 basis and the cell coefficients (for subdivided cells, we can tweak the basis tabulation to be the indicator function) */
1651:         PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);
1652:         PetscFEIntegrateResidual(fe, prob, f, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);
1653:         PetscFEGeomGetChunk(geom,offset,numCells,&chunkGeom);
1654:         PetscFEIntegrateResidual(fe, prob, f, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);
1655:         PetscFEGeomRestoreChunk(geom,offset,numCells,&chunkGeom);
1656:       } else if (id == PETSCFV_CLASSID) {
1657:         PetscFV fv = (PetscFV) obj;

1659:         Ne = numFaces;
1660:         /* Riemann solve over faces (need fields at face centroids) */
1661:         /*   We need to evaluate FE fields at those coordinates */
1662:         PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);
1663:       } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
1664:     }
1665:     /* Loop over domain */
1666:     if (useFEM) {
1667:       /* Add elemVec to locX */
1668:       for (c = cS; c < cE; ++c) {
1669:         const PetscInt cell = cells ? cells[c] : c;
1670:         const PetscInt cind = c - cStart;

1672:         if (mesh->printFEM > 1) {DMPrintCellVector(cell, name, totDim, &elemVec[cind*totDim]);}
1673:         if (ghostLabel) {
1674:           PetscInt ghostVal;

1676:           DMLabelGetValue(ghostLabel,cell,&ghostVal);
1677:           if (ghostVal > 0) continue;
1678:         }
1679:         DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);
1680:       }
1681:     }
1682:     if (useFVM) {
1683:       PetscScalar *fa;
1684:       PetscInt     iface;

1686:       VecGetArray(locF, &fa);
1687:       for (f = 0; f < Nf; ++f) {
1688:         PetscFV      fv;
1689:         PetscObject  obj;
1690:         PetscClassId id;
1691:         PetscInt     foff, pdim;

1693:         PetscDSGetDiscretization(prob, f, &obj);
1694:         PetscDSGetFieldOffset(prob, f, &foff);
1695:         PetscObjectGetClassId(obj, &id);
1696:         if (id != PETSCFV_CLASSID) continue;
1697:         fv   = (PetscFV) obj;
1698:         PetscFVGetNumComponents(fv, &pdim);
1699:         /* Accumulate fluxes to cells */
1700:         for (face = fS, iface = 0; face < fE; ++face) {
1701:           const PetscInt *scells;
1702:           PetscScalar    *fL = NULL, *fR = NULL;
1703:           PetscInt        ghost, d, nsupp, nchild;

1705:           DMLabelGetValue(ghostLabel, face, &ghost);
1706:           DMPlexGetSupportSize(dm, face, &nsupp);
1707:           DMPlexGetTreeChildren(dm, face, &nchild, NULL);
1708:           if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
1709:           DMPlexGetSupport(dm, face, &scells);
1710:           DMLabelGetValue(ghostLabel,scells[0],&ghost);
1711:           if (ghost <= 0) {DMPlexPointLocalFieldRef(dm, scells[0], f, fa, &fL);}
1712:           DMLabelGetValue(ghostLabel,scells[1],&ghost);
1713:           if (ghost <= 0) {DMPlexPointLocalFieldRef(dm, scells[1], f, fa, &fR);}
1714:           for (d = 0; d < pdim; ++d) {
1715:             if (fL) fL[d] -= fluxL[iface*totDim+foff+d];
1716:             if (fR) fR[d] += fluxR[iface*totDim+foff+d];
1717:           }
1718:           ++iface;
1719:         }
1720:       }
1721:       VecRestoreArray(locF, &fa);
1722:     }
1723:     /* Handle time derivative */
1724:     if (locX_t) {
1725:       PetscScalar *x_t, *fa;

1727:       VecGetArray(locF, &fa);
1728:       VecGetArray(locX_t, &x_t);
1729:       for (f = 0; f < Nf; ++f) {
1730:         PetscFV      fv;
1731:         PetscObject  obj;
1732:         PetscClassId id;
1733:         PetscInt     pdim, d;

1735:         PetscDSGetDiscretization(prob, f, &obj);
1736:         PetscObjectGetClassId(obj, &id);
1737:         if (id != PETSCFV_CLASSID) continue;
1738:         fv   = (PetscFV) obj;
1739:         PetscFVGetNumComponents(fv, &pdim);
1740:         for (c = cS; c < cE; ++c) {
1741:           const PetscInt cell = cells ? cells[c] : c;
1742:           PetscScalar   *u_t, *r;

1744:           if (ghostLabel) {
1745:             PetscInt ghostVal;

1747:             DMLabelGetValue(ghostLabel, cell, &ghostVal);
1748:             if (ghostVal > 0) continue;
1749:           }
1750:           DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);
1751:           DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);
1752:           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
1753:         }
1754:       }
1755:       VecRestoreArray(locX_t, &x_t);
1756:       VecRestoreArray(locF, &fa);
1757:     }
1758:     if (useFEM) {
1759:       DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
1760:       DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
1761:     }
1762:     if (useFVM) {
1763:       DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);
1764:       DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);
1765:       DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);
1766:       DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);
1767:       if (dmGrad) {DMRestoreLocalVector(dmGrad, &locGrad);}
1768:     }
1769:   }
1770:   if (useFEM) {ISDestroy(&chunkIS);}
1771:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);

1773:   if (useFEM) {
1774:     DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user);

1776:     if (maxDegree <= 1) {
1777:       DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
1778:       PetscQuadratureDestroy(&affineQuad);
1779:     } else {
1780:       for (f = 0; f < Nf; ++f) {
1781:         DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
1782:         PetscQuadratureDestroy(&quads[f]);
1783:       }
1784:       PetscFree2(quads,geoms);
1785:     }
1786:   }

1788:   /* FEM */
1789:   /* 1: Get sizes from dm and dmAux */
1790:   /* 2: Get geometric data */
1791:   /* 3: Handle boundary values */
1792:   /* 4: Loop over domain */
1793:   /*   Extract coefficients */
1794:   /* Loop over fields */
1795:   /*   Set tiling for FE*/
1796:   /*   Integrate FE residual to get elemVec */
1797:   /*     Loop over subdomain */
1798:   /*       Loop over quad points */
1799:   /*         Transform coords to real space */
1800:   /*         Evaluate field and aux fields at point */
1801:   /*         Evaluate residual at point */
1802:   /*         Transform residual to real space */
1803:   /*       Add residual to elemVec */
1804:   /* Loop over domain */
1805:   /*   Add elemVec to locX */

1807:   /* FVM */
1808:   /* Get geometric data */
1809:   /* If using gradients */
1810:   /*   Compute gradient data */
1811:   /*   Loop over domain faces */
1812:   /*     Count computational faces */
1813:   /*     Reconstruct cell gradient */
1814:   /*   Loop over domain cells */
1815:   /*     Limit cell gradients */
1816:   /* Handle boundary values */
1817:   /* Loop over domain faces */
1818:   /*   Read out field, centroid, normal, volume for each side of face */
1819:   /* Riemann solve over faces */
1820:   /* Loop over domain faces */
1821:   /*   Accumulate fluxes to cells */
1822:   /* TODO Change printFEM to printDisc here */
1823:   if (mesh->printFEM) {
1824:     Vec         locFbc;
1825:     PetscInt    pStart, pEnd, p, maxDof;
1826:     PetscScalar *zeroes;

1828:     VecDuplicate(locF,&locFbc);
1829:     VecCopy(locF,locFbc);
1830:     PetscSectionGetChart(section,&pStart,&pEnd);
1831:     PetscSectionGetMaxDof(section,&maxDof);
1832:     PetscCalloc1(maxDof,&zeroes);
1833:     for (p = pStart; p < pEnd; p++) {
1834:       VecSetValuesSection(locFbc,section,p,zeroes,INSERT_BC_VALUES);
1835:     }
1836:     PetscFree(zeroes);
1837:     DMPrintLocalVec(dm, name, mesh->printTol, locFbc);
1838:     VecDestroy(&locFbc);
1839:   }
1840:   PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);
1841:   return(0);
1842: }

1844: static PetscErrorCode DMPlexComputeResidualFEM_Check_Internal(DM dm, Vec X, Vec X_t, PetscReal t, Vec F, void *user)
1845: {
1846:   DM                dmCh, dmAux;
1847:   Vec               A;
1848:   DMField           coordField = NULL;
1849:   PetscDS           prob, probCh, probAux = NULL;
1850:   PetscSection      section, sectionAux;
1851:   PetscScalar      *elemVec, *elemVecCh, *u, *u_t, *a = NULL;
1852:   PetscInt          Nf, f, numCells, cStart, cEnd, c;
1853:   PetscInt          totDim, totDimAux = 0, diffCell = 0;
1854:   PetscInt          depth;
1855:   PetscInt          maxDegree;
1856:   IS                cellIS;
1857:   DMLabel           depthLabel;
1858:   PetscErrorCode    ierr;

1861:   DMGetSection(dm, &section);
1862:   DMGetDS(dm, &prob);
1863:   PetscDSGetTotalDimension(prob, &totDim);
1864:   PetscSectionGetNumFields(section, &Nf);
1865:   DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
1866:   numCells = cEnd - cStart;
1867:   PetscObjectQuery((PetscObject) dm, "dmCh", (PetscObject *) &dmCh);
1868:   DMGetDS(dmCh, &probCh);
1869:   PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);
1870:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);
1871:   if (dmAux) {
1872:     DMGetSection(dmAux, &sectionAux);
1873:     DMGetDS(dmAux, &probAux);
1874:     PetscDSGetTotalDimension(probAux, &totDimAux);
1875:   }
1876:   VecSet(F, 0.0);
1877:   PetscMalloc3(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,numCells*totDim,&elemVec);
1878:   PetscMalloc1(numCells*totDim,&elemVecCh);
1879:   if (dmAux) {PetscMalloc1(numCells*totDimAux, &a);}
1880:   DMPlexGetDepthLabel(dm, &depthLabel);
1881:   DMPlexGetDepth(dm,&depth);
1882:   DMLabelGetStratumIS(depthLabel,depth,&cellIS);
1883:   DMGetCoordinateField(dm, &coordField);
1884:   for (c = cStart; c < cEnd; ++c) {
1885:     PetscScalar *x = NULL, *x_t = NULL;
1886:     PetscInt     i;

1888:     DMPlexVecGetClosure(dm, section, X, c, NULL, &x);
1889:     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
1890:     DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);
1891:     if (X_t) {
1892:       DMPlexVecGetClosure(dm, section, X_t, c, NULL, &x_t);
1893:       for (i = 0; i < totDim; ++i) u_t[c*totDim+i] = x_t[i];
1894:       DMPlexVecRestoreClosure(dm, section, X_t, c, NULL, &x_t);
1895:     }
1896:     if (dmAux) {
1897:       DM dmAuxPlex;

1899:       DMSNESConvertPlex(dmAux,&dmAuxPlex, PETSC_FALSE);
1900:       DMPlexVecGetClosure(dmAuxPlex, sectionAux, A, c, NULL, &x);
1901:       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
1902:       DMPlexVecRestoreClosure(dmAuxPlex, sectionAux, A, c, NULL, &x);
1903:       DMDestroy(&dmAuxPlex);
1904:     }
1905:   }
1906:   for (f = 0; f < Nf; ++f) {
1907:     PetscFE  fe, feCh;
1908:     PetscInt Nq, Nb;
1909:     /* Conforming batches */
1910:     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1911:     /* Remainder */
1912:     PetscInt Nr, offset;
1913:     PetscQuadrature qGeom = NULL;
1914:     PetscFEGeom *cgeomFEM, *chunkGeom = NULL;

1916:     PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);
1917:     PetscDSGetDiscretization(probCh, f, (PetscObject *) &feCh);
1918:     PetscFEGetDimension(fe, &Nb);
1919:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1920:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
1921:     if (maxDegree <= 1) {
1922:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);
1923:     }
1924:     if (!qGeom) {
1925:       PetscFEGetQuadrature(fe, &qGeom);
1926:       PetscObjectReference((PetscObject)qGeom);
1927:     }
1928:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
1929:     DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
1930:     blockSize = Nb;
1931:     batchSize = numBlocks * blockSize;
1932:      PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1933:     numChunks = numCells / (numBatches*batchSize);
1934:     Ne        = numChunks*numBatches*batchSize;
1935:     Nr        = numCells % (numBatches*batchSize);
1936:     offset    = numCells - Nr;
1937:     PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
1938:     PetscFEIntegrateResidual(fe, prob, f, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);
1939:     PetscFEIntegrateResidual(feCh, prob, f, Ne, chunkGeom, u, u_t, probAux, a, t, elemVecCh);
1940:     PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&chunkGeom);
1941:     PetscFEIntegrateResidual(fe, prob, f, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);
1942:     PetscFEIntegrateResidual(feCh, prob, f, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVecCh[offset*totDim]);
1943:     PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&chunkGeom);
1944:     DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
1945:     PetscQuadratureDestroy(&qGeom);
1946:   }
1947:   ISDestroy(&cellIS);
1948:   for (c = cStart; c < cEnd; ++c) {
1949:     PetscBool diff = PETSC_FALSE;
1950:     PetscInt  d;

1952:     for (d = 0; d < totDim; ++d) if (PetscAbsScalar(elemVec[c*totDim+d] - elemVecCh[c*totDim+d]) > 1.0e-7) {diff = PETSC_TRUE;break;}
1953:     if (diff) {
1954:       PetscPrintf(PetscObjectComm((PetscObject) dm), "Different cell %d\n", c);
1955:       DMPrintCellVector(c, "Residual", totDim, &elemVec[c*totDim]);
1956:       DMPrintCellVector(c, "Check Residual", totDim, &elemVecCh[c*totDim]);
1957:       ++diffCell;
1958:     }
1959:     if (diffCell > 9) break;
1960:     DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*totDim], ADD_ALL_VALUES);
1961:   }
1962:   PetscFree3(u,u_t,elemVec);
1963:   PetscFree(elemVecCh);
1964:   if (dmAux) {PetscFree(a);}
1965:   return(0);
1966: }

1968: /*@
1969:   DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user

1971:   Input Parameters:
1972: + dm - The mesh
1973: . X  - Local solution
1974: - user - The user context

1976:   Output Parameter:
1977: . F  - Local output vector

1979:   Level: developer

1981: .seealso: DMPlexComputeJacobianAction()
1982: @*/
1983: PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1984: {
1985:   PetscObject    check;
1986:   DM             plex;
1987:   IS             cellIS;
1988:   PetscInt       depth;

1992:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1993:   DMPlexGetDepth(plex, &depth);
1994:   DMGetStratumIS(plex, "dim", depth, &cellIS);
1995:   if (!cellIS) {
1996:     DMGetStratumIS(plex, "depth", depth, &cellIS);
1997:   }
1998:   /* The dmCh is used to check two mathematically equivalent discretizations for computational equivalence */
1999:   PetscObjectQuery((PetscObject) plex, "dmCh", &check);
2000:   if (check) {DMPlexComputeResidualFEM_Check_Internal(plex, X, NULL, 0.0, F, user);}
2001:   else       {DMPlexComputeResidual_Internal(plex, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);}
2002:   ISDestroy(&cellIS);
2003:   DMDestroy(&plex);
2004:   return(0);
2005: }

2007: /*@
2008:   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X

2010:   Input Parameters:
2011: + dm - The mesh
2012: - user - The user context

2014:   Output Parameter:
2015: . X  - Local solution

2017:   Level: developer

2019: .seealso: DMPlexComputeJacobianAction()
2020: @*/
2021: PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
2022: {
2023:   DM             plex;

2027:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
2028:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);
2029:   DMDestroy(&plex);
2030:   return(0);
2031: }

2033: PetscErrorCode DMPlexComputeBdJacobian_Single_Internal(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt fieldI, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, DMField coordField, IS facetIS)
2034: {
2035:   DM_Plex       *mesh = (DM_Plex *) dm->data;
2036:   DM             plex = NULL, plexA = NULL;
2037:   PetscDS        prob, probAux = NULL;
2038:   PetscSection   section, sectionAux = NULL;
2039:   PetscSection   globalSection, subSection = NULL;
2040:   Vec            locA = NULL;
2041:   PetscScalar   *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL;
2042:   PetscInt       v;
2043:   PetscInt       Nf, totDim, totDimAux = 0;
2044:   PetscBool      isMatISP;

2048:   DMConvert(dm, DMPLEX, &plex);
2049:   DMGetSection(dm, &section);
2050:   DMGetDS(dm, &prob);
2051:   PetscDSGetNumFields(prob, &Nf);
2052:   PetscDSGetTotalDimension(prob, &totDim);
2053:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);
2054:   if (locA) {
2055:     DM dmAux;

2057:     VecGetDM(locA, &dmAux);
2058:     DMConvert(dmAux, DMPLEX, &plexA);
2059:     DMGetDS(plexA, &probAux);
2060:     PetscDSGetTotalDimension(probAux, &totDimAux);
2061:     DMGetSection(plexA, &sectionAux);
2062:   }

2064:   PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);
2065:   DMGetGlobalSection(dm, &globalSection);
2066:   if (isMatISP) {DMPlexGetSubdomainSection(dm, &subSection);}
2067:   for (v = 0; v < numValues; ++v) {
2068:     PetscFEGeom    *fgeom;
2069:     PetscInt        maxDegree;
2070:     PetscQuadrature qGeom = NULL;
2071:     IS              pointIS;
2072:     const PetscInt *points;
2073:     PetscInt        numFaces, face, Nq;

2075:     DMLabelGetStratumIS(label, values[v], &pointIS);
2076:     if (!pointIS) continue; /* No points with that id on this process */
2077:     {
2078:       IS isectIS;

2080:       /* TODO: Special cases of ISIntersect where it is quick to check a prior if one is a superset of the other */
2081:       ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);
2082:       ISDestroy(&pointIS);
2083:       pointIS = isectIS;
2084:     }
2085:     ISGetLocalSize(pointIS, &numFaces);
2086:     ISGetIndices(pointIS, &points);
2087:     PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim*totDim, &elemMat, locA ? numFaces*totDimAux : 0, &a);
2088:     DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);
2089:     if (maxDegree <= 1) {
2090:       DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);
2091:     }
2092:     if (!qGeom) {
2093:       PetscFE fe;

2095:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
2096:       PetscFEGetFaceQuadrature(fe, &qGeom);
2097:       PetscObjectReference((PetscObject)qGeom);
2098:     }
2099:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
2100:     DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
2101:     for (face = 0; face < numFaces; ++face) {
2102:       const PetscInt point = points[face], *support, *cone;
2103:       PetscScalar   *x     = NULL;
2104:       PetscInt       i, coneSize, faceLoc;

2106:       DMPlexGetSupport(dm, point, &support);
2107:       DMPlexGetConeSize(dm, support[0], &coneSize);
2108:       DMPlexGetCone(dm, support[0], &cone);
2109:       for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break;
2110:       if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %d in cone of support[0] %d", point, support[0]);
2111:       fgeom->face[face][0] = faceLoc;
2112:       DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);
2113:       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
2114:       DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);
2115:       if (locX_t) {
2116:         DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);
2117:         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
2118:         DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);
2119:       }
2120:       if (locA) {
2121:         PetscInt subp;
2122:         DMPlexGetSubpoint(plexA, support[0], &subp);
2123:         DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);
2124:         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
2125:         DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);
2126:       }
2127:     }
2128:     PetscMemzero(elemMat, numFaces*totDim*totDim * sizeof(PetscScalar));
2129:     {
2130:       PetscFE         fe;
2131:       PetscInt        Nb;
2132:       /* Conforming batches */
2133:       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2134:       /* Remainder */
2135:       PetscFEGeom    *chunkGeom = NULL;
2136:       PetscInt        fieldJ, Nr, offset;

2138:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
2139:       PetscFEGetDimension(fe, &Nb);
2140:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
2141:       blockSize = Nb;
2142:       batchSize = numBlocks * blockSize;
2143:       PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
2144:       numChunks = numFaces / (numBatches*batchSize);
2145:       Ne        = numChunks*numBatches*batchSize;
2146:       Nr        = numFaces % (numBatches*batchSize);
2147:       offset    = numFaces - Nr;
2148:       PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);
2149:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2150:         PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
2151:       }
2152:       PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);
2153:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2154:         PetscFEIntegrateBdJacobian(fe, prob, fieldI, fieldJ, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, X_tShift, &elemMat[offset*totDim*totDim]);
2155:       }
2156:       PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);
2157:     }
2158:     for (face = 0; face < numFaces; ++face) {
2159:       const PetscInt point = points[face], *support;

2161:       if (mesh->printFEM > 1) {DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face*totDim*totDim]);}
2162:       DMPlexGetSupport(plex, point, &support);
2163:       if (!isMatISP) {
2164:         DMPlexMatSetClosure(plex, section, globalSection, JacP, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);
2165:       } else {
2166:         Mat lJ;

2168:         MatISGetLocalMat(JacP, &lJ);
2169:         DMPlexMatSetClosure(plex, section, subSection, lJ, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);
2170:       }
2171:     }
2172:     DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
2173:     PetscQuadratureDestroy(&qGeom);
2174:     ISRestoreIndices(pointIS, &points);
2175:     ISDestroy(&pointIS);
2176:     PetscFree4(u, u_t, elemMat, a);
2177:   }
2178:   if (plex)  {DMDestroy(&plex);}
2179:   if (plexA) {DMDestroy(&plexA);}
2180:   return(0);
2181: }

2183: PetscErrorCode DMPlexComputeBdJacobianSingle(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP)
2184: {
2185:   DMField        coordField;
2186:   DMLabel        depthLabel;
2187:   IS             facetIS;
2188:   PetscInt       dim;

2192:   DMGetDimension(dm, &dim);
2193:   DMPlexGetDepthLabel(dm, &depthLabel);
2194:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
2195:   DMGetCoordinateField(dm, &coordField);
2196:   DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);
2197:   return(0);
2198: }

2200: PetscErrorCode DMPlexComputeBdJacobian_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
2201: {
2202:   PetscDS          prob;
2203:   PetscInt         dim, numBd, bd;
2204:   DMLabel          depthLabel;
2205:   DMField          coordField = NULL;
2206:   IS               facetIS;
2207:   PetscErrorCode   ierr;

2210:   DMGetDS(dm, &prob);
2211:   DMPlexGetDepthLabel(dm, &depthLabel);
2212:   DMGetDimension(dm, &dim);
2213:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
2214:   PetscDSGetNumBoundary(prob, &numBd);
2215:   DMGetCoordinateField(dm, &coordField);
2216:   for (bd = 0; bd < numBd; ++bd) {
2217:     DMBoundaryConditionType type;
2218:     const char             *bdLabel;
2219:     DMLabel                 label;
2220:     const PetscInt         *values;
2221:     PetscInt                fieldI, numValues;
2222:     PetscObject             obj;
2223:     PetscClassId            id;

2225:     PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &fieldI, NULL, NULL, NULL, &numValues, &values, NULL);
2226:     PetscDSGetDiscretization(prob, fieldI, &obj);
2227:     PetscObjectGetClassId(obj, &id);
2228:     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
2229:     DMGetLabel(dm, bdLabel, &label);
2230:     DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, fieldI, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);
2231:   }
2232:   ISDestroy(&facetIS);
2233:   return(0);
2234: }

2236: PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user)
2237: {
2238:   DM_Plex        *mesh  = (DM_Plex *) dm->data;
2239:   const char     *name  = "Jacobian";
2240:   DM              dmAux, plex;
2241:   Vec             A;
2242:   DMField         coordField;
2243:   PetscDS         prob, probAux = NULL;
2244:   PetscSection    section, globalSection, subSection, sectionAux;
2245:   PetscScalar    *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
2246:   const PetscInt *cells;
2247:   PetscInt        Nf, fieldI, fieldJ;
2248:   PetscInt        totDim, totDimAux, cStart, cEnd, numCells, c;
2249:   PetscBool       isMatIS, isMatISP, hasJac, hasPrec, hasDyn, hasFV = PETSC_FALSE;
2250:   PetscErrorCode  ierr;

2253:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
2254:   DMGetSection(dm, &section);
2255:   PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);
2256:   DMGetGlobalSection(dm, &globalSection);
2257:   if (isMatISP) {DMPlexGetSubdomainSection(dm, &subSection);}
2258:   DMGetDS(dm, &prob);
2259:   PetscDSGetTotalDimension(prob, &totDim);
2260:   PetscDSHasJacobian(prob, &hasJac);
2261:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
2262:   PetscDSHasDynamicJacobian(prob, &hasDyn);
2263:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
2264:   PetscSectionGetNumFields(section, &Nf);
2265:   ISGetLocalSize(cellIS, &numCells);
2266:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
2267:   PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);
2268:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);
2269:   if (dmAux) {
2270:     DMConvert(dmAux, DMPLEX, &plex);
2271:     DMGetSection(plex, &sectionAux);
2272:     DMGetDS(dmAux, &probAux);
2273:     PetscDSGetTotalDimension(probAux, &totDimAux);
2274:   }
2275:   PetscMalloc5(numCells*totDim,&u,X_t ? numCells*totDim : 0,&u_t,hasJac ? numCells*totDim*totDim : 0,&elemMat,hasPrec ? numCells*totDim*totDim : 0, &elemMatP,hasDyn ? numCells*totDim*totDim : 0, &elemMatD);
2276:   if (dmAux) {PetscMalloc1(numCells*totDimAux, &a);}
2277:   DMGetCoordinateField(dm, &coordField);
2278:   for (c = cStart; c < cEnd; ++c) {
2279:     const PetscInt cell = cells ? cells[c] : c;
2280:     const PetscInt cind = c - cStart;
2281:     PetscScalar   *x = NULL,  *x_t = NULL;
2282:     PetscInt       i;

2284:     DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
2285:     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
2286:     DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
2287:     if (X_t) {
2288:       DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
2289:       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
2290:       DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
2291:     }
2292:     if (dmAux) {
2293:       DMPlexVecGetClosure(plex, sectionAux, A, cell, NULL, &x);
2294:       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
2295:       DMPlexVecRestoreClosure(plex, sectionAux, A, cell, NULL, &x);
2296:     }
2297:   }
2298:   if (hasJac)  {PetscMemzero(elemMat,  numCells*totDim*totDim * sizeof(PetscScalar));}
2299:   if (hasPrec) {PetscMemzero(elemMatP, numCells*totDim*totDim * sizeof(PetscScalar));}
2300:   if (hasDyn)  {PetscMemzero(elemMatD, numCells*totDim*totDim * sizeof(PetscScalar));}
2301:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2302:     PetscClassId    id;
2303:     PetscFE         fe;
2304:     PetscQuadrature qGeom = NULL;
2305:     PetscInt        Nb;
2306:     /* Conforming batches */
2307:     PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2308:     /* Remainder */
2309:     PetscInt        Nr, offset, Nq;
2310:     PetscInt        maxDegree;
2311:     PetscFEGeom     *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

2313:     PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
2314:     PetscObjectGetClassId((PetscObject) fe, &id);
2315:     if (id == PETSCFV_CLASSID) {hasFV = PETSC_TRUE; continue;}
2316:     PetscFEGetDimension(fe, &Nb);
2317:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
2318:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
2319:     if (maxDegree <= 1) {
2320:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);
2321:     }
2322:     if (!qGeom) {
2323:       PetscFEGetQuadrature(fe,&qGeom);
2324:       PetscObjectReference((PetscObject)qGeom);
2325:     }
2326:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
2327:     DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
2328:     blockSize = Nb;
2329:     batchSize = numBlocks * blockSize;
2330:     PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
2331:     numChunks = numCells / (numBatches*batchSize);
2332:     Ne        = numChunks*numBatches*batchSize;
2333:     Nr        = numCells % (numBatches*batchSize);
2334:     offset    = numCells - Nr;
2335:     PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
2336:     PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);
2337:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2338:       if (hasJac) {
2339:         PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
2340:         PetscFEIntegrateJacobian(fe, 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]);
2341:       }
2342:       if (hasPrec) {
2343:         PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP);
2344:         PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatP[offset*totDim*totDim]);
2345:       }
2346:       if (hasDyn) {
2347:         PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);
2348:         PetscFEIntegrateJacobian(fe, 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]);
2349:       }
2350:     }
2351:     PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);
2352:     PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);
2353:     DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
2354:     PetscQuadratureDestroy(&qGeom);
2355:   }
2356:   /*   Add contribution from X_t */
2357:   if (hasDyn) {for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];}
2358:   if (hasFV) {
2359:     PetscClassId id;
2360:     PetscFV      fv;
2361:     PetscInt     offsetI, NcI, NbI = 1, fc, f;

2363:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
2364:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fv);
2365:       PetscDSGetFieldOffset(prob, fieldI, &offsetI);
2366:       PetscObjectGetClassId((PetscObject) fv, &id);
2367:       if (id != PETSCFV_CLASSID) continue;
2368:       /* Put in the identity */
2369:       PetscFVGetNumComponents(fv, &NcI);
2370:       for (c = cStart; c < cEnd; ++c) {
2371:         const PetscInt cind    = c - cStart;
2372:         const PetscInt eOffset = cind*totDim*totDim;
2373:         for (fc = 0; fc < NcI; ++fc) {
2374:           for (f = 0; f < NbI; ++f) {
2375:             const PetscInt i = offsetI + f*NcI+fc;
2376:             if (hasPrec) {
2377:               if (hasJac) {elemMat[eOffset+i*totDim+i] = 1.0;}
2378:               elemMatP[eOffset+i*totDim+i] = 1.0;
2379:             } else {elemMat[eOffset+i*totDim+i] = 1.0;}
2380:           }
2381:         }
2382:       }
2383:     }
2384:     /* No allocated space for FV stuff, so ignore the zero entries */
2385:     MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
2386:   }
2387:   /* Insert values into matrix */
2388:   isMatIS = PETSC_FALSE;
2389:   if (hasPrec && hasJac) {
2390:     PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatIS);
2391:   }
2392:   if (isMatIS && !subSection) {
2393:     DMPlexGetSubdomainSection(dm, &subSection);
2394:   }
2395:   for (c = cStart; c < cEnd; ++c) {
2396:     const PetscInt cell = cells ? cells[c] : c;
2397:     const PetscInt cind = c - cStart;

2399:     if (hasPrec) {
2400:       if (hasJac) {
2401:         if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);}
2402:         if (!isMatIS) {
2403:           DMPlexMatSetClosure(dm, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2404:         } else {
2405:           Mat lJ;

2407:           MatISGetLocalMat(Jac,&lJ);
2408:           DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2409:         }
2410:       }
2411:       if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind*totDim*totDim]);}
2412:       if (!isMatISP) {
2413:         DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);
2414:       } else {
2415:         Mat lJ;

2417:         MatISGetLocalMat(JacP,&lJ);
2418:         DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);
2419:       }
2420:     } else {
2421:       if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);}
2422:       if (!isMatISP) {
2423:         DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2424:       } else {
2425:         Mat lJ;

2427:         MatISGetLocalMat(JacP,&lJ);
2428:         DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2429:       }
2430:     }
2431:   }
2432:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
2433:   if (hasFV) {MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);}
2434:   PetscFree5(u,u_t,elemMat,elemMatP,elemMatD);
2435:   if (dmAux) {
2436:     PetscFree(a);
2437:     DMDestroy(&plex);
2438:   }
2439:   /* Compute boundary integrals */
2440:   DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, user);
2441:   /* Assemble matrix */
2442:   if (hasJac && hasPrec) {
2443:     MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);
2444:     MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);
2445:   }
2446:   MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);
2447:   MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);
2448:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
2449:   return(0);
2450: }

2452: /*@
2453:   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.

2455:   Input Parameters:
2456: + dm - The mesh
2457: . cellIS - 
2458: . t  - The time
2459: . X_tShift - The multiplier for the Jacobian with repsect to X_t
2460: . X  - Local solution vector
2461: . X_t  - Time-derivative of the local solution vector
2462: . Y  - Local input vector
2463: - user - The user context

2465:   Output Parameter:
2466: . Z - Local output vector

2468:   Note:
2469:   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
2470:   like a GPU, or vectorize on a multicore machine.

2472:   Level: developer

2474: .seealso: FormFunctionLocal()
2475: @*/
2476: PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
2477: {
2478:   DM_Plex          *mesh  = (DM_Plex *) dm->data;
2479:   const char       *name  = "Jacobian";
2480:   DM                dmAux, plex, plexAux = NULL;
2481:   Vec               A;
2482:   PetscDS           prob, probAux = NULL;
2483:   PetscQuadrature   quad;
2484:   PetscSection      section, globalSection, sectionAux;
2485:   PetscScalar      *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
2486:   PetscInt          Nf, fieldI, fieldJ;
2487:   PetscInt          totDim, totDimAux = 0;
2488:   const PetscInt   *cells;
2489:   PetscInt          cStart, cEnd, numCells, c;
2490:   PetscBool         hasDyn;
2491:   DMField           coordField;
2492:   PetscErrorCode    ierr;

2495:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
2496:   DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
2497:   if (!cellIS) {
2498:     PetscInt depth;

2500:     DMPlexGetDepth(plex, &depth);
2501:     DMGetStratumIS(plex, "dim", depth, &cellIS);
2502:     if (!cellIS) {DMGetStratumIS(plex, "depth", depth, &cellIS);}
2503:   } else {
2504:     PetscObjectReference((PetscObject) cellIS);
2505:   }
2506:   DMGetSection(dm, &section);
2507:   DMGetGlobalSection(dm, &globalSection);
2508:   DMGetDS(dm, &prob);
2509:   PetscDSGetTotalDimension(prob, &totDim);
2510:   PetscDSHasDynamicJacobian(prob, &hasDyn);
2511:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
2512:   PetscSectionGetNumFields(section, &Nf);
2513:   ISGetLocalSize(cellIS, &numCells);
2514:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
2515:   PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);
2516:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);
2517:   if (dmAux) {
2518:     DMConvert(dmAux, DMPLEX, &plexAux);
2519:     DMGetSection(plexAux, &sectionAux);
2520:     DMGetDS(dmAux, &probAux);
2521:     PetscDSGetTotalDimension(probAux, &totDimAux);
2522:   }
2523:   VecSet(Z, 0.0);
2524:   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);
2525:   if (dmAux) {PetscMalloc1(numCells*totDimAux, &a);}
2526:   DMGetCoordinateField(dm, &coordField);
2527:   for (c = cStart; c < cEnd; ++c) {
2528:     const PetscInt cell = cells ? cells[c] : c;
2529:     const PetscInt cind = c - cStart;
2530:     PetscScalar   *x = NULL,  *x_t = NULL;
2531:     PetscInt       i;

2533:     DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
2534:     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
2535:     DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
2536:     if (X_t) {
2537:       DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
2538:       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
2539:       DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
2540:     }
2541:     if (dmAux) {
2542:       DMPlexVecGetClosure(plexAux, sectionAux, A, cell, NULL, &x);
2543:       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
2544:       DMPlexVecRestoreClosure(plexAux, sectionAux, A, cell, NULL, &x);
2545:     }
2546:     DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);
2547:     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
2548:     DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);
2549:   }
2550:   PetscMemzero(elemMat, numCells*totDim*totDim * sizeof(PetscScalar));
2551:   if (hasDyn)  {PetscMemzero(elemMatD, numCells*totDim*totDim * sizeof(PetscScalar));}
2552:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2553:     PetscFE  fe;
2554:     PetscInt Nb;
2555:     /* Conforming batches */
2556:     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2557:     /* Remainder */
2558:     PetscInt Nr, offset, Nq;
2559:     PetscQuadrature qGeom = NULL;
2560:     PetscInt    maxDegree;
2561:     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

2563:     PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
2564:     PetscFEGetQuadrature(fe, &quad);
2565:     PetscFEGetDimension(fe, &Nb);
2566:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
2567:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
2568:     if (maxDegree <= 1) {DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);}
2569:     if (!qGeom) {
2570:       PetscFEGetQuadrature(fe,&qGeom);
2571:       PetscObjectReference((PetscObject)qGeom);
2572:     }
2573:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
2574:     DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
2575:     blockSize = Nb;
2576:     batchSize = numBlocks * blockSize;
2577:     PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
2578:     numChunks = numCells / (numBatches*batchSize);
2579:     Ne        = numChunks*numBatches*batchSize;
2580:     Nr        = numCells % (numBatches*batchSize);
2581:     offset    = numCells - Nr;
2582:     PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
2583:     PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);
2584:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2585:       PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
2586:       PetscFEIntegrateJacobian(fe, 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]);
2587:       if (hasDyn) {
2588:         PetscFEIntegrateJacobian(fe, prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);
2589:         PetscFEIntegrateJacobian(fe, 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]);
2590:       }
2591:     }
2592:     PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);
2593:     PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);
2594:     DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
2595:     PetscQuadratureDestroy(&qGeom);
2596:   }
2597:   if (hasDyn) {
2598:     for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
2599:   }
2600:   for (c = cStart; c < cEnd; ++c) {
2601:     const PetscInt     cell = cells ? cells[c] : c;
2602:     const PetscInt     cind = c - cStart;
2603:     const PetscBLASInt M = totDim, one = 1;
2604:     const PetscScalar  a = 1.0, b = 0.0;

2606:     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
2607:     if (mesh->printFEM > 1) {
2608:       DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
2609:       DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);
2610:       DMPrintCellVector(c, "Z",  totDim, z);
2611:     }
2612:     DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);
2613:   }
2614:   PetscFree6(u,u_t,elemMat,elemMatD,y,z);
2615:   if (mesh->printFEM) {
2616:     PetscPrintf(PETSC_COMM_WORLD, "Z:\n");
2617:     VecView(Z, PETSC_VIEWER_STDOUT_WORLD);
2618:   }
2619:   PetscFree(a);
2620:   ISDestroy(&cellIS);
2621:   DMDestroy(&plexAux);
2622:   DMDestroy(&plex);
2623:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
2624:   return(0);
2625: }

2627: /*@
2628:   DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.

2630:   Input Parameters:
2631: + dm - The mesh
2632: . X  - Local input vector
2633: - user - The user context

2635:   Output Parameter:
2636: . Jac  - Jacobian matrix

2638:   Note:
2639:   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
2640:   like a GPU, or vectorize on a multicore machine.

2642:   Level: developer

2644: .seealso: FormFunctionLocal()
2645: @*/
2646: PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
2647: {
2648:   DM             plex;
2649:   PetscDS        prob;
2650:   IS             cellIS;
2651:   PetscBool      hasJac, hasPrec;
2652:   PetscInt       depth;

2656:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
2657:   DMPlexGetDepth(plex, &depth);
2658:   DMGetStratumIS(plex, "dim", depth, &cellIS);
2659:   if (!cellIS) {DMGetStratumIS(plex, "depth", depth, &cellIS);}
2660:   DMGetDS(dm, &prob);
2661:   PetscDSHasJacobian(prob, &hasJac);
2662:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
2663:   if (hasJac && hasPrec) {MatZeroEntries(Jac);}
2664:   MatZeroEntries(JacP);
2665:   DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);
2666:   ISDestroy(&cellIS);
2667:   DMDestroy(&plex);
2668:   return(0);
2669: }

2671: /*@
2672:   DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.

2674:   Input Parameters:
2675: + dm - The DM object
2676: . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
2677: . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
2678: - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())

2680:   Level: developer
2681: @*/
2682: PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
2683: {

2687:   DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);
2688:   DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);
2689:   DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);
2690:   return(0);
2691: }

2693: PetscErrorCode DMSNESCheckFromOptions_Internal(SNES snes, DM dm, Vec u, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs)
2694: {
2695:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2696:   PetscDS        prob;
2697:   Mat            J, M;
2698:   Vec            r, b;
2699:   MatNullSpace   nullSpace;
2700:   PetscReal     *error, res = 0.0;
2701:   PetscInt       numFields;
2702:   PetscBool      hasJac, hasPrec;
2703:   PetscInt       Nf, f;

2707:   DMGetDS(dm, &prob);
2708:   PetscDSGetNumFields(prob, &Nf);
2709:   PetscMalloc1(Nf, &exacts);
2710:   for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(prob, f, &exacts[f]);}
2711:   VecDuplicate(u, &r);
2712:   DMCreateMatrix(dm, &J);
2713:   /* TODO Null space for J */
2714:   /* Check discretization error */
2715:   DMGetNumFields(dm, &numFields);
2716:   PetscMalloc1(PetscMax(1, numFields), &error);
2717:   DMProjectFunction(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs, INSERT_ALL_VALUES, u);
2718:   if (numFields > 1) {
2719:     PetscInt f;

2721:     DMComputeL2FieldDiff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs, u, error);
2722:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: [");
2723:     for (f = 0; f < numFields; ++f) {
2724:       if (f) {PetscPrintf(PETSC_COMM_WORLD, ", ");}
2725:       if (error[f] >= 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "%g", (double)error[f]);}
2726:       else                     {PetscPrintf(PETSC_COMM_WORLD, "< 1.0e-11");}
2727:     }
2728:     PetscPrintf(PETSC_COMM_WORLD, "]\n");
2729:   } else {
2730:     DMComputeL2Diff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs, u, &error[0]);
2731:     if (error[0] >= 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)error[0]);}
2732:     else                     {PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");}
2733:   }
2734:   PetscFree(error);
2735:   /* Check residual */
2736:   SNESComputeFunction(snes, u, r);
2737:   VecNorm(r, NORM_2, &res);
2738:   PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);
2739:   VecChop(r, 1.0e-10);
2740:   PetscObjectSetName((PetscObject) r, "Initial Residual");
2741:   PetscObjectSetOptionsPrefix((PetscObject)r,"res_");
2742:   VecViewFromOptions(r, NULL, "-vec_view");
2743:   /* Check Jacobian */
2744:   PetscDSHasJacobian(prob, &hasJac);
2745:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
2746:   if (hasJac && hasPrec) {
2747:     DMCreateMatrix(dm, &M);
2748:     SNESComputeJacobian(snes, u, J, M);
2749:     PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");
2750:     MatViewFromOptions(M, NULL, "-mat_view");
2751:     MatDestroy(&M);
2752:   } else {
2753:     SNESComputeJacobian(snes, u, J, J);
2754:   }
2755:   PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");
2756:   MatViewFromOptions(J, NULL, "-mat_view");
2757:   MatGetNullSpace(J, &nullSpace);
2758:   if (nullSpace) {
2759:     PetscBool isNull;
2760:     MatNullSpaceTest(nullSpace, J, &isNull);
2761:     if (!isNull) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
2762:   }
2763:   VecDuplicate(u, &b);
2764:   VecSet(r, 0.0);
2765:   SNESComputeFunction(snes, r, b);
2766:   MatMult(J, u, r);
2767:   VecAXPY(r, 1.0, b);
2768:   VecDestroy(&b);
2769:   VecNorm(r, NORM_2, &res);
2770:   PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);
2771:   VecChop(r, 1.0e-10);
2772:   PetscObjectSetName((PetscObject) r, "Au - b = Au + F(0)");
2773:   PetscObjectSetOptionsPrefix((PetscObject)r,"linear_res_");
2774:   VecViewFromOptions(r, NULL, "-vec_view");
2775:   VecDestroy(&r);
2776:   MatNullSpaceDestroy(&nullSpace);
2777:   MatDestroy(&J);
2778:   PetscFree(exacts);
2779:   return(0);
2780: }

2782: /*@C
2783:   DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information

2785:   Input Parameters:
2786: + snes - the SNES object
2787: . u    - representative SNES vector
2788: . exactFuncs - pointwise functions of the exact solution for each field
2789: - ctxs - contexts for the functions

2791:   Level: developer
2792: @*/
2793: PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2794: {
2795:   PetscErrorCode (**exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = NULL;
2796:   DM             dm;
2797:   PetscDS        prob;
2798:   Vec            sol;
2799:   PetscBool      check;
2800:   PetscInt       Nf, f;

2804:   PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);
2805:   if (!check) return(0);
2806:   SNESGetDM(snes, &dm);
2807:   DMGetDS(dm, &prob);
2808:   if (!exactFuncs) {
2809:     PetscDSGetNumFields(prob, &Nf);
2810:     PetscMalloc1(Nf, &exact);
2811:     for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(prob, f, &exact[f]);}
2812:   }
2813:   VecDuplicate(u, &sol);
2814:   SNESSetSolution(snes, sol);
2815:   DMSNESCheckFromOptions_Internal(snes, dm, sol, exactFuncs ? exactFuncs : exact, ctxs);
2816:   VecDestroy(&sol);
2817:   PetscFree(exact);
2818:   return(0);
2819: }