Actual source code: dmplexsnes.c

petsc-3.13.6 2020-09-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:   PetscTabulation    T;
548:   const PetscScalar *coords;
549:   PetscScalar    *a;
550:   PetscReal      xir[2];
551:   PetscInt       Nf, p;
552:   const PetscInt dof = ctx->dof;

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

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

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

603:       xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
604:       PetscFEComputeTabulation(fem, 1, xir, 0, T);
605:       for (comp = 0; comp < dof; ++comp) {
606:         a[p*dof+comp] = 0.0;
607:         for (d = 0; d < xSize/dof; ++d) {
608:           a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
609:         }
610:       }
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:   PetscTabulationDestroy(&T);
620:   VecRestoreArray(v, &a);
621:   VecRestoreArrayRead(ctx->coords, &coords);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

869:   Level: beginner

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

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

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

907:   Collective on ctx

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

912:   Level: beginner

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

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

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

933:   Collective on SNES

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

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

944:   Level: intermediate

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

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

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

993: /********************* Residual Computation **************************/


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

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

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

1007:   Level: developer

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

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

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

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

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

1035:   Level: developer

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1520:           if (ghostLabel) {
1521:             PetscInt ghostVal;

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

1549:   if (useFEM) {
1550:     DMPlexComputeBdResidual_Internal(dm, locX, locX_t, t, locF, user);

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

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

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

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

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

1623:   Input Parameters:
1624: + dm - The mesh
1625: . X  - Local solution
1626: - user - The user context

1628:   Output Parameter:
1629: . F  - Local output vector

1631:   Level: developer

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

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

1655: /*@
1656:   DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X

1658:   Input Parameters:
1659: + dm - The mesh
1660: - user - The user context

1662:   Output Parameter:
1663: . X  - Local solution

1665:   Level: developer

1667: .seealso: DMPlexComputeJacobianAction()
1668: @*/
1669: PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1670: {
1671:   DM             plex;

1675:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1676:   DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);
1677:   DMDestroy(&plex);
1678:   return(0);
1679: }

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

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

1709:     VecGetDM(locA, &dmAux);
1710:     DMGetEnclosureRelation(dmAux, dm, &encAux);
1711:     DMConvert(dmAux, DMPLEX, &plexA);
1712:     DMGetDS(plexA, &probAux);
1713:     PetscDSGetTotalDimension(probAux, &totDimAux);
1714:     DMGetLocalSection(plexA, &sectionAux);
1715:   }

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

1728:     DMLabelGetStratumIS(label, values[v], &pointIS);
1729:     if (!pointIS) continue; /* No points with that id on this process */
1730:     {
1731:       IS isectIS;

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

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

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

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

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

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

1838: 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)
1839: {
1840:   DMField        coordField;
1841:   DMLabel        depthLabel;
1842:   IS             facetIS;
1843:   PetscInt       dim;

1847:   DMGetDimension(dm, &dim);
1848:   DMPlexGetDepthLabel(dm, &depthLabel);
1849:   DMLabelGetStratumIS(depthLabel, dim-1, &facetIS);
1850:   DMGetCoordinateField(dm, &coordField);
1851:   DMPlexComputeBdJacobian_Single_Internal(dm, t, label, numValues, values, field, locX, locX_t, X_tShift, Jac, JacP, coordField, facetIS);
1852:   ISDestroy(&facetIS);
1853:   return(0);
1854: }

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

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

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

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

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

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

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

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

