Actual source code: dmplexsnes.c

petsc-3.12.5 2020-03-29
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: /*@C
 44:   DMInterpolationCreate - Creates a DMInterpolationInfo context

 46:   Collective

 48:   Input Parameter:
 49: . comm - the communicator

 51:   Output Parameter:
 52: . ctx - the context

 54:   Level: beginner

 56: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
 57: @*/
 58: PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
 59: {

 64:   PetscNew(ctx);

 66:   (*ctx)->comm   = comm;
 67:   (*ctx)->dim    = -1;
 68:   (*ctx)->nInput = 0;
 69:   (*ctx)->points = NULL;
 70:   (*ctx)->cells  = NULL;
 71:   (*ctx)->n      = -1;
 72:   (*ctx)->coords = NULL;
 73:   return(0);
 74: }

 76: /*@C
 77:   DMInterpolationSetDim - Sets the spatial dimension for the interpolation context

 79:   Not collective

 81:   Input Parameters:
 82: + ctx - the context
 83: - dim - the spatial dimension

 85:   Level: intermediate

 87: .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
 88: @*/
 89: PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
 90: {
 92:   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim);
 93:   ctx->dim = dim;
 94:   return(0);
 95: }

 97: /*@C
 98:   DMInterpolationGetDim - Gets the spatial dimension for the interpolation context

100:   Not collective

102:   Input Parameter:
103: . ctx - the context

105:   Output Parameter:
106: . dim - the spatial dimension

108:   Level: intermediate

110: .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
111: @*/
112: PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
113: {
116:   *dim = ctx->dim;
117:   return(0);
118: }

120: /*@C
121:   DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context

123:   Not collective

125:   Input Parameters:
126: + ctx - the context
127: - dof - the number of fields

129:   Level: intermediate

131: .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
132: @*/
133: PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
134: {
136:   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof);
137:   ctx->dof = dof;
138:   return(0);
139: }

141: /*@C
142:   DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context

144:   Not collective

146:   Input Parameter:
147: . ctx - the context

149:   Output Parameter:
150: . dof - the number of fields

152:   Level: intermediate

154: .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
155: @*/
156: PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
157: {
160:   *dof = ctx->dof;
161:   return(0);
162: }

164: /*@C
165:   DMInterpolationAddPoints - Add points at which we will interpolate the fields

167:   Not collective

169:   Input Parameters:
170: + ctx    - the context
171: . n      - the number of points
172: - points - the coordinates for each point, an array of size n * dim

174:   Note: The coordinate information is copied.

176:   Level: intermediate

178: .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
179: @*/
180: PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
181: {

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

189:   PetscMalloc1(n*ctx->dim, &ctx->points);
190:   PetscArraycpy(ctx->points, points, n*ctx->dim);
191:   return(0);
192: }

194: /*@C
195:   DMInterpolationSetUp - Computea spatial indices that add in point location during interpolation

197:   Collective on ctx

199:   Input Parameters:
200: + ctx - the context
201: . dm  - the DM for the function space used for interpolation
202: - redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.

204:   Level: intermediate

206: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
207: @*/
208: PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
209: {
210:   MPI_Comm          comm = ctx->comm;
211:   PetscScalar       *a;
212:   PetscInt          p, q, i;
213:   PetscMPIInt       rank, size;
214:   PetscErrorCode    ierr;
215:   Vec               pointVec;
216:   PetscSF           cellSF;
217:   PetscLayout       layout;
218:   PetscReal         *globalPoints;
219:   PetscScalar       *globalPointsScalar;
220:   const PetscInt    *ranges;
221:   PetscMPIInt       *counts, *displs;
222:   const PetscSFNode *foundCells;
223:   const PetscInt    *foundPoints;
224:   PetscMPIInt       *foundProcs, *globalProcs;
225:   PetscInt          n, N, numFound;

229:   MPI_Comm_size(comm, &size);
230:   MPI_Comm_rank(comm, &rank);
231:   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
232:   /* Locate points */
233:   n = ctx->nInput;
234:   if (!redundantPoints) {
235:     PetscLayoutCreate(comm, &layout);
236:     PetscLayoutSetBlockSize(layout, 1);
237:     PetscLayoutSetLocalSize(layout, n);
238:     PetscLayoutSetUp(layout);
239:     PetscLayoutGetSize(layout, &N);
240:     /* Communicate all points to all processes */
241:     PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);
242:     PetscLayoutGetRanges(layout, &ranges);
243:     for (p = 0; p < size; ++p) {
244:       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
245:       displs[p] = ranges[p]*ctx->dim;
246:     }
247:     MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);
248:   } else {
249:     N = n;
250:     globalPoints = ctx->points;
251:     counts = displs = NULL;
252:     layout = NULL;
253:   }
254: #if 0
255:   PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);
256:   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
257: #else
258: #if defined(PETSC_USE_COMPLEX)
259:   PetscMalloc1(N*ctx->dim,&globalPointsScalar);
260:   for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
261: #else
262:   globalPointsScalar = globalPoints;
263: #endif
264:   VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);
265:   PetscMalloc2(N,&foundProcs,N,&globalProcs);
266:   for (p = 0; p < N; ++p) {foundProcs[p] = size;}
267:   cellSF = NULL;
268:   DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);
269:   PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);
270: #endif
271:   for (p = 0; p < numFound; ++p) {
272:     if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
273:   }
274:   /* Let the lowest rank process own each point */
275:   MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);
276:   ctx->n = 0;
277:   for (p = 0; p < N; ++p) {
278:     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));
279:     else if (globalProcs[p] == rank) ctx->n++;
280:   }
281:   /* Create coordinates vector and array of owned cells */
282:   PetscMalloc1(ctx->n, &ctx->cells);
283:   VecCreate(comm, &ctx->coords);
284:   VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);
285:   VecSetBlockSize(ctx->coords, ctx->dim);
286:   VecSetType(ctx->coords,VECSTANDARD);
287:   VecGetArray(ctx->coords, &a);
288:   for (p = 0, q = 0, i = 0; p < N; ++p) {
289:     if (globalProcs[p] == rank) {
290:       PetscInt d;

292:       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
293:       ctx->cells[q] = foundCells[q].index;
294:       ++q;
295:     }
296:   }
297:   VecRestoreArray(ctx->coords, &a);
298: #if 0
299:   PetscFree3(foundCells,foundProcs,globalProcs);
300: #else
301:   PetscFree2(foundProcs,globalProcs);
302:   PetscSFDestroy(&cellSF);
303:   VecDestroy(&pointVec);
304: #endif
305:   if ((void*)globalPointsScalar != (void*)globalPoints) {PetscFree(globalPointsScalar);}
306:   if (!redundantPoints) {PetscFree3(globalPoints,counts,displs);}
307:   PetscLayoutDestroy(&layout);
308:   return(0);
309: }

311: /*@C
312:   DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point

314:   Collective on ctx

316:   Input Parameter:
317: . ctx - the context

319:   Output Parameter:
320: . coordinates  - the coordinates of interpolation points

322:   Note: The local vector entries correspond to interpolation points lying on this process, according to the associated DM. This is a borrowed vector that the user should not destroy.

324:   Level: intermediate

326: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
327: @*/
328: PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
329: {
332:   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
333:   *coordinates = ctx->coords;
334:   return(0);
335: }

337: /*@C
338:   DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values

340:   Collective on ctx

342:   Input Parameter:
343: . ctx - the context

345:   Output Parameter:
346: . v  - a vector capable of holding the interpolated field values

348:   Note: This vector should be returned using DMInterpolationRestoreVector().

350:   Level: intermediate

352: .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
353: @*/
354: PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
355: {

360:   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
361:   VecCreate(ctx->comm, v);
362:   VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);
363:   VecSetBlockSize(*v, ctx->dof);
364:   VecSetType(*v,VECSTANDARD);
365:   return(0);
366: }

368: /*@C
369:   DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values

371:   Collective on ctx

373:   Input Parameters:
374: + ctx - the context
375: - v  - a vector capable of holding the interpolated field values

377:   Level: intermediate

379: .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
380: @*/
381: PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
382: {

387:   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
388:   VecDestroy(v);
389:   return(0);
390: }

392: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
393: {
394:   PetscReal      *v0, *J, *invJ, detJ;
395:   const PetscScalar *coords;
396:   PetscScalar    *a;
397:   PetscInt       p;

401:   PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);
402:   VecGetArrayRead(ctx->coords, &coords);
403:   VecGetArray(v, &a);
404:   for (p = 0; p < ctx->n; ++p) {
405:     PetscInt     c = ctx->cells[p];
406:     PetscScalar *x = NULL;
407:     PetscReal    xi[4];
408:     PetscInt     d, f, comp;

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

415:     for (d = 0; d < ctx->dim; ++d) {
416:       xi[d] = 0.0;
417:       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
418:       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];
419:     }
420:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);
421:   }
422:   VecRestoreArray(v, &a);
423:   VecRestoreArrayRead(ctx->coords, &coords);
424:   PetscFree3(v0, J, invJ);
425:   return(0);
426: }

428: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
429: {
430:   PetscReal      *v0, *J, *invJ, detJ;
431:   const PetscScalar *coords;
432:   PetscScalar    *a;
433:   PetscInt       p;

437:   PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);
438:   VecGetArrayRead(ctx->coords, &coords);
439:   VecGetArray(v, &a);
440:   for (p = 0; p < ctx->n; ++p) {
441:     PetscInt       c = ctx->cells[p];
442:     const PetscInt order[3] = {2, 1, 3};
443:     PetscScalar   *x = NULL;
444:     PetscReal      xi[4];
445:     PetscInt       d, f, comp;

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

452:     for (d = 0; d < ctx->dim; ++d) {
453:       xi[d] = 0.0;
454:       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
455:       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];
456:     }
457:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);
458:   }
459:   VecRestoreArray(v, &a);
460:   VecRestoreArrayRead(ctx->coords, &coords);
461:   PetscFree3(v0, J, invJ);
462:   return(0);
463: }

465: PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
466: {
467:   const PetscScalar *vertices = (const PetscScalar*) ctx;
468:   const PetscScalar x0        = vertices[0];
469:   const PetscScalar y0        = vertices[1];
470:   const PetscScalar x1        = vertices[2];
471:   const PetscScalar y1        = vertices[3];
472:   const PetscScalar x2        = vertices[4];
473:   const PetscScalar y2        = vertices[5];
474:   const PetscScalar x3        = vertices[6];
475:   const PetscScalar y3        = vertices[7];
476:   const PetscScalar f_1       = x1 - x0;
477:   const PetscScalar g_1       = y1 - y0;
478:   const PetscScalar f_3       = x3 - x0;
479:   const PetscScalar g_3       = y3 - y0;
480:   const PetscScalar f_01      = x2 - x1 - x3 + x0;
481:   const PetscScalar g_01      = y2 - y1 - y3 + y0;
482:   const PetscScalar *ref;
483:   PetscScalar       *real;
484:   PetscErrorCode    ierr;

487:   VecGetArrayRead(Xref,  &ref);
488:   VecGetArray(Xreal, &real);
489:   {
490:     const PetscScalar p0 = ref[0];
491:     const PetscScalar p1 = ref[1];

493:     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
494:     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
495:   }
496:   PetscLogFlops(28);
497:   VecRestoreArrayRead(Xref,  &ref);
498:   VecRestoreArray(Xreal, &real);
499:   return(0);
500: }

502:  #include <petsc/private/dmimpl.h>
503: PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
504: {
505:   const PetscScalar *vertices = (const PetscScalar*) ctx;
506:   const PetscScalar x0        = vertices[0];
507:   const PetscScalar y0        = vertices[1];
508:   const PetscScalar x1        = vertices[2];
509:   const PetscScalar y1        = vertices[3];
510:   const PetscScalar x2        = vertices[4];
511:   const PetscScalar y2        = vertices[5];
512:   const PetscScalar x3        = vertices[6];
513:   const PetscScalar y3        = vertices[7];
514:   const PetscScalar f_01      = x2 - x1 - x3 + x0;
515:   const PetscScalar g_01      = y2 - y1 - y3 + y0;
516:   const PetscScalar *ref;
517:   PetscErrorCode    ierr;

520:   VecGetArrayRead(Xref,  &ref);
521:   {
522:     const PetscScalar x       = ref[0];
523:     const PetscScalar y       = ref[1];
524:     const PetscInt    rows[2] = {0, 1};
525:     PetscScalar       values[4];

527:     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
528:     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
529:     MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);
530:   }
531:   PetscLogFlops(30);
532:   VecRestoreArrayRead(Xref,  &ref);
533:   MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
534:   MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
535:   return(0);
536: }

