Actual source code: dmplexsnes.c
petsc-3.14.6 2021-03-30
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, NULL, NULL);
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 **************************/
995: /*@
996: DMPlexSNESComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user
998: Input Parameters:
999: + dm - The mesh
1000: . X - Local solution
1001: - user - The user context
1003: Output Parameter:
1004: . F - Local output vector
1006: Level: developer
1008: .seealso: DMPlexComputeJacobianAction()
1009: @*/
1010: PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1011: {
1012: DM plex;
1013: IS allcellIS;
1014: PetscInt Nds, s, depth;
1018: DMGetNumDS(dm, &Nds);
1019: DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
1020: DMPlexGetDepth(plex, &depth);
1021: DMGetStratumIS(plex, "dim", depth, &allcellIS);
1022: if (!allcellIS) {DMGetStratumIS(plex, "depth", depth, &allcellIS);}
1023: for (s = 0; s < Nds; ++s) {
1024: PetscDS ds;
1025: DMLabel label;
1026: IS cellIS;
1028: DMGetRegionNumDS(dm, s, &label, NULL, &ds);
1029: if (!label) {
1030: PetscObjectReference((PetscObject) allcellIS);
1031: cellIS = allcellIS;
1032: } else {
1033: IS pointIS;
1035: DMLabelGetStratumIS(label, 1, &pointIS);
1036: ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);
1037: ISDestroy(&pointIS);
1038: }
1039: DMPlexComputeResidual_Internal(plex, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);
1040: ISDestroy(&cellIS);
1041: }
1042: ISDestroy(&allcellIS);
1043: DMDestroy(&plex);
1044: return(0);
1045: }
1047: /*@
1048: DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1050: Input Parameters:
1051: + dm - The mesh
1052: - user - The user context
1054: Output Parameter:
1055: . X - Local solution
1057: Level: developer
1059: .seealso: DMPlexComputeJacobianAction()
1060: @*/
1061: PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1062: {
1063: DM plex;
1067: DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1068: DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);
1069: DMDestroy(&plex);
1070: return(0);
1071: }
1073: /*@
1074: 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.
1076: Input Parameters:
1077: + dm - The mesh
1078: . cellIS -
1079: . t - The time
1080: . X_tShift - The multiplier for the Jacobian with repsect to X_t
1081: . X - Local solution vector
1082: . X_t - Time-derivative of the local solution vector
1083: . Y - Local input vector
1084: - user - The user context
1086: Output Parameter:
1087: . Z - Local output vector
1089: Note:
1090: We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1091: like a GPU, or vectorize on a multicore machine.
1093: Level: developer
1095: .seealso: FormFunctionLocal()
1096: @*/
1097: PetscErrorCode DMPlexComputeJacobianAction(DM dm, IS cellIS, PetscReal t, PetscReal X_tShift, Vec X, Vec X_t, Vec Y, Vec Z, void *user)
1098: {
1099: DM_Plex *mesh = (DM_Plex *) dm->data;
1100: const char *name = "Jacobian";
1101: DM dmAux, plex, plexAux = NULL;
1102: DMEnclosureType encAux;
1103: Vec A;
1104: PetscDS prob, probAux = NULL;
1105: PetscQuadrature quad;
1106: PetscSection section, globalSection, sectionAux;
1107: PetscScalar *elemMat, *elemMatD, *u, *u_t, *a = NULL, *y, *z;
1108: PetscInt Nf, fieldI, fieldJ;
1109: PetscInt totDim, totDimAux = 0;
1110: const PetscInt *cells;
1111: PetscInt cStart, cEnd, numCells, c;
1112: PetscBool hasDyn;
1113: DMField coordField;
1114: PetscErrorCode ierr;
1117: PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);
1118: DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
1119: if (!cellIS) {
1120: PetscInt depth;
1122: DMPlexGetDepth(plex, &depth);
1123: DMGetStratumIS(plex, "dim", depth, &cellIS);
1124: if (!cellIS) {DMGetStratumIS(plex, "depth", depth, &cellIS);}
1125: } else {
1126: PetscObjectReference((PetscObject) cellIS);
1127: }
1128: DMGetLocalSection(dm, §ion);
1129: DMGetGlobalSection(dm, &globalSection);
1130: ISGetLocalSize(cellIS, &numCells);
1131: ISGetPointRange(cellIS, &cStart, &cEnd, &cells);
1132: DMGetCellDS(dm, cells ? cells[cStart] : cStart, &prob);
1133: PetscDSGetNumFields(prob, &Nf);
1134: PetscDSGetTotalDimension(prob, &totDim);
1135: PetscDSHasDynamicJacobian(prob, &hasDyn);
1136: hasDyn = hasDyn && (X_tShift != 0.0) ? PETSC_TRUE : PETSC_FALSE;
1137: PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);
1138: PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);
1139: if (dmAux) {
1140: DMGetEnclosureRelation(dmAux, dm, &encAux);
1141: DMConvert(dmAux, DMPLEX, &plexAux);
1142: DMGetLocalSection(plexAux, §ionAux);
1143: DMGetDS(dmAux, &probAux);
1144: PetscDSGetTotalDimension(probAux, &totDimAux);
1145: }
1146: VecSet(Z, 0.0);
1147: 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);
1148: if (dmAux) {PetscMalloc1(numCells*totDimAux, &a);}
1149: DMGetCoordinateField(dm, &coordField);
1150: for (c = cStart; c < cEnd; ++c) {
1151: const PetscInt cell = cells ? cells[c] : c;
1152: const PetscInt cind = c - cStart;
1153: PetscScalar *x = NULL, *x_t = NULL;
1154: PetscInt i;
1156: DMPlexVecGetClosure(dm, section, X, cell, NULL, &x);
1157: for (i = 0; i < totDim; ++i) u[cind*totDim+i] = x[i];
1158: DMPlexVecRestoreClosure(dm, section, X, cell, NULL, &x);
1159: if (X_t) {
1160: DMPlexVecGetClosure(dm, section, X_t, cell, NULL, &x_t);
1161: for (i = 0; i < totDim; ++i) u_t[cind*totDim+i] = x_t[i];
1162: DMPlexVecRestoreClosure(dm, section, X_t, cell, NULL, &x_t);
1163: }
1164: if (dmAux) {
1165: PetscInt subcell;
1166: DMGetEnclosurePoint(dmAux, dm, encAux, cell, &subcell);
1167: DMPlexVecGetClosure(plexAux, sectionAux, A, subcell, NULL, &x);
1168: for (i = 0; i < totDimAux; ++i) a[cind*totDimAux+i] = x[i];
1169: DMPlexVecRestoreClosure(plexAux, sectionAux, A, subcell, NULL, &x);
1170: }
1171: DMPlexVecGetClosure(dm, section, Y, cell, NULL, &x);
1172: for (i = 0; i < totDim; ++i) y[cind*totDim+i] = x[i];
1173: DMPlexVecRestoreClosure(dm, section, Y, cell, NULL, &x);
1174: }
1175: PetscArrayzero(elemMat, numCells*totDim*totDim);
1176: if (hasDyn) {PetscArrayzero(elemMatD, numCells*totDim*totDim);}
1177: for (fieldI = 0; fieldI < Nf; ++fieldI) {
1178: PetscFE fe;
1179: PetscInt Nb;
1180: /* Conforming batches */
1181: PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
1182: /* Remainder */
1183: PetscInt Nr, offset, Nq;
1184: PetscQuadrature qGeom = NULL;
1185: PetscInt maxDegree;
1186: PetscFEGeom *cgeomFEM, *chunkGeom = NULL, *remGeom = NULL;
1188: PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);
1189: PetscFEGetQuadrature(fe, &quad);
1190: PetscFEGetDimension(fe, &Nb);
1191: PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);
1192: DMFieldGetDegree(coordField,cellIS,NULL,&maxDegree);
1193: if (maxDegree <= 1) {DMFieldCreateDefaultQuadrature(coordField,cellIS,&qGeom);}
1194: if (!qGeom) {
1195: PetscFEGetQuadrature(fe,&qGeom);
1196: PetscObjectReference((PetscObject)qGeom);
1197: }
1198: PetscQuadratureGetData(qGeom, NULL, NULL, &Nq, NULL, NULL);
1199: DMSNESGetFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
1200: blockSize = Nb;
1201: batchSize = numBlocks * blockSize;
1202: PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);
1203: numChunks = numCells / (numBatches*batchSize);
1204: Ne = numChunks*numBatches*batchSize;
1205: Nr = numCells % (numBatches*batchSize);
1206: offset = numCells - Nr;
1207: PetscFEGeomGetChunk(cgeomFEM,0,offset,&chunkGeom);
1208: PetscFEGeomGetChunk(cgeomFEM,offset,numCells,&remGeom);
1209: for (fieldJ = 0; fieldJ < Nf; ++fieldJ) {
1210: PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMat);
1211: 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]);
1212: if (hasDyn) {
1213: PetscFEIntegrateJacobian(prob, PETSCFE_JACOBIAN_DYN, fieldI, fieldJ, Ne, chunkGeom, u, u_t, probAux, a, t, X_tShift, elemMatD);
1214: 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]);
1215: }
1216: }
1217: PetscFEGeomRestoreChunk(cgeomFEM,offset,numCells,&remGeom);
1218: PetscFEGeomRestoreChunk(cgeomFEM,0,offset,&chunkGeom);
1219: DMSNESRestoreFEGeom(coordField,cellIS,qGeom,PETSC_FALSE,&cgeomFEM);
1220: PetscQuadratureDestroy(&qGeom);
1221: }
1222: if (hasDyn) {
1223: for (c = 0; c < numCells*totDim*totDim; ++c) elemMat[c] += X_tShift*elemMatD[c];
1224: }
1225: for (c = cStart; c < cEnd; ++c) {
1226: const PetscInt cell = cells ? cells[c] : c;
1227: const PetscInt cind = c - cStart;
1228: const PetscBLASInt M = totDim, one = 1;
1229: const PetscScalar a = 1.0, b = 0.0;
1231: PetscStackCallBLAS("BLASgemv", BLASgemv_("N", &M, &M, &a, &elemMat[cind*totDim*totDim], &M, &y[cind*totDim], &one, &b, z, &one));
1232: if (mesh->printFEM > 1) {
1233: DMPrintCellMatrix(c, name, totDim, totDim, &elemMat[cind*totDim*totDim]);
1234: DMPrintCellVector(c, "Y", totDim, &y[cind*totDim]);
1235: DMPrintCellVector(c, "Z", totDim, z);
1236: }
1237: DMPlexVecSetClosure(dm, section, Z, cell, z, ADD_VALUES);
1238: }
1239: PetscFree6(u,u_t,elemMat,elemMatD,y,z);
1240: if (mesh->printFEM) {
1241: PetscPrintf(PetscObjectComm((PetscObject)Z), "Z:\n");
1242: VecView(Z, NULL);
1243: }
1244: PetscFree(a);
1245: ISDestroy(&cellIS);
1246: DMDestroy(&plexAux);
1247: DMDestroy(&plex);
1248: PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);
1249: return(0);
1250: }
1252: /*@
1253: DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1255: Input Parameters:
1256: + dm - The mesh
1257: . X - Local input vector
1258: - user - The user context
1260: Output Parameter:
1261: . Jac - Jacobian matrix
1263: Note:
1264: We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1265: like a GPU, or vectorize on a multicore machine.
1267: Level: developer
1269: .seealso: FormFunctionLocal()
1270: @*/
1271: PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1272: {
1273: DM plex;
1274: IS allcellIS;
1275: PetscBool hasJac, hasPrec;
1276: PetscInt Nds, s, depth;
1280: DMGetNumDS(dm, &Nds);
1281: DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
1282: DMPlexGetDepth(plex, &depth);
1283: DMGetStratumIS(plex, "dim", depth, &allcellIS);
1284: if (!allcellIS) {DMGetStratumIS(plex, "depth", depth, &allcellIS);}
1285: for (s = 0; s < Nds; ++s) {
1286: PetscDS ds;
1287: DMLabel label;
1288: IS cellIS;
1290: DMGetRegionNumDS(dm, s, &label, NULL, &ds);
1291: if (!label) {
1292: PetscObjectReference((PetscObject) allcellIS);
1293: cellIS = allcellIS;
1294: } else {
1295: IS pointIS;
1297: DMLabelGetStratumIS(label, 1, &pointIS);
1298: ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);
1299: ISDestroy(&pointIS);
1300: }
1301: if (!s) {
1302: PetscDSHasJacobian(ds, &hasJac);
1303: PetscDSHasJacobianPreconditioner(ds, &hasPrec);
1304: if (hasJac && hasPrec) {MatZeroEntries(Jac);}
1305: MatZeroEntries(JacP);
1306: }
1307: DMPlexComputeJacobian_Internal(plex, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);
1308: ISDestroy(&cellIS);
1309: }
1310: ISDestroy(&allcellIS);
1311: DMDestroy(&plex);
1312: return(0);
1313: }
1315: /*
1316: MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
1318: Input Parameters:
1319: + X - SNES linearization point
1320: . ovl - index set of overlapping subdomains
1322: Output Parameter:
1323: . J - unassembled (Neumann) local matrix
1325: Level: intermediate
1327: .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
1328: */
1329: static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
1330: {
1331: SNES snes;
1332: Mat pJ;
1333: DM ovldm,origdm;
1334: DMSNES sdm;
1335: PetscErrorCode (*bfun)(DM,Vec,void*);
1336: PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
1337: void *bctx,*jctx;
1341: PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);
1342: if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");
1343: PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);
1344: if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");
1345: MatGetDM(pJ,&ovldm);
1346: DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);
1347: DMSNESSetBoundaryLocal(ovldm,bfun,bctx);
1348: DMSNESGetJacobianLocal(origdm,&jfun,&jctx);
1349: DMSNESSetJacobianLocal(ovldm,jfun,jctx);
1350: PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);
1351: if (!snes) {
1352: SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);
1353: SNESSetDM(snes,ovldm);
1354: PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);
1355: PetscObjectDereference((PetscObject)snes);
1356: }
1357: DMGetDMSNES(ovldm,&sdm);
1358: VecLockReadPush(X);
1359: PetscStackPush("SNES user Jacobian function");
1360: (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);
1361: PetscStackPop;
1362: VecLockReadPop(X);
1363: /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
1364: {
1365: Mat locpJ;
1367: MatISGetLocalMat(pJ,&locpJ);
1368: MatCopy(locpJ,J,SAME_NONZERO_PATTERN);
1369: }
1370: return(0);
1371: }
1373: /*@
1374: DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
1376: Input Parameters:
1377: + dm - The DM object
1378: . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
1379: . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
1380: - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
1382: Level: developer
1383: @*/
1384: PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
1385: {
1389: DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);
1390: DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);
1391: DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);
1392: PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);
1393: return(0);
1394: }
1396: /*@C
1397: DMSNESCheckDiscretization - Check the discretization error of the exact solution
1399: Input Parameters:
1400: + snes - the SNES object
1401: . dm - the DM
1402: . t - the time
1403: . u - a DM vector
1404: - tol - A tolerance for the check, or -1 to print the results instead
1406: Output Parameters:
1407: . error - An array which holds the discretization error in each field, or NULL
1409: Note: The user must call PetscDSSetExactSolution() beforehand
1411: Level: developer
1413: .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian(), PetscDSSetExactSolution()
1414: @*/
1415: PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
1416: {
1417: PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1418: void **ectxs;
1419: PetscReal *err;
1420: MPI_Comm comm;
1421: PetscInt Nf, f;
1422: PetscErrorCode ierr;
1430: DMComputeExactSolution(dm, t, u, NULL);
1431: VecViewFromOptions(u, NULL, "-vec_view");
1433: PetscObjectGetComm((PetscObject) snes, &comm);
1434: DMGetNumFields(dm, &Nf);
1435: PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);
1436: {
1437: PetscInt Nds, s;
1439: DMGetNumDS(dm, &Nds);
1440: for (s = 0; s < Nds; ++s) {
1441: PetscDS ds;
1442: DMLabel label;
1443: IS fieldIS;
1444: const PetscInt *fields;
1445: PetscInt dsNf, f;
1447: DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);
1448: PetscDSGetNumFields(ds, &dsNf);
1449: ISGetIndices(fieldIS, &fields);
1450: for (f = 0; f < dsNf; ++f) {
1451: const PetscInt field = fields[f];
1452: PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]);
1453: }
1454: ISRestoreIndices(fieldIS, &fields);
1455: }
1456: }
1457: if (Nf > 1) {
1458: DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err);
1459: if (tol >= 0.0) {
1460: for (f = 0; f < Nf; ++f) {
1461: 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);
1462: }
1463: } else if (error) {
1464: for (f = 0; f < Nf; ++f) error[f] = err[f];
1465: } else {
1466: PetscPrintf(comm, "L_2 Error: [");
1467: for (f = 0; f < Nf; ++f) {
1468: if (f) {PetscPrintf(comm, ", ");}
1469: PetscPrintf(comm, "%g", (double)err[f]);
1470: }
1471: PetscPrintf(comm, "]\n");
1472: }
1473: } else {
1474: DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]);
1475: if (tol >= 0.0) {
1476: if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1477: } else if (error) {
1478: error[0] = err[0];
1479: } else {
1480: PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]);
1481: }
1482: }
1483: PetscFree3(exacts, ectxs, err);
1484: return(0);
1485: }
1487: /*@C
1488: DMSNESCheckResidual - Check the residual of the exact solution
1490: Input Parameters:
1491: + snes - the SNES object
1492: . dm - the DM
1493: . u - a DM vector
1494: - tol - A tolerance for the check, or -1 to print the results instead
1496: Output Parameters:
1497: . residual - The residual norm of the exact solution, or NULL
1499: Level: developer
1501: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
1502: @*/
1503: PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1504: {
1505: MPI_Comm comm;
1506: Vec r;
1507: PetscReal res;
1515: PetscObjectGetComm((PetscObject) snes, &comm);
1516: DMComputeExactSolution(dm, 0.0, u, NULL);
1517: VecDuplicate(u, &r);
1518: SNESComputeFunction(snes, u, r);
1519: VecNorm(r, NORM_2, &res);
1520: if (tol >= 0.0) {
1521: if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1522: } else if (residual) {
1523: *residual = res;
1524: } else {
1525: PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);
1526: VecChop(r, 1.0e-10);
1527: PetscObjectSetName((PetscObject) r, "Initial Residual");
1528: PetscObjectSetOptionsPrefix((PetscObject)r,"res_");
1529: VecViewFromOptions(r, NULL, "-vec_view");
1530: }
1531: VecDestroy(&r);
1532: return(0);
1533: }
1535: /*@C
1536: DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1538: Input Parameters:
1539: + snes - the SNES object
1540: . dm - the DM
1541: . u - a DM vector
1542: - tol - A tolerance for the check, or -1 to print the results instead
1544: Output Parameters:
1545: + isLinear - Flag indicaing that the function looks linear, or NULL
1546: - convRate - The rate of convergence of the linear model, or NULL
1548: Level: developer
1550: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
1551: @*/
1552: PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1553: {
1554: MPI_Comm comm;
1555: PetscDS ds;
1556: Mat J, M;
1557: MatNullSpace nullspace;
1558: PetscReal slope, intercept;
1559: PetscBool hasJac, hasPrec, isLin = PETSC_FALSE;
1568: PetscObjectGetComm((PetscObject) snes, &comm);
1569: DMComputeExactSolution(dm, 0.0, u, NULL);
1570: /* Create and view matrices */
1571: DMCreateMatrix(dm, &J);
1572: DMGetDS(dm, &ds);
1573: PetscDSHasJacobian(ds, &hasJac);
1574: PetscDSHasJacobianPreconditioner(ds, &hasPrec);
1575: if (hasJac && hasPrec) {
1576: DMCreateMatrix(dm, &M);
1577: SNESComputeJacobian(snes, u, J, M);
1578: PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");
1579: PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");
1580: MatViewFromOptions(M, NULL, "-mat_view");
1581: MatDestroy(&M);
1582: } else {
1583: SNESComputeJacobian(snes, u, J, J);
1584: }
1585: PetscObjectSetName((PetscObject) J, "Jacobian");
1586: PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");
1587: MatViewFromOptions(J, NULL, "-mat_view");
1588: /* Check nullspace */
1589: MatGetNullSpace(J, &nullspace);
1590: if (nullspace) {
1591: PetscBool isNull;
1592: MatNullSpaceTest(nullspace, J, &isNull);
1593: if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
1594: }
1595: MatNullSpaceDestroy(&nullspace);
1596: /* Taylor test */
1597: {
1598: PetscRandom rand;
1599: Vec du, uhat, r, rhat, df;
1600: PetscReal h;
1601: PetscReal *es, *hs, *errors;
1602: PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1603: PetscInt Nv, v;
1605: /* Choose a perturbation direction */
1606: PetscRandomCreate(comm, &rand);
1607: VecDuplicate(u, &du);
1608: VecSetRandom(du, rand);
1609: PetscRandomDestroy(&rand);
1610: VecDuplicate(u, &df);
1611: MatMult(J, du, df);
1612: /* Evaluate residual at u, F(u), save in vector r */
1613: VecDuplicate(u, &r);
1614: SNESComputeFunction(snes, u, r);
1615: /* Look at the convergence of our Taylor approximation as we approach u */
1616: for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
1617: PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);
1618: VecDuplicate(u, &uhat);
1619: VecDuplicate(u, &rhat);
1620: for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
1621: VecWAXPY(uhat, h, du, u);
1622: /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
1623: SNESComputeFunction(snes, uhat, rhat);
1624: VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);
1625: VecNorm(rhat, NORM_2, &errors[Nv]);
1627: es[Nv] = PetscLog10Real(errors[Nv]);
1628: hs[Nv] = PetscLog10Real(h);
1629: }
1630: VecDestroy(&uhat);
1631: VecDestroy(&rhat);
1632: VecDestroy(&df);
1633: VecDestroy(&r);
1634: VecDestroy(&du);
1635: for (v = 0; v < Nv; ++v) {
1636: if ((tol >= 0) && (errors[v] > tol)) break;
1637: else if (errors[v] > PETSC_SMALL) break;
1638: }
1639: if (v == Nv) isLin = PETSC_TRUE;
1640: PetscLinearRegression(Nv, hs, es, &slope, &intercept);
1641: PetscFree3(es, hs, errors);
1642: /* Slope should be about 2 */
1643: if (tol >= 0) {
1644: if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
1645: } else if (isLinear || convRate) {
1646: if (isLinear) *isLinear = isLin;
1647: if (convRate) *convRate = slope;
1648: } else {
1649: if (!isLin) {PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);}
1650: else {PetscPrintf(comm, "Function appears to be linear\n");}
1651: }
1652: }
1653: MatDestroy(&J);
1654: return(0);
1655: }
1657: PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1658: {
1662: DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL);
1663: DMSNESCheckResidual(snes, dm, u, -1.0, NULL);
1664: DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);
1665: return(0);
1666: }
1668: /*@C
1669: DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1671: Input Parameters:
1672: + snes - the SNES object
1673: - u - representative SNES vector
1675: Note: The user must call PetscDSSetExactSolution() beforehand
1677: Level: developer
1678: @*/
1679: PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
1680: {
1681: DM dm;
1682: Vec sol;
1683: PetscBool check;
1687: PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);
1688: if (!check) return(0);
1689: SNESGetDM(snes, &dm);
1690: VecDuplicate(u, &sol);
1691: SNESSetSolution(snes, sol);
1692: DMSNESCheck_Internal(snes, dm, sol);
1693: VecDestroy(&sol);
1694: return(0);
1695: }