Actual source code: dmplexsnes.c
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/snesimpl.h>
3: #include <petscds.h>
4: #include <petsc/private/petscimpl.h>
5: #include <petsc/private/petscfeimpl.h>
7: static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux,
8: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
9: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
10: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[])
11: {
12: p[0] = u[uOff[1]];
13: }
15: /*
16: SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately.
18: Collective on SNES
20: Input Parameters:
21: + snes - The SNES
22: . pfield - The field number for pressure
23: . nullspace - The pressure nullspace
24: . u - The solution vector
25: - ctx - An optional user context
27: Output Parameter:
28: . u - The solution with a continuum pressure integral of zero
30: Notes:
31: If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0. We assume that the nullspace is a single vector given explicitly.
33: Level: developer
35: .seealso: SNESConvergedCorrectPressure()
36: */
37: static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx)
38: {
39: DM dm;
40: PetscDS ds;
41: const Vec *nullvecs;
42: PetscScalar pintd, *intc, *intn;
43: MPI_Comm comm;
44: PetscInt Nf, Nv;
48: PetscObjectGetComm((PetscObject) snes, &comm);
49: SNESGetDM(snes, &dm);
50: if (!dm) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM");
51: if (!nullspace) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace");
52: DMGetDS(dm, &ds);
53: PetscDSSetObjective(ds, pfield, pressure_Private);
54: MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs);
55: if (Nv != 1) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %D", Nv);
56: VecDot(nullvecs[0], u, &pintd);
57: if (PetscAbsScalar(pintd) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));
58: PetscDSGetNumFields(ds, &Nf);
59: PetscMalloc2(Nf, &intc, Nf, &intn);
60: DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx);
61: DMPlexComputeIntegralFEM(dm, u, intc, ctx);
62: VecAXPY(u, -intc[pfield]/intn[pfield], nullvecs[0]);
63: #if defined (PETSC_USE_DEBUG)
64: DMPlexComputeIntegralFEM(dm, u, intc, ctx);
65: if (PetscAbsScalar(intc[pfield]) > PETSC_SMALL) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g\n", (double) PetscRealPart(intc[pfield]));
66: #endif
67: PetscFree2(intc, intn);
68: return(0);
69: }
71: /*@C
72: SNESConvergedCorrectPressure - Convergence test that adds a vector in the nullspace to make the continuum integral of the pressure field equal to zero. This is normally used only to evaluate convergence rates for the pressure accurately. The convergence test itself just mimics SNESConvergedDefault().
74: Logically Collective on SNES
76: Input Parameters:
77: + snes - the SNES context
78: . it - the iteration (0 indicates before any Newton steps)
79: . xnorm - 2-norm of current iterate
80: . snorm - 2-norm of current step
81: . fnorm - 2-norm of function at current iterate
82: - ctx - Optional user context
84: Output Parameter:
85: . reason - SNES_CONVERGED_ITERATING, SNES_CONVERGED_ITS, or SNES_DIVERGED_FNORM_NAN
87: Notes:
88: In order to use this monitor, you must setup several PETSc structures. First fields must be added to the DM, and a PetscDS must be created with discretizations of those fields. We currently assume that the pressure field has index 1. The pressure field must have a nullspace, likely created using the DMSetNullSpaceConstructor() interface. Last we must be able to integrate the pressure over the domain, so the DM attached to the SNES must be a Plex at this time.
90: Level: advanced
92: .seealso: SNESConvergedDefault(), SNESSetConvergenceTest(), DMSetNullSpaceConstructor()
93: @*/
94: PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx)
95: {
96: PetscBool monitorIntegral = PETSC_FALSE;
100: SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx);
101: if (monitorIntegral) {
102: Mat J;
103: Vec u;
104: MatNullSpace nullspace;
105: const Vec *nullvecs;
106: PetscScalar pintd;
108: SNESGetSolution(snes, &u);
109: SNESGetJacobian(snes, &J, NULL, NULL, NULL);
110: MatGetNullSpace(J, &nullspace);
111: MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs);
112: VecDot(nullvecs[0], u, &pintd);
113: PetscInfo1(snes, "SNES: Discrete integral of pressure: %g\n", (double) PetscRealPart(pintd));
114: }
115: if (*reason > 0) {
116: Mat J;
117: Vec u;
118: MatNullSpace nullspace;
119: PetscInt pfield = 1;
121: SNESGetSolution(snes, &u);
122: SNESGetJacobian(snes, &J, NULL, NULL, NULL);
123: MatGetNullSpace(J, &nullspace);
124: SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx);
125: }
126: return(0);
127: }
129: /************************** Interpolation *******************************/
131: static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy)
132: {
133: PetscBool isPlex;
137: PetscObjectTypeCompare((PetscObject) dm, DMPLEX, &isPlex);
138: if (isPlex) {
139: *plex = dm;
140: PetscObjectReference((PetscObject) dm);
141: } else {
142: PetscObjectQuery((PetscObject) dm, "dm_plex", (PetscObject *) plex);
143: if (!*plex) {
144: DMConvert(dm,DMPLEX,plex);
145: PetscObjectCompose((PetscObject) dm, "dm_plex", (PetscObject) *plex);
146: if (copy) {
147: DMCopyDMSNES(dm, *plex);
148: DMCopyAuxiliaryVec(dm, *plex);
149: }
150: } else {
151: PetscObjectReference((PetscObject) *plex);
152: }
153: }
154: return(0);
155: }
157: /*@C
158: DMInterpolationCreate - Creates a DMInterpolationInfo context
160: Collective
162: Input Parameter:
163: . comm - the communicator
165: Output Parameter:
166: . ctx - the context
168: Level: beginner
170: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationDestroy()
171: @*/
172: PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
173: {
178: PetscNew(ctx);
180: (*ctx)->comm = comm;
181: (*ctx)->dim = -1;
182: (*ctx)->nInput = 0;
183: (*ctx)->points = NULL;
184: (*ctx)->cells = NULL;
185: (*ctx)->n = -1;
186: (*ctx)->coords = NULL;
187: return(0);
188: }
190: /*@C
191: DMInterpolationSetDim - Sets the spatial dimension for the interpolation context
193: Not collective
195: Input Parameters:
196: + ctx - the context
197: - dim - the spatial dimension
199: Level: intermediate
201: .seealso: DMInterpolationGetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
202: @*/
203: PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
204: {
206: if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %D", dim);
207: ctx->dim = dim;
208: return(0);
209: }
211: /*@C
212: DMInterpolationGetDim - Gets the spatial dimension for the interpolation context
214: Not collective
216: Input Parameter:
217: . ctx - the context
219: Output Parameter:
220: . dim - the spatial dimension
222: Level: intermediate
224: .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
225: @*/
226: PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
227: {
230: *dim = ctx->dim;
231: return(0);
232: }
234: /*@C
235: DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context
237: Not collective
239: Input Parameters:
240: + ctx - the context
241: - dof - the number of fields
243: Level: intermediate
245: .seealso: DMInterpolationGetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
246: @*/
247: PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
248: {
250: if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %D", dof);
251: ctx->dof = dof;
252: return(0);
253: }
255: /*@C
256: DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context
258: Not collective
260: Input Parameter:
261: . ctx - the context
263: Output Parameter:
264: . dof - the number of fields
266: Level: intermediate
268: .seealso: DMInterpolationSetDof(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
269: @*/
270: PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
271: {
274: *dof = ctx->dof;
275: return(0);
276: }
278: /*@C
279: DMInterpolationAddPoints - Add points at which we will interpolate the fields
281: Not collective
283: Input Parameters:
284: + ctx - the context
285: . n - the number of points
286: - points - the coordinates for each point, an array of size n * dim
288: Note: The coordinate information is copied.
290: Level: intermediate
292: .seealso: DMInterpolationSetDim(), DMInterpolationEvaluate(), DMInterpolationCreate()
293: @*/
294: PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
295: {
299: if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
300: if (ctx->points) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
301: ctx->nInput = n;
303: PetscMalloc1(n*ctx->dim, &ctx->points);
304: PetscArraycpy(ctx->points, points, n*ctx->dim);
305: return(0);
306: }
308: /*@C
309: DMInterpolationSetUp - Compute spatial indices for point location during interpolation
311: Collective on ctx
313: Input Parameters:
314: + ctx - the context
315: . dm - the DM for the function space used for interpolation
316: . redundantPoints - If PETSC_TRUE, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes.
317: - ignoreOutsideDomain - If PETSC_TRUE, ignore points outside the domain, otherwise return an error
319: Level: intermediate
321: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
322: @*/
323: PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain)
324: {
325: MPI_Comm comm = ctx->comm;
326: PetscScalar *a;
327: PetscInt p, q, i;
328: PetscMPIInt rank, size;
329: PetscErrorCode ierr;
330: Vec pointVec;
331: PetscSF cellSF;
332: PetscLayout layout;
333: PetscReal *globalPoints;
334: PetscScalar *globalPointsScalar;
335: const PetscInt *ranges;
336: PetscMPIInt *counts, *displs;
337: const PetscSFNode *foundCells;
338: const PetscInt *foundPoints;
339: PetscMPIInt *foundProcs, *globalProcs;
340: PetscInt n, N, numFound;
344: MPI_Comm_size(comm, &size);
345: MPI_Comm_rank(comm, &rank);
346: if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
347: /* Locate points */
348: n = ctx->nInput;
349: if (!redundantPoints) {
350: PetscLayoutCreate(comm, &layout);
351: PetscLayoutSetBlockSize(layout, 1);
352: PetscLayoutSetLocalSize(layout, n);
353: PetscLayoutSetUp(layout);
354: PetscLayoutGetSize(layout, &N);
355: /* Communicate all points to all processes */
356: PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);
357: PetscLayoutGetRanges(layout, &ranges);
358: for (p = 0; p < size; ++p) {
359: counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
360: displs[p] = ranges[p]*ctx->dim;
361: }
362: MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);
363: } else {
364: N = n;
365: globalPoints = ctx->points;
366: counts = displs = NULL;
367: layout = NULL;
368: }
369: #if 0
370: PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);
371: /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
372: #else
373: #if defined(PETSC_USE_COMPLEX)
374: PetscMalloc1(N*ctx->dim,&globalPointsScalar);
375: for (i=0; i<N*ctx->dim; i++) globalPointsScalar[i] = globalPoints[i];
376: #else
377: globalPointsScalar = globalPoints;
378: #endif
379: VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);
380: PetscMalloc2(N,&foundProcs,N,&globalProcs);
381: for (p = 0; p < N; ++p) {foundProcs[p] = size;}
382: cellSF = NULL;
383: DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF);
384: PetscSFGetGraph(cellSF,NULL,&numFound,&foundPoints,&foundCells);
385: #endif
386: for (p = 0; p < numFound; ++p) {
387: if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank;
388: }
389: /* Let the lowest rank process own each point */
390: MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);
391: ctx->n = 0;
392: for (p = 0; p < N; ++p) {
393: if (globalProcs[p] == size) {
394: if (!ignoreOutsideDomain) 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));
395: else if (rank == 0) ++ctx->n;
396: } else if (globalProcs[p] == rank) ++ctx->n;
397: }
398: /* Create coordinates vector and array of owned cells */
399: PetscMalloc1(ctx->n, &ctx->cells);
400: VecCreate(comm, &ctx->coords);
401: VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);
402: VecSetBlockSize(ctx->coords, ctx->dim);
403: VecSetType(ctx->coords,VECSTANDARD);
404: VecGetArray(ctx->coords, &a);
405: for (p = 0, q = 0, i = 0; p < N; ++p) {
406: if (globalProcs[p] == rank) {
407: PetscInt d;
409: for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
410: ctx->cells[q] = foundCells[q].index;
411: ++q;
412: }
413: if (globalProcs[p] == size && rank == 0) {
414: PetscInt d;
416: for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.;
417: ctx->cells[q] = -1;
418: ++q;
419: }
420: }
421: VecRestoreArray(ctx->coords, &a);
422: #if 0
423: PetscFree3(foundCells,foundProcs,globalProcs);
424: #else
425: PetscFree2(foundProcs,globalProcs);
426: PetscSFDestroy(&cellSF);
427: VecDestroy(&pointVec);
428: #endif
429: if ((void*)globalPointsScalar != (void*)globalPoints) {PetscFree(globalPointsScalar);}
430: if (!redundantPoints) {PetscFree3(globalPoints,counts,displs);}
431: PetscLayoutDestroy(&layout);
432: return(0);
433: }
435: /*@C
436: DMInterpolationGetCoordinates - Gets a Vec with the coordinates of each interpolation point
438: Collective on ctx
440: Input Parameter:
441: . ctx - the context
443: Output Parameter:
444: . coordinates - the coordinates of interpolation points
446: 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.
448: Level: intermediate
450: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
451: @*/
452: PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
453: {
456: if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
457: *coordinates = ctx->coords;
458: return(0);
459: }
461: /*@C
462: DMInterpolationGetVector - Gets a Vec which can hold all the interpolated field values
464: Collective on ctx
466: Input Parameter:
467: . ctx - the context
469: Output Parameter:
470: . v - a vector capable of holding the interpolated field values
472: Note: This vector should be returned using DMInterpolationRestoreVector().
474: Level: intermediate
476: .seealso: DMInterpolationRestoreVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
477: @*/
478: PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
479: {
484: if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
485: VecCreate(ctx->comm, v);
486: VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);
487: VecSetBlockSize(*v, ctx->dof);
488: VecSetType(*v,VECSTANDARD);
489: return(0);
490: }
492: /*@C
493: DMInterpolationRestoreVector - Returns a Vec which can hold all the interpolated field values
495: Collective on ctx
497: Input Parameters:
498: + ctx - the context
499: - v - a vector capable of holding the interpolated field values
501: Level: intermediate
503: .seealso: DMInterpolationGetVector(), DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
504: @*/
505: PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
506: {
511: if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
512: VecDestroy(v);
513: return(0);
514: }
516: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
517: {
518: PetscReal *v0, *J, *invJ, detJ;
519: const PetscScalar *coords;
520: PetscScalar *a;
521: PetscInt p;
525: PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);
526: VecGetArrayRead(ctx->coords, &coords);
527: VecGetArray(v, &a);
528: for (p = 0; p < ctx->n; ++p) {
529: PetscInt c = ctx->cells[p];
530: PetscScalar *x = NULL;
531: PetscReal xi[4];
532: PetscInt d, f, comp;
534: DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);
535: if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
536: DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);
537: for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
539: for (d = 0; d < ctx->dim; ++d) {
540: xi[d] = 0.0;
541: for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
542: 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];
543: }
544: DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);
545: }
546: VecRestoreArray(v, &a);
547: VecRestoreArrayRead(ctx->coords, &coords);
548: PetscFree3(v0, J, invJ);
549: return(0);
550: }
552: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
553: {
554: PetscReal *v0, *J, *invJ, detJ;
555: const PetscScalar *coords;
556: PetscScalar *a;
557: PetscInt p;
561: PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);
562: VecGetArrayRead(ctx->coords, &coords);
563: VecGetArray(v, &a);
564: for (p = 0; p < ctx->n; ++p) {
565: PetscInt c = ctx->cells[p];
566: const PetscInt order[3] = {2, 1, 3};
567: PetscScalar *x = NULL;
568: PetscReal xi[4];
569: PetscInt d, f, comp;
571: DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);
572: if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D", (double)detJ, c);
573: DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);
574: for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
576: for (d = 0; d < ctx->dim; ++d) {
577: xi[d] = 0.0;
578: for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
579: 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];
580: }
581: DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);
582: }
583: VecRestoreArray(v, &a);
584: VecRestoreArrayRead(ctx->coords, &coords);
585: PetscFree3(v0, J, invJ);
586: return(0);
587: }
589: PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
590: {
591: const PetscScalar *vertices = (const PetscScalar*) ctx;
592: const PetscScalar x0 = vertices[0];
593: const PetscScalar y0 = vertices[1];
594: const PetscScalar x1 = vertices[2];
595: const PetscScalar y1 = vertices[3];
596: const PetscScalar x2 = vertices[4];
597: const PetscScalar y2 = vertices[5];
598: const PetscScalar x3 = vertices[6];
599: const PetscScalar y3 = vertices[7];
600: const PetscScalar f_1 = x1 - x0;
601: const PetscScalar g_1 = y1 - y0;
602: const PetscScalar f_3 = x3 - x0;
603: const PetscScalar g_3 = y3 - y0;
604: const PetscScalar f_01 = x2 - x1 - x3 + x0;
605: const PetscScalar g_01 = y2 - y1 - y3 + y0;
606: const PetscScalar *ref;
607: PetscScalar *real;
608: PetscErrorCode ierr;
611: VecGetArrayRead(Xref, &ref);
612: VecGetArray(Xreal, &real);
613: {
614: const PetscScalar p0 = ref[0];
615: const PetscScalar p1 = ref[1];
617: real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
618: real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
619: }
620: PetscLogFlops(28);
621: VecRestoreArrayRead(Xref, &ref);
622: VecRestoreArray(Xreal, &real);
623: return(0);
624: }
626: #include <petsc/private/dmimpl.h>
627: PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
628: {
629: const PetscScalar *vertices = (const PetscScalar*) ctx;
630: const PetscScalar x0 = vertices[0];
631: const PetscScalar y0 = vertices[1];
632: const PetscScalar x1 = vertices[2];
633: const PetscScalar y1 = vertices[3];
634: const PetscScalar x2 = vertices[4];
635: const PetscScalar y2 = vertices[5];
636: const PetscScalar x3 = vertices[6];
637: const PetscScalar y3 = vertices[7];
638: const PetscScalar f_01 = x2 - x1 - x3 + x0;
639: const PetscScalar g_01 = y2 - y1 - y3 + y0;
640: const PetscScalar *ref;
641: PetscErrorCode ierr;
644: VecGetArrayRead(Xref, &ref);
645: {
646: const PetscScalar x = ref[0];
647: const PetscScalar y = ref[1];
648: const PetscInt rows[2] = {0, 1};
649: PetscScalar values[4];
651: values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
652: values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
653: MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);
654: }
655: PetscLogFlops(30);
656: VecRestoreArrayRead(Xref, &ref);
657: MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
658: MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
659: return(0);
660: }
662: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
663: {
664: DM dmCoord;
665: PetscFE fem = NULL;
666: SNES snes;
667: KSP ksp;
668: PC pc;
669: Vec coordsLocal, r, ref, real;
670: Mat J;
671: PetscTabulation T;
672: const PetscScalar *coords;
673: PetscScalar *a;
674: PetscReal xir[2];
675: PetscInt Nf, p;
676: const PetscInt dof = ctx->dof;
680: DMGetNumFields(dm, &Nf);
681: if (Nf) {DMGetField(dm, 0, NULL, (PetscObject *) &fem);}
682: DMGetCoordinatesLocal(dm, &coordsLocal);
683: DMGetCoordinateDM(dm, &dmCoord);
684: SNESCreate(PETSC_COMM_SELF, &snes);
685: SNESSetOptionsPrefix(snes, "quad_interp_");
686: VecCreate(PETSC_COMM_SELF, &r);
687: VecSetSizes(r, 2, 2);
688: VecSetType(r,dm->vectype);
689: VecDuplicate(r, &ref);
690: VecDuplicate(r, &real);
691: MatCreate(PETSC_COMM_SELF, &J);
692: MatSetSizes(J, 2, 2, 2, 2);
693: MatSetType(J, MATSEQDENSE);
694: MatSetUp(J);
695: SNESSetFunction(snes, r, QuadMap_Private, NULL);
696: SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);
697: SNESGetKSP(snes, &ksp);
698: KSPGetPC(ksp, &pc);
699: PCSetType(pc, PCLU);
700: SNESSetFromOptions(snes);
702: VecGetArrayRead(ctx->coords, &coords);
703: VecGetArray(v, &a);
704: PetscFECreateTabulation(fem, 1, 1, xir, 0, &T);
705: for (p = 0; p < ctx->n; ++p) {
706: PetscScalar *x = NULL, *vertices = NULL;
707: PetscScalar *xi;
708: PetscInt c = ctx->cells[p], comp, coordSize, xSize;
710: /* Can make this do all points at once */
711: DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
712: if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 4*2);
713: DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);
714: SNESSetFunction(snes, NULL, NULL, vertices);
715: SNESSetJacobian(snes, NULL, NULL, NULL, vertices);
716: VecGetArray(real, &xi);
717: xi[0] = coords[p*ctx->dim+0];
718: xi[1] = coords[p*ctx->dim+1];
719: VecRestoreArray(real, &xi);
720: SNESSolve(snes, real, ref);
721: VecGetArray(ref, &xi);
722: xir[0] = PetscRealPart(xi[0]);
723: xir[1] = PetscRealPart(xi[1]);
724: if (4*dof != xSize) {
725: PetscInt d;
727: xir[0] = 2.0*xir[0] - 1.0; xir[1] = 2.0*xir[1] - 1.0;
728: PetscFEComputeTabulation(fem, 1, xir, 0, T);
729: for (comp = 0; comp < dof; ++comp) {
730: a[p*dof+comp] = 0.0;
731: for (d = 0; d < xSize/dof; ++d) {
732: a[p*dof+comp] += x[d*dof+comp]*T->T[0][d*dof+comp];
733: }
734: }
735: } else {
736: for (comp = 0; comp < dof; ++comp)
737: 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];
738: }
739: VecRestoreArray(ref, &xi);
740: DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
741: DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);
742: }
743: PetscTabulationDestroy(&T);
744: VecRestoreArray(v, &a);
745: VecRestoreArrayRead(ctx->coords, &coords);
747: SNESDestroy(&snes);
748: VecDestroy(&r);
749: VecDestroy(&ref);
750: VecDestroy(&real);
751: MatDestroy(&J);
752: return(0);
753: }
755: PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
756: {
757: const PetscScalar *vertices = (const PetscScalar*) ctx;
758: const PetscScalar x0 = vertices[0];
759: const PetscScalar y0 = vertices[1];
760: const PetscScalar z0 = vertices[2];
761: const PetscScalar x1 = vertices[9];
762: const PetscScalar y1 = vertices[10];
763: const PetscScalar z1 = vertices[11];
764: const PetscScalar x2 = vertices[6];
765: const PetscScalar y2 = vertices[7];
766: const PetscScalar z2 = vertices[8];
767: const PetscScalar x3 = vertices[3];
768: const PetscScalar y3 = vertices[4];
769: const PetscScalar z3 = vertices[5];
770: const PetscScalar x4 = vertices[12];
771: const PetscScalar y4 = vertices[13];
772: const PetscScalar z4 = vertices[14];
773: const PetscScalar x5 = vertices[15];
774: const PetscScalar y5 = vertices[16];
775: const PetscScalar z5 = vertices[17];
776: const PetscScalar x6 = vertices[18];
777: const PetscScalar y6 = vertices[19];
778: const PetscScalar z6 = vertices[20];
779: const PetscScalar x7 = vertices[21];
780: const PetscScalar y7 = vertices[22];
781: const PetscScalar z7 = vertices[23];
782: const PetscScalar f_1 = x1 - x0;
783: const PetscScalar g_1 = y1 - y0;
784: const PetscScalar h_1 = z1 - z0;
785: const PetscScalar f_3 = x3 - x0;
786: const PetscScalar g_3 = y3 - y0;
787: const PetscScalar h_3 = z3 - z0;
788: const PetscScalar f_4 = x4 - x0;
789: const PetscScalar g_4 = y4 - y0;
790: const PetscScalar h_4 = z4 - z0;
791: const PetscScalar f_01 = x2 - x1 - x3 + x0;
792: const PetscScalar g_01 = y2 - y1 - y3 + y0;
793: const PetscScalar h_01 = z2 - z1 - z3 + z0;
794: const PetscScalar f_12 = x7 - x3 - x4 + x0;
795: const PetscScalar g_12 = y7 - y3 - y4 + y0;
796: const PetscScalar h_12 = z7 - z3 - z4 + z0;
797: const PetscScalar f_02 = x5 - x1 - x4 + x0;
798: const PetscScalar g_02 = y5 - y1 - y4 + y0;
799: const PetscScalar h_02 = z5 - z1 - z4 + z0;
800: const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
801: const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
802: const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
803: const PetscScalar *ref;
804: PetscScalar *real;
805: PetscErrorCode ierr;
808: VecGetArrayRead(Xref, &ref);
809: VecGetArray(Xreal, &real);
810: {
811: const PetscScalar p0 = ref[0];
812: const PetscScalar p1 = ref[1];
813: const PetscScalar p2 = ref[2];
815: 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;
816: 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;
817: 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;
818: }
819: PetscLogFlops(114);
820: VecRestoreArrayRead(Xref, &ref);
821: VecRestoreArray(Xreal, &real);
822: return(0);
823: }
825: PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
826: {
827: const PetscScalar *vertices = (const PetscScalar*) ctx;
828: const PetscScalar x0 = vertices[0];
829: const PetscScalar y0 = vertices[1];
830: const PetscScalar z0 = vertices[2];
831: const PetscScalar x1 = vertices[9];
832: const PetscScalar y1 = vertices[10];
833: const PetscScalar z1 = vertices[11];
834: const PetscScalar x2 = vertices[6];
835: const PetscScalar y2 = vertices[7];
836: const PetscScalar z2 = vertices[8];
837: const PetscScalar x3 = vertices[3];
838: const PetscScalar y3 = vertices[4];
839: const PetscScalar z3 = vertices[5];
840: const PetscScalar x4 = vertices[12];
841: const PetscScalar y4 = vertices[13];
842: const PetscScalar z4 = vertices[14];
843: const PetscScalar x5 = vertices[15];
844: const PetscScalar y5 = vertices[16];
845: const PetscScalar z5 = vertices[17];
846: const PetscScalar x6 = vertices[18];
847: const PetscScalar y6 = vertices[19];
848: const PetscScalar z6 = vertices[20];
849: const PetscScalar x7 = vertices[21];
850: const PetscScalar y7 = vertices[22];
851: const PetscScalar z7 = vertices[23];
852: const PetscScalar f_xy = x2 - x1 - x3 + x0;
853: const PetscScalar g_xy = y2 - y1 - y3 + y0;
854: const PetscScalar h_xy = z2 - z1 - z3 + z0;
855: const PetscScalar f_yz = x7 - x3 - x4 + x0;
856: const PetscScalar g_yz = y7 - y3 - y4 + y0;
857: const PetscScalar h_yz = z7 - z3 - z4 + z0;
858: const PetscScalar f_xz = x5 - x1 - x4 + x0;
859: const PetscScalar g_xz = y5 - y1 - y4 + y0;
860: const PetscScalar h_xz = z5 - z1 - z4 + z0;
861: const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
862: const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
863: const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
864: const PetscScalar *ref;
865: PetscErrorCode ierr;
868: VecGetArrayRead(Xref, &ref);
869: {
870: const PetscScalar x = ref[0];
871: const PetscScalar y = ref[1];
872: const PetscScalar z = ref[2];
873: const PetscInt rows[3] = {0, 1, 2};
874: PetscScalar values[9];
876: values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
877: values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
878: values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
879: values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
880: values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
881: values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
882: values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
883: values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
884: values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
886: MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);
887: }
888: PetscLogFlops(152);
889: VecRestoreArrayRead(Xref, &ref);
890: MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);
891: MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);
892: return(0);
893: }
895: PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
896: {
897: DM dmCoord;
898: SNES snes;
899: KSP ksp;
900: PC pc;
901: Vec coordsLocal, r, ref, real;
902: Mat J;
903: const PetscScalar *coords;
904: PetscScalar *a;
905: PetscInt p;
909: DMGetCoordinatesLocal(dm, &coordsLocal);
910: DMGetCoordinateDM(dm, &dmCoord);
911: SNESCreate(PETSC_COMM_SELF, &snes);
912: SNESSetOptionsPrefix(snes, "hex_interp_");
913: VecCreate(PETSC_COMM_SELF, &r);
914: VecSetSizes(r, 3, 3);
915: VecSetType(r,dm->vectype);
916: VecDuplicate(r, &ref);
917: VecDuplicate(r, &real);
918: MatCreate(PETSC_COMM_SELF, &J);
919: MatSetSizes(J, 3, 3, 3, 3);
920: MatSetType(J, MATSEQDENSE);
921: MatSetUp(J);
922: SNESSetFunction(snes, r, HexMap_Private, NULL);
923: SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);
924: SNESGetKSP(snes, &ksp);
925: KSPGetPC(ksp, &pc);
926: PCSetType(pc, PCLU);
927: SNESSetFromOptions(snes);
929: VecGetArrayRead(ctx->coords, &coords);
930: VecGetArray(v, &a);
931: for (p = 0; p < ctx->n; ++p) {
932: PetscScalar *x = NULL, *vertices = NULL;
933: PetscScalar *xi;
934: PetscReal xir[3];
935: PetscInt c = ctx->cells[p], comp, coordSize, xSize;
937: /* Can make this do all points at once */
938: DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
939: if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %d", coordSize, 8*3);
940: DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);
941: if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %D should be %D", xSize, 8*ctx->dof);
942: SNESSetFunction(snes, NULL, NULL, vertices);
943: SNESSetJacobian(snes, NULL, NULL, NULL, vertices);
944: VecGetArray(real, &xi);
945: xi[0] = coords[p*ctx->dim+0];
946: xi[1] = coords[p*ctx->dim+1];
947: xi[2] = coords[p*ctx->dim+2];
948: VecRestoreArray(real, &xi);
949: SNESSolve(snes, real, ref);
950: VecGetArray(ref, &xi);
951: xir[0] = PetscRealPart(xi[0]);
952: xir[1] = PetscRealPart(xi[1]);
953: xir[2] = PetscRealPart(xi[2]);
954: for (comp = 0; comp < ctx->dof; ++comp) {
955: a[p*ctx->dof+comp] =
956: x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
957: x[3*ctx->dof+comp]* xir[0]*(1-xir[1])*(1-xir[2]) +
958: x[2*ctx->dof+comp]* xir[0]* xir[1]*(1-xir[2]) +
959: x[1*ctx->dof+comp]*(1-xir[0])* xir[1]*(1-xir[2]) +
960: x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])* xir[2] +
961: x[5*ctx->dof+comp]* xir[0]*(1-xir[1])* xir[2] +
962: x[6*ctx->dof+comp]* xir[0]* xir[1]* xir[2] +
963: x[7*ctx->dof+comp]*(1-xir[0])* xir[1]* xir[2];
964: }
965: VecRestoreArray(ref, &xi);
966: DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);
967: DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);
968: }
969: VecRestoreArray(v, &a);
970: VecRestoreArrayRead(ctx->coords, &coords);
972: SNESDestroy(&snes);
973: VecDestroy(&r);
974: VecDestroy(&ref);
975: VecDestroy(&real);
976: MatDestroy(&J);
977: return(0);
978: }
980: /*@C
981: DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points.
983: Input Parameters:
984: + ctx - The DMInterpolationInfo context
985: . dm - The DM
986: - x - The local vector containing the field to be interpolated
988: Output Parameters:
989: . v - The vector containing the interpolated values
991: Note: A suitable v can be obtained using DMInterpolationGetVector().
993: Level: beginner
995: .seealso: DMInterpolationGetVector(), DMInterpolationAddPoints(), DMInterpolationCreate()
996: @*/
997: PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
998: {
999: PetscDS ds;
1000: PetscInt n, p, Nf, field;
1001: PetscBool useDS = PETSC_FALSE;
1008: VecGetLocalSize(v, &n);
1009: 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);
1010: if (!n) return(0);
1011: DMGetDS(dm, &ds);
1012: if (ds) {
1013: useDS = PETSC_TRUE;
1014: PetscDSGetNumFields(ds, &Nf);
1015: for (field = 0; field < Nf; ++field) {
1016: PetscObject obj;
1017: PetscClassId id;
1019: PetscDSGetDiscretization(ds, field, &obj);
1020: PetscObjectGetClassId(obj, &id);
1021: if (id != PETSCFE_CLASSID) {useDS = PETSC_FALSE; break;}
1022: }
1023: }
1024: if (useDS) {
1025: const PetscScalar *coords;
1026: PetscScalar *interpolant;
1027: PetscInt cdim, d;
1029: DMGetCoordinateDim(dm, &cdim);
1030: VecGetArrayRead(ctx->coords, &coords);
1031: VecGetArrayWrite(v, &interpolant);
1032: for (p = 0; p < ctx->n; ++p) {
1033: PetscReal pcoords[3], xi[3];
1034: PetscScalar *xa = NULL;
1035: PetscInt coff = 0, foff = 0, clSize;
1037: if (ctx->cells[p] < 0) continue;
1038: for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p*cdim+d]);
1039: DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi);
1040: DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa);
1041: for (field = 0; field < Nf; ++field) {
1042: PetscTabulation T;
1043: PetscFE fe;
1045: PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);
1046: PetscFECreateTabulation(fe, 1, 1, xi, 0, &T);
1047: {
1048: const PetscReal *basis = T->T[0];
1049: const PetscInt Nb = T->Nb;
1050: const PetscInt Nc = T->Nc;
1051: PetscInt f, fc;
1053: for (fc = 0; fc < Nc; ++fc) {
1054: interpolant[p*ctx->dof+coff+fc] = 0.0;
1055: for (f = 0; f < Nb; ++f) {
1056: interpolant[p*ctx->dof+coff+fc] += xa[foff+f]*basis[(0*Nb + f)*Nc + fc];
1057: }
1058: }
1059: coff += Nc;
1060: foff += Nb;
1061: }
1062: PetscTabulationDestroy(&T);
1063: }
1064: DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa);
1065: if (coff != ctx->dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %D != %D dof specified for interpolation", coff, ctx->dof);
1066: if (foff != clSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE space size %D != %D closure size", foff, clSize);
1067: }
1068: VecRestoreArrayRead(ctx->coords, &coords);
1069: VecRestoreArrayWrite(v, &interpolant);
1070: } else {
1071: DMPolytopeType ct;
1073: /* TODO Check each cell individually */
1074: DMPlexGetCellType(dm, ctx->cells[0], &ct);
1075: switch (ct) {
1076: case DM_POLYTOPE_TRIANGLE: DMInterpolate_Triangle_Private(ctx, dm, x, v);break;
1077: case DM_POLYTOPE_QUADRILATERAL: DMInterpolate_Quad_Private(ctx, dm, x, v);break;
1078: case DM_POLYTOPE_TETRAHEDRON: DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);break;
1079: case DM_POLYTOPE_HEXAHEDRON: DMInterpolate_Hex_Private(ctx, dm, x, v);break;
1080: default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[ct]);
1081: }
1082: }
1083: return(0);
1084: }
1086: /*@C
1087: DMInterpolationDestroy - Destroys a DMInterpolationInfo context
1089: Collective on ctx
1091: Input Parameter:
1092: . ctx - the context
1094: Level: beginner
1096: .seealso: DMInterpolationEvaluate(), DMInterpolationAddPoints(), DMInterpolationCreate()
1097: @*/
1098: PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
1099: {
1104: VecDestroy(&(*ctx)->coords);
1105: PetscFree((*ctx)->points);
1106: PetscFree((*ctx)->cells);
1107: PetscFree(*ctx);
1108: *ctx = NULL;
1109: return(0);
1110: }
1112: /*@C
1113: SNESMonitorFields - Monitors the residual for each field separately
1115: Collective on SNES
1117: Input Parameters:
1118: + snes - the SNES context
1119: . its - iteration number
1120: . fgnorm - 2-norm of residual
1121: - vf - PetscViewerAndFormat of type ASCII
1123: Notes:
1124: This routine prints the residual norm at each iteration.
1126: Level: intermediate
1128: .seealso: SNESMonitorSet(), SNESMonitorDefault()
1129: @*/
1130: PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf)
1131: {
1132: PetscViewer viewer = vf->viewer;
1133: Vec res;
1134: DM dm;
1135: PetscSection s;
1136: const PetscScalar *r;
1137: PetscReal *lnorms, *norms;
1138: PetscInt numFields, f, pStart, pEnd, p;
1139: PetscErrorCode ierr;
1143: SNESGetFunction(snes, &res, NULL, NULL);
1144: SNESGetDM(snes, &dm);
1145: DMGetLocalSection(dm, &s);
1146: PetscSectionGetNumFields(s, &numFields);
1147: PetscSectionGetChart(s, &pStart, &pEnd);
1148: PetscCalloc2(numFields, &lnorms, numFields, &norms);
1149: VecGetArrayRead(res, &r);
1150: for (p = pStart; p < pEnd; ++p) {
1151: for (f = 0; f < numFields; ++f) {
1152: PetscInt fdof, foff, d;
1154: PetscSectionGetFieldDof(s, p, f, &fdof);
1155: PetscSectionGetFieldOffset(s, p, f, &foff);
1156: for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
1157: }
1158: }
1159: VecRestoreArrayRead(res, &r);
1160: MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));
1161: PetscViewerPushFormat(viewer,vf->format);
1162: PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);
1163: PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);
1164: for (f = 0; f < numFields; ++f) {
1165: if (f > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
1166: PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));
1167: }
1168: PetscViewerASCIIPrintf(viewer, "]\n");
1169: PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);
1170: PetscViewerPopFormat(viewer);
1171: PetscFree2(lnorms, norms);
1172: return(0);
1173: }
1175: /********************* Residual Computation **************************/
1177: PetscErrorCode DMPlexGetAllCells_Internal(DM plex, IS *cellIS)
1178: {
1179: PetscInt depth;
1183: DMPlexGetDepth(plex, &depth);
1184: DMGetStratumIS(plex, "dim", depth, cellIS);
1185: if (!*cellIS) {DMGetStratumIS(plex, "depth", depth, cellIS);}
1186: return(0);
1187: }
1189: /*@
1190: DMPlexSNESComputeResidualFEM - Sums the local residual into vector F from the local input X using pointwise functions specified by the user
1192: Input Parameters:
1193: + dm - The mesh
1194: . X - Local solution
1195: - user - The user context
1197: Output Parameter:
1198: . F - Local output vector
1200: Notes:
1201: The residual is summed into F; the caller is responsible for using VecZeroEntries() or otherwise ensuring that any data in F is intentional.
1203: Level: developer
1205: .seealso: DMPlexComputeJacobianAction()
1206: @*/
1207: PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user)
1208: {
1209: DM plex;
1210: IS allcellIS;
1211: PetscInt Nds, s;
1215: DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
1216: DMPlexGetAllCells_Internal(plex, &allcellIS);
1217: DMGetNumDS(dm, &Nds);
1218: for (s = 0; s < Nds; ++s) {
1219: PetscDS ds;
1220: IS cellIS;
1221: PetscFormKey key;
1223: DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);
1224: key.value = 0;
1225: key.field = 0;
1226: key.part = 0;
1227: if (!key.label) {
1228: PetscObjectReference((PetscObject) allcellIS);
1229: cellIS = allcellIS;
1230: } else {
1231: IS pointIS;
1233: key.value = 1;
1234: DMLabelGetStratumIS(key.label, key.value, &pointIS);
1235: ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);
1236: ISDestroy(&pointIS);
1237: }
1238: DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);
1239: ISDestroy(&cellIS);
1240: }
1241: ISDestroy(&allcellIS);
1242: DMDestroy(&plex);
1243: return(0);
1244: }
1246: PetscErrorCode DMSNESComputeResidual(DM dm, Vec X, Vec F, void *user)
1247: {
1248: DM plex;
1249: IS allcellIS;
1250: PetscInt Nds, s;
1254: DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
1255: DMPlexGetAllCells_Internal(plex, &allcellIS);
1256: DMGetNumDS(dm, &Nds);
1257: for (s = 0; s < Nds; ++s) {
1258: PetscDS ds;
1259: DMLabel label;
1260: IS cellIS;
1262: DMGetRegionNumDS(dm, s, &label, NULL, &ds);
1263: {
1264: PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1};
1265: PetscWeakForm wf;
1266: PetscInt Nm = 2, m, Nk = 0, k, kp, off = 0;
1267: PetscFormKey *reskeys;
1269: /* Get unique residual keys */
1270: for (m = 0; m < Nm; ++m) {
1271: PetscInt Nkm;
1272: PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm);
1273: Nk += Nkm;
1274: }
1275: PetscMalloc1(Nk, &reskeys);
1276: for (m = 0; m < Nm; ++m) {
1277: PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys);
1278: }
1279: if (off != Nk) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %D should be %D", off, Nk);
1280: PetscFormKeySort(Nk, reskeys);
1281: for (k = 0, kp = 1; kp < Nk; ++kp) {
1282: if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) {
1283: ++k;
1284: if (kp != k) reskeys[k] = reskeys[kp];
1285: }
1286: }
1287: Nk = k;
1289: PetscDSGetWeakForm(ds, &wf);
1290: for (k = 0; k < Nk; ++k) {
1291: DMLabel label = reskeys[k].label;
1292: PetscInt val = reskeys[k].value;
1294: if (!label) {
1295: PetscObjectReference((PetscObject) allcellIS);
1296: cellIS = allcellIS;
1297: } else {
1298: IS pointIS;
1300: DMLabelGetStratumIS(label, val, &pointIS);
1301: ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);
1302: ISDestroy(&pointIS);
1303: }
1304: DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user);
1305: ISDestroy(&cellIS);
1306: }
1307: PetscFree(reskeys);
1308: }
1309: }
1310: ISDestroy(&allcellIS);
1311: DMDestroy(&plex);
1312: return(0);
1313: }
1315: /*@
1316: DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input X
1318: Input Parameters:
1319: + dm - The mesh
1320: - user - The user context
1322: Output Parameter:
1323: . X - Local solution
1325: Level: developer
1327: .seealso: DMPlexComputeJacobianAction()
1328: @*/
1329: PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user)
1330: {
1331: DM plex;
1335: DMSNESConvertPlex(dm,&plex,PETSC_TRUE);
1336: DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL);
1337: DMDestroy(&plex);
1338: return(0);
1339: }
1341: /*@
1342: DMSNESComputeJacobianAction - Compute the action of the Jacobian J(X) on Y
1344: Input Parameters:
1345: + dm - The DM
1346: . X - Local solution vector
1347: . Y - Local input vector
1348: - user - The user context
1350: Output Parameter:
1351: . F - lcoal output vector
1353: Level: developer
1355: Notes:
1356: Users will typically use DMSNESCreateJacobianMF() followed by MatMult() instead of calling this routine directly.
1358: .seealso: DMSNESCreateJacobianMF(), DMPlexSNESComputeResidualFEM()
1359: @*/
1360: PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user)
1361: {
1362: DM plex;
1363: IS allcellIS;
1364: PetscInt Nds, s;
1368: DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
1369: DMPlexGetAllCells_Internal(plex, &allcellIS);
1370: DMGetNumDS(dm, &Nds);
1371: for (s = 0; s < Nds; ++s) {
1372: PetscDS ds;
1373: DMLabel label;
1374: IS cellIS;
1376: DMGetRegionNumDS(dm, s, &label, NULL, &ds);
1377: {
1378: PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3};
1379: PetscWeakForm wf;
1380: PetscInt Nm = 4, m, Nk = 0, k, kp, off = 0;
1381: PetscFormKey *jackeys;
1383: /* Get unique Jacobian keys */
1384: for (m = 0; m < Nm; ++m) {
1385: PetscInt Nkm;
1386: PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm);
1387: Nk += Nkm;
1388: }
1389: PetscMalloc1(Nk, &jackeys);
1390: for (m = 0; m < Nm; ++m) {
1391: PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys);
1392: }
1393: if (off != Nk) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %D should be %D", off, Nk);
1394: PetscFormKeySort(Nk, jackeys);
1395: for (k = 0, kp = 1; kp < Nk; ++kp) {
1396: if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) {
1397: ++k;
1398: if (kp != k) jackeys[k] = jackeys[kp];
1399: }
1400: }
1401: Nk = k;
1403: PetscDSGetWeakForm(ds, &wf);
1404: for (k = 0; k < Nk; ++k) {
1405: DMLabel label = jackeys[k].label;
1406: PetscInt val = jackeys[k].value;
1408: if (!label) {
1409: PetscObjectReference((PetscObject) allcellIS);
1410: cellIS = allcellIS;
1411: } else {
1412: IS pointIS;
1414: DMLabelGetStratumIS(label, val, &pointIS);
1415: ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);
1416: ISDestroy(&pointIS);
1417: }
1418: DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user);
1419: ISDestroy(&cellIS);
1420: }
1421: PetscFree(jackeys);
1422: }
1423: }
1424: ISDestroy(&allcellIS);
1425: DMDestroy(&plex);
1426: return(0);
1427: }
1429: /*@
1430: DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user.
1432: Input Parameters:
1433: + dm - The mesh
1434: . X - Local input vector
1435: - user - The user context
1437: Output Parameter:
1438: . Jac - Jacobian matrix
1440: Note:
1441: We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1442: like a GPU, or vectorize on a multicore machine.
1444: Level: developer
1446: .seealso: FormFunctionLocal()
1447: @*/
1448: PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user)
1449: {
1450: DM plex;
1451: IS allcellIS;
1452: PetscBool hasJac, hasPrec;
1453: PetscInt Nds, s;
1457: DMSNESConvertPlex(dm, &plex, PETSC_TRUE);
1458: DMPlexGetAllCells_Internal(plex, &allcellIS);
1459: DMGetNumDS(dm, &Nds);
1460: for (s = 0; s < Nds; ++s) {
1461: PetscDS ds;
1462: IS cellIS;
1463: PetscFormKey key;
1465: DMGetRegionNumDS(dm, s, &key.label, NULL, &ds);
1466: key.value = 0;
1467: key.field = 0;
1468: key.part = 0;
1469: if (!key.label) {
1470: PetscObjectReference((PetscObject) allcellIS);
1471: cellIS = allcellIS;
1472: } else {
1473: IS pointIS;
1475: key.value = 1;
1476: DMLabelGetStratumIS(key.label, key.value, &pointIS);
1477: ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS);
1478: ISDestroy(&pointIS);
1479: }
1480: if (!s) {
1481: PetscDSHasJacobian(ds, &hasJac);
1482: PetscDSHasJacobianPreconditioner(ds, &hasPrec);
1483: if (hasJac && hasPrec) {MatZeroEntries(Jac);}
1484: MatZeroEntries(JacP);
1485: }
1486: DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user);
1487: ISDestroy(&cellIS);
1488: }
1489: ISDestroy(&allcellIS);
1490: DMDestroy(&plex);
1491: return(0);
1492: }
1494: struct _DMSNESJacobianMFCtx
1495: {
1496: DM dm;
1497: Vec X;
1498: void *ctx;
1499: };
1501: static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A)
1502: {
1503: struct _DMSNESJacobianMFCtx *ctx;
1504: PetscErrorCode ierr;
1507: MatShellGetContext(A, &ctx);
1508: MatShellSetContext(A, NULL);
1509: DMDestroy(&ctx->dm);
1510: VecDestroy(&ctx->X);
1511: PetscFree(ctx);
1512: return(0);
1513: }
1515: static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z)
1516: {
1517: struct _DMSNESJacobianMFCtx *ctx;
1518: PetscErrorCode ierr;
1521: MatShellGetContext(A, &ctx);
1522: DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx);
1523: return(0);
1524: }
1526: /*@
1527: DMSNESCreateJacobianMF - Create a Mat which computes the action of the Jacobian matrix-free
1529: Collective on dm
1531: Input Parameters:
1532: + dm - The DM
1533: . X - The evaluation point for the Jacobian
1534: - user - A user context, or NULL
1536: Output Parameter:
1537: . J - The Mat
1539: Level: advanced
1541: Notes:
1542: Vec X is kept in Mat J, so updating X then updates the evaluation point.
1544: .seealso: DMSNESComputeJacobianAction()
1545: @*/
1546: PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J)
1547: {
1548: struct _DMSNESJacobianMFCtx *ctx;
1549: PetscInt n, N;
1550: PetscErrorCode ierr;
1553: MatCreate(PetscObjectComm((PetscObject) dm), J);
1554: MatSetType(*J, MATSHELL);
1555: VecGetLocalSize(X, &n);
1556: VecGetSize(X, &N);
1557: MatSetSizes(*J, n, n, N, N);
1558: PetscObjectReference((PetscObject) dm);
1559: PetscObjectReference((PetscObject) X);
1560: PetscMalloc1(1, &ctx);
1561: ctx->dm = dm;
1562: ctx->X = X;
1563: ctx->ctx = user;
1564: MatShellSetContext(*J, ctx);
1565: MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void)) DMSNESJacobianMF_Destroy_Private);
1566: MatShellSetOperation(*J, MATOP_MULT, (void (*)(void)) DMSNESJacobianMF_Mult_Private);
1567: return(0);
1568: }
1570: /*
1571: MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context.
1573: Input Parameters:
1574: + X - SNES linearization point
1575: . ovl - index set of overlapping subdomains
1577: Output Parameter:
1578: . J - unassembled (Neumann) local matrix
1580: Level: intermediate
1582: .seealso: DMCreateNeumannOverlap(), MATIS, PCHPDDMSetAuxiliaryMat()
1583: */
1584: static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx)
1585: {
1586: SNES snes;
1587: Mat pJ;
1588: DM ovldm,origdm;
1589: DMSNES sdm;
1590: PetscErrorCode (*bfun)(DM,Vec,void*);
1591: PetscErrorCode (*jfun)(DM,Vec,Mat,Mat,void*);
1592: void *bctx,*jctx;
1596: PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_MATIS",(PetscObject*)&pJ);
1597: if (!pJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing overlapping Mat");
1598: PetscObjectQuery((PetscObject)ovl,"_DM_Original_HPDDM",(PetscObject*)&origdm);
1599: if (!origdm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing original DM");
1600: MatGetDM(pJ,&ovldm);
1601: DMSNESGetBoundaryLocal(origdm,&bfun,&bctx);
1602: DMSNESSetBoundaryLocal(ovldm,bfun,bctx);
1603: DMSNESGetJacobianLocal(origdm,&jfun,&jctx);
1604: DMSNESSetJacobianLocal(ovldm,jfun,jctx);
1605: PetscObjectQuery((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject*)&snes);
1606: if (!snes) {
1607: SNESCreate(PetscObjectComm((PetscObject)ovl),&snes);
1608: SNESSetDM(snes,ovldm);
1609: PetscObjectCompose((PetscObject)ovl,"_DM_Overlap_HPDDM_SNES",(PetscObject)snes);
1610: PetscObjectDereference((PetscObject)snes);
1611: }
1612: DMGetDMSNES(ovldm,&sdm);
1613: VecLockReadPush(X);
1614: PetscStackPush("SNES user Jacobian function");
1615: (*sdm->ops->computejacobian)(snes,X,pJ,pJ,sdm->jacobianctx);
1616: PetscStackPop;
1617: VecLockReadPop(X);
1618: /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */
1619: {
1620: Mat locpJ;
1622: MatISGetLocalMat(pJ,&locpJ);
1623: MatCopy(locpJ,J,SAME_NONZERO_PATTERN);
1624: }
1625: return(0);
1626: }
1628: /*@
1629: DMPlexSetSNESLocalFEM - Use DMPlex's internal FEM routines to compute SNES boundary values, residual, and Jacobian.
1631: Input Parameters:
1632: + dm - The DM object
1633: . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see PetscDSAddBoundary())
1634: . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see PetscDSSetResidual())
1635: - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see PetscDSSetJacobian())
1637: Level: developer
1638: @*/
1639: PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx)
1640: {
1644: DMSNESSetBoundaryLocal(dm,DMPlexSNESComputeBoundaryFEM,boundaryctx);
1645: DMSNESSetFunctionLocal(dm,DMPlexSNESComputeResidualFEM,residualctx);
1646: DMSNESSetJacobianLocal(dm,DMPlexSNESComputeJacobianFEM,jacobianctx);
1647: PetscObjectComposeFunction((PetscObject)dm,"MatComputeNeumannOverlap_C",MatComputeNeumannOverlap_Plex);
1648: return(0);
1649: }
1651: /*@C
1652: DMSNESCheckDiscretization - Check the discretization error of the exact solution
1654: Input Parameters:
1655: + snes - the SNES object
1656: . dm - the DM
1657: . t - the time
1658: . u - a DM vector
1659: - tol - A tolerance for the check, or -1 to print the results instead
1661: Output Parameters:
1662: . error - An array which holds the discretization error in each field, or NULL
1664: Note: The user must call PetscDSSetExactSolution() beforehand
1666: Level: developer
1668: .seealso: DNSNESCheckFromOptions(), DMSNESCheckResidual(), DMSNESCheckJacobian(), PetscDSSetExactSolution()
1669: @*/
1670: PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[])
1671: {
1672: PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
1673: void **ectxs;
1674: PetscReal *err;
1675: MPI_Comm comm;
1676: PetscInt Nf, f;
1677: PetscErrorCode ierr;
1685: DMComputeExactSolution(dm, t, u, NULL);
1686: VecViewFromOptions(u, NULL, "-vec_view");
1688: PetscObjectGetComm((PetscObject) snes, &comm);
1689: DMGetNumFields(dm, &Nf);
1690: PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err);
1691: {
1692: PetscInt Nds, s;
1694: DMGetNumDS(dm, &Nds);
1695: for (s = 0; s < Nds; ++s) {
1696: PetscDS ds;
1697: DMLabel label;
1698: IS fieldIS;
1699: const PetscInt *fields;
1700: PetscInt dsNf, f;
1702: DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds);
1703: PetscDSGetNumFields(ds, &dsNf);
1704: ISGetIndices(fieldIS, &fields);
1705: for (f = 0; f < dsNf; ++f) {
1706: const PetscInt field = fields[f];
1707: PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]);
1708: }
1709: ISRestoreIndices(fieldIS, &fields);
1710: }
1711: }
1712: if (Nf > 1) {
1713: DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err);
1714: if (tol >= 0.0) {
1715: for (f = 0; f < Nf; ++f) {
1716: 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);
1717: }
1718: } else if (error) {
1719: for (f = 0; f < Nf; ++f) error[f] = err[f];
1720: } else {
1721: PetscPrintf(comm, "L_2 Error: [");
1722: for (f = 0; f < Nf; ++f) {
1723: if (f) {PetscPrintf(comm, ", ");}
1724: PetscPrintf(comm, "%g", (double)err[f]);
1725: }
1726: PetscPrintf(comm, "]\n");
1727: }
1728: } else {
1729: DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0]);
1730: if (tol >= 0.0) {
1731: if (err[0] > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double) err[0], (double) tol);
1732: } else if (error) {
1733: error[0] = err[0];
1734: } else {
1735: PetscPrintf(comm, "L_2 Error: %g\n", (double) err[0]);
1736: }
1737: }
1738: PetscFree3(exacts, ectxs, err);
1739: return(0);
1740: }
1742: /*@C
1743: DMSNESCheckResidual - Check the residual of the exact solution
1745: Input Parameters:
1746: + snes - the SNES object
1747: . dm - the DM
1748: . u - a DM vector
1749: - tol - A tolerance for the check, or -1 to print the results instead
1751: Output Parameters:
1752: . residual - The residual norm of the exact solution, or NULL
1754: Level: developer
1756: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckJacobian()
1757: @*/
1758: PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual)
1759: {
1760: MPI_Comm comm;
1761: Vec r;
1762: PetscReal res;
1770: PetscObjectGetComm((PetscObject) snes, &comm);
1771: DMComputeExactSolution(dm, 0.0, u, NULL);
1772: VecDuplicate(u, &r);
1773: SNESComputeFunction(snes, u, r);
1774: VecNorm(r, NORM_2, &res);
1775: if (tol >= 0.0) {
1776: if (res > tol) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double) res, (double) tol);
1777: } else if (residual) {
1778: *residual = res;
1779: } else {
1780: PetscPrintf(comm, "L_2 Residual: %g\n", (double)res);
1781: VecChop(r, 1.0e-10);
1782: PetscObjectSetName((PetscObject) r, "Initial Residual");
1783: PetscObjectSetOptionsPrefix((PetscObject)r,"res_");
1784: VecViewFromOptions(r, NULL, "-vec_view");
1785: }
1786: VecDestroy(&r);
1787: return(0);
1788: }
1790: /*@C
1791: DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test
1793: Input Parameters:
1794: + snes - the SNES object
1795: . dm - the DM
1796: . u - a DM vector
1797: - tol - A tolerance for the check, or -1 to print the results instead
1799: Output Parameters:
1800: + isLinear - Flag indicaing that the function looks linear, or NULL
1801: - convRate - The rate of convergence of the linear model, or NULL
1803: Level: developer
1805: .seealso: DNSNESCheckFromOptions(), DMSNESCheckDiscretization(), DMSNESCheckResidual()
1806: @*/
1807: PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate)
1808: {
1809: MPI_Comm comm;
1810: PetscDS ds;
1811: Mat J, M;
1812: MatNullSpace nullspace;
1813: PetscReal slope, intercept;
1814: PetscBool hasJac, hasPrec, isLin = PETSC_FALSE;
1823: PetscObjectGetComm((PetscObject) snes, &comm);
1824: DMComputeExactSolution(dm, 0.0, u, NULL);
1825: /* Create and view matrices */
1826: DMCreateMatrix(dm, &J);
1827: DMGetDS(dm, &ds);
1828: PetscDSHasJacobian(ds, &hasJac);
1829: PetscDSHasJacobianPreconditioner(ds, &hasPrec);
1830: if (hasJac && hasPrec) {
1831: DMCreateMatrix(dm, &M);
1832: SNESComputeJacobian(snes, u, J, M);
1833: PetscObjectSetName((PetscObject) M, "Preconditioning Matrix");
1834: PetscObjectSetOptionsPrefix((PetscObject) M, "jacpre_");
1835: MatViewFromOptions(M, NULL, "-mat_view");
1836: MatDestroy(&M);
1837: } else {
1838: SNESComputeJacobian(snes, u, J, J);
1839: }
1840: PetscObjectSetName((PetscObject) J, "Jacobian");
1841: PetscObjectSetOptionsPrefix((PetscObject) J, "jac_");
1842: MatViewFromOptions(J, NULL, "-mat_view");
1843: /* Check nullspace */
1844: MatGetNullSpace(J, &nullspace);
1845: if (nullspace) {
1846: PetscBool isNull;
1847: MatNullSpaceTest(nullspace, J, &isNull);
1848: if (!isNull) SETERRQ(comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid.");
1849: }
1850: /* Taylor test */
1851: {
1852: PetscRandom rand;
1853: Vec du, uhat, r, rhat, df;
1854: PetscReal h;
1855: PetscReal *es, *hs, *errors;
1856: PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1;
1857: PetscInt Nv, v;
1859: /* Choose a perturbation direction */
1860: PetscRandomCreate(comm, &rand);
1861: VecDuplicate(u, &du);
1862: VecSetRandom(du, rand);
1863: PetscRandomDestroy(&rand);
1864: VecDuplicate(u, &df);
1865: MatMult(J, du, df);
1866: /* Evaluate residual at u, F(u), save in vector r */
1867: VecDuplicate(u, &r);
1868: SNESComputeFunction(snes, u, r);
1869: /* Look at the convergence of our Taylor approximation as we approach u */
1870: for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv);
1871: PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors);
1872: VecDuplicate(u, &uhat);
1873: VecDuplicate(u, &rhat);
1874: for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) {
1875: VecWAXPY(uhat, h, du, u);
1876: /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */
1877: SNESComputeFunction(snes, uhat, rhat);
1878: VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df);
1879: VecNorm(rhat, NORM_2, &errors[Nv]);
1881: es[Nv] = PetscLog10Real(errors[Nv]);
1882: hs[Nv] = PetscLog10Real(h);
1883: }
1884: VecDestroy(&uhat);
1885: VecDestroy(&rhat);
1886: VecDestroy(&df);
1887: VecDestroy(&r);
1888: VecDestroy(&du);
1889: for (v = 0; v < Nv; ++v) {
1890: if ((tol >= 0) && (errors[v] > tol)) break;
1891: else if (errors[v] > PETSC_SMALL) break;
1892: }
1893: if (v == Nv) isLin = PETSC_TRUE;
1894: PetscLinearRegression(Nv, hs, es, &slope, &intercept);
1895: PetscFree3(es, hs, errors);
1896: /* Slope should be about 2 */
1897: if (tol >= 0) {
1898: if (!isLin && PetscAbsReal(2 - slope) > tol) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double) slope);
1899: } else if (isLinear || convRate) {
1900: if (isLinear) *isLinear = isLin;
1901: if (convRate) *convRate = slope;
1902: } else {
1903: if (!isLin) {PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double) slope);}
1904: else {PetscPrintf(comm, "Function appears to be linear\n");}
1905: }
1906: }
1907: MatDestroy(&J);
1908: return(0);
1909: }
1911: PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u)
1912: {
1916: DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL);
1917: DMSNESCheckResidual(snes, dm, u, -1.0, NULL);
1918: DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL);
1919: return(0);
1920: }
1922: /*@C
1923: DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information
1925: Input Parameters:
1926: + snes - the SNES object
1927: - u - representative SNES vector
1929: Note: The user must call PetscDSSetExactSolution() beforehand
1931: Level: developer
1932: @*/
1933: PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u)
1934: {
1935: DM dm;
1936: Vec sol;
1937: PetscBool check;
1941: PetscOptionsHasName(((PetscObject)snes)->options,((PetscObject)snes)->prefix, "-dmsnes_check", &check);
1942: if (!check) return(0);
1943: SNESGetDM(snes, &dm);
1944: VecDuplicate(u, &sol);
1945: SNESSetSolution(snes, sol);
1946: DMSNESCheck_Internal(snes, dm, sol);
1947: VecDestroy(&sol);
1948: return(0);
1949: }