xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 0d1cd5e0bf17d8d1a4f7f5396fcab04f005f5dcf)
10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
20bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
30bbe5d31Sbarral #include <pragmatic/cpragmatic.h>
40bbe5d31Sbarral #endif
50bbe5d31Sbarral 
6*0d1cd5e0SMatthew G. Knepley static PetscErrorCode DMPlexLabelToVolumeConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscReal refRatio, PetscReal maxVolumes[])
7*0d1cd5e0SMatthew G. Knepley {
8*0d1cd5e0SMatthew G. Knepley   PetscInt       dim, c;
9*0d1cd5e0SMatthew G. Knepley   PetscErrorCode ierr;
10*0d1cd5e0SMatthew G. Knepley 
11*0d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
12*0d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
13*0d1cd5e0SMatthew G. Knepley   refRatio = refRatio == PETSC_DEFAULT ? (PetscReal) ((PetscInt) 1 << dim) : refRatio;
14*0d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; c++) {
15*0d1cd5e0SMatthew G. Knepley     PetscReal vol;
16*0d1cd5e0SMatthew G. Knepley     PetscInt  closureSize = 0, cl;
17*0d1cd5e0SMatthew G. Knepley     PetscInt *closure     = NULL;
18*0d1cd5e0SMatthew G. Knepley     PetscBool anyRefine   = PETSC_FALSE;
19*0d1cd5e0SMatthew G. Knepley     PetscBool anyCoarsen  = PETSC_FALSE;
20*0d1cd5e0SMatthew G. Knepley     PetscBool anyKeep     = PETSC_FALSE;
21*0d1cd5e0SMatthew G. Knepley 
22*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
23*0d1cd5e0SMatthew G. Knepley     maxVolumes[c - cStart] = vol;
24*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
25*0d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
26*0d1cd5e0SMatthew G. Knepley       const PetscInt point = closure[cl];
27*0d1cd5e0SMatthew G. Knepley       PetscInt       refFlag;
28*0d1cd5e0SMatthew G. Knepley 
29*0d1cd5e0SMatthew G. Knepley       ierr = DMLabelGetValue(adaptLabel, point, &refFlag);CHKERRQ(ierr);
30*0d1cd5e0SMatthew G. Knepley       switch (refFlag) {
31*0d1cd5e0SMatthew G. Knepley       case DM_ADAPT_REFINE:
32*0d1cd5e0SMatthew G. Knepley         anyRefine  = PETSC_TRUE;break;
33*0d1cd5e0SMatthew G. Knepley       case DM_ADAPT_COARSEN:
34*0d1cd5e0SMatthew G. Knepley         anyCoarsen = PETSC_TRUE;break;
35*0d1cd5e0SMatthew G. Knepley       case DM_ADAPT_KEEP:
36*0d1cd5e0SMatthew G. Knepley         anyKeep    = PETSC_TRUE;break;
37*0d1cd5e0SMatthew G. Knepley       case DM_ADAPT_DETERMINE:
38*0d1cd5e0SMatthew G. Knepley         break;
39*0d1cd5e0SMatthew G. Knepley       default:
40*0d1cd5e0SMatthew G. Knepley         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP, "DMPlex does not support refinement flag %D\n", refFlag);break;
41*0d1cd5e0SMatthew G. Knepley       }
42*0d1cd5e0SMatthew G. Knepley       if (anyRefine) break;
43*0d1cd5e0SMatthew G. Knepley     }
44*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
45*0d1cd5e0SMatthew G. Knepley     if (anyRefine) {
46*0d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol / refRatio;
47*0d1cd5e0SMatthew G. Knepley     } else if (anyKeep) {
48*0d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol;
49*0d1cd5e0SMatthew G. Knepley     } else if (anyCoarsen) {
50*0d1cd5e0SMatthew G. Knepley       maxVolumes[c - cStart] = vol * refRatio;
51*0d1cd5e0SMatthew G. Knepley     }
52*0d1cd5e0SMatthew G. Knepley   }
53*0d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
54*0d1cd5e0SMatthew G. Knepley }
55*0d1cd5e0SMatthew G. Knepley 
56*0d1cd5e0SMatthew G. Knepley static PetscErrorCode DMPlexLabelToMetricConstraint(DM dm, DMLabel adaptLabel, PetscInt cStart, PetscInt cEnd, PetscInt vStart, PetscInt vEnd, PetscReal refRatio, Vec *metricVec)
57*0d1cd5e0SMatthew G. Knepley {
58*0d1cd5e0SMatthew G. Knepley   DM              udm, coordDM;
59*0d1cd5e0SMatthew G. Knepley   PetscSection    coordSection;
60*0d1cd5e0SMatthew G. Knepley   Vec             coordinates, mb, mx;
61*0d1cd5e0SMatthew G. Knepley   Mat             A;
62*0d1cd5e0SMatthew G. Knepley   PetscScalar    *metric, *eqns;
63*0d1cd5e0SMatthew G. Knepley   const PetscReal coarseRatio = refRatio == PETSC_DEFAULT ? PetscSqr(0.5) : 1/refRatio;
64*0d1cd5e0SMatthew G. Knepley   PetscInt        dim, Nv, Neq, c, v;
65*0d1cd5e0SMatthew G. Knepley   PetscErrorCode  ierr;
66*0d1cd5e0SMatthew G. Knepley 
67*0d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
68*0d1cd5e0SMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
69*0d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
70*0d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
71*0d1cd5e0SMatthew G. Knepley   ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr);
72*0d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
73*0d1cd5e0SMatthew G. Knepley   Nv   = vEnd - vStart;
74*0d1cd5e0SMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, Nv*PetscSqr(dim), metricVec);CHKERRQ(ierr);
75*0d1cd5e0SMatthew G. Knepley   ierr = VecGetArray(*metricVec, &metric);CHKERRQ(ierr);
76*0d1cd5e0SMatthew G. Knepley   Neq  = (dim*(dim+1))/2;
77*0d1cd5e0SMatthew G. Knepley   ierr = PetscMalloc1(PetscSqr(Neq), &eqns);CHKERRQ(ierr);
78*0d1cd5e0SMatthew G. Knepley   ierr = MatCreateSeqDense(PETSC_COMM_SELF, Neq, Neq, eqns, &A);CHKERRQ(ierr);
79*0d1cd5e0SMatthew G. Knepley   ierr = MatCreateVecs(A, &mx, &mb);CHKERRQ(ierr);
80*0d1cd5e0SMatthew G. Knepley   ierr = VecSet(mb, 1.0);CHKERRQ(ierr);
81*0d1cd5e0SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
82*0d1cd5e0SMatthew G. Knepley     const PetscScalar *sol;
83*0d1cd5e0SMatthew G. Knepley     PetscScalar       *cellCoords = NULL;
84*0d1cd5e0SMatthew G. Knepley     PetscReal          e[3], vol;
85*0d1cd5e0SMatthew G. Knepley     const PetscInt    *cone;
86*0d1cd5e0SMatthew G. Knepley     PetscInt           coneSize, cl, i, j, d, r;
87*0d1cd5e0SMatthew G. Knepley 
88*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
89*0d1cd5e0SMatthew G. Knepley     /* Only works for simplices */
90*0d1cd5e0SMatthew G. Knepley     for (i = 0, r = 0; i < dim+1; ++i) {
91*0d1cd5e0SMatthew G. Knepley       for (j = 0; j < i; ++j, ++r) {
92*0d1cd5e0SMatthew G. Knepley         for (d = 0; d < dim; ++d) e[d] = PetscRealPart(cellCoords[i*dim+d] - cellCoords[j*dim+d]);
93*0d1cd5e0SMatthew G. Knepley         /* FORTRAN ORDERING */
94*0d1cd5e0SMatthew G. Knepley         switch (dim) {
95*0d1cd5e0SMatthew G. Knepley         case 2:
96*0d1cd5e0SMatthew G. Knepley           eqns[0*Neq+r] = PetscSqr(e[0]);
97*0d1cd5e0SMatthew G. Knepley           eqns[1*Neq+r] = 2.0*e[0]*e[1];
98*0d1cd5e0SMatthew G. Knepley           eqns[2*Neq+r] = PetscSqr(e[1]);
99*0d1cd5e0SMatthew G. Knepley           break;
100*0d1cd5e0SMatthew G. Knepley         case 3:
101*0d1cd5e0SMatthew G. Knepley           eqns[0*Neq+r] = PetscSqr(e[0]);
102*0d1cd5e0SMatthew G. Knepley           eqns[1*Neq+r] = 2.0*e[0]*e[1];
103*0d1cd5e0SMatthew G. Knepley           eqns[2*Neq+r] = 2.0*e[0]*e[2];
104*0d1cd5e0SMatthew G. Knepley           eqns[3*Neq+r] = PetscSqr(e[1]);
105*0d1cd5e0SMatthew G. Knepley           eqns[4*Neq+r] = 2.0*e[1]*e[2];
106*0d1cd5e0SMatthew G. Knepley           eqns[5*Neq+r] = PetscSqr(e[2]);
107*0d1cd5e0SMatthew G. Knepley           break;
108*0d1cd5e0SMatthew G. Knepley         }
109*0d1cd5e0SMatthew G. Knepley       }
110*0d1cd5e0SMatthew G. Knepley     }
111*0d1cd5e0SMatthew G. Knepley     ierr = MatSetUnfactored(A);CHKERRQ(ierr);
112*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
113*0d1cd5e0SMatthew G. Knepley     ierr = MatLUFactor(A, NULL, NULL, NULL);CHKERRQ(ierr);
114*0d1cd5e0SMatthew G. Knepley     ierr = MatSolve(A, mb, mx);CHKERRQ(ierr);
115*0d1cd5e0SMatthew G. Knepley     ierr = VecGetArrayRead(mx, &sol);CHKERRQ(ierr);
116*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
117*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
118*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
119*0d1cd5e0SMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) {
120*0d1cd5e0SMatthew G. Knepley       const PetscInt v = cone[cl] - vStart;
121*0d1cd5e0SMatthew G. Knepley 
122*0d1cd5e0SMatthew G. Knepley       if (dim == 2) {
123*0d1cd5e0SMatthew G. Knepley         metric[v*4+0] += vol*coarseRatio*sol[0];
124*0d1cd5e0SMatthew G. Knepley         metric[v*4+1] += vol*coarseRatio*sol[1];
125*0d1cd5e0SMatthew G. Knepley         metric[v*4+2] += vol*coarseRatio*sol[1];
126*0d1cd5e0SMatthew G. Knepley         metric[v*4+3] += vol*coarseRatio*sol[2];
127*0d1cd5e0SMatthew G. Knepley       } else {
128*0d1cd5e0SMatthew G. Knepley         metric[v*9+0] += vol*coarseRatio*sol[0];
129*0d1cd5e0SMatthew G. Knepley         metric[v*9+1] += vol*coarseRatio*sol[1];
130*0d1cd5e0SMatthew G. Knepley         metric[v*9+3] += vol*coarseRatio*sol[1];
131*0d1cd5e0SMatthew G. Knepley         metric[v*9+2] += vol*coarseRatio*sol[2];
132*0d1cd5e0SMatthew G. Knepley         metric[v*9+6] += vol*coarseRatio*sol[2];
133*0d1cd5e0SMatthew G. Knepley         metric[v*9+4] += vol*coarseRatio*sol[3];
134*0d1cd5e0SMatthew G. Knepley         metric[v*9+5] += vol*coarseRatio*sol[4];
135*0d1cd5e0SMatthew G. Knepley         metric[v*9+7] += vol*coarseRatio*sol[4];
136*0d1cd5e0SMatthew G. Knepley         metric[v*9+8] += vol*coarseRatio*sol[5];
137*0d1cd5e0SMatthew G. Knepley       }
138*0d1cd5e0SMatthew G. Knepley     }
139*0d1cd5e0SMatthew G. Knepley     ierr = VecRestoreArrayRead(mx, &sol);CHKERRQ(ierr);
140*0d1cd5e0SMatthew G. Knepley   }
141*0d1cd5e0SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
142*0d1cd5e0SMatthew G. Knepley     const PetscInt *support;
143*0d1cd5e0SMatthew G. Knepley     PetscInt        supportSize, s;
144*0d1cd5e0SMatthew G. Knepley     PetscReal       vol, totVol = 0.0;
145*0d1cd5e0SMatthew G. Knepley 
146*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetSupport(udm, v+vStart, &support);CHKERRQ(ierr);
147*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexGetSupportSize(udm, v+vStart, &supportSize);CHKERRQ(ierr);
148*0d1cd5e0SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {ierr = DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL);CHKERRQ(ierr); totVol += vol;}
149*0d1cd5e0SMatthew G. Knepley     for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol;
150*0d1cd5e0SMatthew G. Knepley   }
151*0d1cd5e0SMatthew G. Knepley   ierr = PetscFree(eqns);CHKERRQ(ierr);
152*0d1cd5e0SMatthew G. Knepley   ierr = VecRestoreArray(*metricVec, &metric);CHKERRQ(ierr);
153*0d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&mx);CHKERRQ(ierr);
154*0d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&mb);CHKERRQ(ierr);
155*0d1cd5e0SMatthew G. Knepley   ierr = MatDestroy(&A);CHKERRQ(ierr);
156*0d1cd5e0SMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
157*0d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
158*0d1cd5e0SMatthew G. Knepley }
159*0d1cd5e0SMatthew G. Knepley 
160*0d1cd5e0SMatthew G. Knepley PetscErrorCode DMPlexRefine_Internal(DM dm, DMLabel adaptLabel, DM *dmRefined)
161*0d1cd5e0SMatthew G. Knepley {
162*0d1cd5e0SMatthew G. Knepley   PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *);
163*0d1cd5e0SMatthew G. Knepley   PetscReal        refinementLimit;
164*0d1cd5e0SMatthew G. Knepley   PetscInt         dim, cStart, cEnd;
165*0d1cd5e0SMatthew G. Knepley   char             genname[1024], *name = NULL;
166*0d1cd5e0SMatthew G. Knepley   PetscBool        isTriangle = PETSC_FALSE, isTetgen = PETSC_FALSE, isCTetgen = PETSC_FALSE, flg, localized;
167*0d1cd5e0SMatthew G. Knepley   PetscErrorCode   ierr;
168*0d1cd5e0SMatthew G. Knepley 
169*0d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
170*0d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
171*0d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetRefinementLimit(dm, &refinementLimit);CHKERRQ(ierr);
172*0d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetRefinementFunction(dm, &refinementFunc);CHKERRQ(ierr);
173*0d1cd5e0SMatthew G. Knepley   if (refinementLimit == 0.0 && !refinementFunc && !adaptLabel) PetscFunctionReturn(0);
174*0d1cd5e0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
175*0d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
176*0d1cd5e0SMatthew G. Knepley   ierr = PetscOptionsGetString(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-dm_plex_generator", genname, 1024, &flg);CHKERRQ(ierr);
177*0d1cd5e0SMatthew G. Knepley   if (flg) name = genname;
178*0d1cd5e0SMatthew G. Knepley   if (name) {
179*0d1cd5e0SMatthew G. Knepley     ierr = PetscStrcmp(name, "triangle", &isTriangle);CHKERRQ(ierr);
180*0d1cd5e0SMatthew G. Knepley     ierr = PetscStrcmp(name, "tetgen",   &isTetgen);CHKERRQ(ierr);
181*0d1cd5e0SMatthew G. Knepley     ierr = PetscStrcmp(name, "ctetgen",  &isCTetgen);CHKERRQ(ierr);
182*0d1cd5e0SMatthew G. Knepley   }
183*0d1cd5e0SMatthew G. Knepley   switch (dim) {
184*0d1cd5e0SMatthew G. Knepley   case 2:
185*0d1cd5e0SMatthew G. Knepley     if (!name || isTriangle) {
186*0d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_TRIANGLE)
187*0d1cd5e0SMatthew G. Knepley       PetscReal *maxVolumes;
188*0d1cd5e0SMatthew G. Knepley       PetscInt  c;
189*0d1cd5e0SMatthew G. Knepley 
190*0d1cd5e0SMatthew G. Knepley       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
191*0d1cd5e0SMatthew G. Knepley       if (adaptLabel) {
192*0d1cd5e0SMatthew G. Knepley         ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
193*0d1cd5e0SMatthew G. Knepley       } else if (refinementFunc) {
194*0d1cd5e0SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
195*0d1cd5e0SMatthew G. Knepley           PetscReal vol, centroid[3];
196*0d1cd5e0SMatthew G. Knepley           PetscReal maxVol;
197*0d1cd5e0SMatthew G. Knepley 
198*0d1cd5e0SMatthew G. Knepley           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
199*0d1cd5e0SMatthew G. Knepley           ierr = (*refinementFunc)(centroid, &maxVol);CHKERRQ(ierr);
200*0d1cd5e0SMatthew G. Knepley           maxVolumes[c - cStart] = (double) maxVol;
201*0d1cd5e0SMatthew G. Knepley         }
202*0d1cd5e0SMatthew G. Knepley       } else {
203*0d1cd5e0SMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
204*0d1cd5e0SMatthew G. Knepley       }
205*0d1cd5e0SMatthew G. Knepley #if !defined(PETSC_USE_REAL_DOUBLE)
206*0d1cd5e0SMatthew G. Knepley       {
207*0d1cd5e0SMatthew G. Knepley         double *mvols;
208*0d1cd5e0SMatthew G. Knepley         ierr = PetscMalloc1(cEnd - cStart,&mvols);CHKERRQ(ierr);
209*0d1cd5e0SMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) mvols[c] = (double)maxVolumes[c];
210*0d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_Triangle(dm, mvols, dmRefined);CHKERRQ(ierr);
211*0d1cd5e0SMatthew G. Knepley         ierr = PetscFree(mvols);CHKERRQ(ierr);
212*0d1cd5e0SMatthew G. Knepley       }
213*0d1cd5e0SMatthew G. Knepley #else
214*0d1cd5e0SMatthew G. Knepley       ierr = DMPlexRefine_Triangle(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
215*0d1cd5e0SMatthew G. Knepley #endif
216*0d1cd5e0SMatthew G. Knepley       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
217*0d1cd5e0SMatthew G. Knepley #else
218*0d1cd5e0SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement needs external package support.\nPlease reconfigure with --download-triangle.");
219*0d1cd5e0SMatthew G. Knepley #endif
220*0d1cd5e0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 2D mesh generation package %s", name);
221*0d1cd5e0SMatthew G. Knepley     break;
222*0d1cd5e0SMatthew G. Knepley   case 3:
223*0d1cd5e0SMatthew G. Knepley     if (!name || isCTetgen || isTetgen) {
224*0d1cd5e0SMatthew G. Knepley       PetscReal *maxVolumes;
225*0d1cd5e0SMatthew G. Knepley       PetscInt   c;
226*0d1cd5e0SMatthew G. Knepley 
227*0d1cd5e0SMatthew G. Knepley       ierr = PetscMalloc1(cEnd - cStart, &maxVolumes);CHKERRQ(ierr);
228*0d1cd5e0SMatthew G. Knepley       if (adaptLabel) {
229*0d1cd5e0SMatthew G. Knepley         ierr = DMPlexLabelToVolumeConstraint(dm, adaptLabel, cStart, cEnd, PETSC_DEFAULT, maxVolumes);CHKERRQ(ierr);
230*0d1cd5e0SMatthew G. Knepley       } else if (refinementFunc) {
231*0d1cd5e0SMatthew G. Knepley         for (c = cStart; c < cEnd; ++c) {
232*0d1cd5e0SMatthew G. Knepley           PetscReal vol, centroid[3];
233*0d1cd5e0SMatthew G. Knepley 
234*0d1cd5e0SMatthew G. Knepley           ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, centroid, NULL);CHKERRQ(ierr);
235*0d1cd5e0SMatthew G. Knepley           ierr = (*refinementFunc)(centroid, &maxVolumes[c-cStart]);CHKERRQ(ierr);
236*0d1cd5e0SMatthew G. Knepley         }
237*0d1cd5e0SMatthew G. Knepley       } else {
238*0d1cd5e0SMatthew G. Knepley         for (c = 0; c < cEnd-cStart; ++c) maxVolumes[c] = refinementLimit;
239*0d1cd5e0SMatthew G. Knepley       }
240*0d1cd5e0SMatthew G. Knepley       if (!name) {
241*0d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
242*0d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
243*0d1cd5e0SMatthew G. Knepley #elif defined(PETSC_HAVE_TETGEN)
244*0d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
245*0d1cd5e0SMatthew G. Knepley #else
246*0d1cd5e0SMatthew G. Knepley         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "External package CTetgen or Tetgen needed.\nPlease reconfigure with '--download-ctetgen' or '--with-clanguage=cxx --download-tetgen'.");
247*0d1cd5e0SMatthew G. Knepley #endif
248*0d1cd5e0SMatthew G. Knepley       } else if (isCTetgen) {
249*0d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_CTETGEN)
250*0d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_CTetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
251*0d1cd5e0SMatthew G. Knepley #else
252*0d1cd5e0SMatthew G. Knepley         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "CTetgen needs external package support.\nPlease reconfigure with --download-ctetgen.");
253*0d1cd5e0SMatthew G. Knepley #endif
254*0d1cd5e0SMatthew G. Knepley       } else {
255*0d1cd5e0SMatthew G. Knepley #if defined(PETSC_HAVE_TETGEN)
256*0d1cd5e0SMatthew G. Knepley         ierr = DMPlexRefine_Tetgen(dm, maxVolumes, dmRefined);CHKERRQ(ierr);
257*0d1cd5e0SMatthew G. Knepley #else
258*0d1cd5e0SMatthew G. Knepley         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Tetgen needs external package support.\nPlease reconfigure with --download-tetgen.");
259*0d1cd5e0SMatthew G. Knepley #endif
260*0d1cd5e0SMatthew G. Knepley       }
261*0d1cd5e0SMatthew G. Knepley       ierr = PetscFree(maxVolumes);CHKERRQ(ierr);
262*0d1cd5e0SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown 3D mesh generation package %s", name);
263*0d1cd5e0SMatthew G. Knepley     break;
264*0d1cd5e0SMatthew G. Knepley   default:
265*0d1cd5e0SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Mesh refinement in dimension %d is not supported.", dim);
266*0d1cd5e0SMatthew G. Knepley   }
267*0d1cd5e0SMatthew G. Knepley   ierr = DMCopyBoundary(dm, *dmRefined);CHKERRQ(ierr);
268*0d1cd5e0SMatthew G. Knepley   if (localized) {ierr = DMLocalizeCoordinates(*dmRefined);CHKERRQ(ierr);}
269*0d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
270*0d1cd5e0SMatthew G. Knepley }
271*0d1cd5e0SMatthew G. Knepley 
272*0d1cd5e0SMatthew G. Knepley PetscErrorCode DMPlexCoarsen_Internal(DM dm, DMLabel adaptLabel, DM *dmCoarsened)
273*0d1cd5e0SMatthew G. Knepley {
274*0d1cd5e0SMatthew G. Knepley   Vec            metricVec;
275*0d1cd5e0SMatthew G. Knepley   PetscInt       cStart, cEnd, vStart, vEnd;
276*0d1cd5e0SMatthew G. Knepley   DMLabel        bdLabel = NULL;
277*0d1cd5e0SMatthew G. Knepley   char           bdLabelName[PETSC_MAX_PATH_LEN];
278*0d1cd5e0SMatthew G. Knepley   PetscBool      localized, flg;
279*0d1cd5e0SMatthew G. Knepley   PetscErrorCode ierr;
280*0d1cd5e0SMatthew G. Knepley 
281*0d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
282*0d1cd5e0SMatthew G. Knepley   ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr);
283*0d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
284*0d1cd5e0SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
285*0d1cd5e0SMatthew G. Knepley   ierr = DMPlexLabelToMetricConstraint(dm, adaptLabel, cStart, cEnd, vStart, vEnd, PETSC_DEFAULT, &metricVec);CHKERRQ(ierr);
286*0d1cd5e0SMatthew G. Knepley   ierr = PetscOptionsGetString(NULL, dm->hdr.prefix, "-dm_plex_coarsen_bd_label", bdLabelName, PETSC_MAX_PATH_LEN-1, &flg);CHKERRQ(ierr);
287*0d1cd5e0SMatthew G. Knepley   if (flg) {ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);}
288*0d1cd5e0SMatthew G. Knepley   ierr = DMAdaptMetric_Plex(dm, metricVec, bdLabel, dmCoarsened);CHKERRQ(ierr);
289*0d1cd5e0SMatthew G. Knepley   ierr = VecDestroy(&metricVec);CHKERRQ(ierr);
290*0d1cd5e0SMatthew G. Knepley   if (localized) {ierr = DMLocalizeCoordinates(*dmCoarsened);CHKERRQ(ierr);}
291*0d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
292*0d1cd5e0SMatthew G. Knepley }
293*0d1cd5e0SMatthew G. Knepley 
294*0d1cd5e0SMatthew G. Knepley PetscErrorCode DMAdaptLabel_Plex(DM dm, DMLabel adaptLabel, DM *dmAdapted)
295*0d1cd5e0SMatthew G. Knepley {
296*0d1cd5e0SMatthew G. Knepley   IS              flagIS;
297*0d1cd5e0SMatthew G. Knepley   const PetscInt *flags;
298*0d1cd5e0SMatthew G. Knepley   PetscInt        defFlag, minFlag, maxFlag, numFlags, f;
299*0d1cd5e0SMatthew G. Knepley   PetscErrorCode  ierr;
300*0d1cd5e0SMatthew G. Knepley 
301*0d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
302*0d1cd5e0SMatthew G. Knepley   ierr = DMLabelGetDefaultValue(adaptLabel, &defFlag);CHKERRQ(ierr);
303*0d1cd5e0SMatthew G. Knepley   minFlag = defFlag;
304*0d1cd5e0SMatthew G. Knepley   maxFlag = defFlag;
305*0d1cd5e0SMatthew G. Knepley   ierr = DMLabelGetValueIS(adaptLabel, &flagIS);CHKERRQ(ierr);
306*0d1cd5e0SMatthew G. Knepley   ierr = ISGetLocalSize(flagIS, &numFlags);CHKERRQ(ierr);
307*0d1cd5e0SMatthew G. Knepley   ierr = ISGetIndices(flagIS, &flags);CHKERRQ(ierr);
308*0d1cd5e0SMatthew G. Knepley   for (f = 0; f < numFlags; ++f) {
309*0d1cd5e0SMatthew G. Knepley     const PetscInt flag = flags[f];
310*0d1cd5e0SMatthew G. Knepley 
311*0d1cd5e0SMatthew G. Knepley     minFlag = PetscMin(minFlag, flag);
312*0d1cd5e0SMatthew G. Knepley     maxFlag = PetscMax(maxFlag, flag);
313*0d1cd5e0SMatthew G. Knepley   }
314*0d1cd5e0SMatthew G. Knepley   ierr = ISRestoreIndices(flagIS, &flags);CHKERRQ(ierr);
315*0d1cd5e0SMatthew G. Knepley   ierr = ISDestroy(&flagIS);CHKERRQ(ierr);
316*0d1cd5e0SMatthew G. Knepley   {
317*0d1cd5e0SMatthew G. Knepley     PetscInt minMaxFlag[2], minMaxFlagGlobal[2];
318*0d1cd5e0SMatthew G. Knepley 
319*0d1cd5e0SMatthew G. Knepley     minMaxFlag[0] =  minFlag;
320*0d1cd5e0SMatthew G. Knepley     minMaxFlag[1] = -maxFlag;
321*0d1cd5e0SMatthew G. Knepley     ierr = MPI_Allreduce(minMaxFlag, minMaxFlagGlobal, 2, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
322*0d1cd5e0SMatthew G. Knepley     minFlag =  minMaxFlagGlobal[0];
323*0d1cd5e0SMatthew G. Knepley     maxFlag = -minMaxFlagGlobal[1];
324*0d1cd5e0SMatthew G. Knepley   }
325*0d1cd5e0SMatthew G. Knepley   if (minFlag == maxFlag) {
326*0d1cd5e0SMatthew G. Knepley     switch (minFlag) {
327*0d1cd5e0SMatthew G. Knepley     case DM_ADAPT_DETERMINE:
328*0d1cd5e0SMatthew G. Knepley       *dmAdapted = NULL;break;
329*0d1cd5e0SMatthew G. Knepley     case DM_ADAPT_REFINE:
330*0d1cd5e0SMatthew G. Knepley       ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE);CHKERRQ(ierr);
331*0d1cd5e0SMatthew G. Knepley       ierr = DMRefine(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
332*0d1cd5e0SMatthew G. Knepley     case DM_ADAPT_COARSEN:
333*0d1cd5e0SMatthew G. Knepley       ierr = DMCoarsen(dm, MPI_COMM_NULL, dmAdapted);CHKERRQ(ierr);break;
334*0d1cd5e0SMatthew G. Knepley     default:
335*0d1cd5e0SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP,"DMPlex does not support refinement flag %D\n", minFlag);break;
336*0d1cd5e0SMatthew G. Knepley     }
337*0d1cd5e0SMatthew G. Knepley   } else {
338*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr);
339*0d1cd5e0SMatthew G. Knepley     ierr = DMPlexRefine_Internal(dm, adaptLabel, dmAdapted);CHKERRQ(ierr);
340*0d1cd5e0SMatthew G. Knepley   }
341*0d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
342*0d1cd5e0SMatthew G. Knepley }
343*0d1cd5e0SMatthew G. Knepley 
344f7bcbcbbSMatthew G. Knepley /*
345*0d1cd5e0SMatthew G. Knepley   DMAdaptMetric_Plex - Generates a new mesh conforming to a metric field.
3460bbe5d31Sbarral 
347f7bcbcbbSMatthew G. Knepley   Input Parameters:
348f7bcbcbbSMatthew G. Knepley + dm - The DM object
34951ff0a1cSMatthew G. Knepley . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise in a LOCAL vector
350*0d1cd5e0SMatthew G. Knepley - bdLabel - Label for boundary tags which are preserved in dmNew, or NULL. Should not be named "_boundary_".
351f7bcbcbbSMatthew G. Knepley 
352f7bcbcbbSMatthew G. Knepley   Output Parameter:
353f7bcbcbbSMatthew G. Knepley . dmNew  - the new DM
354f7bcbcbbSMatthew G. Knepley 
355f7bcbcbbSMatthew G. Knepley   Level: advanced
356f7bcbcbbSMatthew G. Knepley 
357f7bcbcbbSMatthew G. Knepley .seealso: DMCoarsen(), DMRefine()
358f7bcbcbbSMatthew G. Knepley */
359*0d1cd5e0SMatthew G. Knepley PetscErrorCode DMAdaptMetric_Plex(DM dm, Vec vertexMetric, DMLabel bdLabel, DM *dmNew)
360f7bcbcbbSMatthew G. Knepley {
361bc1d7e78SMatthew G. Knepley   MPI_Comm           comm;
362f7bcbcbbSMatthew G. Knepley   const char        *bdName = "_boundary_";
3632b5f8d01SMatthew G. Knepley #if 0
3642b5f8d01SMatthew G. Knepley   DM                 odm = dm;
3652b5f8d01SMatthew G. Knepley #endif
366f7bcbcbbSMatthew G. Knepley   DM                 udm, cdm;
367*0d1cd5e0SMatthew G. Knepley   DMLabel            bdLabelFull;
368*0d1cd5e0SMatthew G. Knepley   const char        *bdLabelName;
3691838b4a6SMatthew G. Knepley   IS                 bdIS, globalVertexNum;
370f7bcbcbbSMatthew G. Knepley   PetscSection       coordSection;
371f7bcbcbbSMatthew G. Knepley   Vec                coordinates;
372f7bcbcbbSMatthew G. Knepley   const PetscScalar *coords, *met;
3731838b4a6SMatthew G. Knepley   const PetscInt    *bdFacesFull, *gV;
3741838b4a6SMatthew G. Knepley   PetscInt          *bdFaces, *bdFaceIds, *l2gv;
375f7bcbcbbSMatthew G. Knepley   PetscReal         *x, *y, *z, *metric;
3762ee120b0SMatthew G. Knepley   PetscInt          *cells;
3771838b4a6SMatthew G. Knepley   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, numLocVertices, v;
3782ee120b0SMatthew G. Knepley   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
379f7bcbcbbSMatthew G. Knepley   PetscBool          flg;
3802ee120b0SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
3812ee120b0SMatthew G. Knepley   DMLabel            bdLabelNew;
382acfe6eadSToby Isaac   double            *coordsNew;
3832ee120b0SMatthew G. Knepley   PetscInt          *bdTags;
3842ee120b0SMatthew G. Knepley   PetscReal         *xNew[3] = {NULL, NULL, NULL};
3852ee120b0SMatthew G. Knepley   PetscInt          *cellsNew;
3862ee120b0SMatthew G. Knepley   PetscInt           d, numCellsNew, numVerticesNew;
3872ee120b0SMatthew G. Knepley   PetscInt           numCornersNew, fStart, fEnd;
3882ee120b0SMatthew G. Knepley #endif
389bc1d7e78SMatthew G. Knepley   PetscMPIInt        numProcs;
390f7bcbcbbSMatthew G. Knepley   PetscErrorCode     ierr;
391f7bcbcbbSMatthew G. Knepley 
392f7bcbcbbSMatthew G. Knepley   PetscFunctionBegin;
393bbc23918SMatthew G. Knepley   /* Check for FEM adjacency flags */
394bc1d7e78SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
395bc1d7e78SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
396*0d1cd5e0SMatthew G. Knepley   if (bdLabel) {
397*0d1cd5e0SMatthew G. Knepley     ierr = DMLabelGetName(bdLabel, &bdLabelName);CHKERRQ(ierr);
398f7bcbcbbSMatthew G. Knepley     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
399bc1d7e78SMatthew G. Knepley     if (flg) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
400f7bcbcbbSMatthew G. Knepley   }
401a2e1bd45SMatthew G. Knepley   /* Add overlap for Pragmatic */
402bbc23918SMatthew G. Knepley #if 0
403bbc23918SMatthew G. Knepley   /* Check for overlap by looking for cell in the SF */
404bbc23918SMatthew G. Knepley   if (!overlapped) {
405a2e1bd45SMatthew G. Knepley     ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr);
406a2e1bd45SMatthew G. Knepley     if (!dm) {dm = odm; ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);}
407bbc23918SMatthew G. Knepley   }
408bbc23918SMatthew G. Knepley #endif
409f7bcbcbbSMatthew G. Knepley   /* Get mesh information */
410f7bcbcbbSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
411f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
412f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
413f7bcbcbbSMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
414f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
415f7bcbcbbSMatthew G. Knepley   numCells    = cEnd - cStart;
416f7bcbcbbSMatthew G. Knepley   numVertices = vEnd - vStart;
417f7bcbcbbSMatthew G. Knepley   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
418f7bcbcbbSMatthew G. Knepley   for (c = 0, coff = 0; c < numCells; ++c) {
419f7bcbcbbSMatthew G. Knepley     const PetscInt *cone;
420f7bcbcbbSMatthew G. Knepley     PetscInt        coneSize, cl;
421f7bcbcbbSMatthew G. Knepley 
422f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
423f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
424f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
425f7bcbcbbSMatthew G. Knepley   }
4261838b4a6SMatthew G. Knepley   ierr = PetscCalloc1(numVertices, &l2gv);CHKERRQ(ierr);
4271838b4a6SMatthew G. Knepley   ierr = DMPlexGetVertexNumbering(udm, &globalVertexNum);CHKERRQ(ierr);
4281838b4a6SMatthew G. Knepley   ierr = ISGetIndices(globalVertexNum, &gV);CHKERRQ(ierr);
4291838b4a6SMatthew G. Knepley   for (v = 0, numLocVertices = 0; v < numVertices; ++v) {
4301838b4a6SMatthew G. Knepley     if (gV[v] >= 0) ++numLocVertices;
4311838b4a6SMatthew G. Knepley     l2gv[v] = gV[v] < 0 ? -(gV[v]+1) : gV[v];
4321838b4a6SMatthew G. Knepley   }
4331838b4a6SMatthew G. Knepley   ierr = ISRestoreIndices(globalVertexNum, &gV);CHKERRQ(ierr);
434f7bcbcbbSMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
435f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
436f7bcbcbbSMatthew G. Knepley   ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr);
437f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
438f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
439f7bcbcbbSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
440f7bcbcbbSMatthew G. Knepley     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
4415f2e139eSMatthew G. Knepley     x[v-vStart] = PetscRealPart(coords[off+0]);
4425f2e139eSMatthew G. Knepley     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
4435f2e139eSMatthew G. Knepley     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
444f7bcbcbbSMatthew G. Knepley   }
445f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
446f7bcbcbbSMatthew G. Knepley   /* Get boundary mesh */
447f7bcbcbbSMatthew G. Knepley   ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr);
448f7bcbcbbSMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr);
449f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
450f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
451f7bcbcbbSMatthew G. Knepley   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
452f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
453f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
454f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
455f7bcbcbbSMatthew G. Knepley 
456f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
457f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
458f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
459f7bcbcbbSMatthew G. Knepley     }
460f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
461f7bcbcbbSMatthew G. Knepley   }
462f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
463f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
464f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
465f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
466f7bcbcbbSMatthew G. Knepley 
467f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
468f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
469f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
470f7bcbcbbSMatthew G. Knepley     }
471f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
472f7bcbcbbSMatthew G. Knepley     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
473f7bcbcbbSMatthew G. Knepley     else         {bdFaceIds[f] = 1;}
474f7bcbcbbSMatthew G. Knepley   }
475f7bcbcbbSMatthew G. Knepley   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
476f7bcbcbbSMatthew G. Knepley   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
477f7bcbcbbSMatthew G. Knepley   /* Get metric */
478f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
4795f2e139eSMatthew G. Knepley   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
480f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
481bbc23918SMatthew G. Knepley #if 0
482a2e1bd45SMatthew G. Knepley   /* Destroy overlap mesh */
483a2e1bd45SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
484bbc23918SMatthew G. Knepley #endif
485f7bcbcbbSMatthew G. Knepley   /* Create new mesh */
486f7bcbcbbSMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
487f7bcbcbbSMatthew G. Knepley   switch (dim) {
488f7bcbcbbSMatthew G. Knepley   case 2:
489cc29c497SMatthew G. Knepley     pragmatic_2d_mpi_init(&numVertices, &numCells, cells, x, y, l2gv, numLocVertices, comm);break;
490f7bcbcbbSMatthew G. Knepley   case 3:
491cc29c497SMatthew G. Knepley     pragmatic_3d_mpi_init(&numVertices, &numCells, cells, x, y, z, l2gv, numLocVertices, comm);break;
4921838b4a6SMatthew G. Knepley   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
493f7bcbcbbSMatthew G. Knepley   }
494f7bcbcbbSMatthew G. Knepley   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
495f7bcbcbbSMatthew G. Knepley   pragmatic_set_metric(metric);
496*0d1cd5e0SMatthew G. Knepley   pragmatic_adapt(((DM_Plex *) dm->data)->remeshBd ? 1 : 0);
4971838b4a6SMatthew G. Knepley   ierr = PetscFree(l2gv);CHKERRQ(ierr);
498f7bcbcbbSMatthew G. Knepley   /* Read out mesh */
4998713909fSMatthew G. Knepley   pragmatic_get_info_mpi(&numVerticesNew, &numCellsNew);
500f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
501f7bcbcbbSMatthew G. Knepley   switch (dim) {
502f7bcbcbbSMatthew G. Knepley   case 2:
503f7bcbcbbSMatthew G. Knepley     numCornersNew = 3;
504f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
5058713909fSMatthew G. Knepley     pragmatic_get_coords_2d_mpi(xNew[0], xNew[1]);
506f7bcbcbbSMatthew G. Knepley     break;
507f7bcbcbbSMatthew G. Knepley   case 3:
508f7bcbcbbSMatthew G. Knepley     numCornersNew = 4;
509f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
5108713909fSMatthew G. Knepley     pragmatic_get_coords_3d_mpi(xNew[0], xNew[1], xNew[2]);
511f7bcbcbbSMatthew G. Knepley     break;
512f7bcbcbbSMatthew G. Knepley   default:
513bc1d7e78SMatthew G. Knepley     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
514f7bcbcbbSMatthew G. Knepley   }
515acfe6eadSToby Isaac   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = (double) xNew[d][v];}
516f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
517f7bcbcbbSMatthew G. Knepley   pragmatic_get_elements(cellsNew);
51851ff0a1cSMatthew G. Knepley   ierr = DMPlexCreateFromCellListParallel(comm, dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, NULL, dmNew);CHKERRQ(ierr);
519f7bcbcbbSMatthew G. Knepley   /* Read out boundary label */
520f7bcbcbbSMatthew G. Knepley   pragmatic_get_boundaryTags(&bdTags);
521f7bcbcbbSMatthew G. Knepley   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
522f7bcbcbbSMatthew G. Knepley   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
523f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
524f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
525f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
526f7bcbcbbSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
527f7bcbcbbSMatthew G. Knepley     /* Only for simplicial meshes */
528f7bcbcbbSMatthew G. Knepley     coff = (c-cStart)*(dim+1);
5298713909fSMatthew G. Knepley     /* d is the local cell number of the vertex opposite to the face we are marking */
530f7bcbcbbSMatthew G. Knepley     for (d = 0; d < dim+1; ++d) {
531f7bcbcbbSMatthew G. Knepley       if (bdTags[coff+d]) {
5328713909fSMatthew G. Knepley         const PetscInt  perm[4][4] = {{-1, -1, -1, -1}, {-1, -1, -1, -1}, {1, 2, 0, -1}, {3, 2, 1, 0}}; /* perm[d] = face opposite */
5338713909fSMatthew G. Knepley         const PetscInt *cone;
534f7bcbcbbSMatthew G. Knepley 
5358713909fSMatthew G. Knepley         /* Mark face opposite to this vertex: This pattern is specified in DMPlexGetRawFaces_Internal() */
5368713909fSMatthew G. Knepley         ierr = DMPlexGetCone(*dmNew, c, &cone);CHKERRQ(ierr);
5378713909fSMatthew G. Knepley         ierr = DMLabelSetValue(bdLabelNew, cone[perm[dim][d]], bdTags[coff+d]);CHKERRQ(ierr);
538f7bcbcbbSMatthew G. Knepley       }
539f7bcbcbbSMatthew G. Knepley     }
540f7bcbcbbSMatthew G. Knepley   }
541f7bcbcbbSMatthew G. Knepley   /* Cleanup */
542f7bcbcbbSMatthew G. Knepley   switch (dim) {
543f7bcbcbbSMatthew G. Knepley   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break;
544f7bcbcbbSMatthew G. Knepley   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break;
545f7bcbcbbSMatthew G. Knepley   }
546f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
547f7bcbcbbSMatthew G. Knepley   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
548f7bcbcbbSMatthew G. Knepley   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
549f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
550f7bcbcbbSMatthew G. Knepley   pragmatic_finalize();
551f7bcbcbbSMatthew G. Knepley #else
552bc1d7e78SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
553f7bcbcbbSMatthew G. Knepley #endif
554f7bcbcbbSMatthew G. Knepley   PetscFunctionReturn(0);
555f7bcbcbbSMatthew G. Knepley }
556