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