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 89371c9d4SSatish Balay typedef enum { 99371c9d4SSatish Balay CENTROID, 109371c9d4SSatish Balay GRID, 119371c9d4SSatish Balay GRID_REPLICATED 129371c9d4SSatish Balay } PointType; 13c4762a1bSJed Brown 14c4762a1bSJed Brown typedef struct { 15c4762a1bSJed Brown PointType pointType; /* Point generation mechanism */ 16c4762a1bSJed Brown } AppCtx; 17c4762a1bSJed Brown 189371c9d4SSatish Balay static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 19c4762a1bSJed Brown PetscInt d, c; 20c4762a1bSJed Brown 2163a3b9bcSJacob Faibussowitsch PetscCheck(Nc == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Something is wrong: %" PetscInt_FMT, Nc); 22c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 23c4762a1bSJed Brown u[c] = 0.0; 24c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[c] += x[d]; 25c4762a1bSJed Brown } 26c4762a1bSJed Brown return 0; 27c4762a1bSJed Brown } 28c4762a1bSJed Brown 299371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 30c4762a1bSJed Brown const char *pointTypes[3] = {"centroid", "grid", "grid_replicated"}; 31c4762a1bSJed Brown PetscInt pt; 32c4762a1bSJed Brown 33c4762a1bSJed Brown PetscFunctionBegin; 34c4762a1bSJed Brown options->pointType = CENTROID; 35d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX"); 36c4762a1bSJed Brown pt = options->pointType; 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL)); 38c4762a1bSJed Brown options->pointType = (PointType)pt; 39d0609cedSBarry Smith PetscOptionsEnd(); 40c4762a1bSJed Brown PetscFunctionReturn(0); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 439371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) { 44c4762a1bSJed Brown PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 469566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 479566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 489566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 49c4762a1bSJed Brown PetscFunctionReturn(0); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 529371c9d4SSatish Balay static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) { 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 } 7463a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", 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 879371c9d4SSatish Balay static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) { 88c4762a1bSJed Brown DM da; 89c4762a1bSJed Brown DMDALocalInfo info; 9030602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 91c4762a1bSJed Brown PetscReal *h; 92c4762a1bSJed Brown PetscMPIInt rank; 93c4762a1bSJed Brown 94362febeeSStefano Zampini PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 969566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 979566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 989566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &ind)); 999566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &h)); 1009371c9d4SSatish Balay h[0] = 1.0 / (N - 1); 1019371c9d4SSatish Balay h[1] = 1.0 / (N - 1); 1029371c9d4SSatish Balay h[2] = 1.0 / (N - 1); 1039566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da)); 1049566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1059566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, N, N, N)); 1069566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 1079566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, 1)); 1089566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1099566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 1109566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 111c4762a1bSJed Brown *Np = info.xm * info.ym * info.zm; 1129566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 113c4762a1bSJed Brown for (k = info.zs; k < info.zs + info.zm; ++k) { 114c4762a1bSJed Brown ind[2] = k; 115c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; ++j) { 116c4762a1bSJed Brown ind[1] = j; 117c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i, ++n) { 118c4762a1bSJed Brown ind[0] = i; 119c4762a1bSJed Brown 120c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d]; 12163a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n)); 122c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1239566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d])); 1249566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 125c4762a1bSJed Brown } 1269566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown } 129c4762a1bSJed Brown } 1309566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1319566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1339566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 134c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 135c4762a1bSJed Brown PetscFunctionReturn(0); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 1389371c9d4SSatish Balay static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) { 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)); 1499371c9d4SSatish Balay h[0] = 1.0 / (N - 1); 1509371c9d4SSatish Balay h[1] = 1.0 / (N - 1); 1519371c9d4SSatish Balay h[2] = 1.0 / (N - 1); 15230602db0SMatthew G. Knepley *Np = N * (dim > 1 ? N : 1) * (dim > 2 ? N : 1); 1539566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 154c4762a1bSJed Brown for (k = 0; k < N; ++k) { 155c4762a1bSJed Brown ind[2] = k; 156c4762a1bSJed Brown for (j = 0; j < N; ++j) { 157c4762a1bSJed Brown ind[1] = j; 158c4762a1bSJed Brown for (i = 0; i < N; ++i, ++n) { 159c4762a1bSJed Brown ind[0] = i; 160c4762a1bSJed Brown 161c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d]; 16263a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n)); 163c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1649566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d])); 1659566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 166c4762a1bSJed Brown } 1679566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown } 170c4762a1bSJed Brown } 1719566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 172c4762a1bSJed Brown *pointsAllProcs = PETSC_TRUE; 1739566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 175c4762a1bSJed Brown PetscFunctionReturn(0); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 1789371c9d4SSatish Balay static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) { 179c4762a1bSJed Brown PetscFunctionBegin; 180c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 181c4762a1bSJed Brown switch (ctx->pointType) { 1829566063dSJacob Faibussowitsch case CENTROID: PetscCall(CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx)); break; 1839566063dSJacob Faibussowitsch case GRID: PetscCall(CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx)); break; 1849566063dSJacob Faibussowitsch case GRID_REPLICATED: PetscCall(CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx)); break; 18598921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int)ctx->pointType); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown PetscFunctionReturn(0); 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 1909371c9d4SSatish Balay int main(int argc, char **argv) { 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 205327415f7SBarry Smith PetscFunctionBeginUser; 2069566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2079566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 2089566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 2099566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2109566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 2119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 2129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 213c4762a1bSJed Brown /* Create points */ 2149566063dSJacob Faibussowitsch PetscCall(CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx)); 215c4762a1bSJed Brown /* Create interpolator */ 2169566063dSJacob Faibussowitsch PetscCall(DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator)); 2179566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDim(interpolator, spaceDim)); 2189566063dSJacob Faibussowitsch PetscCall(DMInterpolationAddPoints(interpolator, Np, pcoords)); 2199566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetUp(interpolator, dm, pointsAllProcs, PETSC_FALSE)); 220c4762a1bSJed Brown /* Check locations */ 221*48a46eb9SPierre Jolivet for (c = 0; c < interpolator->n; ++c) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " is in Cell %" PetscInt_FMT "\n", rank, c, interpolator->cells[c])); 2229566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 2239566063dSJacob Faibussowitsch PetscCall(VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD)); 224c4762a1bSJed Brown /* Setup Discretization */ 22566a0a883SMatthew G. Knepley Nc = dim; 2269566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2279566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, Nc, simplex, NULL, -1, &fe)); 2289566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2299566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2309566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 231c4762a1bSJed Brown /* Create function */ 2329566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nc, &funcs, Nc, &vals)); 233c4762a1bSJed Brown for (c = 0; c < Nc; ++c) funcs[c] = linear; 2349566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lu)); 2359566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu)); 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer)); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank)); 2399566063dSJacob Faibussowitsch PetscCall(VecView(lu, selfviewer)); 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer)); 2419566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 243c4762a1bSJed Brown /* Check interpolant */ 2449566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals)); 2459566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDof(interpolator, Nc)); 2469566063dSJacob Faibussowitsch PetscCall(DMInterpolationEvaluate(interpolator, dm, lu, fieldVals)); 247c4762a1bSJed Brown for (p = 0; p < size; ++p) { 248c4762a1bSJed Brown if (p == rank) { 2499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank)); 2509566063dSJacob Faibussowitsch PetscCall(VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF)); 251c4762a1bSJed Brown } 2529566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)dm)); 253c4762a1bSJed Brown } 2549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolator->coords, &vcoords)); 2559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(fieldVals, &ivals)); 256c4762a1bSJed Brown for (p = 0; p < interpolator->n; ++p) { 257c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 258c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 259c4762a1bSJed Brown PetscReal vcoordsReal[3]; 260c4762a1bSJed Brown PetscInt i; 261c4762a1bSJed Brown 262c4762a1bSJed Brown for (i = 0; i < spaceDim; i++) vcoordsReal[i] = PetscRealPart(vcoords[p * spaceDim + i]); 263c4762a1bSJed Brown #else 264c4762a1bSJed Brown const PetscReal *vcoordsReal = &vcoords[p * spaceDim]; 265c4762a1bSJed Brown #endif 26666a0a883SMatthew G. Knepley (*funcs[c])(dim, 0.0, vcoordsReal, Nc, vals, NULL); 267c4762a1bSJed Brown if (PetscAbsScalar(ivals[p * Nc + c] - vals[c]) > PETSC_SQRT_MACHINE_EPSILON) 26863a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid interpolated value %g != %g (%" PetscInt_FMT ", %" PetscInt_FMT ")", (double)PetscRealPart(ivals[p * Nc + c]), (double)PetscRealPart(vals[c]), p, c); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown } 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolator->coords, &vcoords)); 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(fieldVals, &ivals)); 273c4762a1bSJed Brown /* Cleanup */ 2749566063dSJacob Faibussowitsch PetscCall(PetscFree(pcoords)); 2759566063dSJacob Faibussowitsch PetscCall(PetscFree2(funcs, vals)); 2769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fieldVals)); 2779566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lu)); 2789566063dSJacob Faibussowitsch PetscCall(DMInterpolationDestroy(&interpolator)); 2799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2809566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 281b122ec5aSJacob Faibussowitsch return 0; 282c4762a1bSJed Brown } 283c4762a1bSJed Brown 284c4762a1bSJed Brown /*TEST 285c4762a1bSJed Brown 28630602db0SMatthew G. Knepley testset: 28730602db0SMatthew G. Knepley requires: ctetgen 28830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 28930602db0SMatthew G. Knepley 290c4762a1bSJed Brown test: 291c4762a1bSJed Brown suffix: 0 292c4762a1bSJed Brown test: 293c4762a1bSJed Brown suffix: 1 29430602db0SMatthew G. Knepley args: -dm_refine 2 295c4762a1bSJed Brown test: 296c4762a1bSJed Brown suffix: 2 297c4762a1bSJed Brown nsize: 2 298e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 299c4762a1bSJed Brown test: 300c4762a1bSJed Brown suffix: 3 301c4762a1bSJed Brown nsize: 2 302e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 4 305c4762a1bSJed Brown nsize: 5 306e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 307c4762a1bSJed Brown test: 308c4762a1bSJed Brown suffix: 5 309c4762a1bSJed Brown nsize: 5 310e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 311c4762a1bSJed Brown test: 312c4762a1bSJed Brown suffix: 6 31330602db0SMatthew G. Knepley args: -point_type grid 314c4762a1bSJed Brown test: 315c4762a1bSJed Brown suffix: 7 31630602db0SMatthew G. Knepley args: -dm_refine 2 -point_type grid 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: 8 319c4762a1bSJed Brown nsize: 2 320e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid 321c4762a1bSJed Brown test: 322c4762a1bSJed Brown suffix: 9 32330602db0SMatthew G. Knepley args: -point_type grid_replicated 324c4762a1bSJed Brown test: 325c4762a1bSJed Brown suffix: 10 326c4762a1bSJed Brown nsize: 2 327e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid_replicated 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: 11 330c4762a1bSJed Brown nsize: 2 331e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple -point_type grid_replicated 332c4762a1bSJed Brown 333c4762a1bSJed Brown TEST*/ 334