538: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
539: {
540:   DM             dmCoord;
541:   PetscFE        fem = NULL;
542:   SNES           snes;
543:   KSP            ksp;
544:   PC             pc;
545:   Vec            coordsLocal, r, ref, real;
546:   Mat            J;
547:   const PetscScalar *coords;
548:   PetscScalar    *a;
549:   PetscInt       Nf, p;
550:   const PetscInt dof = ctx->dof;

554:   DMGetNumFields(dm, &Nf);
555:   if (Nf) {DMGetField(dm, 0, NULL, (PetscObject *) &fem);}
556:   DMGetCoordinatesLocal(dm, &coordsLocal);
557:   DMGetCoordinateDM(dm, &dmCoord);
558:   SNESCreate(PETSC_COMM_SELF, &snes);
559:   SNESSetOptionsPrefix(snes, "quad_interp_");
560:   VecCreate(PETSC_COMM_SELF, &r);
561:   VecSetSizes(r, 2, 2);
562:   VecSetType(r,dm->vectype);
563:   VecDuplicate(r, &ref);
564:   VecDuplicate(r, &real);
565:   MatCreate(PETSC_COMM_SELF, &J);
566:   MatSetSizes(J, 2, 2, 2, 2);
567:   MatSetType(J, MATSEQDENSE);
568:   MatSetUp(J);
569:   SNESSetFunction(snes, r, QuadMap_Private, NULL);
570:   SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);
571:   SNESGetKSP(snes, &ksp);
572:   KSPGetPC(ksp, &pc);
573:   PCSetType(pc, PCLU);
574:   SNESSetFromOptions(snes);

576:   VecGetArrayRead(ctx->coords, &coords);
577:   VecGetArray(v, &a);
578:   for (p = 0; p < ctx->n; ++p) {
579:     PetscScalar *x = NULL, *vertices = NULL;
580:     PetscScalar *xi;
581:     PetscReal    xir[2];
582:     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;

584:     /* Can make this do all points at once */
585:     DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
586:     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2);
587:     DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);
588:     SNESSetFunction(snes, NULL, NULL, (void*) vertices);
589:     SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);
590:     VecGetArray(real, &xi);
591:     xi[0]  = coords[p*ctx->dim+0];
592:     xi[1]  = coords[p*ctx->dim+1];
593:     VecRestoreArray(real, &xi);
594:     SNESSolve(snes, real, ref);
595:     VecGetArray(ref, &xi);
596:     xir[0] = PetscRealPart(xi[0]);
597:     xir[1] = PetscRealPart(xi[1]);
598:     if (4*dof != xSize) {
599:       PetscReal *B;
600:       PetscInt   d;

602:       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
603:       PetscFEGetTabulation(fem, 1, xir, &B, NULL, NULL);
604:       for (comp = 0; comp < dof; ++comp) {
605:         a[p*dof+comp] = 0.0;
606:         for (d = 0; d < xSize/dof; ++d) {
607:           a[p*dof+comp] += x[d*dof+comp]*B[d*dof+comp];
608:         }
609:       }
610:       PetscFERestoreTabulation(fem, 1, xir, &B, NULL, NULL);
611:     } else {
612:       for (comp = 0; comp < dof; ++comp)
613:         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];
614:     }
615:     VecRestoreArray(ref, &xi);
616:     DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
617:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);
618:   }
619:   VecRestoreArray(v, &a);
620:   VecRestoreArrayRead(ctx->coords, &coords);

622:   SNESDestroy(&snes);
623:   VecDestroy(&r);
624:   VecDestroy(&ref);
625:   VecDestroy(&real);
626:   MatDestroy(&J);
627:   return(0);
628: }

630: PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
631: {
632:   const PetscScalar *vertices = (const PetscScalar*) ctx;
633:   const PetscScalar x0        = vertices[0];
634:   const PetscScalar y0        = vertices[1];
635:   const PetscScalar z0        = vertices[2];
636:   const PetscScalar x1        = vertices[9];
637:   const PetscScalar y1        = vertices[10];
638:   const PetscScalar z1        = vertices[11];
639:   const PetscScalar x2        = vertices[6];
640:   const PetscScalar y2        = vertices[7];
641:   const PetscScalar z2        = vertices[8];
642:   const PetscScalar x3        = vertices[3];
643:   const PetscScalar y3        = vertices[4];
644:   const PetscScalar z3        = vertices[5];
645:   const PetscScalar x4        = vertices[12];
646:   const PetscScalar y4        = vertices[13];
647:   const PetscScalar z4        = vertices[14];
648:   const PetscScalar x5        = vertices[15];
649:   const PetscScalar y5        = vertices[16];
650:   const PetscScalar z5        = vertices[17];
651:   const PetscScalar x6        = vertices[18];
652:   const PetscScalar y6        = vertices[19];
653:   const PetscScalar z6        = vertices[20];
654:   const PetscScalar x7        = vertices[21];
655:   const PetscScalar y7        = vertices[22];
656:   const PetscScalar z7        = vertices[23];
657:   const PetscScalar f_1       = x1 - x0;
658:   const PetscScalar g_1       = y1 - y0;
659:   const PetscScalar h_1       = z1 - z0;
660:   const PetscScalar f_3       = x3 - x0;
661:   const PetscScalar g_3       = y3 - y0;
662:   const PetscScalar h_3       = z3 - z0;
663:   const PetscScalar f_4       = x4 - x0;
664:   const PetscScalar g_4       = y4 - y0;
665:   const PetscScalar h_4       = z4 - z0;
666:   const PetscScalar f_01      = x2 - x1 - x3 + x0;
667:   const PetscScalar g_01      = y2 - y1 - y3 + y0;
668:   const PetscScalar h_01      = z2 - z1 - z3 + z0;
669:   const PetscScalar f_12      = x7 - x3 - x4 + x0;
670:   const PetscScalar g_12      = y7 - y3 - y4 + y0;
671:   const PetscScalar h_12      = z7 - z3 - z4 + z0;
672:   const PetscScalar f_02      = x5 - x1 - x4 + x0;
673:   const PetscScalar g_02      = y5 - y1 - y4 + y0;
674:   const PetscScalar h_02      = z5 - z1 - z4 + z0;
675:   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
676:   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
677:   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
678:   const PetscScalar *ref;
679:   PetscScalar       *real;
680:   PetscErrorCode    ierr;

683:   VecGetArrayRead(Xref,  &ref);
684:   VecGetArray(Xreal, &real);
685:   {
686:     const PetscScalar p0 = ref[0];
687:     const PetscScalar p1 = ref[1];
688:     const PetscScalar p2 = ref[2];

690:     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;
691:     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;
692:     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;
693:   }
694:   PetscLogFlops(114);
695:   VecRestoreArrayRead(Xref,  &ref);
696:   VecRestoreArray(Xreal, &real);
697:   return(0);
698: }

700: PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
701: {
702:   const PetscScalar *vertices = (const PetscScalar*) ctx;
703:   const PetscScalar x0        = vertices[0];
704:   const PetscScalar y0        = vertices[1];
705:   const PetscScalar z0        = vertices[2];
706:   const PetscScalar x1        = vertices[9];
707:   const PetscScalar y1        = vertices[10];
708:   const PetscScalar z1        = vertices[11];
709:   const PetscScalar x2        = vertices[6];
710:   const PetscScalar y2        = vertices[7];
711:   const PetscScalar z2        = vertices[8];
712:   const PetscScalar x3        = vertices[3];
713:   const PetscScalar y3        = vertices[4];
714:   const PetscScalar z3        = vertices[5];
715:   const PetscScalar x4        = vertices[12];
716:   const PetscScalar y4        = vertices[13];
717:   const PetscScalar z4        = vertices[14];
718:   const PetscScalar x5        = vertices[15];
719:   const PetscScalar y5        = vertices[16];
720:   const PetscScalar z5        = vertices[17];
721:   const PetscScalar x6        = vertices[18];
722:   const PetscScalar y6        = vertices[19];
723:   const PetscScalar z6        = vertices[20];
724:   const PetscScalar x7        = vertices[21];
725:   const PetscScalar y7        = vertices[22];
726:   const PetscScalar z7        = vertices[23];
727:   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
728:   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
729:   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
730:   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
731:   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
732:   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
733:   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
734:   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
735:   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
736:   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
737:   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
738:   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
739:   const PetscScalar *ref;
740:   PetscErrorCode    ierr;

743:   VecGetArrayRead(Xref,  &ref);
744:   {
745:     const PetscScalar x       = ref[0];
746:     const PetscScalar y       = ref[1];
747:     const PetscScalar z       = ref[2];
748:     const PetscInt    rows[3] = {0, 1, 2};
749:     PetscScalar       values[9];

751:     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
752:     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
753:     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
754:     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
755:     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
756:     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
757:     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
758:     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
759:     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;

761:     MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);
762:   }
763:   PetscLogFlops(152);
764:   VecRestoreArrayRead(Xref,  &ref);
765:   MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
766:   MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
767:   return(0);
768: }

770: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
771: {
772:   DM             dmCoord;
773:   SNES           snes;
774:   KSP            ksp;
775:   PC             pc;
776:   Vec            coordsLocal, r, ref, real;
777:   Mat            J;
778:   const PetscScalar *coords;
779:   PetscScalar    *a;
780:   PetscInt       p;

784:   DMGetCoordinatesLocal(dm, &coordsLocal);
785:   DMGetCoordinateDM(dm, &dmCoord);
786:   SNESCreate(PETSC_COMM_SELF, &snes);
787:   SNESSetOptionsPrefix(snes, "hex_interp_");
788:   VecCreate(PETSC_COMM_SELF, &r);
789:   VecSetSizes(r, 3, 3);
790:   VecSetType(r,dm->vectype);
791:   VecDuplicate(r, &ref);
792:   VecDuplicate(r, &real);
793:   MatCreate(PETSC_COMM_SELF, &J);
794:   MatSetSizes(J, 3, 3, 3, 3);
795:   MatSetType(J, MATSEQDENSE);
796:   MatSetUp(J);
797:   SNESSetFunction(snes, r, HexMap_Private, NULL);
798:   SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);
799:   SNESGetKSP(snes, &ksp);
800:   KSPGetPC(ksp, &pc);
801:   PCSetType(pc, PCLU);
802:   SNESSetFromOptions(snes);

804:   VecGetArrayRead(ctx->coords, &coords);
805:   VecGetArray(v, &a);
806:   for (p = 0; p < ctx->n; ++p) {
807:     PetscScalar *x = NULL, *vertices = NULL;
808:     PetscScalar *xi;
809:     PetscReal    xir[3];
810:     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;

812:     /* Can make this do all points at once */
813:     DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
814:     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3);
815:     DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);
816:     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof);
817:     SNESSetFunction(snes, NULL, NULL, (void*) vertices);
818:     SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);
819:     VecGetArray(real, &xi);
820:     xi[0]  = coords[p*ctx->dim+0];
821:     xi[1]  = coords[p*ctx->dim+1];
822:     xi[2]  = coords[p*ctx->dim+2];
823:     VecRestoreArray(real, &xi);
824:     SNESSolve(snes, real, ref);
825:     VecGetArray(ref, &xi);
826:     xir[0] = PetscRealPart(xi[0]);
827:     xir[1] = PetscRealPart(xi[1]);
828:     xir[2] = PetscRealPart(xi[2]);
829:     for (comp = 0; comp < ctx->dof; ++comp) {
830:       a[p*ctx->dof+comp] =
831:         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
832:         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
833:         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
834:         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
835:         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
836:         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
837:         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
838:         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
839:     }
840:     VecRestoreArray(ref, &xi);
841:     DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
842:     DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);
843:   }
844:   VecRestoreArray(v, &a);
845:   VecRestoreArrayRead(ctx->coords, &coords);