2065:     /* Transform to global basis before insertion in Jacobian */
2066:     if (transform) {DMPlexBasisTransformPointTensor_Internal(dm, tdm, tv, cell, PETSC_TRUE, totDim, &elemMat[cind*totDim*totDim]);}
2067:     if (hasPrec) {
2068:       if (hasJac) {
2069:         if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);}
2070:         if (!isMatIS) {
2071:           DMPlexMatSetClosure(dm, section, globalSection, Jac, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2072:         } else {
2073:           Mat lJ;

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

2085:         MatISGetLocalMat(JacP,&lJ);
2086:         DMPlexMatSetClosure(dm, section, subSection, lJ, cell, &elemMatP[cind*totDim*totDim], ADD_VALUES);
2087:       }
2088:     } else {
2089:       if (hasJac) {
2090:         if (mesh->printFEM > 1) {DMPrintCellMatrix(cell, name, totDim, totDim, &elemMat[cind*totDim*totDim]);}
2091:         if (!isMatISP) {
2092:           DMPlexMatSetClosure(dm, section, globalSection, JacP, cell, &elemMat[cind*totDim*totDim], ADD_VALUES);
2093:         } else {
2094:           Mat lJ;

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

2122: /*@
2123:   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.

2125:   Input Parameters:
2126: + dm - The mesh
2127: . cellIS -
2128: . t  - The time
2129: . X_tShift - The multiplier for the Jacobian with repsect to X_t
2130: . X  - Local solution vector
2131: . X_t  - Time-derivative of the local solution vector
2132: . Y  - Local input vector
2133: - user - The user context

2135:   Output Parameter:
2136: . Z - Local output vector

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

2142:   Level: developer

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

2166:   PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
2167:   DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
2168:   if (!cellIS) {
2169:     PetscInt depth;

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

2205:     DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
2206:     for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
2207:     DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
2208:     if (X_t) {
2209:       DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
2210:       for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
2211:       DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
2212:     }
2213:     if (dmAux) {
2214:       PetscInt subcell;
2215:       DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);
2216:       DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);
2217:       for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
2218:       DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);
2219:     }
2220:     DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);
2221:     for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
2222:     DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);
2223:   }
2224:   PetscArrayzero(elemMat, numCells*totDim*totDim);
2225:   if (hasDyn)  {PetscArrayzero(elemMatD, numCells*totDim*totDim);}
2226:   for (fieldI = 0; fieldI < Nf; ++fieldI) {
2227:     PetscFE  fe;
2228:     PetscInt Nb;
2229:     /* Conforming batches */
2230:     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
2231:     /* Remainder */
2232:     PetscInt Nr, offset, Nq;
2233:     PetscQuadrature qGeom = NULL;
2234:     PetscInt    maxDegree;
2235:     PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;

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

2280:     PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
2281:     if (mesh->printFEM > 1) {
2282:       DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
2283:       DMPrintCellVector(c, "Y",  totDim, &y[cind*totDim]);
2284:       DMPrintCellVector(c, "Z",  totDim, z);
2285:     }
2286:     DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);
2287:   }
2288:   PetscFree6(u,u_t,elemMat,elemMatD,y,z);
2289:   if (mesh->printFEM) {
2290:     PetscPrintf(PETSC_COMM_WORLD, "Z:\n");
2291:     VecView(Z, NULL);
2292:   }
2293:   PetscFree(a);
2294:   ISDestroy(&cellIS);
2295:   DMDestroy(&plexAux);
2296:   DMDestroy(&plex);
2297:   PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
2298:   return(0);
2299: }

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

2304:   Input Parameters:
2305: + dm - The mesh
2306: . X  - Local input vector
2307: - user - The user context

2309:   Output Parameter:
2310: . Jac  - Jacobian matrix

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

2316:   Level: developer

2318: .seealso: FormFunctionLocal()
2319: @*/
2320: PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
2321: {
2322:   DM             plex;
2323:   PetscDS        prob;
2324:   IS             cellIS;
2325:   PetscBool      hasJac, hasPrec;
2326:   PetscInt       depth;

2330:   DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
2331:   DMPlexGetDepth(plex, &depth);
2332:   DMGetStratumIS(plex, "dim", depth, &cellIS);
2333:   if (!cellIS) {DMGetStratumIS(plex, "depth", depth, &cellIS);}
2334:   DMGetDS(dm, &prob);
2335:   PetscDSHasJacobian(prob, &hasJac);
2336:   PetscDSHasJacobianPreconditioner(prob, &hasPrec);
2337:   if (hasJac && hasPrec) {MatZeroEntries(Jac);}
2338:   MatZeroEntries(JacP);
2339:   DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);
2340:   ISDestroy(&cellIS);
2341:   DMDestroy(&plex);
2342:   return(0);
2343: }

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

2348:    Input Parameters:
2349: +     X - SNES linearization point
2350: .     ovl - index set of overlapping subdomains

2352:    Output Parameter:
2353: .     J - unassembled (Neumann) local matrix

2355:    Level: intermediate

2357: .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
2358: */
2359: static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
2360: {
2361:   SNES           snes;
2362:   Mat            pJ;
2363:   DM             ovldm,origdm;
2364:   DMSNES         sdm;
2365:   PetscErrorCode (*bfun)(DM,Vec,void*);
2366:   PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
2367:   void           *bctx,*jctx;

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

2397:     MatISGetLocalMat(pJ,&locpJ);
2398:     MatCopy(locpJ,J,SAME_NONZERO_PATTERN);
2399:   }
2400:   return(0);
2401: }

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

