xref: /petsc/src/dm/impls/forest/tests/ex2.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Create a mesh, refine and coarsen simultaneously, and transfer a field\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscds.h>
4*c4762a1bSJed Brown #include <petscdmplex.h>
5*c4762a1bSJed Brown #include <petscdmforest.h>
6*c4762a1bSJed Brown #include <petscoptions.h>
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown static PetscErrorCode AddIdentityLabel(DM dm)
9*c4762a1bSJed Brown {
10*c4762a1bSJed Brown   PetscInt       pStart,pEnd,p;
11*c4762a1bSJed Brown   PetscErrorCode ierr;
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown   PetscFunctionBegin;
14*c4762a1bSJed Brown   ierr = DMCreateLabel(dm, "identity");CHKERRQ(ierr);
15*c4762a1bSJed Brown   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
16*c4762a1bSJed Brown   for (p = pStart; p < pEnd; p++) {ierr = DMSetLabelValue(dm, "identity", p, p);CHKERRQ(ierr);}
17*c4762a1bSJed Brown   PetscFunctionReturn(0);
18*c4762a1bSJed Brown }
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown static PetscErrorCode CreateAdaptivityLabel(DM forest,DMLabel *adaptLabel)
21*c4762a1bSJed Brown {
22*c4762a1bSJed Brown   DMLabel        identLabel;
23*c4762a1bSJed Brown   PetscInt       cStart, cEnd, c;
24*c4762a1bSJed Brown   PetscErrorCode ierr;
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown   PetscFunctionBegin;
27*c4762a1bSJed Brown   ierr = DMLabelCreate(PETSC_COMM_SELF,"adapt",adaptLabel);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = DMLabelSetDefaultValue(*adaptLabel,DM_ADAPT_COARSEN);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = DMGetLabel(forest,"identity",&identLabel);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = DMForestGetCellChart(forest,&cStart,&cEnd);CHKERRQ(ierr);
31*c4762a1bSJed Brown   for (c = cStart; c < cEnd; c++) {
32*c4762a1bSJed Brown     PetscInt basePoint;
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown     ierr = DMLabelGetValue(identLabel,c,&basePoint);CHKERRQ(ierr);
35*c4762a1bSJed Brown     if (!basePoint) {ierr = DMLabelSetValue(*adaptLabel,c,DM_ADAPT_REFINE);CHKERRQ(ierr);}
36*c4762a1bSJed Brown   }
37*c4762a1bSJed Brown   PetscFunctionReturn(0);
38*c4762a1bSJed Brown }
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown static PetscErrorCode LinearFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
41*c4762a1bSJed Brown {
42*c4762a1bSJed Brown   PetscFunctionBeginUser;
43*c4762a1bSJed Brown   u[0] = (x[0] * 2.0 + 1.) + (x[1] * 20.0 + 10.) + ((dim == 3) ? (x[2] * 200.0 + 100.) : 0.);
44*c4762a1bSJed Brown   PetscFunctionReturn(0);
45*c4762a1bSJed Brown }
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown static PetscErrorCode MultiaffineFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
48*c4762a1bSJed Brown {
49*c4762a1bSJed Brown   PetscFunctionBeginUser;
50*c4762a1bSJed Brown   u[0] = (x[0] * 1.0 + 2.0) * (x[1] * 3.0 - 4.0) * ((dim == 3) ? (x[2] * 5.0 + 6.0) : 1.);
51*c4762a1bSJed Brown   PetscFunctionReturn(0);
52*c4762a1bSJed Brown }
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown static PetscErrorCode CoordsFunction(PetscInt dim,PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
55*c4762a1bSJed Brown {
56*c4762a1bSJed Brown   PetscInt f;
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   PetscFunctionBeginUser;
59*c4762a1bSJed Brown   for (f=0;f<Nf;f++) u[f] = x[f];
60*c4762a1bSJed Brown   PetscFunctionReturn(0);
61*c4762a1bSJed Brown }
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown typedef struct _bc_func_ctx
64*c4762a1bSJed Brown {
65*c4762a1bSJed Brown   PetscErrorCode (*func) (PetscInt,PetscReal,const PetscReal [], PetscInt, PetscScalar [], void *);
66*c4762a1bSJed Brown   PetscInt dim;
67*c4762a1bSJed Brown   PetscInt Nf;
68*c4762a1bSJed Brown   void *ctx;
69*c4762a1bSJed Brown }
70*c4762a1bSJed Brown bc_func_ctx;
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown static PetscErrorCode bc_func_fv (PetscReal time, const PetscReal *c, const PetscReal *n, const PetscScalar *xI, PetscScalar *xG, void *ctx)
73*c4762a1bSJed Brown {
74*c4762a1bSJed Brown   bc_func_ctx    *bcCtx;
75*c4762a1bSJed Brown   PetscErrorCode ierr;
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   PetscFunctionBegin;
78*c4762a1bSJed Brown   bcCtx = (bc_func_ctx *) ctx;
79*c4762a1bSJed Brown   ierr = (bcCtx->func)(bcCtx->dim,time,c,bcCtx->Nf,xG,bcCtx->ctx);CHKERRQ(ierr);
80*c4762a1bSJed Brown   PetscFunctionReturn(0);
81*c4762a1bSJed Brown }
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown static PetscErrorCode IdentifyBadPoints (DM dm, Vec vec, PetscReal tol)
84*c4762a1bSJed Brown {
85*c4762a1bSJed Brown   DM             dmplex;
86*c4762a1bSJed Brown   PetscInt       p, pStart, pEnd, maxDof;
87*c4762a1bSJed Brown   Vec            vecLocal;
88*c4762a1bSJed Brown   DMLabel        depthLabel;
89*c4762a1bSJed Brown   PetscSection   section;
90*c4762a1bSJed Brown   PetscErrorCode ierr;
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   PetscFunctionBegin;
93*c4762a1bSJed Brown   ierr = DMCreateLocalVector(dm, &vecLocal);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vecLocal);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vecLocal);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = DMConvert(dm ,DMPLEX, &dmplex);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = DMPlexGetChart(dmplex, &pStart, &pEnd);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = DMPlexGetDepthLabel(dmplex, &depthLabel);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = DMGetLocalSection(dmplex, &section);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = PetscSectionGetMaxDof(section, &maxDof);CHKERRQ(ierr);
101*c4762a1bSJed Brown   for (p = pStart; p < pEnd; p++) {
102*c4762a1bSJed Brown     PetscInt     s, c, cSize, parent, childID, numChildren;
103*c4762a1bSJed Brown     PetscInt     cl, closureSize, *closure = NULL;
104*c4762a1bSJed Brown     PetscScalar *values = NULL;
105*c4762a1bSJed Brown     PetscBool    bad = PETSC_FALSE;
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown     ierr = VecGetValuesSection(vecLocal, section, p, &values);CHKERRQ(ierr);
108*c4762a1bSJed Brown     ierr = PetscSectionGetDof(section, p, &cSize);CHKERRQ(ierr);
109*c4762a1bSJed Brown     for (c = 0; c < cSize; c++) {
110*c4762a1bSJed Brown       PetscReal absDiff = PetscAbsScalar(values[c]);CHKERRQ(ierr);
111*c4762a1bSJed Brown       if (absDiff > tol) {bad = PETSC_TRUE; break;}
112*c4762a1bSJed Brown     }
113*c4762a1bSJed Brown     if (!bad) continue;
114*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "Bad point %D\n", p);CHKERRQ(ierr);
115*c4762a1bSJed Brown     ierr = DMLabelGetValue(depthLabel, p, &s);CHKERRQ(ierr);
116*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "  Depth %D\n", s);CHKERRQ(ierr);
117*c4762a1bSJed Brown     ierr = DMPlexGetTransitiveClosure(dmplex, p, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
118*c4762a1bSJed Brown     for (cl = 0; cl < closureSize; cl++) {
119*c4762a1bSJed Brown       PetscInt cp = closure[2 * cl];
120*c4762a1bSJed Brown       ierr = DMPlexGetTreeParent(dmplex, cp, &parent, &childID);CHKERRQ(ierr);
121*c4762a1bSJed Brown       if (parent != cp) {
122*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "  Closure point %D (%D) child of %D (ID %D)\n", cl, cp, parent, childID);CHKERRQ(ierr);
123*c4762a1bSJed Brown       }
124*c4762a1bSJed Brown       ierr = DMPlexGetTreeChildren(dmplex, cp, &numChildren, NULL);CHKERRQ(ierr);
125*c4762a1bSJed Brown       if (numChildren) {
126*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "  Closure point %D (%D) is parent\n", cl, cp);CHKERRQ(ierr);
127*c4762a1bSJed Brown       }
128*c4762a1bSJed Brown     }
129*c4762a1bSJed Brown     ierr = DMPlexRestoreTransitiveClosure(dmplex, p, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
130*c4762a1bSJed Brown     for (c = 0; c < cSize; c++) {
131*c4762a1bSJed Brown       PetscReal absDiff = PetscAbsScalar(values[c]);CHKERRQ(ierr);
132*c4762a1bSJed Brown       if (absDiff > tol) {
133*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "  Bad dof %D\n", c);CHKERRQ(ierr);
134*c4762a1bSJed Brown       }
135*c4762a1bSJed Brown     }
136*c4762a1bSJed Brown   }
137*c4762a1bSJed Brown   ierr = DMDestroy(&dmplex);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = VecDestroy(&vecLocal);CHKERRQ(ierr);
139*c4762a1bSJed Brown   PetscFunctionReturn(0);
140*c4762a1bSJed Brown }
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown int main(int argc, char **argv)
143*c4762a1bSJed Brown {
144*c4762a1bSJed Brown   MPI_Comm       comm;
145*c4762a1bSJed Brown   DM             base, preForest, postForest;
146*c4762a1bSJed Brown   PetscInt       dim = 2, Nf = 1;
147*c4762a1bSJed Brown   PetscInt       step, adaptSteps = 1;
148*c4762a1bSJed Brown   PetscInt       preCount, postCount;
149*c4762a1bSJed Brown   Vec            preVec, postVecTransfer, postVecExact;
150*c4762a1bSJed Brown   PetscErrorCode (*funcs[1]) (PetscInt,PetscReal,const PetscReal [],PetscInt,PetscScalar [], void *) = {MultiaffineFunction};
151*c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN];
152*c4762a1bSJed Brown   void           *ctxs[1] = {NULL};
153*c4762a1bSJed Brown   const PetscInt cells[] = {3, 3, 3};
154*c4762a1bSJed Brown   PetscReal      diff, tol = PETSC_SMALL;
155*c4762a1bSJed Brown   DMBoundaryType periodicity[] = {DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE};
156*c4762a1bSJed Brown   PetscBool      linear = PETSC_FALSE;
157*c4762a1bSJed Brown   PetscBool      coords = PETSC_FALSE;
158*c4762a1bSJed Brown   PetscBool      useFV = PETSC_FALSE;
159*c4762a1bSJed Brown   PetscBool      conv = PETSC_FALSE;
160*c4762a1bSJed Brown   PetscBool      distribute_base = PETSC_FALSE;
161*c4762a1bSJed Brown   PetscBool      transfer_from_base[2] = {PETSC_TRUE,PETSC_FALSE};
162*c4762a1bSJed Brown   PetscBool      use_bcs = PETSC_TRUE;
163*c4762a1bSJed Brown   PetscBool      periodic = PETSC_FALSE;
164*c4762a1bSJed Brown   bc_func_ctx    bcCtx;
165*c4762a1bSJed Brown   DMLabel        adaptLabel;
166*c4762a1bSJed Brown   size_t         len;
167*c4762a1bSJed Brown   PetscErrorCode ierr;
168*c4762a1bSJed Brown 
169*c4762a1bSJed Brown   filename[0] = '\0';
170*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
171*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
172*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "DMForestTransferVec() Test Options", "DMFOREST");CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The dimension (2 or 3)", "ex2.c", dim, &dim, NULL,2,3);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = PetscOptionsBool("-linear","Transfer a simple linear function", "ex2.c", linear, &linear, NULL);CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = PetscOptionsBool("-coords","Transfer a simple coordinate function", "ex2.c", coords, &coords, NULL);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = PetscOptionsBool("-use_fv","Use a finite volume approximation", "ex2.c", useFV, &useFV, NULL);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = PetscOptionsBool("-test_convert","Test conversion to DMPLEX",NULL,conv,&conv,NULL);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = PetscOptionsString("-filename", "Read the base mesh from file", "ex2.c", filename, filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = PetscOptionsBool("-distribute_base","Distribute base DM", "ex2.c", distribute_base, &distribute_base, NULL);CHKERRQ(ierr);
180*c4762a1bSJed Brown   ierr = PetscOptionsBool("-transfer_from_base","Transfer a vector from base DM to DMForest", "ex2.c", transfer_from_base[0], &transfer_from_base[0], NULL);CHKERRQ(ierr);
181*c4762a1bSJed Brown   transfer_from_base[1] = transfer_from_base[0];
182*c4762a1bSJed Brown   ierr = PetscOptionsBool("-transfer_from_base_steps","Transfer a vector from base DM to the latest DMForest after the adaptivity steps", "ex2.c", transfer_from_base[1], &transfer_from_base[1], NULL);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = PetscOptionsBool("-use_bcs","Use dirichlet boundary conditions", "ex2.c", use_bcs, &use_bcs, NULL);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-adapt_steps","Number of adaptivity steps", "ex2.c", adaptSteps, &adaptSteps, NULL,0);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = PetscOptionsBool("-periodic","Use periodic box mesh", "ex2.c", periodic, &periodic, NULL);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   tol = PetscMax(1.e-10,tol); /* XXX fix for quadruple precision -> why do I need to do this? */
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown   if (periodic) periodicity[0] = periodicity[1] = periodicity[2] = DM_BOUNDARY_PERIODIC;
191*c4762a1bSJed Brown 
192*c4762a1bSJed Brown   if (linear) {
193*c4762a1bSJed Brown     funcs[0] = LinearFunction;
194*c4762a1bSJed Brown   }
195*c4762a1bSJed Brown   if (coords) {
196*c4762a1bSJed Brown     funcs[0] = CoordsFunction;
197*c4762a1bSJed Brown     Nf = dim;
198*c4762a1bSJed Brown   }
199*c4762a1bSJed Brown 
200*c4762a1bSJed Brown   bcCtx.func = funcs[0];
201*c4762a1bSJed Brown   bcCtx.dim  = dim;
202*c4762a1bSJed Brown   bcCtx.Nf   = Nf;
203*c4762a1bSJed Brown   bcCtx.ctx  = NULL;
204*c4762a1bSJed Brown 
205*c4762a1bSJed Brown   /* the base mesh */
206*c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
207*c4762a1bSJed Brown   if (!len) {
208*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(comm,dim,PETSC_FALSE,cells,NULL,NULL,periodicity,PETSC_TRUE,&base);CHKERRQ(ierr);
209*c4762a1bSJed Brown   } else {
210*c4762a1bSJed Brown     DM tdm = NULL;
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm,filename,PETSC_TRUE,&base);CHKERRQ(ierr);
213*c4762a1bSJed Brown     ierr = DMGetDimension(base,&dim);CHKERRQ(ierr);
214*c4762a1bSJed Brown     ierr = DMPlexRefineSimplexToTensor(base,&tdm); CHKERRQ(ierr);
215*c4762a1bSJed Brown     if (tdm) {
216*c4762a1bSJed Brown       ierr = DMDestroy(&base);CHKERRQ(ierr);
217*c4762a1bSJed Brown       base = tdm;
218*c4762a1bSJed Brown     }
219*c4762a1bSJed Brown     use_bcs = PETSC_FALSE;
220*c4762a1bSJed Brown   }
221*c4762a1bSJed Brown   ierr = AddIdentityLabel(base);CHKERRQ(ierr);
222*c4762a1bSJed Brown   if (distribute_base) {
223*c4762a1bSJed Brown     DM               baseParallel;
224*c4762a1bSJed Brown     PetscPartitioner part;
225*c4762a1bSJed Brown 
226*c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(base,&part);CHKERRQ(ierr);
227*c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
228*c4762a1bSJed Brown     ierr = DMPlexDistribute(base,0,NULL,&baseParallel);CHKERRQ(ierr);
229*c4762a1bSJed Brown     if (baseParallel) {
230*c4762a1bSJed Brown       ierr = DMDestroy(&base);CHKERRQ(ierr);
231*c4762a1bSJed Brown       base = baseParallel;
232*c4762a1bSJed Brown     }
233*c4762a1bSJed Brown   }
234*c4762a1bSJed Brown   if (useFV) {
235*c4762a1bSJed Brown     PetscFV      fv;
236*c4762a1bSJed Brown     PetscLimiter limiter;
237*c4762a1bSJed Brown     DM           baseFV;
238*c4762a1bSJed Brown 
239*c4762a1bSJed Brown     ierr = DMPlexConstructGhostCells(base,NULL,NULL,&baseFV);CHKERRQ(ierr);
240*c4762a1bSJed Brown     ierr = DMDestroy(&base);CHKERRQ(ierr);
241*c4762a1bSJed Brown     base = baseFV;
242*c4762a1bSJed Brown     ierr = PetscFVCreate(comm, &fv);CHKERRQ(ierr);
243*c4762a1bSJed Brown     ierr = PetscFVSetSpatialDimension(fv,dim);CHKERRQ(ierr);
244*c4762a1bSJed Brown     ierr = PetscFVSetType(fv,PETSCFVLEASTSQUARES);CHKERRQ(ierr);
245*c4762a1bSJed Brown     ierr = PetscFVSetNumComponents(fv,Nf);CHKERRQ(ierr);
246*c4762a1bSJed Brown     ierr = PetscLimiterCreate(comm,&limiter);CHKERRQ(ierr);
247*c4762a1bSJed Brown     ierr = PetscLimiterSetType(limiter,PETSCLIMITERNONE);CHKERRQ(ierr);
248*c4762a1bSJed Brown     ierr = PetscFVSetLimiter(fv,limiter);CHKERRQ(ierr);
249*c4762a1bSJed Brown     ierr = PetscLimiterDestroy(&limiter);CHKERRQ(ierr);
250*c4762a1bSJed Brown     ierr = PetscFVSetFromOptions(fv);CHKERRQ(ierr);
251*c4762a1bSJed Brown     ierr = DMSetField(base,0,NULL,(PetscObject)fv);CHKERRQ(ierr);
252*c4762a1bSJed Brown     ierr = PetscFVDestroy(&fv);CHKERRQ(ierr);
253*c4762a1bSJed Brown   } else {
254*c4762a1bSJed Brown     PetscFE fe;
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown     ierr = PetscFECreateDefault(comm,dim,Nf,PETSC_FALSE,NULL,PETSC_DEFAULT,&fe);CHKERRQ(ierr);
257*c4762a1bSJed Brown     ierr = DMSetField(base,0,NULL,(PetscObject)fe);CHKERRQ(ierr);
258*c4762a1bSJed Brown     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
259*c4762a1bSJed Brown   }
260*c4762a1bSJed Brown   ierr = DMCreateDS(base);CHKERRQ(ierr);
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   if (use_bcs) {
263*c4762a1bSJed Brown     PetscDS  prob;
264*c4762a1bSJed Brown     PetscInt ids[]   = {1, 2, 3, 4, 5, 6};
265*c4762a1bSJed Brown 
266*c4762a1bSJed Brown     ierr = DMGetDS(base,&prob);CHKERRQ(ierr);
267*c4762a1bSJed Brown     ierr = PetscDSAddBoundary(prob,DM_BC_ESSENTIAL, "bc", "marker", 0, 0, NULL, useFV ? (void(*)(void)) bc_func_fv : (void(*)(void)) funcs[0], 2 * dim, ids, useFV ? (void *) &bcCtx : NULL);CHKERRQ(ierr);
268*c4762a1bSJed Brown   }
269*c4762a1bSJed Brown   ierr = DMViewFromOptions(base,NULL,"-dm_base_view");CHKERRQ(ierr);
270*c4762a1bSJed Brown 
271*c4762a1bSJed Brown   /* the pre adaptivity forest */
272*c4762a1bSJed Brown   ierr = DMCreate(comm,&preForest);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = DMSetType(preForest,(dim == 2) ? DMP4EST : DMP8EST);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = DMCopyDisc(base,preForest);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = DMForestSetBaseDM(preForest,base);CHKERRQ(ierr);
276*c4762a1bSJed Brown   ierr = DMForestSetMinimumRefinement(preForest,0);CHKERRQ(ierr);
277*c4762a1bSJed Brown   ierr = DMForestSetInitialRefinement(preForest,1);CHKERRQ(ierr);
278*c4762a1bSJed Brown   ierr = DMSetFromOptions(preForest);CHKERRQ(ierr);
279*c4762a1bSJed Brown   ierr = DMSetUp(preForest);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = DMViewFromOptions(preForest,NULL,"-dm_pre_view");CHKERRQ(ierr);
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown   /* the pre adaptivity field */
283*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(preForest,&preVec);CHKERRQ(ierr);
284*c4762a1bSJed Brown   ierr = DMProjectFunction(preForest,0.,funcs,ctxs,INSERT_VALUES,preVec);CHKERRQ(ierr);
285*c4762a1bSJed Brown   ierr = VecViewFromOptions(preVec,NULL,"-vec_pre_view");CHKERRQ(ierr);
286*c4762a1bSJed Brown 
287*c4762a1bSJed Brown   /* communicate between base and pre adaptivity forest */
288*c4762a1bSJed Brown   if (transfer_from_base[0]) {
289*c4762a1bSJed Brown     Vec baseVec, baseVecMapped;
290*c4762a1bSJed Brown 
291*c4762a1bSJed Brown     ierr = DMGetGlobalVector(base,&baseVec);CHKERRQ(ierr);
292*c4762a1bSJed Brown     ierr = DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec);CHKERRQ(ierr);
293*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)baseVec,"Function Base");CHKERRQ(ierr);
294*c4762a1bSJed Brown     ierr = VecViewFromOptions(baseVec,NULL,"-vec_base_view");CHKERRQ(ierr);
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown     ierr = DMGetGlobalVector(preForest,&baseVecMapped);CHKERRQ(ierr);
297*c4762a1bSJed Brown     ierr = DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped);CHKERRQ(ierr);
298*c4762a1bSJed Brown     ierr = VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view");CHKERRQ(ierr);
299*c4762a1bSJed Brown 
300*c4762a1bSJed Brown     /* compare */
301*c4762a1bSJed Brown     ierr = VecAXPY(baseVecMapped,-1.,preVec);CHKERRQ(ierr);
302*c4762a1bSJed Brown     ierr = VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view");CHKERRQ(ierr);
303*c4762a1bSJed Brown     ierr = VecNorm(baseVecMapped,NORM_2,&diff);CHKERRQ(ierr);
304*c4762a1bSJed Brown 
305*c4762a1bSJed Brown     /* output */
306*c4762a1bSJed Brown     if (diff < tol) {
307*c4762a1bSJed Brown       ierr = PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n");CHKERRQ(ierr);
308*c4762a1bSJed Brown     } else {
309*c4762a1bSJed Brown       ierr = PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol);CHKERRQ(ierr);
310*c4762a1bSJed Brown     }
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown     ierr = DMRestoreGlobalVector(base,&baseVec);CHKERRQ(ierr);
313*c4762a1bSJed Brown     ierr = DMRestoreGlobalVector(preForest,&baseVecMapped);CHKERRQ(ierr);
314*c4762a1bSJed Brown   }
315*c4762a1bSJed Brown 
316*c4762a1bSJed Brown   for (step = 0; step < adaptSteps; ++step) {
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown     if (!transfer_from_base[1]) {
319*c4762a1bSJed Brown       ierr = PetscObjectGetReference((PetscObject)preForest,&preCount);CHKERRQ(ierr);
320*c4762a1bSJed Brown     }
321*c4762a1bSJed Brown 
322*c4762a1bSJed Brown     /* adapt */
323*c4762a1bSJed Brown     ierr = CreateAdaptivityLabel(preForest,&adaptLabel);CHKERRQ(ierr);
324*c4762a1bSJed Brown     ierr = DMForestTemplate(preForest,comm,&postForest);CHKERRQ(ierr);
325*c4762a1bSJed Brown     if (step) { ierr = DMForestSetAdaptivityLabel(postForest,adaptLabel);CHKERRQ(ierr); }
326*c4762a1bSJed Brown     ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
327*c4762a1bSJed Brown     ierr = DMSetUp(postForest);CHKERRQ(ierr);
328*c4762a1bSJed Brown     ierr = DMViewFromOptions(postForest,NULL,"-dm_post_view");CHKERRQ(ierr);
329*c4762a1bSJed Brown 
330*c4762a1bSJed Brown     /* transfer */
331*c4762a1bSJed Brown     ierr = DMCreateGlobalVector(postForest,&postVecTransfer);CHKERRQ(ierr);
332*c4762a1bSJed Brown     ierr = DMForestTransferVec(preForest,preVec,postForest,postVecTransfer,PETSC_TRUE,0.0);CHKERRQ(ierr);
333*c4762a1bSJed Brown     ierr = VecViewFromOptions(postVecTransfer,NULL,"-vec_post_transfer_view");CHKERRQ(ierr);
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown     /* the exact post adaptivity field */
336*c4762a1bSJed Brown     ierr = DMCreateGlobalVector(postForest,&postVecExact);CHKERRQ(ierr);
337*c4762a1bSJed Brown     ierr = DMProjectFunction(postForest,0.,funcs,ctxs,INSERT_VALUES,postVecExact);CHKERRQ(ierr);
338*c4762a1bSJed Brown     ierr = VecViewFromOptions(postVecExact,NULL,"-vec_post_exact_view");CHKERRQ(ierr);
339*c4762a1bSJed Brown 
340*c4762a1bSJed Brown     /* compare */
341*c4762a1bSJed Brown     ierr = VecAXPY(postVecExact,-1.,postVecTransfer);CHKERRQ(ierr);
342*c4762a1bSJed Brown     ierr = VecViewFromOptions(postVecExact,NULL,"-vec_diff_view");CHKERRQ(ierr);
343*c4762a1bSJed Brown     ierr = VecNorm(postVecExact,NORM_2,&diff);CHKERRQ(ierr);
344*c4762a1bSJed Brown 
345*c4762a1bSJed Brown     /* output */
346*c4762a1bSJed Brown     if (diff < tol) {
347*c4762a1bSJed Brown       ierr = PetscPrintf(comm,"DMForestTransferVec() passes.\n");CHKERRQ(ierr);
348*c4762a1bSJed Brown     } else {
349*c4762a1bSJed Brown       ierr = PetscPrintf(comm,"DMForestTransferVec() fails with error %g and tolerance %g\n",(double)diff,(double)tol);CHKERRQ(ierr);
350*c4762a1bSJed Brown       ierr = IdentifyBadPoints(postForest, postVecExact, tol);CHKERRQ(ierr);
351*c4762a1bSJed Brown     }
352*c4762a1bSJed Brown     ierr = VecDestroy(&postVecExact);CHKERRQ(ierr);
353*c4762a1bSJed Brown 
354*c4762a1bSJed Brown     /* disconnect preForest from postForest if we don't test the transfer throughout the entire refinement process */
355*c4762a1bSJed Brown     if (!transfer_from_base[1]) {
356*c4762a1bSJed Brown       ierr = DMForestSetAdaptivityForest(postForest,NULL);CHKERRQ(ierr);
357*c4762a1bSJed Brown       ierr = PetscObjectGetReference((PetscObject)preForest,&postCount);CHKERRQ(ierr);
358*c4762a1bSJed Brown       if (postCount != preCount) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Adaptation not memory neutral: reference count increase from %d to %d\n",preCount,postCount);
359*c4762a1bSJed Brown     }
360*c4762a1bSJed Brown 
361*c4762a1bSJed Brown     if (conv) {
362*c4762a1bSJed Brown       DM dmConv;
363*c4762a1bSJed Brown 
364*c4762a1bSJed Brown       ierr = DMConvert(postForest,DMPLEX,&dmConv);CHKERRQ(ierr);
365*c4762a1bSJed Brown       ierr = DMViewFromOptions(dmConv,NULL,"-dm_conv_view");CHKERRQ(ierr);
366*c4762a1bSJed Brown       ierr = DMPlexCheckCellShape(dmConv,PETSC_TRUE,PETSC_DETERMINE);CHKERRQ(ierr);
367*c4762a1bSJed Brown       ierr = DMDestroy(&dmConv);CHKERRQ(ierr);
368*c4762a1bSJed Brown     }
369*c4762a1bSJed Brown 
370*c4762a1bSJed Brown     ierr = VecDestroy(&preVec);CHKERRQ(ierr);
371*c4762a1bSJed Brown     ierr = DMDestroy(&preForest);CHKERRQ(ierr);
372*c4762a1bSJed Brown 
373*c4762a1bSJed Brown     preVec    = postVecTransfer;
374*c4762a1bSJed Brown     preForest = postForest;
375*c4762a1bSJed Brown   }
376*c4762a1bSJed Brown 
377*c4762a1bSJed Brown   if (transfer_from_base[1]) {
378*c4762a1bSJed Brown     Vec baseVec, baseVecMapped;
379*c4762a1bSJed Brown 
380*c4762a1bSJed Brown     /* communicate between base and last adapted forest */
381*c4762a1bSJed Brown     ierr = DMGetGlobalVector(base,&baseVec);CHKERRQ(ierr);
382*c4762a1bSJed Brown     ierr = DMProjectFunction(base,0.,funcs,ctxs,INSERT_VALUES,baseVec);CHKERRQ(ierr);
383*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject)baseVec,"Function Base");CHKERRQ(ierr);
384*c4762a1bSJed Brown     ierr = VecViewFromOptions(baseVec,NULL,"-vec_base_view");CHKERRQ(ierr);
385*c4762a1bSJed Brown 
386*c4762a1bSJed Brown     ierr = DMGetGlobalVector(preForest,&baseVecMapped);CHKERRQ(ierr);
387*c4762a1bSJed Brown     ierr = DMForestTransferVecFromBase(preForest,baseVec,baseVecMapped);CHKERRQ(ierr);
388*c4762a1bSJed Brown     ierr = VecViewFromOptions(baseVecMapped,NULL,"-vec_map_base_view");CHKERRQ(ierr);
389*c4762a1bSJed Brown 
390*c4762a1bSJed Brown     /* compare */
391*c4762a1bSJed Brown     ierr = VecAXPY(baseVecMapped,-1.,preVec);CHKERRQ(ierr);
392*c4762a1bSJed Brown     ierr = VecViewFromOptions(baseVecMapped,NULL,"-vec_map_diff_view");CHKERRQ(ierr);
393*c4762a1bSJed Brown     ierr = VecNorm(baseVecMapped,NORM_2,&diff);CHKERRQ(ierr);
394*c4762a1bSJed Brown 
395*c4762a1bSJed Brown     /* output */
396*c4762a1bSJed Brown     if (diff < tol) {
397*c4762a1bSJed Brown       ierr = PetscPrintf(comm,"DMForestTransferVecFromBase() passes.\n");CHKERRQ(ierr);
398*c4762a1bSJed Brown     } else {
399*c4762a1bSJed Brown       ierr = PetscPrintf(comm,"DMForestTransferVecFromBase() fails with error %g and tolerance %g\n",(double)diff,(double)tol);CHKERRQ(ierr);
400*c4762a1bSJed Brown     }
401*c4762a1bSJed Brown 
402*c4762a1bSJed Brown     ierr = DMRestoreGlobalVector(base,&baseVec);CHKERRQ(ierr);
403*c4762a1bSJed Brown     ierr = DMRestoreGlobalVector(preForest,&baseVecMapped);CHKERRQ(ierr);
404*c4762a1bSJed Brown   }
405*c4762a1bSJed Brown 
406*c4762a1bSJed Brown   /* cleanup */
407*c4762a1bSJed Brown   ierr = VecDestroy(&preVec);CHKERRQ(ierr);
408*c4762a1bSJed Brown   ierr = DMDestroy(&preForest);CHKERRQ(ierr);
409*c4762a1bSJed Brown   ierr = DMDestroy(&base);CHKERRQ(ierr);
410*c4762a1bSJed Brown   ierr = PetscFinalize();
411*c4762a1bSJed Brown   return ierr;
412*c4762a1bSJed Brown }
413*c4762a1bSJed Brown 
414*c4762a1bSJed Brown /*TEST
415*c4762a1bSJed Brown 
416*c4762a1bSJed Brown      test:
417*c4762a1bSJed Brown        output_file: output/ex2_2d.out
418*c4762a1bSJed Brown        suffix: p4est_2d
419*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 2 -dim 2
420*c4762a1bSJed Brown        nsize: 3
421*c4762a1bSJed Brown        requires: p4est !single
422*c4762a1bSJed Brown 
423*c4762a1bSJed Brown      test:
424*c4762a1bSJed Brown        output_file: output/ex2_2d.out
425*c4762a1bSJed Brown        suffix: p4est_2d_deg4
426*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 4 -dim 2
427*c4762a1bSJed Brown        requires: p4est !single
428*c4762a1bSJed Brown 
429*c4762a1bSJed Brown      test:
430*c4762a1bSJed Brown        output_file: output/ex2_2d.out
431*c4762a1bSJed Brown        suffix: p4est_2d_deg8
432*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 8 -dim 2
433*c4762a1bSJed Brown        requires: p4est !single
434*c4762a1bSJed Brown 
435*c4762a1bSJed Brown      test:
436*c4762a1bSJed Brown        output_file: output/ex2_steps2.out
437*c4762a1bSJed Brown        suffix: p4est_2d_deg2_steps2
438*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 2 -dim 2 -coords -adapt_steps 2
439*c4762a1bSJed Brown        nsize: 3
440*c4762a1bSJed Brown        requires: p4est !single
441*c4762a1bSJed Brown 
442*c4762a1bSJed Brown      test:
443*c4762a1bSJed Brown        output_file: output/ex2_steps3.out
444*c4762a1bSJed Brown        suffix: p4est_2d_deg3_steps3
445*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 3 -dim 2 -coords -adapt_steps 3
446*c4762a1bSJed Brown        nsize: 3
447*c4762a1bSJed Brown        requires: p4est !single
448*c4762a1bSJed Brown 
449*c4762a1bSJed Brown      test:
450*c4762a1bSJed Brown        output_file: output/ex2_steps3.out
451*c4762a1bSJed Brown        suffix: p4est_2d_deg3_steps3_L2_periodic
452*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 3 -petscdualspace_lagrange_continuity 0 -dim 2 -coords -adapt_steps 3 -periodic -use_bcs 0
453*c4762a1bSJed Brown        nsize: 3
454*c4762a1bSJed Brown        requires: p4est !single
455*c4762a1bSJed Brown 
456*c4762a1bSJed Brown      test:
457*c4762a1bSJed Brown        output_file: output/ex2_steps3.out
458*c4762a1bSJed Brown        suffix: p4est_3d_deg2_steps3_L2_periodic
459*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 2 -petscdualspace_lagrange_continuity 0 -dim 3 -coords -adapt_steps 3 -periodic -use_bcs 0
460*c4762a1bSJed Brown        nsize: 3
461*c4762a1bSJed Brown        requires: p4est !single
462*c4762a1bSJed Brown 
463*c4762a1bSJed Brown      test:
464*c4762a1bSJed Brown        output_file: output/ex2_steps2.out
465*c4762a1bSJed Brown        suffix: p4est_3d_deg2_steps2
466*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 2 -dim 3 -coords -adapt_steps 2
467*c4762a1bSJed Brown        nsize: 3
468*c4762a1bSJed Brown        requires: p4est !single
469*c4762a1bSJed Brown 
470*c4762a1bSJed Brown      test:
471*c4762a1bSJed Brown        output_file: output/ex2_steps3.out
472*c4762a1bSJed Brown        suffix: p4est_3d_deg3_steps3
473*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 3 -dim 3 -coords -adapt_steps 3
474*c4762a1bSJed Brown        nsize: 3
475*c4762a1bSJed Brown        requires: p4est !single
476*c4762a1bSJed Brown 
477*c4762a1bSJed Brown      test:
478*c4762a1bSJed Brown        output_file: output/ex2_2d_fv.out
479*c4762a1bSJed Brown        suffix: p4est_2d_fv
480*c4762a1bSJed Brown        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1
481*c4762a1bSJed Brown        nsize: 3
482*c4762a1bSJed Brown        requires: p4est !single
483*c4762a1bSJed Brown 
484*c4762a1bSJed Brown      test:
485*c4762a1bSJed Brown        TODO: broken (codimension adjacency)
486*c4762a1bSJed Brown        output_file: output/ex2_2d_fv.out
487*c4762a1bSJed Brown        suffix: p4est_2d_fv_adjcodim
488*c4762a1bSJed Brown        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1 -dm_forest_adjacency_codimension 1
489*c4762a1bSJed Brown        nsize: 2
490*c4762a1bSJed Brown        requires: p4est !single
491*c4762a1bSJed Brown 
492*c4762a1bSJed Brown      test:
493*c4762a1bSJed Brown        TODO: broken (dimension adjacency)
494*c4762a1bSJed Brown        output_file: output/ex2_2d_fv.out
495*c4762a1bSJed Brown        suffix: p4est_2d_fv_adjdim
496*c4762a1bSJed Brown        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1 -dm_forest_adjacency_dimension 1
497*c4762a1bSJed Brown        nsize: 2
498*c4762a1bSJed Brown        requires: p4est !single
499*c4762a1bSJed Brown 
500*c4762a1bSJed Brown      test:
501*c4762a1bSJed Brown        output_file: output/ex2_2d_fv.out
502*c4762a1bSJed Brown        suffix: p4est_2d_fv_zerocells
503*c4762a1bSJed Brown        args: -transfer_from_base 0 -use_fv -linear -dim 2 -dm_forest_partition_overlap 1
504*c4762a1bSJed Brown        nsize: 10
505*c4762a1bSJed Brown        requires: p4est !single
506*c4762a1bSJed Brown 
507*c4762a1bSJed Brown      test:
508*c4762a1bSJed Brown        output_file: output/ex2_3d.out
509*c4762a1bSJed Brown        suffix: p4est_3d
510*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 1 -dim 3
511*c4762a1bSJed Brown        nsize: 3
512*c4762a1bSJed Brown        requires: p4est !single
513*c4762a1bSJed Brown 
514*c4762a1bSJed Brown      test:
515*c4762a1bSJed Brown        output_file: output/ex2_3d.out
516*c4762a1bSJed Brown        suffix: p4est_3d_deg3
517*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 3 -dim 3
518*c4762a1bSJed Brown        nsize: 3
519*c4762a1bSJed Brown        requires: p4est !single
520*c4762a1bSJed Brown 
521*c4762a1bSJed Brown      test:
522*c4762a1bSJed Brown        output_file: output/ex2_2d.out
523*c4762a1bSJed Brown        suffix: p4est_2d_deg2_coords
524*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 2 -dim 2 -coords
525*c4762a1bSJed Brown        nsize: 3
526*c4762a1bSJed Brown        requires: p4est !single
527*c4762a1bSJed Brown 
528*c4762a1bSJed Brown      test:
529*c4762a1bSJed Brown        output_file: output/ex2_3d.out
530*c4762a1bSJed Brown        suffix: p4est_3d_deg2_coords
531*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 2 -dim 3 -coords
532*c4762a1bSJed Brown        nsize: 3
533*c4762a1bSJed Brown        requires: p4est !single
534*c4762a1bSJed Brown 
535*c4762a1bSJed Brown      test:
536*c4762a1bSJed Brown        output_file: output/ex2_3d_fv.out
537*c4762a1bSJed Brown        suffix: p4est_3d_fv
538*c4762a1bSJed Brown        args: -transfer_from_base 0 -use_fv -linear -dim 3 -dm_forest_partition_overlap 1
539*c4762a1bSJed Brown        nsize: 3
540*c4762a1bSJed Brown        requires: p4est !single
541*c4762a1bSJed Brown 
542*c4762a1bSJed Brown      test:
543*c4762a1bSJed Brown        suffix: p4est_3d_nans
544*c4762a1bSJed Brown        args: -dim 3 -dm_forest_partition_overlap 1 -test_convert -petscspace_type tensor -petscspace_degree 1
545*c4762a1bSJed Brown        nsize: 2
546*c4762a1bSJed Brown        requires: p4est !single
547*c4762a1bSJed Brown 
548*c4762a1bSJed Brown      test:
549*c4762a1bSJed Brown        TODO: not broken, but the 3D case below is broken, so I do not trust this one
550*c4762a1bSJed Brown        output_file: output/ex2_steps2.out
551*c4762a1bSJed Brown        suffix: p4est_2d_tfb_distributed_nc
552*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 3 -dim 2 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -use_bcs 0 -coords -adapt_steps 2 -distribute_base -petscpartitioner_type shell -petscpartitioner_shell_random
553*c4762a1bSJed Brown        nsize: 3
554*c4762a1bSJed Brown        requires: p4est !single
555*c4762a1bSJed Brown 
556*c4762a1bSJed Brown      test:
557*c4762a1bSJed Brown        TODO: broken
558*c4762a1bSJed Brown        output_file: output/ex2_steps2.out
559*c4762a1bSJed Brown        suffix: p4est_3d_tfb_distributed_nc
560*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 2 -dim 3 -dm_forest_maximum_refinement 2 -dm_p4est_refine_pattern hash -use_bcs 0 -coords -adapt_steps 2 -distribute_base -petscpartitioner_type shell -petscpartitioner_shell_random
561*c4762a1bSJed Brown        nsize: 3
562*c4762a1bSJed Brown        requires: p4est !single
563*c4762a1bSJed Brown 
564*c4762a1bSJed Brown      test:
565*c4762a1bSJed Brown        TODO: broken
566*c4762a1bSJed Brown        output_file: output/ex2_3d.out
567*c4762a1bSJed Brown        suffix: p4est_3d_transfer_fails
568*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 1 -dm_forest_initial_refinement 1
569*c4762a1bSJed Brown        requires: p4est !single
570*c4762a1bSJed Brown 
571*c4762a1bSJed Brown      test:
572*c4762a1bSJed Brown        TODO: broken
573*c4762a1bSJed Brown        output_file: output/ex2_steps2_notfb.out
574*c4762a1bSJed Brown        suffix: p4est_3d_transfer_fails_2
575*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 1 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 2 -dm_forest_initial_refinement 0 -transfer_from_base 0
576*c4762a1bSJed Brown        requires: p4est !single
577*c4762a1bSJed Brown 
578*c4762a1bSJed Brown      test:
579*c4762a1bSJed Brown        output_file: output/ex2_steps2.out
580*c4762a1bSJed Brown        suffix: p4est_3d_multi_transfer_s2t
581*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 3 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 2 -dm_forest_initial_refinement 1 -petscdualspace_lagrange_continuity 0 -use_bcs 0
582*c4762a1bSJed Brown        requires: p4est !single
583*c4762a1bSJed Brown 
584*c4762a1bSJed Brown      test:
585*c4762a1bSJed Brown        output_file: output/ex2_steps2.out
586*c4762a1bSJed Brown        suffix: p4est_3d_coords_transfer_s2t
587*c4762a1bSJed Brown        args: -petscspace_type tensor -petscspace_degree 3 -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -adapt_steps 2 -dm_forest_initial_refinement 1 -petscdualspace_lagrange_continuity 0 -coords -use_bcs 0
588*c4762a1bSJed Brown        requires: p4est !single
589*c4762a1bSJed Brown 
590*c4762a1bSJed Brown TEST*/
591