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 18*d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 19*d71ae5a4SJacob Faibussowitsch { 20c4762a1bSJed Brown PetscInt d, c; 21c4762a1bSJed Brown 2263a3b9bcSJacob Faibussowitsch PetscCheck(Nc == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Something is wrong: %" PetscInt_FMT, Nc); 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 30*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 31*d71ae5a4SJacob Faibussowitsch { 32c4762a1bSJed Brown const char *pointTypes[3] = {"centroid", "grid", "grid_replicated"}; 33c4762a1bSJed Brown PetscInt pt; 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscFunctionBegin; 36c4762a1bSJed Brown options->pointType = CENTROID; 37d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX"); 38c4762a1bSJed Brown pt = options->pointType; 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL)); 40c4762a1bSJed Brown options->pointType = (PointType)pt; 41d0609cedSBarry Smith PetscOptionsEnd(); 42c4762a1bSJed Brown PetscFunctionReturn(0); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 46*d71ae5a4SJacob Faibussowitsch { 47c4762a1bSJed Brown PetscFunctionBegin; 489566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 499566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 52c4762a1bSJed Brown PetscFunctionReturn(0); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 56*d71ae5a4SJacob Faibussowitsch { 57c4762a1bSJed Brown PetscSection coordSection; 58c4762a1bSJed Brown Vec coordsLocal; 59c4762a1bSJed Brown PetscInt spaceDim, p; 60c4762a1bSJed Brown PetscMPIInt rank; 61c4762a1bSJed Brown 62c4762a1bSJed Brown PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 649566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 679566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, Np)); 689566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 69c4762a1bSJed Brown for (p = 0; p < *Np; ++p) { 70c4762a1bSJed Brown PetscScalar *coords = NULL; 71c4762a1bSJed Brown PetscInt size, num, n, d; 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, p, &size, &coords)); 74c4762a1bSJed Brown num = size / spaceDim; 75c4762a1bSJed Brown for (n = 0; n < num; ++n) { 76c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[p * spaceDim + d] += PetscRealPart(coords[n * spaceDim + d]) / num; 77c4762a1bSJed Brown } 7863a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, p)); 79c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 809566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[p * spaceDim + d])); 819566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 82c4762a1bSJed Brown } 839566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 849566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, p, &num, &coords)); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 87c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 88c4762a1bSJed Brown PetscFunctionReturn(0); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 92*d71ae5a4SJacob Faibussowitsch { 93c4762a1bSJed Brown DM da; 94c4762a1bSJed Brown DMDALocalInfo info; 9530602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 96c4762a1bSJed Brown PetscReal *h; 97c4762a1bSJed Brown PetscMPIInt rank; 98c4762a1bSJed Brown 99362febeeSStefano Zampini PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1019566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1029566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 1039566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &ind)); 1049566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &h)); 1059371c9d4SSatish Balay h[0] = 1.0 / (N - 1); 1069371c9d4SSatish Balay h[1] = 1.0 / (N - 1); 1079371c9d4SSatish Balay h[2] = 1.0 / (N - 1); 1089566063dSJacob Faibussowitsch PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da)); 1099566063dSJacob Faibussowitsch PetscCall(DMSetDimension(da, dim)); 1109566063dSJacob Faibussowitsch PetscCall(DMDASetSizes(da, N, N, N)); 1119566063dSJacob Faibussowitsch PetscCall(DMDASetDof(da, 1)); 1129566063dSJacob Faibussowitsch PetscCall(DMDASetStencilWidth(da, 1)); 1139566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1149566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 1159566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 116c4762a1bSJed Brown *Np = info.xm * info.ym * info.zm; 1179566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 118c4762a1bSJed Brown for (k = info.zs; k < info.zs + info.zm; ++k) { 119c4762a1bSJed Brown ind[2] = k; 120c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; ++j) { 121c4762a1bSJed Brown ind[1] = j; 122c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i, ++n) { 123c4762a1bSJed Brown ind[0] = i; 124c4762a1bSJed Brown 125c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d]; 12663a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n)); 127c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1289566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d])); 1299566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 130c4762a1bSJed Brown } 1319566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 1359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1369566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 1379566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1389566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 139c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 140c4762a1bSJed Brown PetscFunctionReturn(0); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 144*d71ae5a4SJacob Faibussowitsch { 14530602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 146c4762a1bSJed Brown PetscReal *h; 147c4762a1bSJed Brown PetscMPIInt rank; 148c4762a1bSJed Brown 14930602db0SMatthew G. Knepley PetscFunctionBeginUser; 1509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1529566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 1539566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &ind)); 1549566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spaceDim, &h)); 1559371c9d4SSatish Balay h[0] = 1.0 / (N - 1); 1569371c9d4SSatish Balay h[1] = 1.0 / (N - 1); 1579371c9d4SSatish Balay h[2] = 1.0 / (N - 1); 15830602db0SMatthew G. Knepley *Np = N * (dim > 1 ? N : 1) * (dim > 2 ? N : 1); 1599566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(*Np * spaceDim, pcoords)); 160c4762a1bSJed Brown for (k = 0; k < N; ++k) { 161c4762a1bSJed Brown ind[2] = k; 162c4762a1bSJed Brown for (j = 0; j < N; ++j) { 163c4762a1bSJed Brown ind[1] = j; 164c4762a1bSJed Brown for (i = 0; i < N; ++i, ++n) { 165c4762a1bSJed Brown ind[0] = i; 166c4762a1bSJed Brown 167c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d]; 16863a3b9bcSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n)); 169c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 1709566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d])); 1719566063dSJacob Faibussowitsch if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 172c4762a1bSJed Brown } 1739566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 1779566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 178c4762a1bSJed Brown *pointsAllProcs = PETSC_TRUE; 1799566063dSJacob Faibussowitsch PetscCall(PetscFree(ind)); 1809566063dSJacob Faibussowitsch PetscCall(PetscFree(h)); 181c4762a1bSJed Brown PetscFunctionReturn(0); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 185*d71ae5a4SJacob Faibussowitsch { 186c4762a1bSJed Brown PetscFunctionBegin; 187c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 188c4762a1bSJed Brown switch (ctx->pointType) { 189*d71ae5a4SJacob Faibussowitsch case CENTROID: 190*d71ae5a4SJacob Faibussowitsch PetscCall(CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx)); 191*d71ae5a4SJacob Faibussowitsch break; 192*d71ae5a4SJacob Faibussowitsch case GRID: 193*d71ae5a4SJacob Faibussowitsch PetscCall(CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx)); 194*d71ae5a4SJacob Faibussowitsch break; 195*d71ae5a4SJacob Faibussowitsch case GRID_REPLICATED: 196*d71ae5a4SJacob Faibussowitsch PetscCall(CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx)); 197*d71ae5a4SJacob Faibussowitsch break; 198*d71ae5a4SJacob Faibussowitsch default: 199*d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int)ctx->pointType); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown PetscFunctionReturn(0); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 205*d71ae5a4SJacob Faibussowitsch { 206c4762a1bSJed Brown AppCtx ctx; 207c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 208c4762a1bSJed Brown DM dm; 209c4762a1bSJed Brown PetscFE fe; 210c4762a1bSJed Brown DMInterpolationInfo interpolator; 211c4762a1bSJed Brown Vec lu, fieldVals; 212c4762a1bSJed Brown PetscScalar *vals; 213c4762a1bSJed Brown const PetscScalar *ivals, *vcoords; 214c4762a1bSJed Brown PetscReal *pcoords; 21530602db0SMatthew G. Knepley PetscBool simplex, pointsAllProcs = PETSC_TRUE; 21666a0a883SMatthew G. Knepley PetscInt dim, spaceDim, Nc, c, Np, p; 217c4762a1bSJed Brown PetscMPIInt rank, size; 218c4762a1bSJed Brown PetscViewer selfviewer; 219c4762a1bSJed Brown 220327415f7SBarry Smith PetscFunctionBeginUser; 2219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2229566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 2239566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 2249566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2259566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &spaceDim)); 2269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 2279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 228c4762a1bSJed Brown /* Create points */ 2299566063dSJacob Faibussowitsch PetscCall(CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx)); 230c4762a1bSJed Brown /* Create interpolator */ 2319566063dSJacob Faibussowitsch PetscCall(DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator)); 2329566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDim(interpolator, spaceDim)); 2339566063dSJacob Faibussowitsch PetscCall(DMInterpolationAddPoints(interpolator, Np, pcoords)); 2349566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetUp(interpolator, dm, pointsAllProcs, PETSC_FALSE)); 235c4762a1bSJed Brown /* Check locations */ 23648a46eb9SPierre 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])); 2379566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 2389566063dSJacob Faibussowitsch PetscCall(VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD)); 239c4762a1bSJed Brown /* Setup Discretization */ 24066a0a883SMatthew G. Knepley Nc = dim; 2419566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2429566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, Nc, simplex, NULL, -1, &fe)); 2439566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2449566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2459566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 246c4762a1bSJed Brown /* Create function */ 2479566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(Nc, &funcs, Nc, &vals)); 248c4762a1bSJed Brown for (c = 0; c < Nc; ++c) funcs[c] = linear; 2499566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lu)); 2509566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu)); 2519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer)); 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank)); 2549566063dSJacob Faibussowitsch PetscCall(VecView(lu, selfviewer)); 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer)); 2569566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 2579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 258c4762a1bSJed Brown /* Check interpolant */ 2599566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals)); 2609566063dSJacob Faibussowitsch PetscCall(DMInterpolationSetDof(interpolator, Nc)); 2619566063dSJacob Faibussowitsch PetscCall(DMInterpolationEvaluate(interpolator, dm, lu, fieldVals)); 262c4762a1bSJed Brown for (p = 0; p < size; ++p) { 263c4762a1bSJed Brown if (p == rank) { 2649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank)); 2659566063dSJacob Faibussowitsch PetscCall(VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF)); 266c4762a1bSJed Brown } 2679566063dSJacob Faibussowitsch PetscCall(PetscBarrier((PetscObject)dm)); 268c4762a1bSJed Brown } 2699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolator->coords, &vcoords)); 2709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(fieldVals, &ivals)); 271c4762a1bSJed Brown for (p = 0; p < interpolator->n; ++p) { 272c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 273c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 274c4762a1bSJed Brown PetscReal vcoordsReal[3]; 275c4762a1bSJed Brown PetscInt i; 276c4762a1bSJed Brown 277c4762a1bSJed Brown for (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 28166a0a883SMatthew G. Knepley (*funcs[c])(dim, 0.0, vcoordsReal, Nc, vals, NULL); 282c4762a1bSJed Brown if (PetscAbsScalar(ivals[p * Nc + c] - vals[c]) > PETSC_SQRT_MACHINE_EPSILON) 28363a3b9bcSJacob 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); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown } 2869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolator->coords, &vcoords)); 2879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(fieldVals, &ivals)); 288c4762a1bSJed Brown /* Cleanup */ 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(pcoords)); 2909566063dSJacob Faibussowitsch PetscCall(PetscFree2(funcs, vals)); 2919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fieldVals)); 2929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lu)); 2939566063dSJacob Faibussowitsch PetscCall(DMInterpolationDestroy(&interpolator)); 2949566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2959566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 296b122ec5aSJacob Faibussowitsch return 0; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown /*TEST 300c4762a1bSJed Brown 30130602db0SMatthew G. Knepley testset: 30230602db0SMatthew G. Knepley requires: ctetgen 30330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 30430602db0SMatthew G. Knepley 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown suffix: 0 307c4762a1bSJed Brown test: 308c4762a1bSJed Brown suffix: 1 30930602db0SMatthew G. Knepley args: -dm_refine 2 310c4762a1bSJed Brown test: 311c4762a1bSJed Brown suffix: 2 312c4762a1bSJed Brown nsize: 2 313e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 314c4762a1bSJed Brown test: 315c4762a1bSJed Brown suffix: 3 316c4762a1bSJed Brown nsize: 2 317e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: 4 320c4762a1bSJed Brown nsize: 5 321e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown suffix: 5 324c4762a1bSJed Brown nsize: 5 325e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 326c4762a1bSJed Brown test: 327c4762a1bSJed Brown suffix: 6 32830602db0SMatthew G. Knepley args: -point_type grid 329c4762a1bSJed Brown test: 330c4762a1bSJed Brown suffix: 7 33130602db0SMatthew G. Knepley args: -dm_refine 2 -point_type grid 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 8 334c4762a1bSJed Brown nsize: 2 335e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid 336c4762a1bSJed Brown test: 337c4762a1bSJed Brown suffix: 9 33830602db0SMatthew G. Knepley args: -point_type grid_replicated 339c4762a1bSJed Brown test: 340c4762a1bSJed Brown suffix: 10 341c4762a1bSJed Brown nsize: 2 342e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid_replicated 343c4762a1bSJed Brown test: 344c4762a1bSJed Brown suffix: 11 345c4762a1bSJed Brown nsize: 2 346e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple -point_type grid_replicated 347c4762a1bSJed Brown 348c4762a1bSJed Brown TEST*/ 349