847:   SNESDestroy(&snes);
848:   VecDestroy(&r);
849:   VecDestroy(&ref);
850:   VecDestroy(&real);
851:   MatDestroy(&J);
852:   return(0);
853: }

855: /*@C
856:   DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.

858:   Input Parameters:
859: + ctx - The DMInterpolationInfo context
860: . dm  - The DM
861: - x   - The local vector containing the field to be interpolated

863:   Output Parameters:
864: . v   - The vector containing the interpolated values

866:   Note: A suitable v can be obtained using DMInterpolationGetVector().

868:   Level: beginner

870: .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
871: @*/
872: PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
873: {
874:   PetscInt       dim, coneSize, n;

881:   VecGetLocalSize(v, &n);
882:   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);
883:   if (n) {
884:     DMGetDimension(dm, &dim);
885:     DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);
886:     if (dim == 2) {
887:       if (coneSize == 3) {
888:         DMInterpolate_Triangle_Private(ctx, dm, x, v);
889:       } else if (coneSize == 4) {
890:         DMInterpolate_Quad_Private(ctx, dm, x, v);
891:       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %D for point interpolation", dim);
892:     } else if (dim == 3) {
893:       if (coneSize == 4) {
894:         DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);
895:       } else {
896:         DMInterpolate_Hex_Private(ctx, dm, x, v);
897:       }
898:     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %D for point interpolation", dim);
899:   }
900:   return(0);
901: }

903: /*@C
904:   DMInterpolationDestroy - Destroys a DMInterpolationInfo context

906:   Collective on ctx

908:   Input Parameter:
909: . ctx - the context

911:   Level: beginner

913: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
914: @*/
915: PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
916: {

921:   VecDestroy(&(*ctx)->coords);
922:   PetscFree((*ctx)->points);
923:   PetscFree((*ctx)->cells);
924:   PetscFree(*ctx);
925:   *ctx = NULL;
926:   return(0);
927: }

929: /*@C
930:   SNESMonitorFields - Monitors the residual for each field separately

932:   Collective on SNES

934:   Input Parameters:
935: + snes   - the SNES context
936: . its    - iteration number
937: . fgnorm - 2-norm of residual
938: - vf  - PetscViewerAndFormat of type ASCII

940:   Notes:
941:   This routine prints the residual norm at each iteration.

943:   Level: intermediate

945: .seealso: SNESMonitorSet(), SNESMonitorDefault()
946: @*/
947: PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
948: {
949:   PetscViewer        viewer = vf->viewer;
950:   Vec                res;
951:   DM                 dm;
952:   PetscSection       s;
953:   const PetscScalar *r;
954:   PetscReal         *lnorms, *norms;
955:   PetscInt           numFields, f, pStart, pEnd, p;
956:   PetscErrorCode     ierr;

960:   SNESGetFunction(snes, &res, 0, 0);
961:   SNESGetDM(snes, &dm);
962:   DMGetLocalSection(dm, &s);
963:   PetscSectionGetNumFields(s, &numFields);
964:   PetscSectionGetChart(s, &pStart, &pEnd);
965:   PetscCalloc2(numFields, &lnorms, numFields, &norms);
966:   VecGetArrayRead(res, &r);
967:   for (p = pStart; p < pEnd; ++p) {
968:     for (f = 0; f < numFields; ++f) {
969:       PetscInt fdof, foff, d;

971:       PetscSectionGetFieldDof(s, p, f, &fdof);
972:       PetscSectionGetFieldOffset(s, p, f, &foff);
973:       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
974:     }
975:   }
976:   VecRestoreArrayRead(res, &r);
977:   MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));
978:   PetscViewerPushFormat(viewer,vf->format);
979:   PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);
980:   PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);
981:   for (f = 0; f < numFields; ++f) {
982:     if (f > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
983:     PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));
984:   }
985:   PetscViewerASCIIPrintf(viewer, "]\n");
986:   PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);
987:   PetscViewerPopFormat(viewer);
988:   PetscFree2(lnorms, norms);
989:   return(0);
990: }

992: /********************* Residual Computation **************************/


995: /*@
996:   DMPlexSNESGetGeometryFVM - Return precomputed geometric data

998:   Input Parameter:
999: . dm - The DM

1001:   Output Parameters:
1002: + facegeom - The values precomputed from face geometry
1003: . cellgeom - The values precomputed from cell geometry
1004: - minRadius - The minimum radius over the mesh of an inscribed sphere in a cell

1006:   Level: developer

1008: .seealso: DMPlexTSSetRHSFunctionLocal()
1009: @*/
1010: PetscErrorCode DMPlexSNESGetGeometryFVM(DM dm, Vec *facegeom, Vec *cellgeom, PetscReal *minRadius)
1011: {
1012:   DM             plex;

1017:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1018:   DMPlexGetDataFVM(plex, NULL, cellgeom, facegeom, NULL);
1019:   if (minRadius) {DMPlexGetMinRadius(plex, minRadius);}
1020:   DMDestroy(&plex);
1021:   return(0);
1022: }

1024: /*@
1025:   DMPlexSNESGetGradientDM - Return gradient data layout

1027:   Input Parameters:
1028: + dm - The DM
1029: - fv - The PetscFV

1031:   Output Parameter:
1032: . dmGrad - The layout for gradient values

1034:   Level: developer

1036: .seealso: DMPlexSNESGetGeometryFVM()
1037: @*/
1038: PetscErrorCode DMPlexSNESGetGradientDM(DM dm, PetscFV fv, DM *dmGrad)
1039: {
1040:   DM             plex;
1041:   PetscBool      computeGradients;

1048:   PetscFVGetComputeGradients(fv, &computeGradients);
1049:   if (!computeGradients) {*dmGrad = NULL; return(0);}
1050:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1051:   DMPlexGetDataFVM(plex, fv, NULL, NULL, dmGrad);
1052:   DMDestroy(&plex);
1053:   return(0);
1054: }

1056: 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)
1057: {
1058:   DM_Plex         *mesh = (DM_Plex *) dm->data;
1059:   DM               plex = NULL, plexA = NULL;
1060:   PetscDS          prob, probAux = NULL;
1061:   PetscSection     section, sectionAux = NULL;
1062:   Vec              locA = NULL;
1063:   PetscScalar     *u = NULL, *u_t = NULL, *a = NULL, *elemVec = NULL;
1064:   PetscInt         v;
1065:   PetscInt         totDim, totDimAux = 0;
1066:   PetscErrorCode   ierr;

1069:   DMConvert(dm, DMPLEX, &plex);
1070:   DMGetLocalSection(dm, &section);
1071:   DMGetDS(dm, &prob);
1072:   PetscDSGetTotalDimension(prob, &totDim);
1073:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);
1074:   if (locA) {
1075:     DM dmAux;

1077:     VecGetDM(locA, &dmAux);
1078:     DMConvert(dmAux, DMPLEX, &plexA);
1079:     DMGetDS(plexA, &probAux);
1080:     PetscDSGetTotalDimension(probAux, &totDimAux);
1081:     DMGetLocalSection(plexA, &sectionAux);
1082:   }
1083:   for (v = 0; v < numValues; ++v) {
1084:     PetscFEGeom    *fgeom;
1085:     PetscInt        maxDegree;
1086:     PetscQuadrature qGeom = NULL;
1087:     IS              pointIS;
1088:     const PetscInt *points;
1089:     PetscInt        numFaces, face, Nq;

1091:     DMLabelGetStratumIS(label, values[v], &pointIS);
1092:     if (!pointIS) continue; /* No points with that id on this process */
1093:     {
1094:       IS isectIS;

1096:       /* TODO: Special cases of ISIntersect where it is quick to check a priori if one is a superset of the other */
1097:       ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);
1098:       ISDestroy(&pointIS);
1099:       pointIS = isectIS;
1100:     }
1101:     ISGetLocalSize(pointIS,&numFaces);
1102:     ISGetIndices(pointIS,&points);
1103:     PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim, &elemVec, locA ? numFaces*totDimAux : 0, &a);
1104:     DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);
1105:     if (maxDegree <= 1) {
1106:       DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);
1107:     }
1108:     if (!qGeom) {
1109:       PetscFE fe;

1111:       PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1112:       PetscFEGetFaceQuadrature(fe, &qGeom);
1113:       PetscObjectReference((PetscObject)qGeom);
1114:     }
1115:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
1116:     DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
1117:     for (face = 0; face < numFaces; ++face) {
1118:       const PetscInt point = points[face], *support, *cone;
1119:       PetscScalar   *x     = NULL;
1120:       PetscInt       i, coneSize, faceLoc;

1122:       DMPlexGetSupport(dm, point, &support);
1123:       DMPlexGetConeSize(dm, support[0], &coneSize);
1124:       DMPlexGetCone(dm, support[0], &cone);
1125:       for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break;
1126:       if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %D in cone of support[0] %D", point, support[0]);
1127:       fgeom->face[face][0] = faceLoc;
1128:       DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);
1129:       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
1130:       DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);
1131:       if (locX_t) {
1132:         DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);
1133:         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
1134:         DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);
1135:       }
1136:       if (locA) {
1137:         PetscInt subp;

1139:         DMPlexGetAuxiliaryPoint(plex, plexA, support[0], &subp);
1140:         DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);
1141:         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
1142:         DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);
1143:       }
1144:     }
1145:     PetscArrayzero(elemVec, numFaces*totDim);
1146:     {
1147:       PetscFE         fe;
1148:       PetscInt        Nb;
1149:       PetscFEGeom     *chunkGeom = NULL;
1150:       /* Conforming batches */
1151:       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1152:       /* Remainder */
1153:       PetscInt        Nr, offset;

1155:       PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);
1156:       PetscFEGetDimension(fe, &Nb);
1157:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1158:       /* TODO: documentation is unclear about what is going on with these numbers: how should Nb / Nq factor in ? */
1159:       blockSize = Nb;
1160:       batchSize = numBlocks * blockSize;
1161:        PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1162:       numChunks = numFaces / (numBatches*batchSize);
1163:       Ne        = numChunks*numBatches*batchSize;
1164:       Nr        = numFaces % (numBatches*batchSize);
1165:       offset    = numFaces - Nr;
1166:       PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);
1167:       PetscFEIntegrateBdResidual(prob, field, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);
1168:       PetscFEGeomRestoreChunk(fgeom, 0, offset, &chunkGeom);
1169:       PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);
1170:       PetscFEIntegrateBdResidual(prob, field, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, a ? &a[offset*totDimAux] : NULL, t, &elemVec[offset*totDim]);
1171:       PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);
1172:     }
1173:     for (face = 0; face < numFaces; ++face) {
1174:       const PetscInt point = points[face], *support;

1176:       if (mesh->printFEM > 1) {DMPrintCellVector(point, "BdResidual", totDim, &elemVec[face*totDim]);}
1177:       DMPlexGetSupport(plex, point, &support);
1178:       DMPlexVecSetClosure(plex, NULL, locF, support[0], &elemVec[face*totDim], ADD_ALL_VALUES);
1179:     }
1180:     DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
1181:     PetscQuadratureDestroy(&qGeom);
1182:     ISRestoreIndices(pointIS, &points);
1183:     ISDestroy(&pointIS);
1184:     PetscFree4(u, u_t, elemVec, a);
1185:   }
1186:   if (plex)  {DMDestroy(&plex);}
1187:   if (plexA) {DMDestroy(&plexA);}
1188:   return(0);
1189: }

1191: PetscErrorCode DMPlexComputeBdResidualSingle(DM dm, PetscReal t, DMLabel label, PetscInt numValues, const PetscInt values[], PetscInt field, Vec locX, Vec locX_t, Vec locF)
1192: {
1193:   DMField        coordField;
1194:   DMLabel        depthLabel;
1195:   IS             facetIS;
1196:   PetscInt       dim;

1200:   DMGetDimension(dm, &dim);
1201:   DMPlexGetDepthLabel(dm, &depthLabel);
1202:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
1203:   DMGetCoordinateField(dm, &coordField);
1204:   DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);
1205:   ISDestroy(&facetIS);
1206:   return(0);
1207: }

