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 PointType pointType; /* Point generation mechanism */ 12c4762a1bSJed Brown } AppCtx; 13c4762a1bSJed Brown 1466a0a883SMatthew G. Knepley static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscInt d, c; 17c4762a1bSJed Brown 18e00437b9SBarry Smith PetscCheck(Nc == 3,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Something is wrong: %D", Nc); 19c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 20c4762a1bSJed Brown u[c] = 0.0; 21c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[c] += x[d]; 22c4762a1bSJed Brown } 23c4762a1bSJed Brown return 0; 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 26c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 27c4762a1bSJed Brown { 28c4762a1bSJed Brown const char *pointTypes[3] = {"centroid", "grid", "grid_replicated"}; 29c4762a1bSJed Brown PetscInt pt; 30c4762a1bSJed Brown 31c4762a1bSJed Brown PetscFunctionBegin; 32c4762a1bSJed Brown options->pointType = CENTROID; 33*d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX"); 34c4762a1bSJed Brown pt = options->pointType; 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL)); 36c4762a1bSJed Brown options->pointType = (PointType) pt; 37*d0609cedSBarry Smith PetscOptionsEnd(); 38c4762a1bSJed Brown PetscFunctionReturn(0); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 42c4762a1bSJed Brown { 43c4762a1bSJed Brown PetscFunctionBegin; 449566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 459566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 469566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 479566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 48c4762a1bSJed Brown PetscFunctionReturn(0); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 52c4762a1bSJed Brown { 53c4762a1bSJed Brown PetscSection coordSection; 54c4762a1bSJed Brown Vec coordsLocal; 55c4762a1bSJed Brown PetscInt spaceDim, p; 56c4762a1bSJed Brown PetscMPIInt rank; 57c4762a1bSJed Brown 58c4762a1bSJed Brown PetscFunctionBegin; 599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 629566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 639566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, Np)); 649566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 65c4762a1bSJed Brown for (p = 0; p < *Np; ++p) { 66c4762a1bSJed Brown PetscScalar *coords = NULL; 67c4762a1bSJed Brown PetscInt size, num, n, d; 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, p, &size, &coords)); 70c4762a1bSJed Brown num = size/spaceDim; 71c4762a1bSJed Brown for (n = 0; n < num; ++n) { 72c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[p*spaceDim+d] += PetscRealPart(coords[n*spaceDim+d]) / num; 73c4762a1bSJed Brown } 749566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, p)); 75c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 769566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[p*spaceDim+d])); 779566063dSJacob Faibussowitsch if (d < spaceDim-1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 78c4762a1bSJed Brown } 799566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 809566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, p, &num, &coords)); 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 83c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 84c4762a1bSJed Brown PetscFunctionReturn(0); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87c4762a1bSJed Brown static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown DM da; 90c4762a1bSJed Brown DMDALocalInfo info; 9130602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 92c4762a1bSJed Brown PetscReal *h; 93c4762a1bSJed Brown PetscMPIInt rank; 94c4762a1bSJed Brown 95362febeeSStefano Zampini PetscFunctionBegin; 969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 979566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 989566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 999566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim,&ind)); 1009566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim,&h)); 101c4762a1bSJed Brown h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1); 1029566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject) dm), &da)); 1039566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1049566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, N, N, N)); 1059566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 1069566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, 1)); 1079566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1089566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 1099566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 110c4762a1bSJed Brown *Np = info.xm * info.ym * info.zm; 1119566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 112c4762a1bSJed Brown for (k = info.zs; k < info.zs + info.zm; ++k) { 113c4762a1bSJed Brown ind[2] = k; 114c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; ++j) { 115c4762a1bSJed Brown ind[1] = j; 116c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i, ++n) { 117c4762a1bSJed Brown ind[0] = i; 118c4762a1bSJed Brown 119c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n*spaceDim+d] = ind[d]*h[d]; 1209566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n)); 121c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1229566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d])); 1239566063dSJacob Faibussowitsch if (d < spaceDim-1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 124c4762a1bSJed Brown } 1259566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown } 128c4762a1bSJed Brown } 1299566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1309566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 133c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 134c4762a1bSJed Brown PetscFunctionReturn(0); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 138c4762a1bSJed Brown { 13930602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 140c4762a1bSJed Brown PetscReal *h; 141c4762a1bSJed Brown PetscMPIInt rank; 142c4762a1bSJed Brown 14330602db0SMatthew G. Knepley PetscFunctionBeginUser; 1449566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1469566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 1479566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim,&ind)); 1489566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim,&h)); 149c4762a1bSJed Brown h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1); 15030602db0SMatthew G. Knepley *Np = N * (dim > 1 ? N : 1) * (dim > 2 ? N : 1); 1519566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 152c4762a1bSJed Brown for (k = 0; k < N; ++k) { 153c4762a1bSJed Brown ind[2] = k; 154c4762a1bSJed Brown for (j = 0; j < N; ++j) { 155c4762a1bSJed Brown ind[1] = j; 156c4762a1bSJed Brown for (i = 0; i < N; ++i, ++n) { 157c4762a1bSJed Brown ind[0] = i; 158c4762a1bSJed Brown 159c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n*spaceDim+d] = ind[d]*h[d]; 1609566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n)); 161c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1629566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d])); 1639566063dSJacob Faibussowitsch if (d < spaceDim-1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 164c4762a1bSJed Brown } 1659566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown } 168c4762a1bSJed Brown } 1699566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 170c4762a1bSJed Brown *pointsAllProcs = PETSC_TRUE; 1719566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1729566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 173c4762a1bSJed Brown PetscFunctionReturn(0); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown PetscFunctionBegin; 179c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 180c4762a1bSJed Brown switch (ctx->pointType) { 1819566063dSJacob Faibussowitsch case CENTROID: PetscCall(CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx));break; 1829566063dSJacob Faibussowitsch case GRID: PetscCall(CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx));break; 1839566063dSJacob Faibussowitsch case GRID_REPLICATED: PetscCall(CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx));break; 18498921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int) ctx->pointType); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown PetscFunctionReturn(0); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown int main(int argc, char **argv) 190c4762a1bSJed Brown { 191c4762a1bSJed Brown AppCtx ctx; 192c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 193c4762a1bSJed Brown DM dm; 194c4762a1bSJed Brown PetscFE fe; 195c4762a1bSJed Brown DMInterpolationInfo interpolator; 196c4762a1bSJed Brown Vec lu, fieldVals; 197c4762a1bSJed Brown PetscScalar *vals; 198c4762a1bSJed Brown const PetscScalar *ivals, *vcoords; 199c4762a1bSJed Brown PetscReal *pcoords; 20030602db0SMatthew G. Knepley PetscBool simplex, pointsAllProcs=PETSC_TRUE; 20166a0a883SMatthew G. Knepley PetscInt dim, spaceDim, Nc, c, Np, p; 202c4762a1bSJed Brown PetscMPIInt rank, size; 203c4762a1bSJed Brown PetscViewer selfviewer; 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 2069566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 2079566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 2089566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2099566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 2109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 2119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 212c4762a1bSJed Brown /* Create points */ 2139566063dSJacob Faibussowitsch PetscCall(CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx)); 214c4762a1bSJed Brown /* Create interpolator */ 2159566063dSJacob Faibussowitsch PetscCall(DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator)); 2169566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDim(interpolator, spaceDim)); 2179566063dSJacob Faibussowitsch PetscCall(DMInterpolationAddPoints(interpolator, Np, pcoords)); 2189566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetUp(interpolator, dm, pointsAllProcs, PETSC_FALSE)); 219c4762a1bSJed Brown /* Check locations */ 220c4762a1bSJed Brown for (c = 0; c < interpolator->n; ++c) { 2219566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D is in Cell %D\n", rank, c, interpolator->cells[c])); 222c4762a1bSJed Brown } 2239566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 2249566063dSJacob Faibussowitsch PetscCall(VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD)); 225c4762a1bSJed Brown /* Setup Discretization */ 22666a0a883SMatthew G. Knepley Nc = dim; 2279566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2289566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, Nc, simplex, NULL, -1, &fe)); 2299566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 2309566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2319566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 232c4762a1bSJed Brown /* Create function */ 2339566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nc, &funcs, Nc, &vals)); 234c4762a1bSJed Brown for (c = 0; c < Nc; ++c) funcs[c] = linear; 2359566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lu)); 2369566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu)); 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer)); 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank)); 2409566063dSJacob Faibussowitsch PetscCall(VecView(lu,selfviewer)); 2419566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer)); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 244c4762a1bSJed Brown /* Check interpolant */ 2459566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals)); 2469566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDof(interpolator, Nc)); 2479566063dSJacob Faibussowitsch PetscCall(DMInterpolationEvaluate(interpolator, dm, lu, fieldVals)); 248c4762a1bSJed Brown for (p = 0; p < size; ++p) { 249c4762a1bSJed Brown if (p == rank) { 2509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank)); 2519566063dSJacob Faibussowitsch PetscCall(VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF)); 252c4762a1bSJed Brown } 2539566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject) dm)); 254c4762a1bSJed Brown } 2559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolator->coords, &vcoords)); 2569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(fieldVals, &ivals)); 257c4762a1bSJed Brown for (p = 0; p < interpolator->n; ++p) { 258c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 259c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 260c4762a1bSJed Brown PetscReal vcoordsReal[3]; 261c4762a1bSJed Brown PetscInt i; 262c4762a1bSJed Brown 263c4762a1bSJed Brown for (i = 0; i < spaceDim; i++) vcoordsReal[i] = PetscRealPart(vcoords[p * spaceDim + i]); 264c4762a1bSJed Brown #else 265c4762a1bSJed Brown const PetscReal *vcoordsReal = &vcoords[p*spaceDim]; 266c4762a1bSJed Brown #endif 26766a0a883SMatthew G. Knepley (*funcs[c])(dim, 0.0, vcoordsReal, Nc, vals, NULL); 268c4762a1bSJed Brown if (PetscAbsScalar(ivals[p*Nc+c] - vals[c]) > PETSC_SQRT_MACHINE_EPSILON) 26998921bdaSJacob Faibussowitsch SETERRQ(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); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolator->coords, &vcoords)); 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(fieldVals, &ivals)); 274c4762a1bSJed Brown /* Cleanup */ 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(pcoords)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree2(funcs, vals)); 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fieldVals)); 2789566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lu)); 2799566063dSJacob Faibussowitsch PetscCall(DMInterpolationDestroy(&interpolator)); 2809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2819566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 282b122ec5aSJacob Faibussowitsch return 0; 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285c4762a1bSJed Brown /*TEST 286c4762a1bSJed Brown 28730602db0SMatthew G. Knepley testset: 28830602db0SMatthew G. Knepley requires: ctetgen 28930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 29030602db0SMatthew G. Knepley 291c4762a1bSJed Brown test: 292c4762a1bSJed Brown suffix: 0 293c4762a1bSJed Brown test: 294c4762a1bSJed Brown suffix: 1 29530602db0SMatthew G. Knepley args: -dm_refine 2 296c4762a1bSJed Brown test: 297c4762a1bSJed Brown suffix: 2 298c4762a1bSJed Brown nsize: 2 299e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown suffix: 3 302c4762a1bSJed Brown nsize: 2 303e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown suffix: 4 306c4762a1bSJed Brown nsize: 5 307e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: 5 310c4762a1bSJed Brown nsize: 5 311e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 312c4762a1bSJed Brown test: 313c4762a1bSJed Brown suffix: 6 31430602db0SMatthew G. Knepley args: -point_type grid 315c4762a1bSJed Brown test: 316c4762a1bSJed Brown suffix: 7 31730602db0SMatthew G. Knepley args: -dm_refine 2 -point_type grid 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: 8 320c4762a1bSJed Brown nsize: 2 321e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown suffix: 9 32430602db0SMatthew G. Knepley args: -point_type grid_replicated 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 10 327c4762a1bSJed Brown nsize: 2 328e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid_replicated 329c4762a1bSJed Brown test: 330c4762a1bSJed Brown suffix: 11 331c4762a1bSJed Brown nsize: 2 332e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple -point_type grid_replicated 333c4762a1bSJed Brown 334c4762a1bSJed Brown TEST*/ 335