2406:   Input Parameters:
2407: + dm - The DM object
2408: . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
2409: . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
2410: - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())

2412:   Level: developer
2413: @*/
2414: PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
2415: {

2419:   DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);
2420:   DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);
2421:   DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);
2422:   PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);
2423:   return(0);
2424: }

2426: /*@C
2427:   DMSNESCheckDiscretization - Check the discretization error of the exact solution

2429:   Input Parameters:
2430: + snes - the SNES object
2431: . dm   - the DM
2432: . u    - a DM vector
2433: . exactFuncs - pointwise functions of the exact solution for each field
2434: . ctxs - contexts for the functions
2435: - tol  - A tolerance for the check, or -1 to print the results instead

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

2440:   Level: developer

2442: .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian()
2443: @*/
2444: 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[])
2445: {
2446:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2447:   void            **ectxs;
2448:   MPI_Comm          comm;
2449:   PetscDS           ds;
2450:   PetscReal        *err;
2451:   PetscInt          Nf, f;
2452:   PetscErrorCode    ierr;

2459:   PetscObjectGetComm((PetscObject) snes, &comm);
2460:   DMGetDS(dm, &ds);
2461:   DMGetNumFields(dm, &Nf);
2462:   PetscMalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);
2463:   for (f = 0; f < Nf; ++f) {PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);}
2464:   DMProjectFunction(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, INSERT_ALL_VALUES, u);
2465:   PetscObjectSetName((PetscObject) u, "Exact Solution");
2466:   PetscObjectSetOptionsPrefix((PetscObject) u, "exact_");
2467:   VecViewFromOptions(u, NULL, "-vec_view");
2468:   if (Nf > 1) {
2469:     DMComputeL2FieldDiff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs, u, err);
2470:     if (tol >= 0.0) {
2471:       for (f = 0; f < Nf; ++f) {
2472:         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);
2473:       }
2474:     } else if (error) {
2475:       for (f = 0; f < Nf; ++f) error[f] = err[f];
2476:     } else {
2477:       PetscPrintf(comm, "L_2 Error: [");
2478:       for (f = 0; f < Nf; ++f) {
2479:         if (f) {PetscPrintf(comm, ", ");}
2480:         PetscPrintf(comm, "%g", (double)err[f]);
2481:       }
2482:       PetscPrintf(comm, "]\n");
2483:     }
2484:   } else {
2485:     DMComputeL2Diff(dm, 0.0, exactFuncs ? exactFuncs : exacts, ctxs ? ctxs : ectxs , u, &err[0]);
2486:     if (tol >= 0.0) {
2487:       if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
2488:     } else if (error) {
2489:       error[0] = err[0];
2490:     } else {
2491:       PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0]);
2492:     }
2493:   }
2494:   PetscFree3(exacts, ectxs, err);
2495:   return(0);
2496: }

2498: /*@C
2499:   DMSNESCheckResidual - Check the residual of the exact solution

2501:   Input Parameters:
2502: + snes - the SNES object
2503: . dm   - the DM
2504: . u    - a DM vector
2505: - tol  - A tolerance for the check, or -1 to print the results instead

2507:   Output Parameters:
2508: . residual - The residual norm of the exact solution, or NULL

2510:   Level: developer

2512: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
2513: @*/
2514: PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
2515: {
2516:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2517:   void            **ectxs;
2518:   MPI_Comm          comm;
2519:   PetscDS           ds;
2520:   Vec               r;
2521:   PetscReal         res;
2522:   PetscInt          Nf, f;
2523:   PetscBool         computeSol = PETSC_FALSE;
2524:   PetscErrorCode    ierr;

2531:   PetscObjectGetComm((PetscObject) snes, &comm);
2532:   DMGetDS(dm, &ds);
2533:   DMGetNumFields(dm, &Nf);
2534:   PetscMalloc2(Nf, &exacts, Nf, &ectxs);
2535:   for (f = 0; f < Nf; ++f) {
2536:     PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);
2537:     if (exacts[f]) computeSol = PETSC_TRUE;
2538:   }
2539:   if (computeSol) {DMProjectFunction(dm, 0.0, exacts, ectxs, INSERT_ALL_VALUES, u);}
2540:   PetscFree2(exacts, ectxs);

