xref: /petsc/src/snes/tests/ex2.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Interpolation Tests for Plex\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscsnes.h>
4*c4762a1bSJed Brown #include <petscdmplex.h>
5*c4762a1bSJed Brown #include <petscdmda.h>
6*c4762a1bSJed Brown #include <petscds.h>
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown typedef enum {CENTROID, GRID, GRID_REPLICATED} PointType;
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown typedef struct {
11*c4762a1bSJed Brown   PetscInt      dim;                          /* The topological mesh dimension */
12*c4762a1bSJed Brown   PetscBool     cellSimplex;                  /* Use simplices or hexes */
13*c4762a1bSJed Brown   char          filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
14*c4762a1bSJed Brown   PointType     pointType;                    /* Point generation mechanism */
15*c4762a1bSJed Brown } AppCtx;
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown static PetscInt Nc = 3;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
20*c4762a1bSJed Brown {
21*c4762a1bSJed Brown   PetscInt d, c;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
24*c4762a1bSJed Brown     u[c] = 0.0;
25*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) u[c] += x[d];
26*c4762a1bSJed Brown   }
27*c4762a1bSJed Brown   return 0;
28*c4762a1bSJed Brown }
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
31*c4762a1bSJed Brown {
32*c4762a1bSJed Brown   const char    *pointTypes[3] = {"centroid", "grid", "grid_replicated"};
33*c4762a1bSJed Brown   PetscInt       pt;
34*c4762a1bSJed Brown   PetscErrorCode ierr;
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown   PetscFunctionBegin;
37*c4762a1bSJed Brown   options->dim           = 3;
38*c4762a1bSJed Brown   options->cellSimplex   = PETSC_TRUE;
39*c4762a1bSJed Brown   options->filename[0]   = '\0';
40*c4762a1bSJed Brown   options->pointType     = CENTROID;
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX");CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex2.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex2.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = PetscOptionsString("-filename", "The mesh file", "ex2.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
46*c4762a1bSJed Brown   pt   = options->pointType;
47*c4762a1bSJed Brown   ierr = PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL);CHKERRQ(ierr);
48*c4762a1bSJed Brown   options->pointType = (PointType) pt;
49*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   PetscFunctionReturn(0);
52*c4762a1bSJed Brown }
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
55*c4762a1bSJed Brown {
56*c4762a1bSJed Brown   PetscInt       dim         = ctx->dim;
57*c4762a1bSJed Brown   PetscBool      cellSimplex = ctx->cellSimplex;
58*c4762a1bSJed Brown   const char    *filename    = ctx->filename;
59*c4762a1bSJed Brown   const PetscInt cells[3]    = {1, 1, 1};
60*c4762a1bSJed Brown   size_t         len;
61*c4762a1bSJed Brown   PetscMPIInt    rank, size;
62*c4762a1bSJed Brown   PetscErrorCode ierr;
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown   PetscFunctionBegin;
65*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
68*c4762a1bSJed Brown   if (len) {ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr);}
69*c4762a1bSJed Brown   else     {ierr = DMPlexCreateBoxMesh(comm, dim, cellSimplex, cells, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);}
70*c4762a1bSJed Brown   {
71*c4762a1bSJed Brown     DM               distributedMesh = NULL;
72*c4762a1bSJed Brown     PetscPartitioner part;
73*c4762a1bSJed Brown 
74*c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
75*c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown     /* Distribute mesh over processes */
78*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
79*c4762a1bSJed Brown     if (distributedMesh) {
80*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
81*c4762a1bSJed Brown       *dm  = distributedMesh;
82*c4762a1bSJed Brown     }
83*c4762a1bSJed Brown   }
84*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
87*c4762a1bSJed Brown   PetscFunctionReturn(0);
88*c4762a1bSJed Brown }
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
91*c4762a1bSJed Brown {
92*c4762a1bSJed Brown   PetscSection   coordSection;
93*c4762a1bSJed Brown   Vec            coordsLocal;
94*c4762a1bSJed Brown   PetscInt       spaceDim, p;
95*c4762a1bSJed Brown   PetscMPIInt    rank;
96*c4762a1bSJed Brown   PetscErrorCode ierr;
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   PetscFunctionBegin;
99*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, NULL, Np);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = PetscCalloc1(*Np * spaceDim, pcoords);CHKERRQ(ierr);
105*c4762a1bSJed Brown   for (p = 0; p < *Np; ++p) {
106*c4762a1bSJed Brown     PetscScalar *coords = NULL;
107*c4762a1bSJed Brown     PetscInt     size, num, n, d;
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown     ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, p, &size, &coords);CHKERRQ(ierr);
110*c4762a1bSJed Brown     num  = size/spaceDim;
111*c4762a1bSJed Brown     for (n = 0; n < num; ++n) {
112*c4762a1bSJed Brown       for (d = 0; d < spaceDim; ++d) (*pcoords)[p*spaceDim+d] += PetscRealPart(coords[n*spaceDim+d]) / num;
113*c4762a1bSJed Brown     }
114*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, p);CHKERRQ(ierr);
115*c4762a1bSJed Brown     for (d = 0; d < spaceDim; ++d) {
116*c4762a1bSJed Brown       ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[p*spaceDim+d]);CHKERRQ(ierr);
117*c4762a1bSJed Brown       if (d < spaceDim-1) {ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);}
118*c4762a1bSJed Brown     }
119*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n");CHKERRQ(ierr);
120*c4762a1bSJed Brown     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, p, &num, &coords);CHKERRQ(ierr);
121*c4762a1bSJed Brown   }
122*c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr);
123*c4762a1bSJed Brown   *pointsAllProcs = PETSC_FALSE;
124*c4762a1bSJed Brown   PetscFunctionReturn(0);
125*c4762a1bSJed Brown }
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
128*c4762a1bSJed Brown {
129*c4762a1bSJed Brown   DM             da;
130*c4762a1bSJed Brown   DMDALocalInfo  info;
131*c4762a1bSJed Brown   PetscInt       N = 3, n = 0, spaceDim, i, j, k, *ind, d;
132*c4762a1bSJed Brown   PetscReal      *h;
133*c4762a1bSJed Brown   PetscMPIInt    rank;
134*c4762a1bSJed Brown   PetscErrorCode ierr;
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = PetscCalloc1(spaceDim,&ind);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = PetscCalloc1(spaceDim,&h);CHKERRQ(ierr);
140*c4762a1bSJed Brown   h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1);
141*c4762a1bSJed Brown   ierr = DMDACreate(PetscObjectComm((PetscObject) dm), &da);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = DMSetDimension(da, ctx->dim);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = DMDASetSizes(da, N, N, N);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = DMDASetDof(da, 1);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = DMDASetStencilWidth(da, 1);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr);
149*c4762a1bSJed Brown   *Np  = info.xm * info.ym * info.zm;
150*c4762a1bSJed Brown   ierr = PetscCalloc1(*Np * spaceDim, pcoords);CHKERRQ(ierr);
151*c4762a1bSJed Brown   for (k = info.zs; k < info.zs + info.zm; ++k) {
152*c4762a1bSJed Brown     ind[2] = k;
153*c4762a1bSJed Brown     for (j = info.ys; j < info.ys + info.ym; ++j) {
154*c4762a1bSJed Brown       ind[1] = j;
155*c4762a1bSJed Brown       for (i = info.xs; i < info.xs + info.xm; ++i, ++n) {
156*c4762a1bSJed Brown         ind[0] = i;
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) (*pcoords)[n*spaceDim+d] = ind[d]*h[d];
159*c4762a1bSJed Brown         ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n);CHKERRQ(ierr);
160*c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) {
161*c4762a1bSJed Brown           ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d]);CHKERRQ(ierr);
162*c4762a1bSJed Brown           if (d < spaceDim-1) {ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);}
163*c4762a1bSJed Brown         }
164*c4762a1bSJed Brown         ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n");CHKERRQ(ierr);
165*c4762a1bSJed Brown       }
166*c4762a1bSJed Brown     }
167*c4762a1bSJed Brown   }
168*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = PetscFree(ind);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = PetscFree(h);CHKERRQ(ierr);
172*c4762a1bSJed Brown   *pointsAllProcs = PETSC_FALSE;
173*c4762a1bSJed Brown   PetscFunctionReturn(0);
174*c4762a1bSJed Brown }
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
177*c4762a1bSJed Brown {
178*c4762a1bSJed Brown   PetscInt       N = 3, n = 0, spaceDim, i, j, k, *ind, d;
179*c4762a1bSJed Brown   PetscReal      *h;
180*c4762a1bSJed Brown   PetscMPIInt    rank;
181*c4762a1bSJed Brown   PetscErrorCode ierr;
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = PetscCalloc1(spaceDim,&ind);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = PetscCalloc1(spaceDim,&h);CHKERRQ(ierr);
187*c4762a1bSJed Brown   h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1);
188*c4762a1bSJed Brown   *Np  = N * (ctx->dim > 1 ? N : 1) * (ctx->dim > 2 ? N : 1);
189*c4762a1bSJed Brown   ierr = PetscCalloc1(*Np * spaceDim, pcoords);CHKERRQ(ierr);
190*c4762a1bSJed Brown   for (k = 0; k < N; ++k) {
191*c4762a1bSJed Brown     ind[2] = k;
192*c4762a1bSJed Brown     for (j = 0; j < N; ++j) {
193*c4762a1bSJed Brown       ind[1] = j;
194*c4762a1bSJed Brown       for (i = 0; i < N; ++i, ++n) {
195*c4762a1bSJed Brown         ind[0] = i;
196*c4762a1bSJed Brown 
197*c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) (*pcoords)[n*spaceDim+d] = ind[d]*h[d];
198*c4762a1bSJed Brown         ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n);CHKERRQ(ierr);
199*c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) {
200*c4762a1bSJed Brown           ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d]);CHKERRQ(ierr);
201*c4762a1bSJed Brown           if (d < spaceDim-1) {ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);}
202*c4762a1bSJed Brown         }
203*c4762a1bSJed Brown         ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n");CHKERRQ(ierr);
204*c4762a1bSJed Brown       }
205*c4762a1bSJed Brown     }
206*c4762a1bSJed Brown   }
207*c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr);
208*c4762a1bSJed Brown   *pointsAllProcs = PETSC_TRUE;
209*c4762a1bSJed Brown   ierr = PetscFree(ind);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = PetscFree(h);CHKERRQ(ierr);
211*c4762a1bSJed Brown   PetscFunctionReturn(0);
212*c4762a1bSJed Brown }
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
215*c4762a1bSJed Brown {
216*c4762a1bSJed Brown   PetscErrorCode ierr;
217*c4762a1bSJed Brown 
218*c4762a1bSJed Brown   PetscFunctionBegin;
219*c4762a1bSJed Brown   *pointsAllProcs = PETSC_FALSE;
220*c4762a1bSJed Brown   switch (ctx->pointType) {
221*c4762a1bSJed Brown   case CENTROID:        ierr = CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx);CHKERRQ(ierr);break;
222*c4762a1bSJed Brown   case GRID:            ierr = CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx);CHKERRQ(ierr);break;
223*c4762a1bSJed Brown   case GRID_REPLICATED: ierr = CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx);CHKERRQ(ierr);break;
224*c4762a1bSJed Brown   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int) ctx->pointType);
225*c4762a1bSJed Brown   }
226*c4762a1bSJed Brown   PetscFunctionReturn(0);
227*c4762a1bSJed Brown }
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown int main(int argc, char **argv)
230*c4762a1bSJed Brown {
231*c4762a1bSJed Brown   AppCtx              ctx;
232*c4762a1bSJed Brown   PetscErrorCode   (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
233*c4762a1bSJed Brown   DM                  dm;
234*c4762a1bSJed Brown   PetscFE             fe;
235*c4762a1bSJed Brown   DMInterpolationInfo interpolator;
236*c4762a1bSJed Brown   Vec                 lu, fieldVals;
237*c4762a1bSJed Brown   PetscScalar        *vals;
238*c4762a1bSJed Brown   const PetscScalar  *ivals, *vcoords;
239*c4762a1bSJed Brown   PetscReal          *pcoords;
240*c4762a1bSJed Brown   PetscBool           pointsAllProcs=PETSC_TRUE;
241*c4762a1bSJed Brown   PetscInt            spaceDim, c, Np, p;
242*c4762a1bSJed Brown   PetscMPIInt         rank, size;
243*c4762a1bSJed Brown   PetscViewer         selfviewer;
244*c4762a1bSJed Brown   PetscErrorCode      ierr;
245*c4762a1bSJed Brown 
246*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
247*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
248*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr);
249*c4762a1bSJed Brown   ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr);
250*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
251*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
252*c4762a1bSJed Brown   /* Create points */
253*c4762a1bSJed Brown   ierr = CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx);CHKERRQ(ierr);
254*c4762a1bSJed Brown   /* Create interpolator */
255*c4762a1bSJed Brown   ierr = DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator);CHKERRQ(ierr);
256*c4762a1bSJed Brown   ierr = DMInterpolationSetDim(interpolator, spaceDim);CHKERRQ(ierr);
257*c4762a1bSJed Brown   ierr = DMInterpolationAddPoints(interpolator, Np, pcoords);CHKERRQ(ierr);
258*c4762a1bSJed Brown   ierr = DMInterpolationSetUp(interpolator, dm, pointsAllProcs);CHKERRQ(ierr);
259*c4762a1bSJed Brown   /* Check locations */
260*c4762a1bSJed Brown   for (c = 0; c < interpolator->n; ++c) {
261*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D is in Cell %D\n", rank, c, interpolator->cells[c]);CHKERRQ(ierr);
262*c4762a1bSJed Brown   }
263*c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr);
264*c4762a1bSJed Brown   ierr = VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
265*c4762a1bSJed Brown   /* Setup Discretization */
266*c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), ctx.dim, Nc, ctx.cellSimplex, NULL, -1, &fe);CHKERRQ(ierr);
267*c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
270*c4762a1bSJed Brown   /* Create function */
271*c4762a1bSJed Brown   ierr = PetscCalloc2(Nc, &funcs, Nc, &vals);CHKERRQ(ierr);
272*c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) funcs[c] = linear;
273*c4762a1bSJed Brown   ierr = DMGetLocalVector(dm, &lu);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
276*c4762a1bSJed Brown   ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr);
277*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank);CHKERRQ(ierr);
278*c4762a1bSJed Brown   ierr = VecView(lu,selfviewer);CHKERRQ(ierr);
279*c4762a1bSJed Brown   ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
282*c4762a1bSJed Brown   /* Check interpolant */
283*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals);CHKERRQ(ierr);
284*c4762a1bSJed Brown   ierr = DMInterpolationSetDof(interpolator, Nc);CHKERRQ(ierr);
285*c4762a1bSJed Brown   ierr = DMInterpolationEvaluate(interpolator, dm, lu, fieldVals);CHKERRQ(ierr);
286*c4762a1bSJed Brown   for (p = 0; p < size; ++p) {
287*c4762a1bSJed Brown     if (p == rank) {
288*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank);CHKERRQ(ierr);
289*c4762a1bSJed Brown       ierr = VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
290*c4762a1bSJed Brown     }
291*c4762a1bSJed Brown     ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
292*c4762a1bSJed Brown   }
293*c4762a1bSJed Brown   ierr = VecGetArrayRead(interpolator->coords, &vcoords);CHKERRQ(ierr);
294*c4762a1bSJed Brown   ierr = VecGetArrayRead(fieldVals, &ivals);CHKERRQ(ierr);
295*c4762a1bSJed Brown   for (p = 0; p < interpolator->n; ++p) {
296*c4762a1bSJed Brown     for (c = 0; c < Nc; ++c) {
297*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
298*c4762a1bSJed Brown       PetscReal vcoordsReal[3];
299*c4762a1bSJed Brown       PetscInt  i;
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown       for (i = 0; i < spaceDim; i++) vcoordsReal[i] = PetscRealPart(vcoords[p * spaceDim + i]);
302*c4762a1bSJed Brown #else
303*c4762a1bSJed Brown       const PetscReal *vcoordsReal = &vcoords[p*spaceDim];
304*c4762a1bSJed Brown #endif
305*c4762a1bSJed Brown       (*funcs[c])(ctx.dim, 0.0, vcoordsReal, 1, vals, NULL);
306*c4762a1bSJed Brown       if (PetscAbsScalar(ivals[p*Nc+c] - vals[c]) > PETSC_SQRT_MACHINE_EPSILON)
307*c4762a1bSJed Brown         SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid interpolated value %g != %g (%D, %D)", (double) PetscRealPart(ivals[p*Nc+c]), (double) PetscRealPart(vals[c]), p, c);
308*c4762a1bSJed Brown     }
309*c4762a1bSJed Brown   }
310*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(interpolator->coords, &vcoords);CHKERRQ(ierr);
311*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(fieldVals, &ivals);CHKERRQ(ierr);
312*c4762a1bSJed Brown   /* Cleanup */
313*c4762a1bSJed Brown   ierr = PetscFree(pcoords);CHKERRQ(ierr);
314*c4762a1bSJed Brown   ierr = PetscFree2(funcs, vals);CHKERRQ(ierr);
315*c4762a1bSJed Brown   ierr = VecDestroy(&fieldVals);CHKERRQ(ierr);
316*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(dm, &lu);CHKERRQ(ierr);
317*c4762a1bSJed Brown   ierr = DMInterpolationDestroy(&interpolator);CHKERRQ(ierr);
318*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
319*c4762a1bSJed Brown   ierr = PetscFinalize();
320*c4762a1bSJed Brown   return ierr;
321*c4762a1bSJed Brown }
322*c4762a1bSJed Brown 
323*c4762a1bSJed Brown /*TEST
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown   test:
326*c4762a1bSJed Brown     suffix: 0
327*c4762a1bSJed Brown     requires: ctetgen
328*c4762a1bSJed Brown     args: -petscspace_degree 1
329*c4762a1bSJed Brown   test:
330*c4762a1bSJed Brown     suffix: 1
331*c4762a1bSJed Brown     requires: ctetgen
332*c4762a1bSJed Brown     args: -petscspace_degree 1 -dm_refine 2
333*c4762a1bSJed Brown   test:
334*c4762a1bSJed Brown     suffix: 2
335*c4762a1bSJed Brown     requires: ctetgen
336*c4762a1bSJed Brown     nsize: 2
337*c4762a1bSJed Brown     args: -petscspace_degree 1 -petscpartitioner_type simple
338*c4762a1bSJed Brown   test:
339*c4762a1bSJed Brown     suffix: 3
340*c4762a1bSJed Brown     requires: ctetgen
341*c4762a1bSJed Brown     nsize: 2
342*c4762a1bSJed Brown     args: -petscspace_degree 1 -dm_refine 2 -petscpartitioner_type simple
343*c4762a1bSJed Brown   test:
344*c4762a1bSJed Brown     suffix: 4
345*c4762a1bSJed Brown     requires: ctetgen
346*c4762a1bSJed Brown     nsize: 5
347*c4762a1bSJed Brown     args: -petscspace_degree 1 -petscpartitioner_type simple
348*c4762a1bSJed Brown   test:
349*c4762a1bSJed Brown     suffix: 5
350*c4762a1bSJed Brown     requires: ctetgen
351*c4762a1bSJed Brown     nsize: 5
352*c4762a1bSJed Brown     args: -petscspace_degree 1 -dm_refine 2 -petscpartitioner_type simple
353*c4762a1bSJed Brown   test:
354*c4762a1bSJed Brown     suffix: 6
355*c4762a1bSJed Brown     requires: ctetgen
356*c4762a1bSJed Brown     args: -petscspace_degree 1 -point_type grid
357*c4762a1bSJed Brown   test:
358*c4762a1bSJed Brown     suffix: 7
359*c4762a1bSJed Brown     requires: ctetgen
360*c4762a1bSJed Brown     args: -petscspace_degree 1 -dm_refine 2 -point_type grid
361*c4762a1bSJed Brown   test:
362*c4762a1bSJed Brown     suffix: 8
363*c4762a1bSJed Brown     requires: ctetgen
364*c4762a1bSJed Brown     nsize: 2
365*c4762a1bSJed Brown     args: -petscspace_degree 1 -point_type grid -petscpartitioner_type simple
366*c4762a1bSJed Brown   test:
367*c4762a1bSJed Brown     suffix: 9
368*c4762a1bSJed Brown     requires: ctetgen
369*c4762a1bSJed Brown     args: -petscspace_degree 1 -point_type grid_replicated
370*c4762a1bSJed Brown   test:
371*c4762a1bSJed Brown     suffix: 10
372*c4762a1bSJed Brown     requires: ctetgen
373*c4762a1bSJed Brown     nsize: 2
374*c4762a1bSJed Brown     args: -petscspace_degree 1 -point_type grid_replicated -petscpartitioner_type simple
375*c4762a1bSJed Brown   test:
376*c4762a1bSJed Brown     suffix: 11
377*c4762a1bSJed Brown     requires: ctetgen
378*c4762a1bSJed Brown     nsize: 2
379*c4762a1bSJed Brown     args: -petscspace_degree 1 -dm_refine 2 -point_type grid_replicated -petscpartitioner_type simple
380*c4762a1bSJed Brown 
381*c4762a1bSJed Brown TEST*/
382