1209: PetscErrorCode DMPlexComputeBdResidual_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1210: {
1211:   PetscDS        prob;
1212:   PetscInt       numBd, bd;
1213:   DMField        coordField = NULL;
1214:   IS             facetIS    = NULL;
1215:   DMLabel        depthLabel;
1216:   PetscInt       dim;

1220:   DMGetDS(dm, &prob);
1221:   DMPlexGetDepthLabel(dm, &depthLabel);
1222:   DMGetDimension(dm, &dim);
1223:   DMLabelGetStratumIS(depthLabel,dim - 1,&facetIS);
1224:   PetscDSGetNumBoundary(prob, &numBd);
1225:   for (bd = 0; bd < numBd; ++bd) {
1226:     DMBoundaryConditionType type;
1227:     const char             *bdLabel;
1228:     DMLabel                 label;
1229:     const PetscInt         *values;
1230:     PetscInt                field, numValues;
1231:     PetscObject             obj;
1232:     PetscClassId            id;

1234:     PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &field, NULL, NULL, NULL, &numValues, &values, NULL);
1235:     PetscDSGetDiscretization(prob, field, &obj);
1236:     PetscObjectGetClassId(obj, &id);
1237:     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
1238:     if (!facetIS) {
1239:       DMLabel  depthLabel;
1240:       PetscInt dim;

1242:       DMPlexGetDepthLabel(dm, &depthLabel);
1243:       DMGetDimension(dm, &dim);
1244:       DMLabelGetStratumIS(depthLabel, dim - 1, &facetIS);
1245:     }
1246:     DMGetCoordinateField(dm, &coordField);
1247:     DMGetLabel(dm, bdLabel, &label);
1248:     DMPlexComputeBdResidual_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, locF, coordField, facetIS);
1249:   }
1250:   ISDestroy(&facetIS);
1251:   return(0);
1252: }

1254: PetscErrorCode DMPlexComputeResidual_Internal(DM dm, IS cellIS, PetscReal time, Vec locX, Vec locX_t, PetscReal t, Vec locF, void *user)
1255: {
1256:   DM_Plex         *mesh       = (DM_Plex *) dm->data;
1257:   const char      *name       = "Residual";
1258:   DM               dmAux      = NULL;
1259:   DM               dmGrad     = NULL;
1260:   DMLabel          ghostLabel = NULL;
1261:   PetscDS          prob       = NULL;
1262:   PetscDS          probAux    = NULL;
1263:   PetscSection     section    = NULL;
1264:   PetscBool        useFEM     = PETSC_FALSE;
1265:   PetscBool        useFVM     = PETSC_FALSE;
1266:   PetscBool        isImplicit = (locX_t || time == PETSC_MIN_REAL) ? PETSC_TRUE : PETSC_FALSE;
1267:   PetscFV          fvm        = NULL;
1268:   PetscFVCellGeom *cgeomFVM   = NULL;
1269:   PetscFVFaceGeom *fgeomFVM   = NULL;
1270:   DMField          coordField = NULL;
1271:   Vec              locA, cellGeometryFVM = NULL, faceGeometryFVM = NULL, grad, locGrad = NULL;
1272:   PetscScalar     *u = NULL, *u_t, *a, *uL, *uR;
1273:   IS               chunkIS;
1274:   const PetscInt  *cells;
1275:   PetscInt         cStart, cEnd, numCells;
1276:   PetscInt         Nf, f, totDim, totDimAux, numChunks, cellChunkSize, faceChunkSize, chunk, fStart, fEnd;
1277:   PetscInt         maxDegree = PETSC_MAX_INT;
1278:   PetscQuadrature  affineQuad = NULL, *quads = NULL;
1279:   PetscFEGeom     *affineGeom = NULL, **geoms = NULL;
1280:   PetscErrorCode   ierr;

1283:   PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);
1284:   /* TODO The places where we have to use isFE are probably the member functions for the PetscDisc class */
1285:   /* TODO The FVM geometry is over-manipulated. Make the precalc functions return exactly what we need */
1286:   /* FEM+FVM */
1287:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
1288:   DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
1289:   /* 1: Get sizes from dm and dmAux */
1290:   DMGetLocalSection(dm, &section);
1291:   DMGetLabel(dm, "ghost", &ghostLabel);
1292:   DMGetCellDS(dm, cStart, &prob);
1293:   PetscDSGetNumFields(prob, &Nf);
1294:   PetscDSGetTotalDimension(prob, &totDim);
1295:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);
1296:   if (locA) {
1297:     PetscInt subcell;
1298:     DMPlexGetAuxiliaryPoint(dm, dmAux, cStart, &subcell);
1299:     VecGetDM(locA, &dmAux);
1300:     DMGetCellDS(dmAux, subcell, &probAux);
1301:     PetscDSGetTotalDimension(probAux, &totDimAux);
1302:   }
1303:   /* 2: Get geometric data */
1304:   for (f = 0; f < Nf; ++f) {
1305:     PetscObject  obj;
1306:     PetscClassId id;
1307:     PetscBool    fimp;

1309:     PetscDSGetImplicit(prob, f, &fimp);
1310:     if (isImplicit != fimp) continue;
1311:     PetscDSGetDiscretization(prob, f, &obj);
1312:     PetscObjectGetClassId(obj, &id);
1313:     if (id == PETSCFE_CLASSID) {useFEM = PETSC_TRUE;}
1314:     if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;}
1315:   }
1316:   if (useFEM) {
1317:     DMGetCoordinateField(dm, &coordField);
1318:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
1319:     if (maxDegree <= 1) {
1320:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&affineQuad);
1321:       if (affineQuad) {
1322:         DMSNESGetFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
1323:       }
1324:     } else {
1325:       PetscCalloc2(Nf,&quads,Nf,&geoms);
1326:       for (f = 0; f < Nf; ++f) {
1327:         PetscObject  obj;
1328:         PetscClassId id;
1329:         PetscBool    fimp;

1331:         PetscDSGetImplicit(prob, f, &fimp);
1332:         if (isImplicit != fimp) continue;
1333:         PetscDSGetDiscretization(prob, f, &obj);
1334:         PetscObjectGetClassId(obj, &id);
1335:         if (id == PETSCFE_CLASSID) {
1336:           PetscFE fe = (PetscFE) obj;

1338:           PetscFEGetQuadrature(fe, &quads[f]);
1339:           PetscObjectReference((PetscObject)quads[f]);
1340:           DMSNESGetFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
1341:         }
1342:       }
1343:     }
1344:   }
1345:   if (useFVM) {
1346:     DMPlexSNESGetGeometryFVM(dm, &faceGeometryFVM, &cellGeometryFVM, NULL);
1347:     VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);
1348:     VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);
1349:     /* Reconstruct and limit cell gradients */
1350:     DMPlexSNESGetGradientDM(dm, fvm, &dmGrad);
1351:     if (dmGrad) {
1352:       DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);
1353:       DMGetGlobalVector(dmGrad, &grad);
1354:       DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, locX, grad);
1355:       /* Communicate gradient values */
1356:       DMGetLocalVector(dmGrad, &locGrad);
1357:       DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);
1358:       DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);
1359:       DMRestoreGlobalVector(dmGrad, &grad);
1360:     }
1361:     /* Handle non-essential (e.g. outflow) boundary values */
1362:     DMPlexInsertBoundaryValues(dm, PETSC_FALSE, locX, time, faceGeometryFVM, cellGeometryFVM, locGrad);
1363:   }
1364:   /* Loop over chunks */
1365:   if (useFEM) {ISCreate(PETSC_COMM_SELF, &chunkIS);}
1366:   numCells      = cEnd - cStart;
1367:   numChunks     = 1;
1368:   cellChunkSize = numCells/numChunks;
1369:   faceChunkSize = (fEnd - fStart)/numChunks;
1370:   numChunks     = PetscMin(1,numCells);
1371:   for (chunk = 0; chunk < numChunks; ++chunk) {
1372:     PetscScalar     *elemVec, *fluxL, *fluxR;
1373:     PetscReal       *vol;
1374:     PetscFVFaceGeom *fgeom;
1375:     PetscInt         cS = cStart+chunk*cellChunkSize, cE = PetscMin(cS+cellChunkSize, cEnd), numCells = cE - cS, c;
1376:     PetscInt         fS = fStart+chunk*faceChunkSize, fE = PetscMin(fS+faceChunkSize, fEnd), numFaces = 0, face;

1378:     /* Extract field coefficients */
1379:     if (useFEM) {
1380:       ISGetPointSubrange(chunkIS, cS, cE, cells);
1381:       DMPlexGetCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
1382:       DMGetWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
1383:       PetscArrayzero(elemVec, numCells*totDim);
1384:     }
1385:     if (useFVM) {
1386:       DMPlexGetFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);
1387:       DMPlexGetFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);
1388:       DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);
1389:       DMGetWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);
1390:       PetscArrayzero(fluxL, numFaces*totDim);
1391:       PetscArrayzero(fluxR, numFaces*totDim);
1392:     }
1393:     /* TODO We will interlace both our field coefficients (u, u_t, uL, uR, etc.) and our output (elemVec, fL, fR). I think this works */
1394:     /* Loop over fields */
1395:     for (f = 0; f < Nf; ++f) {
1396:       PetscObject  obj;
1397:       PetscClassId id;
1398:       PetscBool    fimp;
1399:       PetscInt     numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset;

1401:       PetscDSGetImplicit(prob, f, &fimp);
1402:       if (isImplicit != fimp) continue;
1403:       PetscDSGetDiscretization(prob, f, &obj);
1404:       PetscObjectGetClassId(obj, &id);
1405:       if (id == PETSCFE_CLASSID) {
1406:         PetscFE         fe = (PetscFE) obj;
1407:         PetscFEGeom    *geom = affineGeom ? affineGeom : geoms[f];
1408:         PetscFEGeom    *chunkGeom = NULL;
1409:         PetscQuadrature quad = affineQuad ? affineQuad : quads[f];
1410:         PetscInt        Nq, Nb;

1412:         PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1413:         PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);
1414:         PetscFEGetDimension(fe, &Nb);
1415:         blockSize = Nb;
1416:         batchSize = numBlocks * blockSize;
1417:         PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1418:         numChunks = numCells / (numBatches*batchSize);
1419:         Ne        = numChunks*numBatches*batchSize;
1420:         Nr        = numCells % (numBatches*batchSize);
1421:         offset    = numCells - Nr;
1422:         /* Integrate FE residual to get elemVec (need fields at quadrature points) */
1423:         /*   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) */
1424:         PetscFEGeomGetChunk(geom,0,offset,&chunkGeom);
1425:         PetscFEIntegrateResidual(prob, f, Ne, chunkGeom, u, u_t, probAux, a, t, elemVec);
1426:         PetscFEGeomGetChunk(geom,offset,numCells,&chunkGeom);
1427:         PetscFEIntegrateResidual(prob, f, Nr, chunkGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, &elemVec[offset*totDim]);
1428:         PetscFEGeomRestoreChunk(geom,offset,numCells,&chunkGeom);
1429:       } else if (id == PETSCFV_CLASSID) {
1430:         PetscFV fv = (PetscFV) obj;

1432:         Ne = numFaces;
1433:         /* Riemann solve over faces (need fields at face centroids) */
1434:         /*   We need to evaluate FE fields at those coordinates */
1435:         PetscFVIntegrateRHSFunction(fv, prob, f, Ne, fgeom, vol, uL, uR, fluxL, fluxR);
1436:       } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %D", f);
1437:     }
1438:     /* Loop over domain */
1439:     if (useFEM) {
1440:       /* Add elemVec to locX */
1441:       for (c = cS; c < cE; ++c) {
1442:         const PetscInt cell = cells ? cells[c] : c;
1443:         const PetscInt cind = c - cStart;

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

1449:           DMLabelGetValue(ghostLabel,cell,&ghostVal);
1450:           if (ghostVal > 0) continue;
1451:         }
1452:         DMPlexVecSetClosure(dm, section, locF, cell, &elemVec[cind*totDim], ADD_ALL_VALUES);
1453:       }
1454:     }
1455:     if (useFVM) {
1456:       PetscScalar *fa;
1457:       PetscInt     iface;

1459:       VecGetArray(locF, &fa);
1460:       for (f = 0; f < Nf; ++f) {
1461:         PetscFV      fv;
1462:         PetscObject  obj;
1463:         PetscClassId id;
1464:         PetscInt     foff, pdim;

1466:         PetscDSGetDiscretization(prob, f, &obj);
1467:         PetscDSGetFieldOffset(prob, f, &foff);
1468:         PetscObjectGetClassId(obj, &id);
1469:         if (id != PETSCFV_CLASSID) continue;
1470:         fv   = (PetscFV) obj;
1471:         PetscFVGetNumComponents(fv, &pdim);
1472:         /* Accumulate fluxes to cells */
1473:         for (face = fS, iface = 0; face < fE; ++face) {
1474:           const PetscInt *scells;
1475:           PetscScalar    *fL = NULL, *fR = NULL;
1476:           PetscInt        ghost, d, nsupp, nchild;

1478:           DMLabelGetValue(ghostLabel, face, &ghost);
1479:           DMPlexGetSupportSize(dm, face, &nsupp);
1480:           DMPlexGetTreeChildren(dm, face, &nchild, NULL);
1481:           if (ghost >= 0 || nsupp > 2 || nchild > 0) continue;
1482:           DMPlexGetSupport(dm, face, &scells);
1483:           DMLabelGetValue(ghostLabel,scells[0],&ghost);
1484:           if (ghost <= 0) {DMPlexPointLocalFieldRef(dm, scells[0], f, fa, &fL);}
1485:           DMLabelGetValue(ghostLabel,scells[1],&ghost);
1486:           if (ghost <= 0) {DMPlexPointLocalFieldRef(dm, scells[1], f, fa, &fR);}
1487:           for (d = 0; d < pdim; ++d) {
1488:             if (fL) fL[d] -= fluxL[iface*totDim+foff+d];
1489:             if (fR) fR[d] += fluxR[iface*totDim+foff+d];
1490:           }
1491:           ++iface;
1492:         }
1493:       }
1494:       VecRestoreArray(locF, &fa);
1495:     }
1496:     /* Handle time derivative */
1497:     if (locX_t) {
1498:       PetscScalar *x_t, *fa;

1500:       VecGetArray(locF, &fa);
1501:       VecGetArray(locX_t, &x_t);
1502:       for (f = 0; f < Nf; ++f) {
1503:         PetscFV      fv;
1504:         PetscObject  obj;
1505:         PetscClassId id;
1506:         PetscInt     pdim, d;

1508:         PetscDSGetDiscretization(prob, f, &obj);
1509:         PetscObjectGetClassId(obj, &id);
1510:         if (id != PETSCFV_CLASSID) continue;
1511:         fv   = (PetscFV) obj;
1512:         PetscFVGetNumComponents(fv, &pdim);
1513:         for (c = cS; c < cE; ++c) {
1514:           const PetscInt cell = cells ? cells[c] : c;
1515:           PetscScalar   *u_t, *r;

1517:           if (ghostLabel) {
1518:             PetscInt ghostVal;

1520:             DMLabelGetValue(ghostLabel, cell, &ghostVal);
1521:             if (ghostVal > 0) continue;
1522:           }
1523:           DMPlexPointLocalFieldRead(dm, cell, f, x_t, &u_t);
1524:           DMPlexPointLocalFieldRef(dm, cell, f, fa, &r);
1525:           for (d = 0; d < pdim; ++d) r[d] += u_t[d];
1526:         }
1527:       }
1528:       VecRestoreArray(locX_t, &x_t);
1529:       VecRestoreArray(locF, &fa);
1530:     }
1531:     if (useFEM) {
1532:       DMPlexRestoreCellFields(dm, chunkIS, locX, locX_t, locA, &u, &u_t, &a);
1533:       DMRestoreWorkArray(dm, numCells*totDim, MPIU_SCALAR, &elemVec);
1534:     }
1535:     if (useFVM) {
1536:       DMPlexRestoreFaceFields(dm, fS, fE, locX, locX_t, faceGeometryFVM, cellGeometryFVM, locGrad, &numFaces, &uL, &uR);
1537:       DMPlexRestoreFaceGeometry(dm, fS, fE, faceGeometryFVM, cellGeometryFVM, &numFaces, &fgeom, &vol);
1538:       DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxL);
1539:       DMRestoreWorkArray(dm, numFaces*totDim, MPIU_SCALAR, &fluxR);
1540:       if (dmGrad) {DMRestoreLocalVector(dmGrad, &locGrad);}
1541:     }
1542:   }
1543:   if (useFEM) {ISDestroy(&chunkIS);}
1544:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);

