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 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 19d71ae5a4SJacob Faibussowitsch { 20c4762a1bSJed Brown PetscInt d, c; 21c4762a1bSJed Brown 22*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2363a3b9bcSJacob Faibussowitsch PetscCheck(Nc == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Something is wrong: %" PetscInt_FMT, Nc); 24c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 25c4762a1bSJed Brown u[c] = 0.0; 26c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[c] += x[d]; 27c4762a1bSJed Brown } 28*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 31d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 32d71ae5a4SJacob Faibussowitsch { 33c4762a1bSJed Brown const char *pointTypes[3] = {"centroid", "grid", "grid_replicated"}; 34c4762a1bSJed Brown PetscInt pt; 35c4762a1bSJed Brown 36c4762a1bSJed Brown PetscFunctionBegin; 37c4762a1bSJed Brown options->pointType = CENTROID; 38d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX"); 39c4762a1bSJed Brown pt = options->pointType; 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL)); 41c4762a1bSJed Brown options->pointType = (PointType)pt; 42d0609cedSBarry Smith PetscOptionsEnd(); 43*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 47d71ae5a4SJacob Faibussowitsch { 48c4762a1bSJed Brown PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 509566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 519566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 53*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 57d71ae5a4SJacob Faibussowitsch { 58c4762a1bSJed Brown PetscSection coordSection; 59c4762a1bSJed Brown Vec coordsLocal; 60c4762a1bSJed Brown PetscInt spaceDim, p; 61c4762a1bSJed Brown PetscMPIInt rank; 62c4762a1bSJed Brown 63c4762a1bSJed Brown PetscFunctionBegin; 649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 679566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 689566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, Np)); 699566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 70c4762a1bSJed Brown for (p = 0; p < *Np; ++p) { 71c4762a1bSJed Brown PetscScalar *coords = NULL; 72c4762a1bSJed Brown PetscInt size, num, n, d; 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, p, &size, &coords)); 75c4762a1bSJed Brown num = size / spaceDim; 76c4762a1bSJed Brown for (n = 0; n < num; ++n) { 77c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[p * spaceDim + d] += PetscRealPart(coords[n * spaceDim + d]) / num; 78c4762a1bSJed Brown } 7963a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, p)); 80c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 819566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[p * spaceDim + d])); 829566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 83c4762a1bSJed Brown } 849566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 859566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, p, &num, &coords)); 86c4762a1bSJed Brown } 879566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 88c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 89*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 93d71ae5a4SJacob Faibussowitsch { 94c4762a1bSJed Brown DM da; 95c4762a1bSJed Brown DMDALocalInfo info; 9630602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 97c4762a1bSJed Brown PetscReal *h; 98c4762a1bSJed Brown PetscMPIInt rank; 99c4762a1bSJed Brown 100362febeeSStefano Zampini PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1029566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1039566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 1049566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &ind)); 1059566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &h)); 1069371c9d4SSatish Balay h[0] = 1.0 / (N - 1); 1079371c9d4SSatish Balay h[1] = 1.0 / (N - 1); 1089371c9d4SSatish Balay h[2] = 1.0 / (N - 1); 1099566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da)); 1109566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1119566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, N, N, N)); 1129566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 1139566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, 1)); 1149566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1159566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 1169566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 117c4762a1bSJed Brown *Np = info.xm * info.ym * info.zm; 1189566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 119c4762a1bSJed Brown for (k = info.zs; k < info.zs + info.zm; ++k) { 120c4762a1bSJed Brown ind[2] = k; 121c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; ++j) { 122c4762a1bSJed Brown ind[1] = j; 123c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i, ++n) { 124c4762a1bSJed Brown ind[0] = i; 125c4762a1bSJed Brown 126c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d]; 12763a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n)); 128c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1299566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d])); 1309566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 131c4762a1bSJed Brown } 1329566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 1369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1379566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 1389566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1399566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 140c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 141*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 145d71ae5a4SJacob Faibussowitsch { 14630602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 147c4762a1bSJed Brown PetscReal *h; 148c4762a1bSJed Brown PetscMPIInt rank; 149c4762a1bSJed Brown 15030602db0SMatthew G. Knepley PetscFunctionBeginUser; 1519566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1539566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 1549566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &ind)); 1559566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &h)); 1569371c9d4SSatish Balay h[0] = 1.0 / (N - 1); 1579371c9d4SSatish Balay h[1] = 1.0 / (N - 1); 1589371c9d4SSatish Balay h[2] = 1.0 / (N - 1); 15930602db0SMatthew G. Knepley *Np = N * (dim > 1 ? N : 1) * (dim > 2 ? N : 1); 1609566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 161c4762a1bSJed Brown for (k = 0; k < N; ++k) { 162c4762a1bSJed Brown ind[2] = k; 163c4762a1bSJed Brown for (j = 0; j < N; ++j) { 164c4762a1bSJed Brown ind[1] = j; 165c4762a1bSJed Brown for (i = 0; i < N; ++i, ++n) { 166c4762a1bSJed Brown ind[0] = i; 167c4762a1bSJed Brown 168c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d]; 16963a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n)); 170c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1719566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d])); 1729566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 173c4762a1bSJed Brown } 1749566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 177c4762a1bSJed Brown } 1789566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 179c4762a1bSJed Brown *pointsAllProcs = PETSC_TRUE; 1809566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1819566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 182*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 186d71ae5a4SJacob Faibussowitsch { 187c4762a1bSJed Brown PetscFunctionBegin; 188c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 189c4762a1bSJed Brown switch (ctx->pointType) { 190d71ae5a4SJacob Faibussowitsch case CENTROID: 191d71ae5a4SJacob Faibussowitsch PetscCall(CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx)); 192d71ae5a4SJacob Faibussowitsch break; 193d71ae5a4SJacob Faibussowitsch case GRID: 194d71ae5a4SJacob Faibussowitsch PetscCall(CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx)); 195d71ae5a4SJacob Faibussowitsch break; 196d71ae5a4SJacob Faibussowitsch case GRID_REPLICATED: 197d71ae5a4SJacob Faibussowitsch PetscCall(CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx)); 198d71ae5a4SJacob Faibussowitsch break; 199d71ae5a4SJacob Faibussowitsch default: 200d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int)ctx->pointType); 201c4762a1bSJed Brown } 202*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 206d71ae5a4SJacob Faibussowitsch { 207c4762a1bSJed Brown AppCtx ctx; 208c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 209c4762a1bSJed Brown DM dm; 210c4762a1bSJed Brown PetscFE fe; 211c4762a1bSJed Brown DMInterpolationInfo interpolator; 212c4762a1bSJed Brown Vec lu, fieldVals; 213c4762a1bSJed Brown PetscScalar *vals; 214c4762a1bSJed Brown const PetscScalar *ivals, *vcoords; 215c4762a1bSJed Brown PetscReal *pcoords; 21630602db0SMatthew G. Knepley PetscBool simplex, pointsAllProcs = PETSC_TRUE; 21766a0a883SMatthew G. Knepley PetscInt dim, spaceDim, Nc, c, Np, p; 218c4762a1bSJed Brown PetscMPIInt rank, size; 219c4762a1bSJed Brown PetscViewer selfviewer; 220c4762a1bSJed Brown 221327415f7SBarry Smith PetscFunctionBeginUser; 2229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2239566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 2249566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 2259566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2269566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 2279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 2289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 229c4762a1bSJed Brown /* Create points */ 2309566063dSJacob Faibussowitsch PetscCall(CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx)); 231c4762a1bSJed Brown /* Create interpolator */ 2329566063dSJacob Faibussowitsch PetscCall(DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator)); 2339566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDim(interpolator, spaceDim)); 2349566063dSJacob Faibussowitsch PetscCall(DMInterpolationAddPoints(interpolator, Np, pcoords)); 2359566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetUp(interpolator, dm, pointsAllProcs, PETSC_FALSE)); 236c4762a1bSJed Brown /* Check locations */ 23748a46eb9SPierre 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])); 2389566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 2399566063dSJacob Faibussowitsch PetscCall(VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD)); 240c4762a1bSJed Brown /* Setup Discretization */ 24166a0a883SMatthew G. Knepley Nc = dim; 2429566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2439566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, Nc, simplex, NULL, -1, &fe)); 2449566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2459566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2469566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 247c4762a1bSJed Brown /* Create function */ 2489566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nc, &funcs, Nc, &vals)); 249c4762a1bSJed Brown for (c = 0; c < Nc; ++c) funcs[c] = linear; 2509566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lu)); 2519566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer)); 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank)); 2559566063dSJacob Faibussowitsch PetscCall(VecView(lu, selfviewer)); 2569566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer)); 2579566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 2589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 259c4762a1bSJed Brown /* Check interpolant */ 2609566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals)); 2619566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDof(interpolator, Nc)); 2629566063dSJacob Faibussowitsch PetscCall(DMInterpolationEvaluate(interpolator, dm, lu, fieldVals)); 263c4762a1bSJed Brown for (p = 0; p < size; ++p) { 264c4762a1bSJed Brown if (p == rank) { 2659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank)); 2669566063dSJacob Faibussowitsch PetscCall(VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF)); 267c4762a1bSJed Brown } 2689566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)dm)); 269c4762a1bSJed Brown } 2709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolator->coords, &vcoords)); 2719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(fieldVals, &ivals)); 272c4762a1bSJed Brown for (p = 0; p < interpolator->n; ++p) { 273c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 274c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 275c4762a1bSJed Brown PetscReal vcoordsReal[3]; 276c4762a1bSJed Brown 277*3ba16761SJacob Faibussowitsch for (PetscInt i = 0; i < spaceDim; i++) vcoordsReal[i] = PetscRealPart(vcoords[p * spaceDim + i]); 278c4762a1bSJed Brown #else 279c4762a1bSJed Brown const PetscReal *vcoordsReal = &vcoords[p * spaceDim]; 280c4762a1bSJed Brown #endif 281*3ba16761SJacob Faibussowitsch PetscCall((*funcs[c])(dim, 0.0, vcoordsReal, Nc, vals, NULL)); 282*3ba16761SJacob Faibussowitsch PetscCheck(PetscAbsScalar(ivals[p * Nc + c] - vals[c]) <= PETSC_SQRT_MACHINE_EPSILON, 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); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown } 2859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolator->coords, &vcoords)); 2869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(fieldVals, &ivals)); 287c4762a1bSJed Brown /* Cleanup */ 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(pcoords)); 2899566063dSJacob Faibussowitsch PetscCall(PetscFree2(funcs, vals)); 2909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fieldVals)); 2919566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lu)); 2929566063dSJacob Faibussowitsch PetscCall(DMInterpolationDestroy(&interpolator)); 2939566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2949566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 295b122ec5aSJacob Faibussowitsch return 0; 296c4762a1bSJed Brown } 297c4762a1bSJed Brown 298c4762a1bSJed Brown /*TEST 299c4762a1bSJed Brown 30030602db0SMatthew G. Knepley testset: 30130602db0SMatthew G. Knepley requires: ctetgen 30230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 30330602db0SMatthew G. Knepley 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown suffix: 0 306c4762a1bSJed Brown test: 307c4762a1bSJed Brown suffix: 1 30830602db0SMatthew G. Knepley args: -dm_refine 2 309c4762a1bSJed Brown test: 310c4762a1bSJed Brown suffix: 2 311c4762a1bSJed Brown nsize: 2 312e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown suffix: 3 315c4762a1bSJed Brown nsize: 2 316e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: 4 319c4762a1bSJed Brown nsize: 5 320e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 321c4762a1bSJed Brown test: 322c4762a1bSJed Brown suffix: 5 323c4762a1bSJed Brown nsize: 5 324e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 6 32730602db0SMatthew G. Knepley args: -point_type grid 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: 7 33030602db0SMatthew G. Knepley args: -dm_refine 2 -point_type grid 331c4762a1bSJed Brown test: 332c4762a1bSJed Brown suffix: 8 333c4762a1bSJed Brown nsize: 2 334e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid 335c4762a1bSJed Brown test: 336c4762a1bSJed Brown suffix: 9 33730602db0SMatthew G. Knepley args: -point_type grid_replicated 338c4762a1bSJed Brown test: 339c4762a1bSJed Brown suffix: 10 340c4762a1bSJed Brown nsize: 2 341e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid_replicated 342c4762a1bSJed Brown test: 343c4762a1bSJed Brown suffix: 11 344c4762a1bSJed Brown nsize: 2 345e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple -point_type grid_replicated 346c4762a1bSJed Brown 347c4762a1bSJed Brown TEST*/ 348