Actual source code: dmmeshsnes.c
petsc-3.3-p5 2012-12-01
1: #include <petsc-private/meshimpl.h> /*I "petscdmmesh.h" I*/
2: #include <petscsnes.h> /*I "petscsnes.h" I*/
6: PetscErrorCode DMMeshInterpolationCreate(DM dm, DMMeshInterpolationInfo *ctx) {
12: PetscMalloc(sizeof(struct _DMMeshInterpolationInfo), ctx);
13: (*ctx)->dim = -1;
14: (*ctx)->nInput = 0;
15: (*ctx)->points = PETSC_NULL;
16: (*ctx)->cells = PETSC_NULL;
17: (*ctx)->n = -1;
18: (*ctx)->coords = PETSC_NULL;
19: return(0);
20: }
24: PetscErrorCode DMMeshInterpolationSetDim(DM dm, PetscInt dim, DMMeshInterpolationInfo ctx) {
27: if ((dim < 1) || (dim > 3)) {SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);}
28: ctx->dim = dim;
29: return(0);
30: }
34: PetscErrorCode DMMeshInterpolationGetDim(DM dm, PetscInt *dim, DMMeshInterpolationInfo ctx) {
38: *dim = ctx->dim;
39: return(0);
40: }
44: PetscErrorCode DMMeshInterpolationSetDof(DM dm, PetscInt dof, DMMeshInterpolationInfo ctx) {
47: if (dof < 1) {SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);}
48: ctx->dof = dof;
49: return(0);
50: }
54: PetscErrorCode DMMeshInterpolationGetDof(DM dm, PetscInt *dof, DMMeshInterpolationInfo ctx) {
58: *dof = ctx->dof;
59: return(0);
60: }
64: PetscErrorCode DMMeshInterpolationAddPoints(DM dm, PetscInt n, PetscReal points[], DMMeshInterpolationInfo ctx) {
69: if (ctx->dim < 0) {
70: SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
71: }
72: if (ctx->points) {
73: SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
74: }
75: ctx->nInput = n;
76: PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);
77: PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));
78: return(0);
79: }
83: PetscErrorCode DMMeshInterpolationSetUp(DM dm, DMMeshInterpolationInfo ctx, PetscBool redundantPoints) {
84: ALE::Obj<PETSC_MESH_TYPE> m;
85: MPI_Comm comm = ((PetscObject) dm)->comm;
86: PetscScalar *a;
87: PetscInt p, q, i;
88: PetscMPIInt rank, size;
93: MPI_Comm_size(comm, &size);
94: MPI_Comm_rank(comm, &rank);
95: DMMeshGetMesh(dm, m);
96: if (ctx->dim < 0) {
97: SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
98: }
99: // Locate points
100: PetscLayout layout;
101: PetscReal *globalPoints;
102: const PetscInt *ranges;
103: PetscMPIInt *counts, *displs;
104: PetscInt *foundCells;
105: PetscMPIInt *foundProcs, *globalProcs;
106: PetscInt n = ctx->nInput, N;
108: if (!redundantPoints) {
109: PetscLayoutCreate(comm, &layout);
110: PetscLayoutSetBlockSize(layout, 1);
111: PetscLayoutSetLocalSize(layout, n);
112: PetscLayoutSetUp(layout);
113: PetscLayoutGetSize(layout, &N);
114: /* Communicate all points to all processes */
115: PetscMalloc3(N*ctx->dim,PetscReal,&globalPoints,size,PetscMPIInt,&counts,size,PetscMPIInt,&displs);
116: PetscLayoutGetRanges(layout, &ranges);
117: for(p = 0; p < size; ++p) {
118: counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
119: displs[p] = ranges[p]*ctx->dim;
120: }
121: MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);
122: } else {
123: N = n;
124: globalPoints = ctx->points;
125: }
126: PetscMalloc3(N,PetscInt,&foundCells,N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);
127: for(p = 0; p < N; ++p) {
128: foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]);
129: if (foundCells[p] >= 0) {
130: foundProcs[p] = rank;
131: } else {
132: foundProcs[p] = size;
133: }
134: }
135: /* Let the lowest rank process own each point */
136: MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);
137: ctx->n = 0;
138: for(p = 0; p < N; ++p) {
139: if (globalProcs[p] == size) {
140: SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %g %g %g not located in mesh", p, globalPoints[p*ctx->dim+0], ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0, ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0);
141: } else if (globalProcs[p] == rank) {
142: ctx->n++;
143: }
144: }
145: /* Create coordinates vector and array of owned cells */
146: PetscMalloc(ctx->n * sizeof(PetscInt), &ctx->cells);
147: VecCreate(comm, &ctx->coords);
148: VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);
149: VecSetBlockSize(ctx->coords, ctx->dim);
150: VecSetFromOptions(ctx->coords);
151: VecGetArray(ctx->coords, &a);
152: for(p = 0, q = 0, i = 0; p < N; ++p) {
153: if (globalProcs[p] == rank) {
154: PetscInt d;
156: for(d = 0; d < ctx->dim; ++d, ++i) {
157: a[i] = globalPoints[p*ctx->dim+d];
158: }
159: ctx->cells[q++] = foundCells[p];
160: }
161: }
162: VecRestoreArray(ctx->coords, &a);
163: PetscFree3(foundCells,foundProcs,globalProcs);
164: if (!redundantPoints) {
165: PetscFree3(globalPoints,counts,displs);
166: PetscLayoutDestroy(&layout);
167: }
168: return(0);
169: }
173: PetscErrorCode DMMeshInterpolationGetCoordinates(DM dm, Vec *coordinates, DMMeshInterpolationInfo ctx) {
177: if (!ctx->coords) {SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
178: *coordinates = ctx->coords;
179: return(0);
180: }
184: PetscErrorCode DMMeshInterpolationGetVector(DM dm, Vec *v, DMMeshInterpolationInfo ctx) {
190: if (!ctx->coords) {SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
191: VecCreate(((PetscObject) dm)->comm, v);
192: VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);
193: VecSetBlockSize(*v, ctx->dof);
194: VecSetFromOptions(*v);
195: return(0);
196: }
200: PetscErrorCode DMMeshInterpolationRestoreVector(DM dm, Vec *v, DMMeshInterpolationInfo ctx) {
206: if (!ctx->coords) {SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");}
207: VecDestroy(v);
208: return(0);
209: }
213: PetscErrorCode DMMeshInterpolate_Simplex_Private(DM dm, SectionReal x, Vec v, DMMeshInterpolationInfo ctx) {
214: #ifdef PETSC_HAVE_SIEVE
215: ALE::Obj<PETSC_MESH_TYPE> m;
216: ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
217: PetscInt p;
221: DMMeshGetMesh(dm, m);
222: SectionRealGetSection(x, s);
223: const ALE::Obj<PETSC_MESH_TYPE::real_section_type>& coordinates = m->getRealSection("coordinates");
224: PetscReal *v0, *J, *invJ, detJ;
225: PetscScalar *a, *coords;
227: PetscMalloc3(ctx->dim,PetscReal,&v0,ctx->dim*ctx->dim,PetscReal,&J,ctx->dim*ctx->dim,PetscReal,&invJ);
228: VecGetArray(ctx->coords, &coords);
229: VecGetArray(v, &a);
230: for(p = 0; p < ctx->n; ++p) {
231: PetscInt e = ctx->cells[p];
232: PetscReal xi[4];
233: PetscInt d, f, comp;
235: if ((ctx->dim+1)*ctx->dof != m->sizeWithBC(s, e)) {SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_ARG_SIZ, "Invalid restrict size %d should be %d", m->sizeWithBC(s, e), (ctx->dim+1)*ctx->dof);}
236: m->computeElementGeometry(coordinates, e, v0, J, invJ, detJ);
237: const PetscScalar *c = m->restrictClosure(s, e); /* Must come after geom, since it uses closure temp space*/
238: for(comp = 0; comp < ctx->dof; ++comp) {
239: a[p*ctx->dof+comp] = c[0*ctx->dof+comp];
240: }
241: for(d = 0; d < ctx->dim; ++d) {
242: xi[d] = 0.0;
243: for(f = 0; f < ctx->dim; ++f) {
244: xi[d] += invJ[d*ctx->dim+f]*0.5*(coords[p*ctx->dim+f] - v0[f]);
245: }
246: for(comp = 0; comp < ctx->dof; ++comp) {
247: a[p*ctx->dof+comp] += (c[(d+1)*ctx->dof+comp] - c[0*ctx->dof+comp])*xi[d];
248: }
249: }
250: }
251: VecRestoreArray(v, &a);
252: VecRestoreArray(ctx->coords, &coords);
253: PetscFree3(v0, J, invJ);
254: return(0);
255: #else
256: SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Interpolation only work with DMMesh currently.");
257: #endif
258: }
262: PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
263: {
264: const PetscScalar*vertices = (const PetscScalar *) ctx;
265: const PetscScalar x0 = vertices[0];
266: const PetscScalar y0 = vertices[1];
267: const PetscScalar x1 = vertices[2];
268: const PetscScalar y1 = vertices[3];
269: const PetscScalar x2 = vertices[4];
270: const PetscScalar y2 = vertices[5];
271: const PetscScalar x3 = vertices[6];
272: const PetscScalar y3 = vertices[7];
273: const PetscScalar f_1 = x1 - x0;
274: const PetscScalar g_1 = y1 - y0;
275: const PetscScalar f_3 = x3 - x0;
276: const PetscScalar g_3 = y3 - y0;
277: const PetscScalar f_01 = x2 - x1 - x3 + x0;
278: const PetscScalar g_01 = y2 - y1 - y3 + y0;
279: PetscScalar *ref, *real;
280: PetscErrorCode ierr;
283: VecGetArray(Xref, &ref);
284: VecGetArray(Xreal, &real);
285: {
286: const PetscScalar p0 = ref[0];
287: const PetscScalar p1 = ref[1];
289: real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
290: real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
291: }
292: PetscLogFlops(28);
293: VecRestoreArray(Xref, &ref);
294: VecRestoreArray(Xreal, &real);
295: return(0);
296: }
300: PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
301: {
302: const PetscScalar*vertices = (const PetscScalar *) ctx;
303: const PetscScalar x0 = vertices[0];
304: const PetscScalar y0 = vertices[1];
305: const PetscScalar x1 = vertices[2];
306: const PetscScalar y1 = vertices[3];
307: const PetscScalar x2 = vertices[4];
308: const PetscScalar y2 = vertices[5];
309: const PetscScalar x3 = vertices[6];
310: const PetscScalar y3 = vertices[7];
311: const PetscScalar f_01 = x2 - x1 - x3 + x0;
312: const PetscScalar g_01 = y2 - y1 - y3 + y0;
313: PetscScalar *ref;
314: PetscErrorCode ierr;
317: VecGetArray(Xref, &ref);
318: {
319: const PetscScalar x = ref[0];
320: const PetscScalar y = ref[1];
321: const PetscInt rows[2] = {0, 1};
322: const PetscScalar values[4] = {(x1 - x0 + f_01*y) * 0.5, (x3 - x0 + f_01*x) * 0.5,
323: (y1 - y0 + g_01*y) * 0.5, (y3 - y0 + g_01*x) * 0.5};
324: MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);
325: }
326: PetscLogFlops(30);
327: VecRestoreArray(Xref, &ref);
328: MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);
329: MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);
330: return(0);
331: }
335: PetscErrorCode DMMeshInterpolate_Quad_Private(DM dm, SectionReal x, Vec v, DMMeshInterpolationInfo ctx) {
336: #ifdef PETSC_HAVE_SIEVE
337: SNES snes;
338: KSP ksp;
339: PC pc;
340: Vec r, ref, real;
341: Mat J;
342: PetscScalar vertices[8];
344: ALE::Obj<PETSC_MESH_TYPE> m;
345: ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
346: PetscInt p;
350: SNESCreate(PETSC_COMM_SELF, &snes);
351: SNESSetOptionsPrefix(snes, "quad_interp_");
352: VecCreate(PETSC_COMM_SELF, &r);
353: VecSetSizes(r, 2, 2);
354: VecSetFromOptions(r);
355: VecDuplicate(r, &ref);
356: VecDuplicate(r, &real);
357: MatCreate(PETSC_COMM_SELF, &J);
358: MatSetSizes(J, 2, 2, 2, 2);
359: MatSetType(J, MATSEQDENSE);
360: MatSetUp(J);
361: SNESSetFunction(snes, r, QuadMap_Private, vertices);
362: SNESSetJacobian(snes, J, J, QuadJacobian_Private, vertices);
363: SNESGetKSP(snes, &ksp);
364: KSPGetPC(ksp, &pc);
365: PCSetType(pc, PCLU);
366: SNESSetFromOptions(snes);
368: DMMeshGetMesh(dm, m);
369: SectionRealGetSection(x, s);
370: const ALE::Obj<PETSC_MESH_TYPE::real_section_type>& coordinates = m->getRealSection("coordinates");
371: PetscScalar *a, *coords;
373: VecGetArray(ctx->coords, &coords);
374: VecGetArray(v, &a);
375: for(p = 0; p < ctx->n; ++p) {
376: PetscScalar *xi;
377: PetscInt e = ctx->cells[p], comp;
379: if (4*ctx->dof != m->sizeWithBC(s, e)) {SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_ARG_SIZ, "Invalid restrict size %d should be %d", m->sizeWithBC(s, e), 4*ctx->dof);}
380: /* Can make this do all points at once */
381: {
382: const PetscReal *v = m->restrictClosure(coordinates, e);
383: for(PetscInt i = 0; i < 8; ++i) vertices[i] = v[i];
384: }
385: const PetscScalar *c = m->restrictClosure(s, e); /* Must come after geom, since it uses closure temp space*/
386: VecGetArray(real, &xi);
387: xi[0] = coords[p*ctx->dim+0];
388: xi[1] = coords[p*ctx->dim+1];
389: VecRestoreArray(real, &xi);
390: SNESSolve(snes, real, ref);
391: VecGetArray(ref, &xi);
392: for(comp = 0; comp < ctx->dof; ++comp) {
393: a[p*ctx->dof+comp] = c[0*ctx->dof+comp]*(1 - xi[0])*(1 - xi[1]) + c[1*ctx->dof+comp]*xi[0]*(1 - xi[1]) + c[2*ctx->dof+comp]*xi[0]*xi[1] + c[3*ctx->dof+comp]*(1 - xi[0])*xi[1];
394: }
395: VecRestoreArray(ref, &xi);
396: }
397: VecRestoreArray(v, &a);
398: VecRestoreArray(ctx->coords, &coords);
400: SNESDestroy(&snes);
401: VecDestroy(&r);
402: VecDestroy(&ref);
403: VecDestroy(&real);
404: MatDestroy(&J);
405: return(0);
406: #else
407: SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Interpolation only work with DMMesh currently.");
408: #endif
409: }
413: PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
414: {
415: const PetscScalar*vertices = (const PetscScalar *) ctx;
416: const PetscScalar x0 = vertices[0];
417: const PetscScalar y0 = vertices[1];
418: const PetscScalar z0 = vertices[2];
419: const PetscScalar x1 = vertices[3];
420: const PetscScalar y1 = vertices[4];
421: const PetscScalar z1 = vertices[5];
422: const PetscScalar x2 = vertices[6];
423: const PetscScalar y2 = vertices[7];
424: const PetscScalar z2 = vertices[8];
425: const PetscScalar x3 = vertices[9];
426: const PetscScalar y3 = vertices[10];
427: const PetscScalar z3 = vertices[11];
428: const PetscScalar x4 = vertices[12];
429: const PetscScalar y4 = vertices[13];
430: const PetscScalar z4 = vertices[14];
431: const PetscScalar x5 = vertices[15];
432: const PetscScalar y5 = vertices[16];
433: const PetscScalar z5 = vertices[17];
434: const PetscScalar x6 = vertices[18];
435: const PetscScalar y6 = vertices[19];
436: const PetscScalar z6 = vertices[20];
437: const PetscScalar x7 = vertices[21];
438: const PetscScalar y7 = vertices[22];
439: const PetscScalar z7 = vertices[23];
440: const PetscScalar f_1 = x1 - x0;
441: const PetscScalar g_1 = y1 - y0;
442: const PetscScalar h_1 = z1 - z0;
443: const PetscScalar f_3 = x3 - x0;
444: const PetscScalar g_3 = y3 - y0;
445: const PetscScalar h_3 = z3 - z0;
446: const PetscScalar f_4 = x4 - x0;
447: const PetscScalar g_4 = y4 - y0;
448: const PetscScalar h_4 = z4 - z0;
449: const PetscScalar f_01 = x2 - x1 - x3 + x0;
450: const PetscScalar g_01 = y2 - y1 - y3 + y0;
451: const PetscScalar h_01 = z2 - z1 - z3 + z0;
452: const PetscScalar f_12 = x7 - x3 - x4 + x0;
453: const PetscScalar g_12 = y7 - y3 - y4 + y0;
454: const PetscScalar h_12 = z7 - z3 - z4 + z0;
455: const PetscScalar f_02 = x5 - x1 - x4 + x0;
456: const PetscScalar g_02 = y5 - y1 - y4 + y0;
457: const PetscScalar h_02 = z5 - z1 - z4 + z0;
458: const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
459: const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
460: const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
461: PetscScalar *ref, *real;
462: PetscErrorCode ierr;
465: VecGetArray(Xref, &ref);
466: VecGetArray(Xreal, &real);
467: {
468: const PetscScalar p0 = ref[0];
469: const PetscScalar p1 = ref[1];
470: const PetscScalar p2 = ref[2];
472: 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;
473: 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;
474: 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;
475: }
476: PetscLogFlops(114);
477: VecRestoreArray(Xref, &ref);
478: VecRestoreArray(Xreal, &real);
479: return(0);
480: }
484: PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx)
485: {
486: const PetscScalar*vertices = (const PetscScalar *) ctx;
487: const PetscScalar x0 = vertices[0];
488: const PetscScalar y0 = vertices[1];
489: const PetscScalar z0 = vertices[2];
490: const PetscScalar x1 = vertices[3];
491: const PetscScalar y1 = vertices[4];
492: const PetscScalar z1 = vertices[5];
493: const PetscScalar x2 = vertices[6];
494: const PetscScalar y2 = vertices[7];
495: const PetscScalar z2 = vertices[8];
496: const PetscScalar x3 = vertices[9];
497: const PetscScalar y3 = vertices[10];
498: const PetscScalar z3 = vertices[11];
499: const PetscScalar x4 = vertices[12];
500: const PetscScalar y4 = vertices[13];
501: const PetscScalar z4 = vertices[14];
502: const PetscScalar x5 = vertices[15];
503: const PetscScalar y5 = vertices[16];
504: const PetscScalar z5 = vertices[17];
505: const PetscScalar x6 = vertices[18];
506: const PetscScalar y6 = vertices[19];
507: const PetscScalar z6 = vertices[20];
508: const PetscScalar x7 = vertices[21];
509: const PetscScalar y7 = vertices[22];
510: const PetscScalar z7 = vertices[23];
511: const PetscScalar f_xy = x2 - x1 - x3 + x0;
512: const PetscScalar g_xy = y2 - y1 - y3 + y0;
513: const PetscScalar h_xy = z2 - z1 - z3 + z0;
514: const PetscScalar f_yz = x7 - x3 - x4 + x0;
515: const PetscScalar g_yz = y7 - y3 - y4 + y0;
516: const PetscScalar h_yz = z7 - z3 - z4 + z0;
517: const PetscScalar f_xz = x5 - x1 - x4 + x0;
518: const PetscScalar g_xz = y5 - y1 - y4 + y0;
519: const PetscScalar h_xz = z5 - z1 - z4 + z0;
520: const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
521: const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
522: const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
523: PetscScalar *ref;
524: PetscErrorCode ierr;
527: VecGetArray(Xref, &ref);
528: {
529: const PetscScalar x = ref[0];
530: const PetscScalar y = ref[1];
531: const PetscScalar z = ref[2];
532: const PetscInt rows[3] = {0, 1, 2};
533: const PetscScalar values[9] = {
534: (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0,
535: (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0,
536: (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0,
537: (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0,
538: (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0,
539: (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0,
540: (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0,
541: (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0,
542: (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0};
543: MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);
544: }
545: PetscLogFlops(152);
546: VecRestoreArray(Xref, &ref);
547: MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);
548: MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);
549: return(0);
550: }
554: PetscErrorCode DMMeshInterpolate_Hex_Private(DM dm, SectionReal x, Vec v, DMMeshInterpolationInfo ctx) {
555: #ifdef PETSC_HAVE_SIEVE
556: SNES snes;
557: KSP ksp;
558: PC pc;
559: Vec r, ref, real;
560: Mat J;
561: PetscScalar vertices[24];
563: ALE::Obj<PETSC_MESH_TYPE> m;
564: ALE::Obj<PETSC_MESH_TYPE::real_section_type> s;
565: PetscInt p;
569: SNESCreate(PETSC_COMM_SELF, &snes);
570: SNESSetOptionsPrefix(snes, "hex_interp_");
571: VecCreate(PETSC_COMM_SELF, &r);
572: VecSetSizes(r, 3, 3);
573: VecSetFromOptions(r);
574: VecDuplicate(r, &ref);
575: VecDuplicate(r, &real);
576: MatCreate(PETSC_COMM_SELF, &J);
577: MatSetSizes(J, 3, 3, 3, 3);
578: MatSetType(J, MATSEQDENSE);
579: MatSetUp(J);
580: SNESSetFunction(snes, r, HexMap_Private, vertices);
581: SNESSetJacobian(snes, J, J, HexJacobian_Private, vertices);
582: SNESGetKSP(snes, &ksp);
583: KSPGetPC(ksp, &pc);
584: PCSetType(pc, PCLU);
585: SNESSetFromOptions(snes);
587: DMMeshGetMesh(dm, m);
588: SectionRealGetSection(x, s);
589: const ALE::Obj<PETSC_MESH_TYPE::real_section_type>& coordinates = m->getRealSection("coordinates");
590: PetscScalar *a, *coords;
592: VecGetArray(ctx->coords, &coords);
593: VecGetArray(v, &a);
594: for(p = 0; p < ctx->n; ++p) {
595: PetscScalar *xi;
596: PetscInt e = ctx->cells[p], comp;
598: if (8*ctx->dof != m->sizeWithBC(s, e)) {SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_ARG_SIZ, "Invalid restrict size %d should be %d", m->sizeWithBC(s, e), 8*ctx->dof);}
599: /* Can make this do all points at once */
600: {
601: const PetscReal *v = m->restrictClosure(coordinates, e);
602: for(PetscInt i = 0; i < 24; ++i) vertices[i] = v[i];
603: }
604: const PetscScalar *c = m->restrictClosure(s, e); /* Must come after geom, since it uses closure temp space*/
605: VecGetArray(real, &xi);
606: xi[0] = coords[p*ctx->dim+0];
607: xi[1] = coords[p*ctx->dim+1];
608: xi[2] = coords[p*ctx->dim+2];
609: VecRestoreArray(real, &xi);
610: SNESSolve(snes, real, ref);
611: VecGetArray(ref, &xi);
612: for(comp = 0; comp < ctx->dof; ++comp) {
613: a[p*ctx->dof+comp] =
614: c[0*ctx->dof+comp]*(1-xi[0])*(1-xi[1])*(1-xi[2]) +
615: c[1*ctx->dof+comp]* xi[0]*(1-xi[1])*(1-xi[2]) +
616: c[2*ctx->dof+comp]* xi[0]* xi[1]*(1-xi[2]) +
617: c[3*ctx->dof+comp]*(1-xi[0])* xi[1]*(1-xi[2]) +
618: c[4*ctx->dof+comp]*(1-xi[0])*(1-xi[1])* xi[2] +
619: c[5*ctx->dof+comp]* xi[0]*(1-xi[1])* xi[2] +
620: c[6*ctx->dof+comp]* xi[0]* xi[1]* xi[2] +
621: c[7*ctx->dof+comp]*(1-xi[0])* xi[1]* xi[2];
622: }
623: VecRestoreArray(ref, &xi);
624: }
625: VecRestoreArray(v, &a);
626: VecRestoreArray(ctx->coords, &coords);
628: SNESDestroy(&snes);
629: VecDestroy(&r);
630: VecDestroy(&ref);
631: VecDestroy(&real);
632: MatDestroy(&J);
633: return(0);
634: #else
635: SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Interpolation only work with DMMesh currently.");
636: #endif
637: }
641: PetscErrorCode DMMeshInterpolationEvaluate(DM dm, SectionReal x, Vec v, DMMeshInterpolationInfo ctx) {
642: PetscInt dim, coneSize, n;
649: VecGetLocalSize(v, &n);
650: if (n != ctx->n*ctx->dof) {SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %d should be %d", n, ctx->n*ctx->dof);}
651: if (n) {
652: DMMeshGetDimension(dm, &dim);
653: DMMeshGetConeSize(dm, ctx->cells[0], &coneSize);
654: if (dim == 2) {
655: if (coneSize == 3) {
656: DMMeshInterpolate_Simplex_Private(dm, x, v, ctx);
657: } else if (coneSize == 4) {
658: DMMeshInterpolate_Quad_Private(dm, x, v, ctx);
659: } else {
660: SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
661: }
662: } else if (dim == 3) {
663: if (coneSize == 4) {
664: DMMeshInterpolate_Simplex_Private(dm, x, v, ctx);
665: } else {
666: DMMeshInterpolate_Hex_Private(dm, x, v, ctx);
667: }
668: } else {
669: SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
670: }
671: }
672: return(0);
673: }
677: PetscErrorCode DMMeshInterpolationDestroy(DM dm, DMMeshInterpolationInfo *ctx) {
683: VecDestroy(&(*ctx)->coords);
684: PetscFree((*ctx)->points);
685: PetscFree((*ctx)->cells);
686: PetscFree(*ctx);
687: *ctx = PETSC_NULL;
688: return(0);
689: }