1546:   if (useFEM) {
1547:     DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user);

1549:     if (maxDegree <= 1) {
1550:       DMSNESRestoreFEGeom(coordField,cellIS,affineQuad,PETSC_FALSE,&affineGeom);
1551:       PetscQuadratureDestroy(&affineQuad);
1552:     } else {
1553:       for (f = 0; f < Nf; ++f) {
1554:         DMSNESRestoreFEGeom(coordField,cellIS,quads[f],PETSC_FALSE,&geoms[f]);
1555:         PetscQuadratureDestroy(&quads[f]);
1556:       }
1557:       PetscFree2(quads,geoms);
1558:     }
1559:   }

1561:   /* FEM */
1562:   /* 1: Get sizes from dm and dmAux */
1563:   /* 2: Get geometric data */
1564:   /* 3: Handle boundary values */
1565:   /* 4: Loop over domain */
1566:   /*   Extract coefficients */
1567:   /* Loop over fields */
1568:   /*   Set tiling for FE*/
1569:   /*   Integrate FE residual to get elemVec */
1570:   /*     Loop over subdomain */
1571:   /*       Loop over quad points */
1572:   /*         Transform coords to real space */
1573:   /*         Evaluate field and aux fields at point */
1574:   /*         Evaluate residual at point */
1575:   /*         Transform residual to real space */
1576:   /*       Add residual to elemVec */
1577:   /* Loop over domain */
1578:   /*   Add elemVec to locX */

1580:   /* FVM */
1581:   /* Get geometric data */
1582:   /* If using gradients */
1583:   /*   Compute gradient data */
1584:   /*   Loop over domain faces */
1585:   /*     Count computational faces */
1586:   /*     Reconstruct cell gradient */
1587:   /*   Loop over domain cells */
1588:   /*     Limit cell gradients */
1589:   /* Handle boundary values */
1590:   /* Loop over domain faces */
1591:   /*   Read out field, centroid, normal, volume for each side of face */
1592:   /* Riemann solve over faces */
1593:   /* Loop over domain faces */
1594:   /*   Accumulate fluxes to cells */
1595:   /* TODO Change printFEM to printDisc here */
1596:   if (mesh->printFEM) {
1597:     Vec         locFbc;
1598:     PetscInt    pStart, pEnd, p, maxDof;
1599:     PetscScalar *zeroes;

1601:     VecDuplicate(locF,&locFbc);
1602:     VecCopy(locF,locFbc);
1603:     PetscSectionGetChart(section,&pStart,&pEnd);
1604:     PetscSectionGetMaxDof(section,&maxDof);
1605:     PetscCalloc1(maxDof,&zeroes);
1606:     for (p = pStart; p < pEnd; p++) {
1607:       VecSetValuesSection(locFbc,section,p,zeroes,INSERT_BC_VALUES);
1608:     }
1609:     PetscFree(zeroes);
1610:     DMPrintLocalVec(dm, name, mesh->printTol, locFbc);
1611:     VecDestroy(&locFbc);
1612:   }
1613:   PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);
1614:   return(0);
1615: }

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

1620:   Input Parameters:
1621: + dm - The mesh
1622: . X  - Local solution
1623: - user - The user context

1625:   Output Parameter:
1626: . F  - Local output vector

1628:   Level: developer

1630: .seealso: DMPlexComputeJacobianAction()
1631: @*/
1632: PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1633: {
1634:   DM             plex;
1635:   IS             cellIS;
1636:   PetscInt       depth;

1640:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1641:   DMPlexGetDepth(plex, &depth);
1642:   DMGetStratumIS(plex, "dim", depth, &cellIS);
1643:   if (!cellIS) {
1644:     DMGetStratumIS(plex, "depth", depth, &cellIS);
1645:   }
1646:   DMPlexComputeResidual_Internal(plex, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);
1647:   ISDestroy(&cellIS);
1648:   DMDestroy(&plex);
1649:   return(0);
1650: }

1652: /*@
1653:   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X

1655:   Input Parameters:
1656: + dm - The mesh
1657: - user - The user context

1659:   Output Parameter:
1660: . X  - Local solution

1662:   Level: developer

1664: .seealso: DMPlexComputeJacobianAction()
1665: @*/
1666: PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1667: {
1668:   DM             plex;

1672:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1673:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);
1674:   DMDestroy(&plex);
1675:   return(0);
1676: }

1678: 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)
1679: {
1680:   DM_Plex       *mesh = (DM_Plex *) dm->data;
1681:   DM             plex = NULL, plexA = NULL, tdm;
1682:   PetscDS        prob, probAux = NULL;
1683:   PetscSection   section, sectionAux = NULL;
1684:   PetscSection   globalSection, subSection = NULL;
1685:   Vec            locA = NULL, tv;
1686:   PetscScalar   *u = NULL, *u_t = NULL, *a = NULL, *elemMat = NULL;
1687:   PetscInt       v;
1688:   PetscInt       Nf, totDim, totDimAux = 0;
1689:   PetscBool      isMatISP, transform;

1693:   DMConvert(dm, DMPLEX, &plex);
1694:   DMHasBasisTransform(dm, &transform);
1695:   DMGetBasisTransformDM_Internal(dm, &tdm);
1696:   DMGetBasisTransformVec_Internal(dm, &tv);
1697:   DMGetLocalSection(dm, &section);
1698:   DMGetDS(dm, &prob);
1699:   PetscDSGetNumFields(prob, &Nf);
1700:   PetscDSGetTotalDimension(prob, &totDim);
1701:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &locA);
1702:   if (locA) {
1703:     DM dmAux;

1705:     VecGetDM(locA, &dmAux);
1706:     DMConvert(dmAux, DMPLEX, &plexA);
1707:     DMGetDS(plexA, &probAux);
1708:     PetscDSGetTotalDimension(probAux, &totDimAux);
1709:     DMGetLocalSection(plexA, &sectionAux);
1710:   }