2542:   VecDuplicate(u, &r);
2543:   SNESComputeFunction(snes, u, r);
2544:   VecNorm(r, NORM_2, &res);
2545:   if (tol >= 0.0) {
2546:     if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
2547:   } else if (residual) {
2548:     *residual = res;
2549:   } else {
2550:     PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);
2551:     VecChop(r, 1.0e-10);
2552:     PetscObjectSetName((PetscObject) r, "Initial Residual");
2553:     PetscObjectSetOptionsPrefix((PetscObject)r,"res_");
2554:     VecViewFromOptions(r, NULL, "-vec_view");
2555:   }
2556:   VecDestroy(&r);
2557:   return(0);
2558: }

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

2563:   Input Parameters:
2564: + snes - the SNES object
2565: . dm   - the DM
2566: . u    - a DM vector
2567: - tol  - A tolerance for the check, or -1 to print the results instead

2569:   Output Parameters:
2570: + isLinear - Flag indicaing that the function looks linear, or NULL
2571: - convRate - The rate of convergence of the linear model, or NULL

2573:   Level: developer

2575: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
2576: @*/
2577: PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
2578: {
2579:   PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
2580:   void            **ectxs;
2581:   MPI_Comm          comm;
2582:   PetscDS           ds;
2583:   Mat               J, M;
2584:   MatNullSpace      nullspace;
2585:   PetscReal         slope, intercept;
2586:   PetscInt          Nf, f;
2587:   PetscBool         hasJac, hasPrec, isLin = PETSC_FALSE, computeSol = PETSC_FALSE;
2588:   PetscErrorCode    ierr;

2596:   PetscObjectGetComm((PetscObject) snes, &comm);
2597:   DMGetDS(dm, &ds);
2598:   DMGetNumFields(dm, &Nf);
2599:   PetscMalloc2(Nf, &exacts, Nf, &ectxs);
2600:   for (f = 0; f < Nf; ++f) {
2601:     PetscDSGetExactSolution(ds, f, &exacts[f], &ectxs[f]);
2602:     if (exacts[f]) computeSol = PETSC_TRUE;
2603:   }
2604:   if (computeSol) {DMProjectFunction(dm, 0.0, exacts, ectxs, INSERT_ALL_VALUES, u);}
2605:   PetscFree2(exacts, ectxs);

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

2641:     /* Choose a perturbation direction */
2642:     PetscRandomCreate(comm, &rand);
2643:     VecDuplicate(u, &du);
2644:     VecSetRandom(du, rand); 
2645:     PetscRandomDestroy(&rand);
2646:     VecDuplicate(u, &df);
2647:     MatMult(J, du, df);
2648:     /* Evaluate residual at u, F(u), save in vector r */
2649:     VecDuplicate(u, &r);
2650:     SNESComputeFunction(snes, u, r);
2651:     /* Look at the convergence of our Taylor approximation as we approach u */
2652:     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
2653:     PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);
2654:     VecDuplicate(u, &uhat);
2655:     VecDuplicate(u, &rhat);
2656:     for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
2657:       VecWAXPY(uhat, h, du, u);
2658:       /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
2659:       SNESComputeFunction(snes, uhat, rhat);
2660:       VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);
2661:       VecNorm(rhat, NORM_2, &errors[Nv]);

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

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

2698:   DMSNESCheckDiscretization(snes, dm, u, exactFuncs, ctxs, -1.0, NULL);
2699:   DMSNESCheckResidual(snes, dm, u, -1.0, NULL);
2700:   DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);
2701:   return(0);
2702: }

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

2707:   Input Parameters:
2708: + snes - the SNES object
2709: . u    - representative SNES vector
2710: . exactFuncs - pointwise functions of the exact solution for each field
2711: - ctxs - contexts for the functions

2713:   Level: developer
2714: @*/
2715: PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u, PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx), void **ctxs)
2716: {
2717:   DM             dm;
2718:   Vec            sol;
2719:   PetscBool      check;

2723:   PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);
2724:   if (!check) return(0);
2725:   SNESGetDM(snes, &dm);
2726:   VecDuplicate(u, &sol);
2727:   SNESSetSolution(snes, sol);
2728:   DMSNESCheck_Internal(snes, dm, sol, exactFuncs, ctxs);
2729:   VecDestroy(&sol);
2730:   return(0);
2731: }