xref: /petsc/src/snes/tests/ex2.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 
14c3af198dSMatthew G. Knepley typedef enum {
15c3af198dSMatthew G. Knepley   CONSTANT,
16c3af198dSMatthew G. Knepley   LINEAR
17c3af198dSMatthew G. Knepley } FuncType;
18c3af198dSMatthew G. Knepley 
19c4762a1bSJed Brown typedef struct {
20c3af198dSMatthew G. Knepley   PointType pointType; // Point generation mechanism
21c3af198dSMatthew G. Knepley   FuncType  funcType;  // Type of interpolated function
22c3af198dSMatthew G. Knepley   PetscBool useFV;     // Use finite volume, instead of finite element
23c4762a1bSJed Brown } AppCtx;
24c4762a1bSJed Brown 
constant(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)25*2a8381b2SBarry Smith static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
26c3af198dSMatthew G. Knepley {
27c3af198dSMatthew G. Knepley   PetscFunctionBeginUser;
28c3af198dSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) u[c] = c + 1.;
29c3af198dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
30c3af198dSMatthew G. Knepley }
31c3af198dSMatthew G. Knepley 
linear(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)32*2a8381b2SBarry Smith static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
33d71ae5a4SJacob Faibussowitsch {
343ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
3563a3b9bcSJacob Faibussowitsch   PetscCheck(Nc == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Something is wrong: %" PetscInt_FMT, Nc);
36c3af198dSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) {
37c4762a1bSJed Brown     u[c] = 0.0;
38c3af198dSMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) u[c] += x[d];
39c4762a1bSJed Brown   }
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41c4762a1bSJed Brown }
42c4762a1bSJed Brown 
ProcessOptions(MPI_Comm comm,AppCtx * options)43d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown   const char *pointTypes[3] = {"centroid", "grid", "grid_replicated"};
46c3af198dSMatthew G. Knepley   const char *funcTypes[2]  = {"constant", "linear"};
47c3af198dSMatthew G. Knepley   PetscInt    pt, fn;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   PetscFunctionBegin;
50c4762a1bSJed Brown   options->pointType = CENTROID;
51c3af198dSMatthew G. Knepley   options->funcType  = LINEAR;
52c3af198dSMatthew G. Knepley   options->useFV     = PETSC_FALSE;
53c3af198dSMatthew G. Knepley 
54d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Interpolation Options", "DMPLEX");
55c4762a1bSJed Brown   pt = options->pointType;
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-point_type", "The point type", "ex2.c", pointTypes, 3, pointTypes[options->pointType], &pt, NULL));
57c4762a1bSJed Brown   options->pointType = (PointType)pt;
58c3af198dSMatthew G. Knepley   fn                 = options->funcType;
59c3af198dSMatthew G. Knepley   PetscCall(PetscOptionsEList("-func_type", "The function type", "ex2.c", funcTypes, 2, funcTypes[options->funcType], &fn, NULL));
60c3af198dSMatthew G. Knepley   options->funcType = (FuncType)fn;
61c3af198dSMatthew G. Knepley   PetscCall(PetscOptionsBool("-use_fv", "Use finite volumes, instead of finite elements", "ex2.c", options->useFV, &options->useFV, NULL));
62d0609cedSBarry Smith   PetscOptionsEnd();
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,AppCtx * ctx,DM * dm)66d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
67d71ae5a4SJacob Faibussowitsch {
68c4762a1bSJed Brown   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
709566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
719566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
729566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
CreatePoints_Centroid(DM dm,PetscInt * Np,PetscReal ** pcoords,PetscBool * pointsAllProcs,AppCtx * ctx)76d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_Centroid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
77d71ae5a4SJacob Faibussowitsch {
78c4762a1bSJed Brown   PetscSection coordSection;
79c4762a1bSJed Brown   Vec          coordsLocal;
80c4762a1bSJed Brown   PetscInt     spaceDim, p;
81c4762a1bSJed Brown   PetscMPIInt  rank;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   PetscFunctionBegin;
849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
859566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal));
869566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordSection));
879566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &spaceDim));
889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, Np));
899566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(*Np * spaceDim, pcoords));
90c4762a1bSJed Brown   for (p = 0; p < *Np; ++p) {
91c4762a1bSJed Brown     PetscScalar *coords = NULL;
92c4762a1bSJed Brown     PetscInt     size, num, n, d;
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(dm, coordSection, coordsLocal, p, &size, &coords));
95c4762a1bSJed Brown     num = size / spaceDim;
96c4762a1bSJed Brown     for (n = 0; n < num; ++n) {
97c4762a1bSJed Brown       for (d = 0; d < spaceDim; ++d) (*pcoords)[p * spaceDim + d] += PetscRealPart(coords[n * spaceDim + d]) / num;
98c4762a1bSJed Brown     }
9963a3b9bcSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, p));
100c4762a1bSJed Brown     for (d = 0; d < spaceDim; ++d) {
1019566063dSJacob Faibussowitsch       PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[p * spaceDim + d]));
1029566063dSJacob Faibussowitsch       if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", "));
103c4762a1bSJed Brown     }
1049566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n"));
1059566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(dm, coordSection, coordsLocal, p, &num, &coords));
106c4762a1bSJed Brown   }
1079566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
108c4762a1bSJed Brown   *pointsAllProcs = PETSC_FALSE;
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
CreatePoints_Grid(DM dm,PetscInt * Np,PetscReal ** pcoords,PetscBool * pointsAllProcs,AppCtx * ctx)112d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_Grid(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown   DM            da;
115c4762a1bSJed Brown   DMDALocalInfo info;
11630602db0SMatthew G. Knepley   PetscInt      N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d;
117c4762a1bSJed Brown   PetscReal    *h;
118c4762a1bSJed Brown   PetscMPIInt   rank;
119c4762a1bSJed Brown 
120362febeeSStefano Zampini   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1229566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &spaceDim));
1249566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(spaceDim, &ind));
1259566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(spaceDim, &h));
1269371c9d4SSatish Balay   h[0] = 1.0 / (N - 1);
1279371c9d4SSatish Balay   h[1] = 1.0 / (N - 1);
1289371c9d4SSatish Balay   h[2] = 1.0 / (N - 1);
1299566063dSJacob Faibussowitsch   PetscCall(DMDACreate(PetscObjectComm((PetscObject)dm), &da));
1309566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(da, dim));
1319566063dSJacob Faibussowitsch   PetscCall(DMDASetSizes(da, N, N, N));
1329566063dSJacob Faibussowitsch   PetscCall(DMDASetDof(da, 1));
1339566063dSJacob Faibussowitsch   PetscCall(DMDASetStencilWidth(da, 1));
1349566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1359566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1369566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
137c4762a1bSJed Brown   *Np = info.xm * info.ym * info.zm;
1389566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(*Np * spaceDim, pcoords));
139c4762a1bSJed Brown   for (k = info.zs; k < info.zs + info.zm; ++k) {
140c4762a1bSJed Brown     ind[2] = k;
141c4762a1bSJed Brown     for (j = info.ys; j < info.ys + info.ym; ++j) {
142c4762a1bSJed Brown       ind[1] = j;
143c4762a1bSJed Brown       for (i = info.xs; i < info.xs + info.xm; ++i, ++n) {
144c4762a1bSJed Brown         ind[0] = i;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d];
14763a3b9bcSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n));
148c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) {
1499566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d]));
1509566063dSJacob Faibussowitsch           if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", "));
151c4762a1bSJed Brown         }
1529566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n"));
153c4762a1bSJed Brown       }
154c4762a1bSJed Brown     }
155c4762a1bSJed Brown   }
1569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1579566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
1589566063dSJacob Faibussowitsch   PetscCall(PetscFree(ind));
1599566063dSJacob Faibussowitsch   PetscCall(PetscFree(h));
160c4762a1bSJed Brown   *pointsAllProcs = PETSC_FALSE;
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162c4762a1bSJed Brown }
163c4762a1bSJed Brown 
CreatePoints_GridReplicated(DM dm,PetscInt * Np,PetscReal ** pcoords,PetscBool * pointsAllProcs,AppCtx * ctx)164d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints_GridReplicated(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
165d71ae5a4SJacob Faibussowitsch {
16630602db0SMatthew G. Knepley   PetscInt    N = 3, n = 0, dim, spaceDim, i, j, k, *ind, d;
167c4762a1bSJed Brown   PetscReal  *h;
168c4762a1bSJed Brown   PetscMPIInt rank;
169c4762a1bSJed Brown 
17030602db0SMatthew G. Knepley   PetscFunctionBeginUser;
1719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1729566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1739566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &spaceDim));
1749566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(spaceDim, &ind));
1759566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(spaceDim, &h));
1769371c9d4SSatish Balay   h[0] = 1.0 / (N - 1);
1779371c9d4SSatish Balay   h[1] = 1.0 / (N - 1);
1789371c9d4SSatish Balay   h[2] = 1.0 / (N - 1);
17930602db0SMatthew G. Knepley   *Np  = N * (dim > 1 ? N : 1) * (dim > 2 ? N : 1);
1809566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(*Np * spaceDim, pcoords));
181c4762a1bSJed Brown   for (k = 0; k < N; ++k) {
182c4762a1bSJed Brown     ind[2] = k;
183c4762a1bSJed Brown     for (j = 0; j < N; ++j) {
184c4762a1bSJed Brown       ind[1] = j;
185c4762a1bSJed Brown       for (i = 0; i < N; ++i, ++n) {
186c4762a1bSJed Brown         ind[0] = i;
187c4762a1bSJed Brown 
188c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) (*pcoords)[n * spaceDim + d] = ind[d] * h[d];
18963a3b9bcSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d]Point %" PetscInt_FMT " (", rank, n));
190c4762a1bSJed Brown         for (d = 0; d < spaceDim; ++d) {
1919566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%g", (double)(*pcoords)[n * spaceDim + d]));
1929566063dSJacob Faibussowitsch           if (d < spaceDim - 1) PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ", "));
193c4762a1bSJed Brown         }
1949566063dSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, ")\n"));
195c4762a1bSJed Brown       }
196c4762a1bSJed Brown     }
197c4762a1bSJed Brown   }
1989566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
199c4762a1bSJed Brown   *pointsAllProcs = PETSC_TRUE;
2009566063dSJacob Faibussowitsch   PetscCall(PetscFree(ind));
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(h));
2023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
203c4762a1bSJed Brown }
204c4762a1bSJed Brown 
CreatePoints(DM dm,PetscInt * Np,PetscReal ** pcoords,PetscBool * pointsAllProcs,AppCtx * ctx)205d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreatePoints(DM dm, PetscInt *Np, PetscReal **pcoords, PetscBool *pointsAllProcs, AppCtx *ctx)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown   PetscFunctionBegin;
208c4762a1bSJed Brown   *pointsAllProcs = PETSC_FALSE;
209c4762a1bSJed Brown   switch (ctx->pointType) {
210d71ae5a4SJacob Faibussowitsch   case CENTROID:
211d71ae5a4SJacob Faibussowitsch     PetscCall(CreatePoints_Centroid(dm, Np, pcoords, pointsAllProcs, ctx));
212d71ae5a4SJacob Faibussowitsch     break;
213d71ae5a4SJacob Faibussowitsch   case GRID:
214d71ae5a4SJacob Faibussowitsch     PetscCall(CreatePoints_Grid(dm, Np, pcoords, pointsAllProcs, ctx));
215d71ae5a4SJacob Faibussowitsch     break;
216d71ae5a4SJacob Faibussowitsch   case GRID_REPLICATED:
217d71ae5a4SJacob Faibussowitsch     PetscCall(CreatePoints_GridReplicated(dm, Np, pcoords, pointsAllProcs, ctx));
218d71ae5a4SJacob Faibussowitsch     break;
219d71ae5a4SJacob Faibussowitsch   default:
220d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid point generation type %d", (int)ctx->pointType);
221c4762a1bSJed Brown   }
2223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
223c4762a1bSJed Brown }
224c4762a1bSJed Brown 
CreateDiscretization(DM dm,PetscInt Nc,AppCtx * ctx)225c3af198dSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, PetscInt Nc, AppCtx *ctx)
226c3af198dSMatthew G. Knepley {
227c3af198dSMatthew G. Knepley   PetscFunctionBegin;
228c3af198dSMatthew G. Knepley   if (ctx->useFV) {
229c3af198dSMatthew G. Knepley     PetscFV  fv;
230c3af198dSMatthew G. Knepley     PetscInt cdim;
231c3af198dSMatthew G. Knepley 
232c3af198dSMatthew G. Knepley     PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)dm), &fv));
233c3af198dSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)fv, "phi"));
234c3af198dSMatthew G. Knepley     PetscCall(PetscFVSetFromOptions(fv));
235c3af198dSMatthew G. Knepley     PetscCall(PetscFVSetNumComponents(fv, Nc));
236c3af198dSMatthew G. Knepley     PetscCall(DMGetCoordinateDim(dm, &cdim));
237c3af198dSMatthew G. Knepley     PetscCall(PetscFVSetSpatialDimension(fv, cdim));
238c3af198dSMatthew G. Knepley     PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fv));
239c3af198dSMatthew G. Knepley     PetscCall(PetscFVDestroy(&fv));
240c3af198dSMatthew G. Knepley     PetscCall(DMCreateDS(dm));
241c3af198dSMatthew G. Knepley   } else {
242c3af198dSMatthew G. Knepley     PetscFE        fe;
243c3af198dSMatthew G. Knepley     DMPolytopeType ct;
244c3af198dSMatthew G. Knepley     PetscInt       dim, cStart;
245c3af198dSMatthew G. Knepley 
246c3af198dSMatthew G. Knepley     PetscCall(DMGetDimension(dm, &dim));
247c3af198dSMatthew G. Knepley     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
248c3af198dSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dm, cStart, &ct));
249c3af198dSMatthew G. Knepley     PetscCall(PetscFECreateByCell(PetscObjectComm((PetscObject)dm), dim, Nc, ct, NULL, -1, &fe));
250c3af198dSMatthew G. Knepley     PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
251c3af198dSMatthew G. Knepley     PetscCall(PetscFEDestroy(&fe));
252c3af198dSMatthew G. Knepley     PetscCall(DMCreateDS(dm));
253c3af198dSMatthew G. Knepley   }
254c3af198dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
255c3af198dSMatthew G. Knepley }
256c3af198dSMatthew G. Knepley 
main(int argc,char ** argv)257d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
258d71ae5a4SJacob Faibussowitsch {
259c4762a1bSJed Brown   AppCtx ctx;
260*2a8381b2SBarry Smith   PetscErrorCode (**funcs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
261c4762a1bSJed Brown   DM                  dm;
262c4762a1bSJed Brown   DMInterpolationInfo interpolator;
263c4762a1bSJed Brown   Vec                 lu, fieldVals;
264c4762a1bSJed Brown   PetscScalar        *vals;
265c4762a1bSJed Brown   const PetscScalar  *ivals, *vcoords;
266c4762a1bSJed Brown   PetscReal          *pcoords;
267c3af198dSMatthew G. Knepley   PetscBool           pointsAllProcs = PETSC_TRUE;
26866a0a883SMatthew G. Knepley   PetscInt            dim, spaceDim, Nc, c, Np, p;
269c4762a1bSJed Brown   PetscMPIInt         rank, size;
270c4762a1bSJed Brown   PetscViewer         selfviewer;
271c4762a1bSJed Brown 
272327415f7SBarry Smith   PetscFunctionBeginUser;
2739566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2749566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
2759566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
2769566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2779566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &spaceDim));
2789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
280c4762a1bSJed Brown   /* Create points */
2819566063dSJacob Faibussowitsch   PetscCall(CreatePoints(dm, &Np, &pcoords, &pointsAllProcs, &ctx));
282c4762a1bSJed Brown   /* Create interpolator */
2839566063dSJacob Faibussowitsch   PetscCall(DMInterpolationCreate(PETSC_COMM_WORLD, &interpolator));
2849566063dSJacob Faibussowitsch   PetscCall(DMInterpolationSetDim(interpolator, spaceDim));
2859566063dSJacob Faibussowitsch   PetscCall(DMInterpolationAddPoints(interpolator, Np, pcoords));
2869566063dSJacob Faibussowitsch   PetscCall(DMInterpolationSetUp(interpolator, dm, pointsAllProcs, PETSC_FALSE));
287c4762a1bSJed Brown   /* Check locations */
28848a46eb9SPierre 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]));
2899566063dSJacob Faibussowitsch   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
2909566063dSJacob Faibussowitsch   PetscCall(VecView(interpolator->coords, PETSC_VIEWER_STDOUT_WORLD));
29166a0a883SMatthew G. Knepley   Nc = dim;
292c3af198dSMatthew G. Knepley   PetscCall(CreateDiscretization(dm, Nc, &ctx));
293c4762a1bSJed Brown   /* Create function */
2949566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(Nc, &funcs, Nc, &vals));
295c3af198dSMatthew G. Knepley   switch (ctx.funcType) {
296c3af198dSMatthew G. Knepley   case CONSTANT:
297c3af198dSMatthew G. Knepley     for (c = 0; c < Nc; ++c) funcs[c] = constant;
298c3af198dSMatthew G. Knepley     break;
299c3af198dSMatthew G. Knepley   case LINEAR:
300c4762a1bSJed Brown     for (c = 0; c < Nc; ++c) funcs[c] = linear;
301c3af198dSMatthew G. Knepley     break;
302c3af198dSMatthew G. Knepley   default:
303c3af198dSMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid function type: %d", (int)ctx.funcType);
304c3af198dSMatthew G. Knepley   }
3059566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lu));
3069566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, lu));
3079566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
3089566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer));
3099566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(selfviewer, "[%d]solution\n", rank));
3109566063dSJacob Faibussowitsch   PetscCall(VecView(lu, selfviewer));
3119566063dSJacob Faibussowitsch   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &selfviewer));
3129566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
313c4762a1bSJed Brown   /* Check interpolant */
3149566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, interpolator->n * Nc, &fieldVals));
3159566063dSJacob Faibussowitsch   PetscCall(DMInterpolationSetDof(interpolator, Nc));
3169566063dSJacob Faibussowitsch   PetscCall(DMInterpolationEvaluate(interpolator, dm, lu, fieldVals));
317c4762a1bSJed Brown   for (p = 0; p < size; ++p) {
318c4762a1bSJed Brown     if (p == rank) {
3199566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d]Field values\n", rank));
3209566063dSJacob Faibussowitsch       PetscCall(VecView(fieldVals, PETSC_VIEWER_STDOUT_SELF));
321c4762a1bSJed Brown     }
3229566063dSJacob Faibussowitsch     PetscCall(PetscBarrier((PetscObject)dm));
323c4762a1bSJed Brown   }
3249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(interpolator->coords, &vcoords));
3259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(fieldVals, &ivals));
326c4762a1bSJed Brown   for (p = 0; p < interpolator->n; ++p) {
327c4762a1bSJed Brown     for (c = 0; c < Nc; ++c) {
328c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
329c4762a1bSJed Brown       PetscReal vcoordsReal[3];
330c4762a1bSJed Brown 
3313ba16761SJacob Faibussowitsch       for (PetscInt i = 0; i < spaceDim; i++) vcoordsReal[i] = PetscRealPart(vcoords[p * spaceDim + i]);
332c4762a1bSJed Brown #else
333c4762a1bSJed Brown       const PetscReal *vcoordsReal = &vcoords[p * spaceDim];
334c4762a1bSJed Brown #endif
3353ba16761SJacob Faibussowitsch       PetscCall((*funcs[c])(dim, 0.0, vcoordsReal, Nc, vals, NULL));
3363ba16761SJacob 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);
337c4762a1bSJed Brown     }
338c4762a1bSJed Brown   }
3399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(interpolator->coords, &vcoords));
3409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(fieldVals, &ivals));
341c4762a1bSJed Brown   /* Cleanup */
3429566063dSJacob Faibussowitsch   PetscCall(PetscFree(pcoords));
3439566063dSJacob Faibussowitsch   PetscCall(PetscFree2(funcs, vals));
3449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fieldVals));
3459566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lu));
3469566063dSJacob Faibussowitsch   PetscCall(DMInterpolationDestroy(&interpolator));
3479566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3489566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
349b122ec5aSJacob Faibussowitsch   return 0;
350c4762a1bSJed Brown }
351c4762a1bSJed Brown 
352c4762a1bSJed Brown /*TEST
353c4762a1bSJed Brown 
35430602db0SMatthew G. Knepley   testset:
35530602db0SMatthew G. Knepley     requires: ctetgen
35630602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -petscspace_degree 1
35730602db0SMatthew G. Knepley 
358c4762a1bSJed Brown     test:
359c4762a1bSJed Brown       suffix: 0
360c4762a1bSJed Brown     test:
361c4762a1bSJed Brown       suffix: 1
36230602db0SMatthew G. Knepley       args: -dm_refine 2
363c4762a1bSJed Brown     test:
364c4762a1bSJed Brown       suffix: 2
365c4762a1bSJed Brown       nsize: 2
366e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple
367c4762a1bSJed Brown     test:
368c4762a1bSJed Brown       suffix: 3
369c4762a1bSJed Brown       nsize: 2
370e600fa54SMatthew G. Knepley       args: -dm_refine 2 -petscpartitioner_type simple
371c4762a1bSJed Brown     test:
372c4762a1bSJed Brown       suffix: 4
373c4762a1bSJed Brown       nsize: 5
374e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple
375c4762a1bSJed Brown     test:
376c4762a1bSJed Brown       suffix: 5
377c4762a1bSJed Brown       nsize: 5
378e600fa54SMatthew G. Knepley       args: -dm_refine 2 -petscpartitioner_type simple
379c4762a1bSJed Brown     test:
380c4762a1bSJed Brown       suffix: 6
38130602db0SMatthew G. Knepley       args: -point_type grid
382c4762a1bSJed Brown     test:
383c4762a1bSJed Brown       suffix: 7
38430602db0SMatthew G. Knepley       args: -dm_refine 2 -point_type grid
385c4762a1bSJed Brown     test:
386c4762a1bSJed Brown       suffix: 8
387c4762a1bSJed Brown       nsize: 2
388e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -point_type grid
389c4762a1bSJed Brown     test:
390c4762a1bSJed Brown       suffix: 9
39130602db0SMatthew G. Knepley       args: -point_type grid_replicated
392c4762a1bSJed Brown     test:
393c4762a1bSJed Brown       suffix: 10
394c4762a1bSJed Brown       nsize: 2
395e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -point_type grid_replicated
396c4762a1bSJed Brown     test:
397c4762a1bSJed Brown       suffix: 11
398c4762a1bSJed Brown       nsize: 2
399e600fa54SMatthew G. Knepley       args: -dm_refine 2 -petscpartitioner_type simple -point_type grid_replicated
400c4762a1bSJed Brown 
401c3af198dSMatthew G. Knepley   test:
402c3af198dSMatthew G. Knepley     suffix: fv_0
403c3af198dSMatthew G. Knepley     requires: triangle
404c3af198dSMatthew G. Knepley     args: -use_fv -func_type constant
405c3af198dSMatthew G. Knepley 
406c4762a1bSJed Brown TEST*/
407