1712:   PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);
1713:   DMGetGlobalSection(dm, &globalSection);
1714:   if (isMatISP) {DMPlexGetSubdomainSection(dm, &subSection);}
1715:   for (v = 0; v < numValues; ++v) {
1716:     PetscFEGeom    *fgeom;
1717:     PetscInt        maxDegree;
1718:     PetscQuadrature qGeom = NULL;
1719:     IS              pointIS;
1720:     const PetscInt *points;
1721:     PetscInt        numFaces, face, Nq;

1723:     DMLabelGetStratumIS(label, values[v], &pointIS);
1724:     if (!pointIS) continue; /* No points with that id on this process */
1725:     {
1726:       IS isectIS;

1728:       /* TODO: Special cases of ISIntersect where it is quick to check a prior if one is a superset of the other */
1729:       ISIntersect_Caching_Internal(facetIS,pointIS,&isectIS);
1730:       ISDestroy(&pointIS);
1731:       pointIS = isectIS;
1732:     }
1733:     ISGetLocalSize(pointIS, &numFaces);
1734:     ISGetIndices(pointIS, &points);
1735:     PetscMalloc4(numFaces*totDim, &u, locX_t ? numFaces*totDim : 0, &u_t, numFaces*totDim*totDim, &elemMat, locA ? numFaces*totDimAux : 0, &a);
1736:     DMFieldGetDegree(coordField,pointIS,NULL,&maxDegree);
1737:     if (maxDegree <= 1) {
1738:       DMFieldCreateDefaultQuadrature(coordField,pointIS,&qGeom);
1739:     }
1740:     if (!qGeom) {
1741:       PetscFE fe;

1743:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1744:       PetscFEGetFaceQuadrature(fe, &qGeom);
1745:       PetscObjectReference((PetscObject)qGeom);
1746:     }
1747:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
1748:     DMSNESGetFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
1749:     for (face = 0; face < numFaces; ++face) {
1750:       const PetscInt point = points[face], *support, *cone;
1751:       PetscScalar   *x     = NULL;
1752:       PetscInt       i, coneSize, faceLoc;

1754:       DMPlexGetSupport(dm, point, &support);
1755:       DMPlexGetConeSize(dm, support[0], &coneSize);
1756:       DMPlexGetCone(dm, support[0], &cone);
1757:       for (faceLoc = 0; faceLoc < coneSize; ++faceLoc) if (cone[faceLoc] == point) break;
1758:       if (faceLoc == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not find face %D in cone of support[0] %D", point, support[0]);
1759:       fgeom->face[face][0] = faceLoc;
1760:       DMPlexVecGetClosure(plex, section, locX, support[0], NULL, &x);
1761:       for (i = 0; i < totDim; ++i) u[face*totDim+i] = x[i];
1762:       DMPlexVecRestoreClosure(plex, section, locX, support[0], NULL, &x);
1763:       if (locX_t) {
1764:         DMPlexVecGetClosure(plex, section, locX_t, support[0], NULL, &x);
1765:         for (i = 0; i < totDim; ++i) u_t[face*totDim+i] = x[i];
1766:         DMPlexVecRestoreClosure(plex, section, locX_t, support[0], NULL, &x);
1767:       }
1768:       if (locA) {
1769:         PetscInt subp;
1770:         DMPlexGetSubpoint(plexA, support[0], &subp);
1771:         DMPlexVecGetClosure(plexA, sectionAux, locA, subp, NULL, &x);
1772:         for (i = 0; i < totDimAux; ++i) a[face*totDimAux+i] = x[i];
1773:         DMPlexVecRestoreClosure(plexA, sectionAux, locA, subp, NULL, &x);
1774:       }
1775:     }
1776:     PetscArrayzero(elemMat, numFaces*totDim*totDim);
1777:     {
1778:       PetscFE         fe;
1779:       PetscInt        Nb;
1780:       /* Conforming batches */
1781:       PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1782:       /* Remainder */
1783:       PetscFEGeom    *chunkGeom = NULL;
1784:       PetscInt        fieldJ, Nr, offset;

1786:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1787:       PetscFEGetDimension(fe, &Nb);
1788:       PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1789:       blockSize = Nb;
1790:       batchSize = numBlocks * blockSize;
1791:       PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1792:       numChunks = numFaces / (numBatches*batchSize);
1793:       Ne        = numChunks*numBatches*batchSize;
1794:       Nr        = numFaces % (numBatches*batchSize);
1795:       offset    = numFaces - Nr;
1796:       PetscFEGeomGetChunk(fgeom,0,offset,&chunkGeom);
1797:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1798:         PetscFEIntegrateBdJacobian(prob, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
1799:       }
1800:       PetscFEGeomGetChunk(fgeom,offset,numFaces,&chunkGeom);
1801:       for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1802:         PetscFEIntegrateBdJacobian(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]);
1803:       }
1804:       PetscFEGeomRestoreChunk(fgeom,offset,numFaces,&chunkGeom);
1805:     }
1806:     for (face = 0; face < numFaces; ++face) {
1807:       const PetscInt point = points[face], *support;

1809:       /* Transform to global basis before insertion in Jacobian */
1810:       DMPlexGetSupport(plex, point, &support);
1811:       if (transform) {DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, support[0], PETSC_TRUE, totDim, &elemMat[face*totDim*totDim]);}
1812:       if (mesh->printFEM > 1) {DMPrintCellMatrix(point, "BdJacobian", totDim, totDim, &elemMat[face*totDim*totDim]);}
1813:       if (!isMatISP) {
1814:         DMPlexMatSetClosure(plex, section, globalSection, JacP, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);
1815:       } else {
1816:         Mat lJ;

1818:         MatISGetLocalMat(JacP, &lJ);
1819:         DMPlexMatSetClosure(plex, section, subSection, lJ, support[0], &elemMat[face*totDim*totDim], ADD_VALUES);
1820:       }
1821:     }
1822:     DMSNESRestoreFEGeom(coordField,pointIS,qGeom,PETSC_TRUE,&fgeom);
1823:     PetscQuadratureDestroy(&qGeom);
1824:     ISRestoreIndices(pointIS, &points);
1825:     ISDestroy(&pointIS);
1826:     PetscFree4(u, u_t, elemMat, a);
1827:   }
1828:   if (plex)  {DMDestroy(&plex);}
1829:   if (plexA) {DMDestroy(&plexA);}
1830:   return(0);
1831: }

1833: 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)
1834: {
1835:   DMField        coordField;
1836:   DMLabel        depthLabel;
1837:   IS             facetIS;
1838:   PetscInt       dim;

1842:   DMGetDimension(dm, &dim);
1843:   DMPlexGetDepthLabel(dm, &depthLabel);
1844:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
1845:   DMGetCoordinateField(dm, &coordField);
1846:   DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);
1847:   ISDestroy(&facetIS);
1848:   return(0);
1849: }

1851: PetscErrorCode DMPlexComputeBdJacobian_Internal(DM dm, Vec locX, Vec locX_t, PetscReal t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
1852: {
1853:   PetscDS          prob;
1854:   PetscInt         dim, numBd, bd;
1855:   DMLabel          depthLabel;
1856:   DMField          coordField = NULL;
1857:   IS               facetIS;
1858:   PetscErrorCode   ierr;

1861:   DMGetDS(dm, &prob);
1862:   DMPlexGetDepthLabel(dm, &depthLabel);
1863:   DMGetDimension(dm, &dim);
1864:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
1865:   PetscDSGetNumBoundary(prob, &numBd);
1866:   DMGetCoordinateField(dm, &coordField);
1867:   for (bd = 0; bd < numBd; ++bd) {
1868:     DMBoundaryConditionType type;
1869:     const char             *bdLabel;
1870:     DMLabel                 label;
1871:     const PetscInt         *values;
1872:     PetscInt                fieldI, numValues;
1873:     PetscObject             obj;
1874:     PetscClassId            id;

1876:     PetscDSGetBoundary(prob, bd, &type, NULL, &bdLabel, &fieldI, NULL, NULL, NULL, &numValues, &values, NULL);
1877:     PetscDSGetDiscretization(prob, fieldI, &obj);
1878:     PetscObjectGetClassId(obj, &id);
1879:     if ((id != PETSCFE_CLASSID) || (type & DM_BC_ESSENTIAL)) continue;
1880:     DMGetLabel(dm, bdLabel, &label);
1881:     DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, fieldI, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);
1882:   }
1883:   ISDestroy(&facetIS);
1884:   return(0);
1885: }

1887: PetscErrorCode DMPlexComputeJacobian_Internal(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Mat Jac, Mat JacP,void *user)
1888: {
1889:   DM_Plex        *mesh  = (DM_Plex *) dm->data;
1890:   const char     *name  = "Jacobian";
1891:   DM              dmAux, plex, tdm;
1892:   Vec             A, tv;
1893:   DMField         coordField;
1894:   PetscDS         prob, probAux = NULL;
1895:   PetscSection    section, globalSection, subSection, sectionAux;
1896:   PetscScalar    *elemMat, *elemMatP, *elemMatD, *u, *u_t, *a = NULL;
1897:   const PetscInt *cells;
1898:   PetscInt        Nf, fieldI, fieldJ;
1899:   PetscInt        totDim, totDimAux, cStart, cEnd, numCells, c;
1900:   PetscBool       isMatIS, isMatISP, hasJac, hasPrec, hasDyn, hasFV = PETSC_FALSE, transform;
1901:   PetscErrorCode  ierr;

1904:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
1905:   DMHasBasisTransform(dm, &transform);
1906:   DMGetBasisTransformDM_Internal(dm, &tdm);
1907:   DMGetBasisTransformVec_Internal(dm, &tv);
1908:   DMGetLocalSection(dm, &section);
1909:   PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatISP);
1910:   DMGetGlobalSection(dm, &globalSection);
1911:   if (isMatISP) {DMPlexGetSubdomainSection(dm, &subSection);}
1912:   DMGetDS(dm, &prob);
1913:   PetscDSGetTotalDimension(prob, &totDim);
1914:   PetscDSHasJacobian(prob, &hasJac);
1915:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
1916:   /* user passed in the same matrix, avoid double contributions and
1917:      only assemble the Jacobian */
1918:   if (hasJac && Jac == JacP) hasPrec = PETSC_FALSE;
1919:   PetscDSHasDynamicJacobian(prob, &hasDyn);
1920:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
1921:   PetscSectionGetNumFields(section, &Nf);
1922:   ISGetLocalSize(cellIS, &numCells);
1923:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
1924:   PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);
1925:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);
1926:   if (dmAux) {
1927:     DMConvert(dmAux, DMPLEX, &plex);
1928:     DMGetLocalSection(plex, &sectionAux);
1929:     DMGetDS(dmAux, &probAux);
1930:     PetscDSGetTotalDimension(probAux, &totDimAux);
1931:   }
1932:   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);
1933:   if (dmAux) {PetscMalloc1(numCells*totDimAux, &a);}
1934:   DMGetCoordinateField(dm, &coordField);
1935:   for (c = cStart; c < cEnd; ++c) {
1936:     const PetscInt cell = cells ? cells[c] : c;
1937:     const PetscInt cind = c - cStart;
1938:     PetscScalar   *x = NULL,  *x_t = NULL;
1939:     PetscInt       i;

1941:     DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
1942:     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
1943:     DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
1944:     if (X_t) {
1945:       DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
1946:       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
1947:       DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
1948:     }
1949:     if (dmAux) {
1950:       PetscInt subcell;
1951:       DMPlexGetAuxiliaryPoint(dm, dmAux, cell, &subcell);
1952:       DMPlexVecGetClosure(plex, sectionAux, A, subcell, NULL, &x);
1953:       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
1954:       DMPlexVecRestoreClosure(plex, sectionAux, A, subcell, NULL, &x);
1955:     }
1956:   }
1957:   if (hasJac)  {PetscArrayzero(elemMat,  numCells*totDim*totDim);}
1958:   if (hasPrec) {PetscArrayzero(elemMatP, numCells*totDim*totDim);}
1959:   if (hasDyn)  {PetscArrayzero(elemMatD, numCells*totDim*totDim);}
1960:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
1961:     PetscClassId    id;
1962:     PetscFE         fe;
1963:     PetscQuadrature qGeom = NULL;
1964:     PetscInt        Nb;
1965:     /* Conforming batches */
1966:     PetscInt        numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1967:     /* Remainder */
1968:     PetscInt        Nr, offset, Nq;
1969:     PetscInt        maxDegree;
1970:     PetscFEGeom     *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

