xref: /petsc/src/dm/impls/plex/tests/ex19.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests mesh adaptation with DMPlex and pragmatic.\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown #include <petscksp.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown typedef struct {
8*c4762a1bSJed Brown   DM        dm;
9*c4762a1bSJed Brown   /* Definition of the test case (mesh and metric field) */
10*c4762a1bSJed Brown   PetscInt  dim;                         /* The topological mesh dimension */
11*c4762a1bSJed Brown   char      mshNam[PETSC_MAX_PATH_LEN];  /* Name of the mesh filename if any */
12*c4762a1bSJed Brown   PetscInt  nbrVerEdge;                  /* Number of vertices per edge if unit square/cube generated */
13*c4762a1bSJed Brown   char      bdLabel[PETSC_MAX_PATH_LEN]; /* Name of the label marking boundary facets */
14*c4762a1bSJed Brown   PetscInt  metOpt;                      /* Different choices of metric */
15*c4762a1bSJed Brown   PetscReal hmax, hmin;                  /* Max and min sizes prescribed by the metric */
16*c4762a1bSJed Brown   PetscBool doL2;                        /* Test L2 projection */
17*c4762a1bSJed Brown } AppCtx;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
20*c4762a1bSJed Brown {
21*c4762a1bSJed Brown   PetscErrorCode ierr;
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown   PetscFunctionBegin;
24*c4762a1bSJed Brown   options->dim        = 2;
25*c4762a1bSJed Brown   ierr = PetscStrcpy(options->mshNam, "");CHKERRQ(ierr);
26*c4762a1bSJed Brown   options->nbrVerEdge = 5;
27*c4762a1bSJed Brown   ierr = PetscStrcpy(options->bdLabel, "");CHKERRQ(ierr);
28*c4762a1bSJed Brown   options->metOpt     = 1;
29*c4762a1bSJed Brown   options->hmin       = 0.05;
30*c4762a1bSJed Brown   options->hmax       = 0.5;
31*c4762a1bSJed Brown   options->doL2       = PETSC_FALSE;
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Adaptation Options", "DMPLEX");CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex19.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = PetscOptionsString("-msh", "Name of the mesh filename if any", "ex19.c", options->mshNam, options->mshNam, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-nbrVerEdge", "Number of vertices per edge if unit square/cube generated", "ex19.c", options->nbrVerEdge, &options->nbrVerEdge, NULL,0);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = PetscOptionsString("-bdLabel", "Name of the label marking boundary facets", "ex19.c", options->bdLabel, options->bdLabel, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-met", "Different choices of metric", "ex19.c", options->metOpt, &options->metOpt, NULL,0);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscOptionsReal("-hmax", "Max size prescribed by the metric", "ex19.c", options->hmax, &options->hmax, NULL);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = PetscOptionsReal("-hmin", "Min size prescribed by the metric", "ex19.c", options->hmin, &options->hmin, NULL);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = PetscOptionsBool("-do_L2", "Test L2 projection", "ex19.c", options->doL2, &options->doL2, NULL);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   PetscFunctionReturn(0);
45*c4762a1bSJed Brown }
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user)
48*c4762a1bSJed Brown {
49*c4762a1bSJed Brown   PetscBool      flag;
50*c4762a1bSJed Brown   PetscErrorCode ierr;
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   PetscFunctionBegin;
53*c4762a1bSJed Brown   ierr = PetscStrcmp(user->mshNam, "", &flag);CHKERRQ(ierr);
54*c4762a1bSJed Brown   if (flag) {
55*c4762a1bSJed Brown     PetscInt faces[3];
56*c4762a1bSJed Brown     faces[0] = user->nbrVerEdge-1;
57*c4762a1bSJed Brown     faces[1] = user->nbrVerEdge-1;
58*c4762a1bSJed Brown     faces[2] = user->nbrVerEdge-1;
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(comm, user->dim, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, &user->dm);CHKERRQ(ierr);
61*c4762a1bSJed Brown   } else {
62*c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm, user->mshNam, PETSC_TRUE, &user->dm);CHKERRQ(ierr);
63*c4762a1bSJed Brown     ierr = DMGetDimension(user->dm, &user->dim);CHKERRQ(ierr);
64*c4762a1bSJed Brown   }
65*c4762a1bSJed Brown   {
66*c4762a1bSJed Brown     DM distributedMesh = NULL;
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown     /* Distribute mesh over processes */
69*c4762a1bSJed Brown     ierr = DMPlexDistribute(user->dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
70*c4762a1bSJed Brown     if (distributedMesh) {
71*c4762a1bSJed Brown       ierr = DMDestroy(&user->dm);CHKERRQ(ierr);
72*c4762a1bSJed Brown       user->dm  = distributedMesh;
73*c4762a1bSJed Brown     }
74*c4762a1bSJed Brown   }
75*c4762a1bSJed Brown   ierr = DMSetFromOptions(user->dm);CHKERRQ(ierr);
76*c4762a1bSJed Brown   PetscFunctionReturn(0);
77*c4762a1bSJed Brown }
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown static PetscErrorCode ComputeMetric(DM dm, AppCtx *user, Vec *metric)
80*c4762a1bSJed Brown {
81*c4762a1bSJed Brown   DM                 cdm, mdm;
82*c4762a1bSJed Brown   PetscSection       csec, msec;
83*c4762a1bSJed Brown   Vec                coordinates;
84*c4762a1bSJed Brown   const PetscScalar *coords;
85*c4762a1bSJed Brown   PetscScalar       *met;
86*c4762a1bSJed Brown   PetscReal          h, *lambda, lbd, lmax;
87*c4762a1bSJed Brown   PetscInt           pStart, pEnd, p, d;
88*c4762a1bSJed Brown   const PetscInt     dim = user->dim, Nd = dim*dim;
89*c4762a1bSJed Brown   PetscErrorCode     ierr;
90*c4762a1bSJed Brown 
91*c4762a1bSJed Brown   PetscFunctionBeginUser;
92*c4762a1bSJed Brown   ierr = PetscCalloc1(PetscMax(3, dim),&lambda);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = DMClone(cdm, &mdm);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = DMGetLocalSection(cdm, &csec);CHKERRQ(ierr);
96*c4762a1bSJed Brown 
97*c4762a1bSJed Brown   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &msec);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = PetscSectionSetNumFields(msec, 1);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = PetscSectionSetFieldComponents(msec, 0, Nd);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = PetscSectionGetChart(csec, &pStart, &pEnd);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = PetscSectionSetChart(msec, pStart, pEnd);CHKERRQ(ierr);
102*c4762a1bSJed Brown   for (p = pStart; p < pEnd; ++p) {
103*c4762a1bSJed Brown     ierr = PetscSectionSetDof(msec, p, Nd);CHKERRQ(ierr);
104*c4762a1bSJed Brown     ierr = PetscSectionSetFieldDof(msec, p, 0, Nd);CHKERRQ(ierr);
105*c4762a1bSJed Brown   }
106*c4762a1bSJed Brown   ierr = PetscSectionSetUp(msec);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = DMSetLocalSection(mdm, msec);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = PetscSectionDestroy(&msec);CHKERRQ(ierr);
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = DMCreateLocalVector(mdm, metric);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = VecGetArray(*metric, &met);CHKERRQ(ierr);
114*c4762a1bSJed Brown   for (p = pStart; p < pEnd; ++p) {
115*c4762a1bSJed Brown     PetscScalar       *pcoords;
116*c4762a1bSJed Brown     PetscScalar       *pmet;
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown     ierr = DMPlexPointLocalRead(cdm, p, coords, &pcoords);CHKERRQ(ierr);
119*c4762a1bSJed Brown     switch (user->metOpt) {
120*c4762a1bSJed Brown     case 0:
121*c4762a1bSJed Brown       lbd = 1/(user->hmax*user->hmax);
122*c4762a1bSJed Brown       lambda[0] = lambda[1] = lambda[2] = lbd;
123*c4762a1bSJed Brown       break;
124*c4762a1bSJed Brown     case 1:
125*c4762a1bSJed Brown       h = user->hmax - (user->hmax-user->hmin)*PetscRealPart(pcoords[0]);
126*c4762a1bSJed Brown       h = h*h;
127*c4762a1bSJed Brown       lmax = 1/(user->hmax*user->hmax);
128*c4762a1bSJed Brown       lambda[0] = 1/h;
129*c4762a1bSJed Brown       lambda[1] = lmax;
130*c4762a1bSJed Brown       lambda[2] = lmax;
131*c4762a1bSJed Brown       break;
132*c4762a1bSJed Brown     case 2:
133*c4762a1bSJed Brown       h = user->hmax*PetscAbsReal(((PetscReal) 1.0)-PetscExpReal(-PetscAbsScalar(pcoords[0]-(PetscReal)0.5))) + user->hmin;
134*c4762a1bSJed Brown       lbd = 1/(h*h);
135*c4762a1bSJed Brown       lmax = 1/(user->hmax*user->hmax);
136*c4762a1bSJed Brown       lambda[0] = lbd;
137*c4762a1bSJed Brown       lambda[1] = lmax;
138*c4762a1bSJed Brown       lambda[2] = lmax;
139*c4762a1bSJed Brown       break;
140*c4762a1bSJed Brown     default:
141*c4762a1bSJed Brown       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1 or 2, cannot be %d", user->metOpt);
142*c4762a1bSJed Brown     }
143*c4762a1bSJed Brown     /* Only set the diagonal */
144*c4762a1bSJed Brown     ierr = DMPlexPointLocalRef(mdm, p, met, &pmet);CHKERRQ(ierr);
145*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) pmet[d*(dim+1)] = lambda[d];
146*c4762a1bSJed Brown   }
147*c4762a1bSJed Brown   ierr = VecRestoreArray(*metric, &met);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = DMDestroy(&mdm);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = PetscFree(lambda);CHKERRQ(ierr);
151*c4762a1bSJed Brown   PetscFunctionReturn(0);
152*c4762a1bSJed Brown }
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
155*c4762a1bSJed Brown {
156*c4762a1bSJed Brown   u[0] = x[0] + x[1];
157*c4762a1bSJed Brown   return 0;
158*c4762a1bSJed Brown }
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
161*c4762a1bSJed Brown                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
162*c4762a1bSJed Brown                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
163*c4762a1bSJed Brown                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
164*c4762a1bSJed Brown {
165*c4762a1bSJed Brown   g0[0] = 1.0;
166*c4762a1bSJed Brown }
167*c4762a1bSJed Brown 
168*c4762a1bSJed Brown static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user)
169*c4762a1bSJed Brown {
170*c4762a1bSJed Brown   PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
171*c4762a1bSJed Brown   KSP              ksp;
172*c4762a1bSJed Brown   PetscDS          prob;
173*c4762a1bSJed Brown   PetscFE          fe;
174*c4762a1bSJed Brown   Mat              Interp, mass;
175*c4762a1bSJed Brown   Vec              u, ua, scaling, ones, massLumped, rhs, uproj;
176*c4762a1bSJed Brown   PetscReal        error;
177*c4762a1bSJed Brown   PetscInt         dim;
178*c4762a1bSJed Brown   MPI_Comm         comm;
179*c4762a1bSJed Brown   PetscErrorCode   ierr;
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown   PetscFunctionBeginUser;
182*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, NULL, -1, &fe);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = PetscDSSetDiscretization(prob, 0, (PetscObject) fe);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
188*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, NULL, -1, &fe);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = DMGetDS(dma, &prob);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = PetscDSSetDiscretization(prob, 0, (PetscObject) fe);CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
192*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, identity, NULL, NULL, NULL);CHKERRQ(ierr);
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown   funcs[0] = linear;
195*c4762a1bSJed Brown   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
196*c4762a1bSJed Brown   ierr = DMGetGlobalVector(dma, &ua);CHKERRQ(ierr);
197*c4762a1bSJed Brown   ierr = DMGetGlobalVector(dma, &ones);CHKERRQ(ierr);
198*c4762a1bSJed Brown   ierr = DMGetGlobalVector(dma, &massLumped);CHKERRQ(ierr);
199*c4762a1bSJed Brown   ierr = DMGetGlobalVector(dma, &rhs);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr = DMGetGlobalVector(dma, &uproj);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "Original");CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-orig_vec_view");CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = DMComputeL2Diff(dm, 0.0, funcs, NULL, u, &error);CHKERRQ(ierr);
205*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double) error);CHKERRQ(ierr);
206*c4762a1bSJed Brown   ierr = DMCreateInterpolation(dm, dma, &Interp, &scaling);CHKERRQ(ierr);
207*c4762a1bSJed Brown   ierr = MatInterpolate(Interp, u, ua);CHKERRQ(ierr);
208*c4762a1bSJed Brown   ierr = MatDestroy(&Interp);CHKERRQ(ierr);
209*c4762a1bSJed Brown   ierr = VecDestroy(&scaling);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) ua, "Interpolation");CHKERRQ(ierr);
211*c4762a1bSJed Brown   ierr = VecViewFromOptions(ua, NULL, "-interp_vec_view");CHKERRQ(ierr);
212*c4762a1bSJed Brown   ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, ua, &error);CHKERRQ(ierr);
213*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double) error);CHKERRQ(ierr);
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown   ierr = VecSet(ones, 1.0);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = DMPlexComputeJacobianAction(dma, NULL, 0, 0, ua, NULL, ones, massLumped, user);CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) massLumped, "Lumped mass");CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = VecViewFromOptions(massLumped, NULL, "-mass_vec_view");CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = DMCreateMassMatrix(dm, dma, &mass);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = MatMult(mass, u, rhs);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = MatDestroy(&mass);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = VecViewFromOptions(rhs, NULL, "-lumped_rhs_view");CHKERRQ(ierr);
223*c4762a1bSJed Brown   ierr = VecPointwiseDivide(uproj, rhs, massLumped);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) uproj, "Different Lumped Projection");CHKERRQ(ierr);
225*c4762a1bSJed Brown   ierr = VecViewFromOptions(uproj, NULL, "-lumped_rhs_vec_view");CHKERRQ(ierr);
226*c4762a1bSJed Brown   ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr);
227*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Lumped (rhs) L2 Error: %g\n", (double) error);CHKERRQ(ierr);
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   ierr = DMCreateMatrix(dma, &mass);CHKERRQ(ierr);
230*c4762a1bSJed Brown   ierr = DMPlexSNESComputeJacobianFEM(dma, ua, mass, mass, user);CHKERRQ(ierr);
231*c4762a1bSJed Brown   ierr = MatViewFromOptions(mass, NULL, "-mass_mat_view");CHKERRQ(ierr);
232*c4762a1bSJed Brown   ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr);
233*c4762a1bSJed Brown   ierr = KSPSetOperators(ksp, mass, mass);CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = KSPSolve(ksp, rhs, uproj);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
237*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) uproj, "Full Projection");CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = VecViewFromOptions(uproj, NULL, "-proj_vec_view");CHKERRQ(ierr);
239*c4762a1bSJed Brown   ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double) error);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = MatDestroy(&mass);CHKERRQ(ierr);
242*c4762a1bSJed Brown 
243*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dma, &ua);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dma, &ones);CHKERRQ(ierr);
246*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dma, &massLumped);CHKERRQ(ierr);
247*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dma, &rhs);CHKERRQ(ierr);
248*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dma, &uproj);CHKERRQ(ierr);
249*c4762a1bSJed Brown   PetscFunctionReturn(0);
250*c4762a1bSJed Brown }
251*c4762a1bSJed Brown 
252*c4762a1bSJed Brown int main (int argc, char * argv[]) {
253*c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
254*c4762a1bSJed Brown   DMLabel        bdLabel = NULL;
255*c4762a1bSJed Brown   MPI_Comm       comm;
256*c4762a1bSJed Brown   DM             dma, odm;
257*c4762a1bSJed Brown   Vec            metric;
258*c4762a1bSJed Brown   size_t         len;
259*c4762a1bSJed Brown   PetscErrorCode ierr;
260*c4762a1bSJed Brown 
261*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
262*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
263*c4762a1bSJed Brown   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown   ierr = CreateMesh(comm, &user);CHKERRQ(ierr);
266*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) user.dm, "DMinit");CHKERRQ(ierr);
267*c4762a1bSJed Brown   ierr = DMViewFromOptions(user.dm, NULL, "-init_dm_view");CHKERRQ(ierr);
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown   odm  = user.dm;
270*c4762a1bSJed Brown   ierr = DMPlexDistributeOverlap(odm, 1, NULL, &user.dm);CHKERRQ(ierr);
271*c4762a1bSJed Brown   if (!user.dm) {user.dm = odm;}
272*c4762a1bSJed Brown   else          {ierr = DMDestroy(&odm);CHKERRQ(ierr);}
273*c4762a1bSJed Brown   ierr = ComputeMetric(user.dm, &user, &metric);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = PetscStrlen(user.bdLabel, &len);CHKERRQ(ierr);
275*c4762a1bSJed Brown   if (len) {
276*c4762a1bSJed Brown     ierr = DMCreateLabel(user.dm, user.bdLabel);CHKERRQ(ierr);
277*c4762a1bSJed Brown     ierr = DMGetLabel(user.dm, user.bdLabel, &bdLabel);CHKERRQ(ierr);
278*c4762a1bSJed Brown   }
279*c4762a1bSJed Brown   ierr = DMAdaptMetric(user.dm, metric, bdLabel, &dma);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) dma, "DMadapt");CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject) dma, "adapt_");CHKERRQ(ierr);
282*c4762a1bSJed Brown   ierr = DMViewFromOptions(dma, NULL, "-dm_view");CHKERRQ(ierr);
283*c4762a1bSJed Brown   if (user.doL2) {ierr = TestL2Projection(user.dm, dma, &user);CHKERRQ(ierr);}
284*c4762a1bSJed Brown   ierr = DMDestroy(&dma);CHKERRQ(ierr);
285*c4762a1bSJed Brown   ierr = VecDestroy(&metric);CHKERRQ(ierr);
286*c4762a1bSJed Brown   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
287*c4762a1bSJed Brown   PetscFinalize();
288*c4762a1bSJed Brown   return 0;
289*c4762a1bSJed Brown }
290*c4762a1bSJed Brown 
291*c4762a1bSJed Brown /*TEST
292*c4762a1bSJed Brown 
293*c4762a1bSJed Brown   test:
294*c4762a1bSJed Brown     suffix: 0
295*c4762a1bSJed Brown     requires: pragmatic
296*c4762a1bSJed Brown     TODO: broken
297*c4762a1bSJed Brown     args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view
298*c4762a1bSJed Brown   test:
299*c4762a1bSJed Brown     suffix: 1
300*c4762a1bSJed Brown     requires: pragmatic
301*c4762a1bSJed Brown     TODO: broken
302*c4762a1bSJed Brown     args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 1 -bdLabel marker -met 2 -init_dm_view -adapt_dm_view
303*c4762a1bSJed Brown   test:
304*c4762a1bSJed Brown     suffix: 2
305*c4762a1bSJed Brown     requires: pragmatic
306*c4762a1bSJed Brown     TODO: broken
307*c4762a1bSJed Brown     args: -dim 3 -nbrVerEdge 5 -met 2 -init_dm_view -adapt_dm_view
308*c4762a1bSJed Brown   test:
309*c4762a1bSJed Brown     suffix: 3
310*c4762a1bSJed Brown     requires: pragmatic
311*c4762a1bSJed Brown     TODO: broken
312*c4762a1bSJed Brown     args: -dim 3 -nbrVerEdge 5 -bdLabel marker -met 2 -init_dm_view -adapt_dm_view
313*c4762a1bSJed Brown   test:
314*c4762a1bSJed Brown     suffix: 4
315*c4762a1bSJed Brown     requires: pragmatic
316*c4762a1bSJed Brown     TODO: broken
317*c4762a1bSJed Brown     nsize: 2
318*c4762a1bSJed Brown     args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view
319*c4762a1bSJed Brown   test:
320*c4762a1bSJed Brown     suffix: 5
321*c4762a1bSJed Brown     requires: pragmatic
322*c4762a1bSJed Brown     TODO: broken
323*c4762a1bSJed Brown     nsize: 4
324*c4762a1bSJed Brown     args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view
325*c4762a1bSJed Brown   test:
326*c4762a1bSJed Brown     suffix: 6
327*c4762a1bSJed Brown     requires: pragmatic
328*c4762a1bSJed Brown     TODO: broken
329*c4762a1bSJed Brown     nsize: 2
330*c4762a1bSJed Brown     args: -dim 3 -nbrVerEdge 10 -dm_plex_separate_marker 0 -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view
331*c4762a1bSJed Brown   test:
332*c4762a1bSJed Brown     suffix: 7
333*c4762a1bSJed Brown     requires: pragmatic
334*c4762a1bSJed Brown     TODO: broken
335*c4762a1bSJed Brown     nsize: 5
336*c4762a1bSJed Brown     args: -dim 2 -nbrVerEdge 20 -dm_plex_separate_marker 0 -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view
337*c4762a1bSJed Brown   test:
338*c4762a1bSJed Brown     suffix: proj_0
339*c4762a1bSJed Brown     requires: pragmatic
340*c4762a1bSJed Brown     TODO: broken
341*c4762a1bSJed Brown     args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -init_dm_view -adapt_dm_view -do_L2 -petscspace_degree 1 -petscfe_default_quadrature_order 1 -dm_plex_hash_location -pc_type lu
342*c4762a1bSJed Brown   test:
343*c4762a1bSJed Brown     suffix: proj_1
344*c4762a1bSJed Brown     requires: pragmatic
345*c4762a1bSJed Brown     TODO: broken
346*c4762a1bSJed Brown     args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 0 -init_dm_view -adapt_dm_view -do_L2 -petscspace_degree 2 -petscfe_default_quadrature_order 4 -dm_plex_hash_location -pc_type lu
347*c4762a1bSJed Brown 
348*c4762a1bSJed Brown TEST*/
349