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 182c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 PetscErrorCode ierr; 31c4762a1bSJed Brown 32c4762a1bSJed Brown PetscFunctionBegin; 33c4762a1bSJed Brown options->pointType = CENTROID; 34c4762a1bSJed Brown 35c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX");CHKERRQ(ierr); 36c4762a1bSJed Brown pt = options->pointType; 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL)); 38c4762a1bSJed Brown options->pointType = (PointType) pt; 391e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 40c4762a1bSJed Brown PetscFunctionReturn(0); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown PetscFunctionBegin; 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 50c4762a1bSJed Brown PetscFunctionReturn(0); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53c4762a1bSJed Brown static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 54c4762a1bSJed Brown { 55c4762a1bSJed Brown PetscSection coordSection; 56c4762a1bSJed Brown Vec coordsLocal; 57c4762a1bSJed Brown PetscInt spaceDim, p; 58c4762a1bSJed Brown PetscMPIInt rank; 59c4762a1bSJed Brown 60c4762a1bSJed Brown PetscFunctionBegin; 61*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordsLocal)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(dm, &coordSection)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &spaceDim)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, NULL, Np)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(*Np * spaceDim, pcoords)); 67c4762a1bSJed Brown for (p = 0; p < *Np; ++p) { 68c4762a1bSJed Brown PetscScalar *coords = NULL; 69c4762a1bSJed Brown PetscInt size, num, n, d; 70c4762a1bSJed Brown 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(dm, coordSection, coordsLocal, p, &size, &coords)); 72c4762a1bSJed Brown num = size/spaceDim; 73c4762a1bSJed Brown for (n = 0; n < num; ++n) { 74c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[p*spaceDim+d] += PetscRealPart(coords[n*spaceDim+d]) / num; 75c4762a1bSJed Brown } 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, p)); 77c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[p*spaceDim+d])); 79*5f80ce2aSJacob Faibussowitsch if (d < spaceDim-1) CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 80c4762a1bSJed Brown } 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, p, &num, &coords)); 83c4762a1bSJed Brown } 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 85c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 86c4762a1bSJed Brown PetscFunctionReturn(0); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 90c4762a1bSJed Brown { 91c4762a1bSJed Brown DM da; 92c4762a1bSJed Brown DMDALocalInfo info; 9330602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 94c4762a1bSJed Brown PetscReal *h; 95c4762a1bSJed Brown PetscMPIInt rank; 96c4762a1bSJed Brown 97362febeeSStefano Zampini PetscFunctionBegin; 98*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &spaceDim)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(spaceDim,&ind)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(spaceDim,&h)); 103c4762a1bSJed Brown h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate(PetscObjectComm((PetscObject) dm), &da)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(da, dim)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetSizes(da, N, N, N)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetDof(da, 1)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetStencilWidth(da, 1)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da, &info)); 112c4762a1bSJed Brown *Np = info.xm * info.ym * info.zm; 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(*Np * spaceDim, pcoords)); 114c4762a1bSJed Brown for (k = info.zs; k < info.zs + info.zm; ++k) { 115c4762a1bSJed Brown ind[2] = k; 116c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; ++j) { 117c4762a1bSJed Brown ind[1] = j; 118c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; ++i, ++n) { 119c4762a1bSJed Brown ind[0] = i; 120c4762a1bSJed Brown 121c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[n*spaceDim+d] = ind[d]*h[d]; 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n)); 123c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d])); 125*5f80ce2aSJacob Faibussowitsch if (d < spaceDim-1) CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 126c4762a1bSJed Brown } 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown } 130c4762a1bSJed Brown } 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ind)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(h)); 135c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 136c4762a1bSJed Brown PetscFunctionReturn(0); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 140c4762a1bSJed Brown { 14130602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 142c4762a1bSJed Brown PetscReal *h; 143c4762a1bSJed Brown PetscMPIInt rank; 144c4762a1bSJed Brown 14530602db0SMatthew G. Knepley PetscFunctionBeginUser; 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 147*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &spaceDim)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(spaceDim,&ind)); 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(spaceDim,&h)); 151c4762a1bSJed Brown h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1); 15230602db0SMatthew G. Knepley *Np = N * (dim > 1 ? N : 1) * (dim > 2 ? N : 1); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(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]; 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n)); 163c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d])); 165*5f80ce2aSJacob Faibussowitsch if (d < spaceDim-1) CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ")); 166c4762a1bSJed Brown } 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n")); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown } 170c4762a1bSJed Brown } 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 172c4762a1bSJed Brown *pointsAllProcs = PETSC_TRUE; 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ind)); 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(h)); 175c4762a1bSJed Brown PetscFunctionReturn(0); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown PetscFunctionBegin; 181c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 182c4762a1bSJed Brown switch (ctx->pointType) { 183*5f80ce2aSJacob Faibussowitsch case CENTROID: CHKERRQ(CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx));break; 184*5f80ce2aSJacob Faibussowitsch case GRID: CHKERRQ(CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx));break; 185*5f80ce2aSJacob Faibussowitsch case GRID_REPLICATED: CHKERRQ(CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx));break; 18698921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int) ctx->pointType); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown PetscFunctionReturn(0); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown int main(int argc, char **argv) 192c4762a1bSJed Brown { 193c4762a1bSJed Brown AppCtx ctx; 194c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 195c4762a1bSJed Brown DM dm; 196c4762a1bSJed Brown PetscFE fe; 197c4762a1bSJed Brown DMInterpolationInfo interpolator; 198c4762a1bSJed Brown Vec lu, fieldVals; 199c4762a1bSJed Brown PetscScalar *vals; 200c4762a1bSJed Brown const PetscScalar *ivals, *vcoords; 201c4762a1bSJed Brown PetscReal *pcoords; 20230602db0SMatthew G. Knepley PetscBool simplex, pointsAllProcs=PETSC_TRUE; 20366a0a883SMatthew G. Knepley PetscInt dim, spaceDim, Nc, c, Np, p; 204c4762a1bSJed Brown PetscMPIInt rank, size; 205c4762a1bSJed Brown PetscViewer selfviewer; 206c4762a1bSJed Brown PetscErrorCode ierr; 207c4762a1bSJed Brown 208c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 212*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &spaceDim)); 213*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 214*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 215c4762a1bSJed Brown /* Create points */ 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx)); 217c4762a1bSJed Brown /* Create interpolator */ 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator)); 219*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationSetDim(interpolator, spaceDim)); 220*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationAddPoints(interpolator, Np, pcoords)); 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationSetUp(interpolator, dm, pointsAllProcs, PETSC_FALSE)); 222c4762a1bSJed Brown /* Check locations */ 223c4762a1bSJed Brown for (c = 0; c < interpolator->n; ++c) { 224*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D is in Cell %D\n", rank, c, interpolator->cells[c])); 225c4762a1bSJed Brown } 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD)); 228c4762a1bSJed Brown /* Setup Discretization */ 22966a0a883SMatthew G. Knepley Nc = dim; 230*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 231*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, Nc, simplex, NULL, -1, &fe)); 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 233*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 234*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 235c4762a1bSJed Brown /* Create function */ 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc2(Nc, &funcs, Nc, &vals)); 237c4762a1bSJed Brown for (c = 0; c < Nc; ++c) funcs[c] = linear; 238*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lu)); 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu)); 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer)); 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank)); 243*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(lu,selfviewer)); 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer)); 245*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 247c4762a1bSJed Brown /* Check interpolant */ 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals)); 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationSetDof(interpolator, Nc)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationEvaluate(interpolator, dm, lu, fieldVals)); 251c4762a1bSJed Brown for (p = 0; p < size; ++p) { 252c4762a1bSJed Brown if (p == rank) { 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank)); 254*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF)); 255c4762a1bSJed Brown } 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBarrier((PetscObject) dm)); 257c4762a1bSJed Brown } 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(interpolator->coords, &vcoords)); 259*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(fieldVals, &ivals)); 260c4762a1bSJed Brown for (p = 0; p < interpolator->n; ++p) { 261c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 262c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 263c4762a1bSJed Brown PetscReal vcoordsReal[3]; 264c4762a1bSJed Brown PetscInt i; 265c4762a1bSJed Brown 266c4762a1bSJed Brown for (i = 0; i < spaceDim; i++) vcoordsReal[i] = PetscRealPart(vcoords[p * spaceDim + i]); 267c4762a1bSJed Brown #else 268c4762a1bSJed Brown const PetscReal *vcoordsReal = &vcoords[p*spaceDim]; 269c4762a1bSJed Brown #endif 27066a0a883SMatthew G. Knepley (*funcs[c])(dim, 0.0, vcoordsReal, Nc, vals, NULL); 271c4762a1bSJed Brown if (PetscAbsScalar(ivals[p*Nc+c] - vals[c]) > PETSC_SQRT_MACHINE_EPSILON) 27298921bdaSJacob 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); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown } 275*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(interpolator->coords, &vcoords)); 276*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(fieldVals, &ivals)); 277c4762a1bSJed Brown /* Cleanup */ 278*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pcoords)); 279*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(funcs, vals)); 280*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&fieldVals)); 281*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lu)); 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMInterpolationDestroy(&interpolator)); 283*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 284c4762a1bSJed Brown ierr = PetscFinalize(); 285c4762a1bSJed Brown return ierr; 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /*TEST 289c4762a1bSJed Brown 29030602db0SMatthew G. Knepley testset: 29130602db0SMatthew G. Knepley requires: ctetgen 29230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 29330602db0SMatthew G. Knepley 294c4762a1bSJed Brown test: 295c4762a1bSJed Brown suffix: 0 296c4762a1bSJed Brown test: 297c4762a1bSJed Brown suffix: 1 29830602db0SMatthew G. Knepley args: -dm_refine 2 299c4762a1bSJed Brown test: 300c4762a1bSJed Brown suffix: 2 301c4762a1bSJed Brown nsize: 2 302e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 3 305c4762a1bSJed Brown nsize: 2 306e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 307c4762a1bSJed Brown test: 308c4762a1bSJed Brown suffix: 4 309c4762a1bSJed Brown nsize: 5 310e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple 311c4762a1bSJed Brown test: 312c4762a1bSJed Brown suffix: 5 313c4762a1bSJed Brown nsize: 5 314e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple 315c4762a1bSJed Brown test: 316c4762a1bSJed Brown suffix: 6 31730602db0SMatthew G. Knepley args: -point_type grid 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: 7 32030602db0SMatthew G. Knepley args: -dm_refine 2 -point_type grid 321c4762a1bSJed Brown test: 322c4762a1bSJed Brown suffix: 8 323c4762a1bSJed Brown nsize: 2 324e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 9 32730602db0SMatthew G. Knepley args: -point_type grid_replicated 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: 10 330c4762a1bSJed Brown nsize: 2 331e600fa54SMatthew G. Knepley args: -petscpartitioner_type simple -point_type grid_replicated 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 11 334c4762a1bSJed Brown nsize: 2 335e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple -point_type grid_replicated 336c4762a1bSJed Brown 337c4762a1bSJed Brown TEST*/ 338