1972:     PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1973:     PetscObjectGetClassId((PetscObject) fe, &id);
1974:     if (id == PETSCFV_CLASSID) {hasFV = PETSC_TRUE; continue;}
1975:     PetscFEGetDimension(fe, &Nb);
1976:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1977:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
1978:     if (maxDegree <= 1) {
1979:       DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);
1980:     }
1981:     if (!qGeom) {
1982:       PetscFEGetQuadrature(fe,&qGeom);
1983:       PetscObjectReference((PetscObject)qGeom);
1984:     }
1985:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
1986:     DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
1987:     blockSize = Nb;
1988:     batchSize = numBlocks * blockSize;
1989:     PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1990:     numChunks = numCells / (numBatches*batchSize);
1991:     Ne        = numChunks*numBatches*batchSize;
1992:     Nr        = numCells % (numBatches*batchSize);
1993:     offset    = numCells - Nr;
1994:     PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
1995:     PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);
1996:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1997:       if (hasJac) {
1998:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
1999:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);
2000:       }
2001:       if (hasPrec) {
2002:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_PRE, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatP);
2003:         PetscFEIntegrateJacobian(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]);
2004:       }
2005:       if (hasDyn) {
2006:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);
2007:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);
2008:       }
2009:     }
2010:     PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);
2011:     PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);
2012:     DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
2013:     PetscQuadratureDestroy(&qGeom);
2014:   }
2015:   /*   Add contribution from X_t */
2016:   if (hasDyn) {for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];}
2017:   if (hasFV) {
2018:     PetscClassId id;
2019:     PetscFV      fv;
2020:     PetscInt     offsetI, NcI, NbI = 1, fc, f;

2022:     for (fieldI = 0; fieldI < Nf; ++fieldI) {
2023:       PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fv);
2024:       PetscDSGetFieldOffset(prob, fieldI, &offsetI);
2025:       PetscObjectGetClassId((PetscObject) fv, &id);
2026:       if (id != PETSCFV_CLASSID) continue;
2027:       /* Put in the identity */
2028:       PetscFVGetNumComponents(fv, &NcI);
2029:       for (c = cStart; c < cEnd; ++c) {
2030:         const PetscInt cind    = c - cStart;
2031:         const PetscInt eOffset = cind*totDim*totDim;
2032:         for (fc = 0; fc < NcI; ++fc) {
2033:           for (f = 0; f < NbI; ++f) {
2034:             const PetscInt i = offsetI + f*NcI+fc;
2035:             if (hasPrec) {
2036:               if (hasJac) {elemMat[eOffset+i*totDim+i] = 1.0;}
2037:               elemMatP[eOffset+i*totDim+i] = 1.0;
2038:             } else {elemMat[eOffset+i*totDim+i] = 1.0;}
2039:           }
2040:         }
2041:       }
2042:     }
2043:     /* No allocated space for FV stuff, so ignore the zero entries */
2044:     MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);
2045:   }
2046:   /* Insert values into matrix */
2047:   isMatIS = PETSC_FALSE;
2048:   if (hasPrec && hasJac) {
2049:     PetscObjectTypeCompare((PetscObject) JacP, MATIS, &isMatIS);
2050:   }
2051:   if (isMatIS && !subSection) {
2052:     DMPlexGetSubdomainSection(dm, &subSection);
2053:   }
2054:   for (c = cStart; c < cEnd; ++c) {
2055:     const PetscInt cell = cells ? cells[c] : c;
2056:     const PetscInt cind = c - cStart;

2058:     /* Transform to global basis before insertion in Jacobian */
2059:     if (transform) {DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind*totDim*totDim]);}
2060:     if (hasPrec) {
2061:       if (hasJac) {
2062:         if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);}
2063:         if (!isMatIS) {
2064:           DMPlexMatSetClosure(dm, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2065:         } else {
2066:           Mat lJ;

2068:           MatISGetLocalMat(Jac,&lJ);
2069:           DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2070:         }
2071:       }
2072:       if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMatP[cind*totDim*totDim]);}
2073:       if (!isMatISP) {
2074:         DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);
2075:       } else {
2076:         Mat lJ;

2078:         MatISGetLocalMat(JacP,&lJ);
2079:         DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);
2080:       }
2081:     } else {
2082:       if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);}
2083:       if (!isMatISP) {
2084:         DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2085:       } else {
2086:         Mat lJ;

2088:         MatISGetLocalMat(JacP,&lJ);
2089:         DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2090:       }
2091:     }
2092:   }
2093:   ISRestorePointRange(cellIS, &cStart, &cEnd, &cells);
2094:   if (hasFV) {MatSetOption(JacP, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);}
2095:   PetscFree5(u,u_t,elemMat,elemMatP,elemMatD);
2096:   if (dmAux) {
2097:     PetscFree(a);
2098:     DMDestroy(&plex);
2099:   }
2100:   /* Compute boundary integrals */
2101:   DMPlexComputeBdJacobian_Internal(dm, X, X_t, t, X_tShift, Jac, JacP, user);
2102:   /* Assemble matrix */
2103:   if (hasJac && hasPrec) {
2104:     MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);
2105:     MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);
2106:   }
2107:   MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);
2108:   MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);
2109:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
2110:   return(0);
2111: }

2113: /*@
2114:   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.

2116:   Input Parameters:
2117: + dm - The mesh
2118: . cellIS -
2119: . t  - The time
2120: . X_tShift - The multiplier for the Jacobian with repsect to X_t
2121: . X  - Local solution vector
2122: . X_t  - Time-derivative of the local solution vector
2123: . Y  - Local input vector
2124: - user - The user context

2126:   Output Parameter:
2127: . Z - Local output vector

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

2133:   Level: developer

2135: .seealso: FormFunctionLocal()
2136: @*/
2137: PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
2138: {
2139:   DM_Plex          *mesh  = (DM_Plex *) dm->data;
2140:   const char       *name  = "Jacobian";
2141:   DM                dmAux, plex, plexAux = NULL;
2142:   Vec               A;
2143:   PetscDS           prob, probAux = NULL;
2144:   PetscQuadrature   quad;
2145:   PetscSection      section, globalSection, sectionAux;
2146:   PetscScalar      *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
2147:   PetscInt          Nf, fieldI, fieldJ;
2148:   PetscInt          totDim, totDimAux = 0;
2149:   const PetscInt   *cells;
2150:   PetscInt          cStart, cEnd, numCells, c;
2151:   PetscBool         hasDyn;
2152:   DMField           coordField;
2153:   PetscErrorCode    ierr;

2156:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
2157:   DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
2158:   if (!cellIS) {
2159:     PetscInt depth;

2161:     DMPlexGetDepth(plex, &depth);
2162:     DMGetStratumIS(plex, "dim", depth, &cellIS);
2163:     if (!cellIS) {DMGetStratumIS(plex, "depth", depth, &cellIS);}
2164:   } else {
2165:     PetscObjectReference((PetscObject) cellIS);
2166:   }
2167:   DMGetLocalSection(dm, &section);
2168:   DMGetGlobalSection(dm, &globalSection);
2169:   DMGetDS(dm, &prob);
2170:   PetscDSGetTotalDimension(prob, &totDim);
2171:   PetscDSHasDynamicJacobian(prob, &hasDyn);
2172:   hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
2173:   PetscSectionGetNumFields(section, &Nf);
2174:   ISGetLocalSize(cellIS, &numCells);
2175:   ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
2176:   PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);
2177:   PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);
2178:   if (dmAux) {
2179:     DMConvert(dmAux, DMPLEX, &plexAux);
2180:     DMGetLocalSection(plexAux, &sectionAux);
2181:     DMGetDS(dmAux, &probAux);
2182:     PetscDSGetTotalDimension(probAux, &totDimAux);
2183:   }
2184:   VecSet(Z, 0.0);
2185:   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);
2186:   if (dmAux) {PetscMalloc1(numCells*totDimAux, &a);}
2187:   DMGetCoordinateField(dm, &coordField);
2188:   for (c = cStart; c < cEnd; ++c) {
2189:     const PetscInt cell = cells ? cells[c] : c;
2190:     const PetscInt cind = c - cStart;
2191:     PetscScalar   *x = NULL,  *x_t = NULL;
2192:     PetscInt       i;

2194:     DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
2195:     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
2196:     DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
2197:     if (X_t) {
2198:       DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
2199:       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
2200:       DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
2201:     }
2202:     if (dmAux) {
2203:       PetscInt subcell;
2204:       DMPlexGetAuxiliaryPoint(dm, dmAux, cell, &subcell);
2205:       DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);
2206:       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
2207:       DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);
2208:     }
2209:     DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);
2210:     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
2211:     DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);
2212:   }
2213:   PetscArrayzero(elemMat, numCells*totDim*totDim);
2214:   if (hasDyn)  {PetscArrayzero(elemMatD, numCells*totDim*totDim);}
2215:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2216:     PetscFE  fe;
2217:     PetscInt Nb;
2218:     /* Conforming batches */
2219:     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2220:     /* Remainder */
2221:     PetscInt Nr, offset, Nq;
2222:     PetscQuadrature qGeom = NULL;
2223:     PetscInt    maxDegree;
2224:     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

2226:     PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
2227:     PetscFEGetQuadrature(fe, &quad);
2228:     PetscFEGetDimension(fe, &Nb);
2229:     PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
2230:     DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
2231:     if (maxDegree <= 1) {DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);}
2232:     if (!qGeom) {
2233:       PetscFEGetQuadrature(fe,&qGeom);
2234:       PetscObjectReference((PetscObject)qGeom);
2235:     }
2236:     PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
2237:     DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
2238:     blockSize = Nb;
2239:     batchSize = numBlocks * blockSize;
2240:     PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
2241:     numChunks = numCells / (numBatches*batchSize);
2242:     Ne        = numChunks*numBatches*batchSize;
2243:     Nr        = numCells % (numBatches*batchSize);
2244:     offset    = numCells - Nr;
2245:     PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
2246:     PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);
2247:     for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
2248:       PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
2249:       PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMat[offset*totDim*totDim]);
2250:       if (hasDyn) {
2251:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);
2252:         PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Nr, remGeom, &u[offset*totDim], u_t ? &u_t[offset*totDim] : NULL, probAux, &a[offset*totDimAux], t, X_tShift, &elemMatD[offset*totDim*totDim]);
2253:       }
2254:     }
2255:     PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);
2256:     PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);
2257:     DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
2258:     PetscQuadratureDestroy(&qGeom);
2259:   }
2260:   if (hasDyn) {
2261:     for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
2262:   }
2263:   for (c = cStart; c < cEnd; ++c) {
2264:     const PetscInt     cell = cells ? cells[c] : c;
2265:     const PetscInt     cind = c - cStart;
2266:     const PetscBLASInt M = totDim, one = 1;
2267:     const PetscScalar  a = 1.0, b = 0.0;

2269:     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
2270:     if (mesh->printFEM > 1) {
2271:       DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
2272:       DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);
2273:       DMPrintCellVector(c, "Z",  totDim, z);
2274:     }
2275:     DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);
2276:   }
2277:   PetscFree6(u,u_t,elemMat,elemMatD,y,z);
2278:   if (mesh->printFEM) {
2279:     PetscPrintf(PETSC_COMM_WORLD, "Z:\n");
2280:     VecView(Z, NULL);
2281:   }
2282:   PetscFree(a);
2283:   ISDestroy(&cellIS);
2284:   DMDestroy(&plexAux);
2285:   DMDestroy(&plex);
2286:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
2287:   return(0);
2288: }

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

2293:   Input Parameters:
2294: + dm - The mesh
2295: . X  - Local input vector
2296: - user - The user context

2298:   Output Parameter:
2299: . Jac  - Jacobian matrix

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

2305:   Level: developer

2307: .seealso: FormFunctionLocal()
2308: @*/
2309: PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
2310: {
2311:   DM             plex;
2312:   PetscDS        prob;
2313:   IS             cellIS;
2314:   PetscBool      hasJac, hasPrec;
2315:   PetscInt       depth;

2319:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
2320:   DMPlexGetDepth(plex, &depth);
2321:   DMGetStratumIS(plex, "dim", depth, &cellIS);
2322:   if (!cellIS) {DMGetStratumIS(plex, "depth", depth, &cellIS);}
2323:   DMGetDS(dm, &prob);
2324:   PetscDSHasJacobian(prob, &hasJac);
2325:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
2326:   if (hasJac && hasPrec) {MatZeroEntries(Jac);}
2327:   MatZeroEntries(JacP);
2328:   DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);
2329:   ISDestroy(&cellIS);
2330:   DMDestroy(&plex);
2331:   return(0);
2332: }

