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 14c4762a1bSJed Brown static PetscInt Nc = 3; 15c4762a1bSJed Brown 16c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown PetscInt d, c; 19c4762a1bSJed Brown 20c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 21c4762a1bSJed Brown u[c] = 0.0; 22c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[c] += x[d]; 23c4762a1bSJed Brown } 24c4762a1bSJed Brown return 0; 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown const char *pointTypes[3] = {"centroid", "grid", "grid_replicated"}; 30c4762a1bSJed Brown PetscInt pt; 31c4762a1bSJed Brown PetscErrorCode ierr; 32c4762a1bSJed Brown 33c4762a1bSJed Brown PetscFunctionBegin; 34c4762a1bSJed Brown options->pointType = CENTROID; 35c4762a1bSJed Brown 36c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX");CHKERRQ(ierr); 37c4762a1bSJed Brown pt = options->pointType; 38c4762a1bSJed Brown ierr = PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL);CHKERRQ(ierr); 39c4762a1bSJed Brown options->pointType = (PointType) pt; 40c4762a1bSJed Brown ierr = PetscOptionsEnd(); 41c4762a1bSJed Brown PetscFunctionReturn(0); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown PetscErrorCode ierr; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscFunctionBegin; 4930602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 5030602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 53c4762a1bSJed Brown PetscFunctionReturn(0); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown PetscSection coordSection; 59c4762a1bSJed Brown Vec coordsLocal; 60c4762a1bSJed Brown PetscInt spaceDim, p; 61c4762a1bSJed Brown PetscMPIInt rank; 62c4762a1bSJed Brown PetscErrorCode ierr; 63c4762a1bSJed Brown 64c4762a1bSJed Brown PetscFunctionBegin; 65ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 66c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, NULL, Np);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = PetscCalloc1(*Np * spaceDim, pcoords);CHKERRQ(ierr); 71c4762a1bSJed Brown for (p = 0; p < *Np; ++p) { 72c4762a1bSJed Brown PetscScalar *coords = NULL; 73c4762a1bSJed Brown PetscInt size, num, n, d; 74c4762a1bSJed Brown 75c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dm, coordSection, coordsLocal, p, &size, &coords);CHKERRQ(ierr); 76c4762a1bSJed Brown num = size/spaceDim; 77c4762a1bSJed Brown for (n = 0; n < num; ++n) { 78c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) (*pcoords)[p*spaceDim+d] += PetscRealPart(coords[n*spaceDim+d]) / num; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, p);CHKERRQ(ierr); 81c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 82c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[p*spaceDim+d]);CHKERRQ(ierr); 83c4762a1bSJed Brown if (d < spaceDim-1) {ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);} 84c4762a1bSJed Brown } 85c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n");CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, p, &num, &coords);CHKERRQ(ierr); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr); 89c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 90c4762a1bSJed Brown PetscFunctionReturn(0); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 94c4762a1bSJed Brown { 95c4762a1bSJed Brown DM da; 96c4762a1bSJed Brown DMDALocalInfo info; 9730602db0SMatthew G. Knepley PetscInt N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d; 98c4762a1bSJed Brown PetscReal *h; 99c4762a1bSJed Brown PetscMPIInt rank; 100c4762a1bSJed Brown PetscErrorCode ierr; 101c4762a1bSJed Brown 102*362febeeSStefano Zampini PetscFunctionBegin; 103ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 10430602db0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = PetscCalloc1(spaceDim,&ind);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = PetscCalloc1(spaceDim,&h);CHKERRQ(ierr); 108c4762a1bSJed Brown h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1); 109c4762a1bSJed Brown ierr = DMDACreate(PetscObjectComm((PetscObject) dm), &da);CHKERRQ(ierr); 11030602db0SMatthew G. Knepley ierr = DMSetDimension(da, dim);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = DMDASetSizes(da, N, N, N);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = DMDASetDof(da, 1);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = DMDASetStencilWidth(da, 1);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = DMSetUp(da);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = DMDAGetLocalInfo(da, &info);CHKERRQ(ierr); 117c4762a1bSJed Brown *Np = info.xm * info.ym * info.zm; 118c4762a1bSJed Brown ierr = PetscCalloc1(*Np * spaceDim, pcoords);CHKERRQ(ierr); 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]; 127c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n);CHKERRQ(ierr); 128c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 129c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d]);CHKERRQ(ierr); 130c4762a1bSJed Brown if (d < spaceDim-1) {ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);} 131c4762a1bSJed Brown } 132c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n");CHKERRQ(ierr); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 136c4762a1bSJed Brown ierr = DMDestroy(&da);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = PetscFree(ind);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = PetscFree(h);CHKERRQ(ierr); 140c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 145c4762a1bSJed Brown { 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 PetscErrorCode ierr; 150c4762a1bSJed Brown 15130602db0SMatthew G. Knepley PetscFunctionBeginUser; 15230602db0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 153ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 154c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = PetscCalloc1(spaceDim,&ind);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = PetscCalloc1(spaceDim,&h);CHKERRQ(ierr); 157c4762a1bSJed Brown h[0] = 1.0/(N-1); h[1] = 1.0/(N-1); h[2] = 1.0/(N-1); 15830602db0SMatthew G. Knepley *Np = N * (dim > 1 ? N : 1) * (dim > 2 ? N : 1); 159c4762a1bSJed Brown ierr = PetscCalloc1(*Np * spaceDim, pcoords);CHKERRQ(ierr); 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]; 168c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D (", rank, n);CHKERRQ(ierr); 169c4762a1bSJed Brown for (d = 0; d < spaceDim; ++d) { 170c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n*spaceDim+d]);CHKERRQ(ierr); 171c4762a1bSJed Brown if (d < spaceDim-1) {ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", ");CHKERRQ(ierr);} 172c4762a1bSJed Brown } 173c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n");CHKERRQ(ierr); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 177c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr); 178c4762a1bSJed Brown *pointsAllProcs = PETSC_TRUE; 179c4762a1bSJed Brown ierr = PetscFree(ind);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = PetscFree(h);CHKERRQ(ierr); 181c4762a1bSJed Brown PetscFunctionReturn(0); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown PetscErrorCode ierr; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBegin; 189c4762a1bSJed Brown *pointsAllProcs = PETSC_FALSE; 190c4762a1bSJed Brown switch (ctx->pointType) { 191c4762a1bSJed Brown case CENTROID: ierr = CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx);CHKERRQ(ierr);break; 192c4762a1bSJed Brown case GRID: ierr = CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx);CHKERRQ(ierr);break; 193c4762a1bSJed Brown case GRID_REPLICATED: ierr = CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx);CHKERRQ(ierr);break; 194c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int) ctx->pointType); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown PetscFunctionReturn(0); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown int main(int argc, char **argv) 200c4762a1bSJed Brown { 201c4762a1bSJed Brown AppCtx ctx; 202c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 203c4762a1bSJed Brown DM dm; 204c4762a1bSJed Brown PetscFE fe; 205c4762a1bSJed Brown DMInterpolationInfo interpolator; 206c4762a1bSJed Brown Vec lu, fieldVals; 207c4762a1bSJed Brown PetscScalar *vals; 208c4762a1bSJed Brown const PetscScalar *ivals, *vcoords; 209c4762a1bSJed Brown PetscReal *pcoords; 21030602db0SMatthew G. Knepley PetscBool simplex, pointsAllProcs=PETSC_TRUE; 21130602db0SMatthew G. Knepley PetscInt dim, spaceDim, c, Np, p; 212c4762a1bSJed Brown PetscMPIInt rank, size; 213c4762a1bSJed Brown PetscViewer selfviewer; 214c4762a1bSJed Brown PetscErrorCode ierr; 215c4762a1bSJed Brown 216c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 217c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr); 21930602db0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &spaceDim);CHKERRQ(ierr); 221ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 222ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 223c4762a1bSJed Brown /* Create points */ 224c4762a1bSJed Brown ierr = CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx);CHKERRQ(ierr); 225c4762a1bSJed Brown /* Create interpolator */ 226c4762a1bSJed Brown ierr = DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = DMInterpolationSetDim(interpolator, spaceDim);CHKERRQ(ierr); 228c4762a1bSJed Brown ierr = DMInterpolationAddPoints(interpolator, Np, pcoords);CHKERRQ(ierr); 22952aa1562SMatthew G. Knepley ierr = DMInterpolationSetUp(interpolator, dm, pointsAllProcs, PETSC_FALSE);CHKERRQ(ierr); 230c4762a1bSJed Brown /* Check locations */ 231c4762a1bSJed Brown for (c = 0; c < interpolator->n; ++c) { 232c4762a1bSJed Brown ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %D is in Cell %D\n", rank, c, interpolator->cells[c]);CHKERRQ(ierr); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 236c4762a1bSJed Brown /* Setup Discretization */ 23730602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 23830602db0SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, Nc, simplex, NULL, -1, &fe);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 242c4762a1bSJed Brown /* Create function */ 243c4762a1bSJed Brown ierr = PetscCalloc2(Nc, &funcs, Nc, &vals);CHKERRQ(ierr); 244c4762a1bSJed Brown for (c = 0; c < Nc; ++c) funcs[c] = linear; 245c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lu);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = VecView(lu,selfviewer);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 254c4762a1bSJed Brown /* Check interpolant */ 255c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = DMInterpolationSetDof(interpolator, Nc);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = DMInterpolationEvaluate(interpolator, dm, lu, fieldVals);CHKERRQ(ierr); 258c4762a1bSJed Brown for (p = 0; p < size; ++p) { 259c4762a1bSJed Brown if (p == rank) { 260c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank);CHKERRQ(ierr); 261c4762a1bSJed Brown ierr = VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 262c4762a1bSJed Brown } 263c4762a1bSJed Brown ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown ierr = VecGetArrayRead(interpolator->coords, &vcoords);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = VecGetArrayRead(fieldVals, &ivals);CHKERRQ(ierr); 267c4762a1bSJed Brown for (p = 0; p < interpolator->n; ++p) { 268c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 269c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 270c4762a1bSJed Brown PetscReal vcoordsReal[3]; 271c4762a1bSJed Brown PetscInt i; 272c4762a1bSJed Brown 273c4762a1bSJed Brown for (i = 0; i < spaceDim; i++) vcoordsReal[i] = PetscRealPart(vcoords[p * spaceDim + i]); 274c4762a1bSJed Brown #else 275c4762a1bSJed Brown const PetscReal *vcoordsReal = &vcoords[p*spaceDim]; 276c4762a1bSJed Brown #endif 27730602db0SMatthew G. Knepley (*funcs[c])(dim, 0.0, vcoordsReal, 1, vals, NULL); 278c4762a1bSJed Brown if (PetscAbsScalar(ivals[p*Nc+c] - vals[c]) > PETSC_SQRT_MACHINE_EPSILON) 279c4762a1bSJed Brown SETERRQ4(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); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown } 282c4762a1bSJed Brown ierr = VecRestoreArrayRead(interpolator->coords, &vcoords);CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = VecRestoreArrayRead(fieldVals, &ivals);CHKERRQ(ierr); 284c4762a1bSJed Brown /* Cleanup */ 285c4762a1bSJed Brown ierr = PetscFree(pcoords);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = PetscFree2(funcs, vals);CHKERRQ(ierr); 287c4762a1bSJed Brown ierr = VecDestroy(&fieldVals);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lu);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = DMInterpolationDestroy(&interpolator);CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = PetscFinalize(); 292c4762a1bSJed Brown return ierr; 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown /*TEST 296c4762a1bSJed Brown 29730602db0SMatthew G. Knepley testset: 29830602db0SMatthew G. Knepley requires: ctetgen 29930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 30030602db0SMatthew G. Knepley 301c4762a1bSJed Brown test: 302c4762a1bSJed Brown suffix: 0 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 1 30530602db0SMatthew G. Knepley args: -dm_refine 2 306c4762a1bSJed Brown test: 307c4762a1bSJed Brown suffix: 2 308c4762a1bSJed Brown nsize: 2 30930602db0SMatthew G. Knepley args: -dm_distribute -petscpartitioner_type simple 310c4762a1bSJed Brown test: 311c4762a1bSJed Brown suffix: 3 312c4762a1bSJed Brown nsize: 2 31330602db0SMatthew G. Knepley args: -dm_refine 2 -dm_distribute -petscpartitioner_type simple 314c4762a1bSJed Brown test: 315c4762a1bSJed Brown suffix: 4 316c4762a1bSJed Brown nsize: 5 31730602db0SMatthew G. Knepley args: -dm_distribute -petscpartitioner_type simple 318c4762a1bSJed Brown test: 319c4762a1bSJed Brown suffix: 5 320c4762a1bSJed Brown nsize: 5 32130602db0SMatthew G. Knepley args: -dm_refine 2 -dm_distribute -petscpartitioner_type simple 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown suffix: 6 32430602db0SMatthew G. Knepley args: -point_type grid 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 7 32730602db0SMatthew G. Knepley args: -dm_refine 2 -point_type grid 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: 8 330c4762a1bSJed Brown nsize: 2 33130602db0SMatthew G. Knepley args: -dm_distribute -petscpartitioner_type simple -point_type grid 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 9 33430602db0SMatthew G. Knepley args: -point_type grid_replicated 335c4762a1bSJed Brown test: 336c4762a1bSJed Brown suffix: 10 337c4762a1bSJed Brown nsize: 2 33830602db0SMatthew G. Knepley args: -dm_distribute -petscpartitioner_type simple -point_type grid_replicated 339c4762a1bSJed Brown test: 340c4762a1bSJed Brown suffix: 11 341c4762a1bSJed Brown nsize: 2 34230602db0SMatthew G. Knepley args: -dm_refine 2 -dm_distribute -petscpartitioner_type simple -point_type grid_replicated 343c4762a1bSJed Brown 344c4762a1bSJed Brown TEST*/ 345