2334: /*
2335:      MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.

2337:    Input Parameters:
2338: +     X - SNES linearization point
2339: .     ovl - index set of overlapping subdomains

2341:    Output Parameter:
2342: .     J - unassembled (Neumann) local matrix

2344:    Level: intermediate

2346: .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
2347: */
2348: static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
2349: {
2350:   SNES           snes;
2351:   Mat            pJ;
2352:   DM             ovldm,origdm;
2353:   DMSNES         sdm;
2354:   PetscErrorCode (*bfun)(DM,Vec,void*);
2355:   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
2356:   void           *bctx,*jctx;

2360:   PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);
2361:   if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");
2362:   PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);
2363:   if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");
2364:   MatGetDM(pJ,&ovldm);
2365:   DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);
2366:   DMSNESSetBoundaryLocal(ovldm,bfun,bctx);
2367:   DMSNESGetJacobianLocal(origdm,&jfun,&jctx);
2368:   DMSNESSetJacobianLocal(ovldm,jfun,jctx);
2369:   PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);
2370:   if (!snes) {
2371:     SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);
2372:     SNESSetDM(snes,ovldm);
2373:     PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);
2374:     PetscObjectDereference((PetscObject)snes);
2375:   }
2376:   DMGetDMSNES(ovldm,&sdm);
2377:   VecLockReadPush(X);
2378:   PetscStackPush("SNES user Jacobian function");
2379:   (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);
2380:   PetscStackPop;
2381:   VecLockReadPop(X);
2382:   /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
2383:   {
2384:     Mat locpJ;

2386:     MatISGetLocalMat(pJ,&locpJ);
2387:     MatCopy(locpJ,J,SAME_NONZERO_PATTERN);
2388:   }
2389:   return(0);
2390: }

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

2395:   Input Parameters:
2396: + dm - The DM object
2397: . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
2398: . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
2399: - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())

2401:   Level: developer
2402: @*/
2403: PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
2404: {

2408:   DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);
2409:   DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);
2410:   DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);
2411:   PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);
2412:   return(0);
2413: }

2415: /*@C
2416:   DMSNESCheckDiscretization - Check the discretization error of the exact solution

2418:   Input Parameters:
2419: + snes - the SNES object
2420: . dm   - the DM
2421: . u    - a DM vector
2422: . exactFuncs - pointwise functions of the exact solution for each field
2423: . ctxs - contexts for the functions
2424: . tol  - A tolerance for the check, or -1 to print the results instead

2426:   Output Parameters:
2427: . error - An array which holds the discretization error in each field, or NULL

2429:   Level: developer

2431: .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian()
2432: @*/
2433: PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, Vec u, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs, PetscReal tol, PetscReal error[])
2434: {
2435:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2436:   void            **ectxs;
2437:   MPI_Comm          comm;
2438:   PetscDS           ds;
2439:   PetscReal        *err;
2440:   PetscInt          Nf, f;
2441:   PetscErrorCode    ierr;

2448:   PetscObjectGetComm((PetscObject) snes, &comm);
2449:   DMGetDS(dm, &ds);
2450:   DMGetNumFields(dm, &Nf);
2451:   PetscMalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);
2452:   for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);}
2453:   DMProjectFunction(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, INSERT_ALL_VALUES, u);
2454:   PetscObjectSetName((PetscObject) u, "Exact Solution");
2455:   PetscObjectSetOptionsPrefix((PetscObject) u, "exact_");
2456:   VecViewFromOptions(u, NULL, "-vec_view");
2457:   if (Nf > 1) {
2458:     DMComputeL2FieldDiff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, u, err);
2459:     if (tol >= 0.0) {
2460:       for (f = 0; f < Nf; ++f) {
2461:         if (err[f] > tol) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %D exceeds tolerance %g", (double) err[f], f, (double) tol);
2462:       }
2463:     } else if (error) {
2464:       for (f = 0; f < Nf; ++f) error[f] = err[f];
2465:     } else {
2466:       PetscPrintf(comm, "L_2 Error: [");
2467:       for (f = 0; f < Nf; ++f) {
2468:         if (f) {PetscPrintf(comm, ", ");}
2469:         PetscPrintf(comm, "%g", (double)err[f]);
2470:       }
2471:       PetscPrintf(comm, "]\n");
2472:     }
2473:   } else {
2474:     DMComputeL2Diff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs , u, &err[0]);
2475:     if (tol >= 0.0) {
2476:       if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
2477:     } else if (error) {
2478:       error[0] = err[0];
2479:     } else {
2480:       PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0]);
2481:     }
2482:   }
2483:   PetscFree3(exacts, ectxs, err);
2484:   return(0);
2485: }

2487: /*@C
2488:   DMSNESCheckResidual - Check the residual of the exact solution

2490:   Input Parameters:
2491: + snes - the SNES object
2492: . dm   - the DM
2493: . u    - a DM vector
2494: . tol  - A tolerance for the check, or -1 to print the results instead

2496:   Output Parameters:
2497: . residual - The residual norm of the exact solution, or NULL

2499:   Level: developer

2501: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
2502: @*/
2503: PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
2504: {
2505:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2506:   void            **ectxs;
2507:   MPI_Comm          comm;
2508:   PetscDS           ds;
2509:   Vec               r;
2510:   PetscReal         res;
2511:   PetscInt          Nf, f;
2512:   PetscBool         computeSol = PETSC_FALSE;
2513:   PetscErrorCode    ierr;

2520:   PetscObjectGetComm((PetscObject) snes, &comm);
2521:   DMGetDS(dm, &ds);
2522:   DMGetNumFields(dm, &Nf);
2523:   PetscMalloc2(Nf, &exacts, Nf, &ectxs);
2524:   for (f = 0; f < Nf; ++f) {
2525:     PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);
2526:     if (exacts[f]) computeSol = PETSC_TRUE;
2527:   }
2528:   if (computeSol) {DMProjectFunction(dm, 0.0, exacts, ectxs, INSERT_ALL_VALUES, u);}
2529:   PetscFree2(exacts, ectxs);

2531:   VecDuplicate(u, &r);
2532:   SNESComputeFunction(snes, u, r);
2533:   VecNorm(r, NORM_2, &res);
2534:   if (tol >= 0.0) {
2535:     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
2536:   } else if (residual) {
2537:     *residual = res;
2538:   } else {
2539:     PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);
2540:     VecChop(r, 1.0e-10);
2541:     PetscObjectSetName((PetscObject) r, "Initial Residual");
2542:     PetscObjectSetOptionsPrefix((PetscObject)r,"res_");
2543:     VecViewFromOptions(r, NULL, "-vec_view");
2544:   }
2545:   VecDestroy(&r);
2546:   return(0);
2547: }

2549: /*@C
2550:   DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test

2552:   Input Parameters:
2553: + snes - the SNES object
2554: . dm   - the DM
2555: . u    - a DM vector
2556: . tol  - A tolerance for the check, or -1 to print the results instead

2558:   Output Parameters:
2559: + isLinear - Flag indicaing that the function looks linear, or NULL
2560: - convRate - The rate of convergence of the linear model, or NULL

2562:   Level: developer

2564: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
2565: @*/
2566: PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
2567: {
2568:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2569:   void            **ectxs;
2570:   MPI_Comm          comm;
2571:   PetscDS           ds;
2572:   Mat               J, M;
2573:   MatNullSpace      nullspace;
2574:   PetscReal         slope, intercept;
2575:   PetscInt          Nf, f;
2576:   PetscBool         hasJac, hasPrec, isLin = PETSC_FALSE, computeSol = PETSC_FALSE;
2577:   PetscErrorCode    ierr;

2585:   PetscObjectGetComm((PetscObject) snes, &comm);
2586:   DMGetDS(dm, &ds);
2587:   DMGetNumFields(dm, &Nf);
2588:   PetscMalloc2(Nf, &exacts, Nf, &ectxs);
2589:   for (f = 0; f < Nf; ++f) {
2590:     PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);
2591:     if (exacts[f]) computeSol = PETSC_TRUE;
2592:   }
2593:   if (computeSol) {DMProjectFunction(dm, 0.0, exacts, ectxs, INSERT_ALL_VALUES, u);}
2594:   PetscFree2(exacts, ectxs);

2596:   /* Create and view matrices */
2597:   DMCreateMatrix(dm, &J);
2598:   PetscDSHasJacobian(ds, &hasJac);
2599:   PetscDSHasJacobianPreconditioner(ds, &hasPrec);
2600:   if (hasJac && hasPrec) {
2601:     DMCreateMatrix(dm, &M);
2602:     SNESComputeJacobian(snes, u, J, M);
2603:     PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");
2604:     PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");
2605:     MatViewFromOptions(M, NULL, "-mat_view");
2606:     MatDestroy(&M);
2607:   } else {
2608:     SNESComputeJacobian(snes, u, J, J);
2609:   }
2610:   PetscObjectSetName((PetscObject) J, "Jacobian");
2611:   PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");
2612:   MatViewFromOptions(J, NULL, "-mat_view");
2613:   /* Check nullspace */
2614:   MatGetNullSpace(J, &nullspace);
2615:   if (nullspace) {
2616:     PetscBool isNull;
2617:     MatNullSpaceTest(nullspace, J, &isNull);
2618:     if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
2619:   }
2620:   MatNullSpaceDestroy(&nullspace);
2621:   /* Taylor test */
2622:   {
2623:     PetscRandom rand;
2624:     Vec         du, uhat, r, rhat, df;
2625:     PetscReal   h;
2626:     PetscReal  *es, *hs, *errors;
2627:     PetscReal   hMax = 1.0, hMin = 1e-6, hMult = 0.1;
2628:     PetscInt    Nv, v;

2630:     /* Choose a perturbation direction */
2631:     PetscRandomCreate(comm, &rand);
2632:     VecDuplicate(u, &du);
2633:     VecSetRandom(du, rand);
2634:     PetscRandomDestroy(&rand);
2635:     VecDuplicate(u, &df);
2636:     MatMult(J, du, df);
2637:     /* Evaluate residual at u, F(u), save in vector r */
2638:     VecDuplicate(u, &r);
2639:     SNESComputeFunction(snes, u, r);
2640:     /* Look at the convergence of our Taylor approximation as we approach u */
2641:     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
2642:     PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);
2643:     VecDuplicate(u, &uhat);
2644:     VecDuplicate(u, &rhat);
2645:     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
2646:       VecWAXPY(uhat, h, du, u);
2647:       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
2648:       SNESComputeFunction(snes, uhat, rhat);
2649:       VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);
2650:       VecNorm(rhat, NORM_2, &errors[Nv]);

2652:       es[Nv] = PetscLog10Real(errors[Nv]);
2653:       hs[Nv] = PetscLog10Real(h);
2654:     }
2655:     VecDestroy(&uhat);
2656:     VecDestroy(&rhat);
2657:     VecDestroy(&df);
2658:     VecDestroy(&r);
2659:     VecDestroy(&du);
2660:     for (v = 0; v < Nv; ++v) {
2661:       if ((tol >= 0) && (errors[v] > tol)) break;
2662:       else if (errors[v] > PETSC_SMALL)    break;
2663:     }
2664:     if (v == Nv) isLin = PETSC_TRUE;
2665:     PetscLinearRegression(Nv, hs, es, &slope, &intercept);
2666:     PetscFree3(es, hs, errors);
2667:     /* Slope should be about 2 */
2668:     if (tol >= 0) {
2669:       if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
2670:     } else if (isLinear || convRate) {
2671:       if (isLinear) *isLinear = isLin;
2672:       if (convRate) *convRate = slope;
2673:     } else {
2674:       if (!isLin) {PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);}
2675:       else        {PetscPrintf(comm, "Function appears to be linear\n");}
2676:     }
2677:   }
2678:   MatDestroy(&J);
2679:   return(0);
2680: }

2682: PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx), void **ctxs)
2683: {

2687:   DMSNESCheckDiscretization(snes, dm, u, exactFuncs, ctxs, -1.0, NULL);
2688:   DMSNESCheckResidual(snes, dm, u, -1.0, NULL);
2689:   DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);
2690:   return(0);
2691: }

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

2696:   Input Parameters:
2697: + snes - the SNES object
2698: . u    - representative SNES vector
2699: . exactFuncs - pointwise functions of the exact solution for each field
2700: - ctxs - contexts for the functions

2702:   Level: developer
2703: @*/
2704: PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2705: {
2706:   DM             dm;
2707:   Vec            sol;
2708:   PetscBool      check;

2712:   PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);
2713:   if (!check) return(0);
2714:   SNESGetDM(snes, &dm);
2715:   VecDuplicate(u, &sol);
2716:   SNESSetSolution(snes, sol);
2717:   DMSNESCheck_Internal(snes, dm, sol, exactFuncs, ctxs);
2718:   VecDestroy(&sol);
2719:   return(0);
2720: }