xref: /petsc/src/dm/impls/plex/tests/ex46.c (revision a12d352d5eb7e163a98d06e375dd80b20fcd305b)
1*a12d352dSMatthew G. Knepley static char help[] = "Tests 1D nested mesh refinement.\n\n";
2*a12d352dSMatthew G. Knepley 
3*a12d352dSMatthew G. Knepley #include <petscdmplex.h>
4*a12d352dSMatthew G. Knepley #include <petscds.h>
5*a12d352dSMatthew G. Knepley 
6*a12d352dSMatthew G. Knepley typedef struct {
7*a12d352dSMatthew G. Knepley   PetscInt             Nr;       /* Number of refinements */
8*a12d352dSMatthew G. Knepley   PetscSimplePointFunc funcs[1]; /* Functions to test */
9*a12d352dSMatthew G. Knepley } AppCtx;
10*a12d352dSMatthew G. Knepley 
11*a12d352dSMatthew G. Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
12*a12d352dSMatthew G. Knepley {
13*a12d352dSMatthew G. Knepley   u[0] = 1.;
14*a12d352dSMatthew G. Knepley   return 0;
15*a12d352dSMatthew G. Knepley }
16*a12d352dSMatthew G. Knepley 
17*a12d352dSMatthew G. Knepley static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
18*a12d352dSMatthew G. Knepley {
19*a12d352dSMatthew G. Knepley   u[0] = 2.*x[0] + 1.;
20*a12d352dSMatthew G. Knepley   return 0;
21*a12d352dSMatthew G. Knepley }
22*a12d352dSMatthew G. Knepley 
23*a12d352dSMatthew G. Knepley static PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
24*a12d352dSMatthew G. Knepley {
25*a12d352dSMatthew G. Knepley   u[0] = 3.*x[0]*x[0] + 2.*x[0] + 1.;
26*a12d352dSMatthew G. Knepley   return 0;
27*a12d352dSMatthew G. Knepley }
28*a12d352dSMatthew G. Knepley 
29*a12d352dSMatthew G. Knepley static PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
30*a12d352dSMatthew G. Knepley {
31*a12d352dSMatthew G. Knepley   u[0] = 4.*x[0]*x[0]*x[0] + 3.*x[0]*x[0] + 2.*x[0] + 1.;
32*a12d352dSMatthew G. Knepley   return 0;
33*a12d352dSMatthew G. Knepley }
34*a12d352dSMatthew G. Knepley 
35*a12d352dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
36*a12d352dSMatthew G. Knepley {
37*a12d352dSMatthew G. Knepley   PetscErrorCode ierr;
38*a12d352dSMatthew G. Knepley 
39*a12d352dSMatthew G. Knepley   PetscFunctionBeginUser;
40*a12d352dSMatthew G. Knepley   options->Nr = 1;
41*a12d352dSMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "1D Refinement Options", "DMPLEX");CHKERRQ(ierr);
42*a12d352dSMatthew G. Knepley   ierr = PetscOptionsInt("-num_refine", "Refine cycles", "ex46.c", options->Nr, &options->Nr, NULL);CHKERRQ(ierr);
43*a12d352dSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
44*a12d352dSMatthew G. Knepley   PetscFunctionReturn(0);
45*a12d352dSMatthew G. Knepley }
46*a12d352dSMatthew G. Knepley 
47*a12d352dSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
48*a12d352dSMatthew G. Knepley {
49*a12d352dSMatthew G. Knepley   PetscErrorCode ierr;
50*a12d352dSMatthew G. Knepley 
51*a12d352dSMatthew G. Knepley   PetscFunctionBeginUser;
52*a12d352dSMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
53*a12d352dSMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
54*a12d352dSMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
55*a12d352dSMatthew G. Knepley   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
56*a12d352dSMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
57*a12d352dSMatthew G. Knepley   PetscFunctionReturn(0);
58*a12d352dSMatthew G. Knepley }
59*a12d352dSMatthew G. Knepley 
60*a12d352dSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
61*a12d352dSMatthew G. Knepley {
62*a12d352dSMatthew G. Knepley   DM             cdm = dm;
63*a12d352dSMatthew G. Knepley   PetscFE        fe;
64*a12d352dSMatthew G. Knepley   PetscSpace     sp;
65*a12d352dSMatthew G. Knepley   PetscInt       dim, deg;
66*a12d352dSMatthew G. Knepley   PetscErrorCode ierr;
67*a12d352dSMatthew G. Knepley 
68*a12d352dSMatthew G. Knepley   PetscFunctionBeginUser;
69*a12d352dSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
70*a12d352dSMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, -1, &fe);CHKERRQ(ierr);
71*a12d352dSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) fe, "scalar");CHKERRQ(ierr);
72*a12d352dSMatthew G. Knepley   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
73*a12d352dSMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
74*a12d352dSMatthew G. Knepley   while (cdm) {
75*a12d352dSMatthew G. Knepley     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
76*a12d352dSMatthew G. Knepley     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
77*a12d352dSMatthew G. Knepley   }
78*a12d352dSMatthew G. Knepley   ierr = PetscFEGetBasisSpace(fe, &sp);CHKERRQ(ierr);
79*a12d352dSMatthew G. Knepley   ierr = PetscSpaceGetDegree(sp, &deg, NULL);CHKERRQ(ierr);
80*a12d352dSMatthew G. Knepley   switch (deg) {
81*a12d352dSMatthew G. Knepley   case 0: user->funcs[0] = constant;break;
82*a12d352dSMatthew G. Knepley   case 1: user->funcs[0] = linear;break;
83*a12d352dSMatthew G. Knepley   case 2: user->funcs[0] = quadratic;break;
84*a12d352dSMatthew G. Knepley   case 3: user->funcs[0] = cubic;break;
85*a12d352dSMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Could not determine function to test for degree %D", deg);
86*a12d352dSMatthew G. Knepley   }
87*a12d352dSMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
88*a12d352dSMatthew G. Knepley   PetscFunctionReturn(0);
89*a12d352dSMatthew G. Knepley }
90*a12d352dSMatthew G. Knepley 
91*a12d352dSMatthew G. Knepley static PetscErrorCode CheckError(DM dm, Vec u, PetscSimplePointFunc funcs[])
92*a12d352dSMatthew G. Knepley {
93*a12d352dSMatthew G. Knepley   PetscReal      error, tol = PETSC_SMALL;
94*a12d352dSMatthew G. Knepley   MPI_Comm       comm;
95*a12d352dSMatthew G. Knepley   PetscErrorCode ierr;
96*a12d352dSMatthew G. Knepley 
97*a12d352dSMatthew G. Knepley   PetscFunctionBeginUser;
98*a12d352dSMatthew G. Knepley   ierr = DMComputeL2Diff(dm, 0.0, funcs, NULL, u, &error);CHKERRQ(ierr);
99*a12d352dSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
100*a12d352dSMatthew G. Knepley   if (error > tol) {ierr = PetscPrintf(comm, "Function tests FAIL at tolerance %g error %g\n", (double)tol,(double) error);CHKERRQ(ierr);}
101*a12d352dSMatthew G. Knepley   else             {ierr = PetscPrintf(comm, "Function tests pass at tolerance %g\n", (double)tol);CHKERRQ(ierr);}
102*a12d352dSMatthew G. Knepley   PetscFunctionReturn(0);
103*a12d352dSMatthew G. Knepley }
104*a12d352dSMatthew G. Knepley 
105*a12d352dSMatthew G. Knepley int main(int argc, char **argv)
106*a12d352dSMatthew G. Knepley {
107*a12d352dSMatthew G. Knepley   DM             dm;
108*a12d352dSMatthew G. Knepley   Vec            u;
109*a12d352dSMatthew G. Knepley   AppCtx         user;
110*a12d352dSMatthew G. Knepley   PetscInt       cStart, cEnd, c, r;
111*a12d352dSMatthew G. Knepley   PetscErrorCode ierr;
112*a12d352dSMatthew G. Knepley 
113*a12d352dSMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
114*a12d352dSMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
115*a12d352dSMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
116*a12d352dSMatthew G. Knepley   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
117*a12d352dSMatthew G. Knepley   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
118*a12d352dSMatthew G. Knepley   ierr = DMProjectFunction(dm, 0.0, user.funcs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
119*a12d352dSMatthew G. Knepley   ierr = CheckError(dm, u, user.funcs);CHKERRQ(ierr);
120*a12d352dSMatthew G. Knepley   for (r = 0; r < user.Nr; ++r) {
121*a12d352dSMatthew G. Knepley     DM      adm;
122*a12d352dSMatthew G. Knepley     DMLabel adapt;
123*a12d352dSMatthew G. Knepley     Vec     au;
124*a12d352dSMatthew G. Knepley     Mat     Interp;
125*a12d352dSMatthew G. Knepley 
126*a12d352dSMatthew G. Knepley     ierr = DMLabelCreate(PETSC_COMM_SELF, "adapt", &adapt);CHKERRQ(ierr);
127*a12d352dSMatthew G. Knepley     ierr = DMLabelSetDefaultValue(adapt, DM_ADAPT_COARSEN);CHKERRQ(ierr);
128*a12d352dSMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
129*a12d352dSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
130*a12d352dSMatthew G. Knepley       if (c % 2) {ierr = DMLabelSetValue(adapt, c, DM_ADAPT_REFINE);CHKERRQ(ierr);}
131*a12d352dSMatthew G. Knepley     }
132*a12d352dSMatthew G. Knepley     ierr = DMAdaptLabel(dm, adapt, &adm);CHKERRQ(ierr);
133*a12d352dSMatthew G. Knepley     ierr = DMLabelDestroy(&adapt);CHKERRQ(ierr);
134*a12d352dSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) adm, "Adapted Mesh");CHKERRQ(ierr);
135*a12d352dSMatthew G. Knepley     ierr = DMViewFromOptions(adm, NULL, "-dm_view");CHKERRQ(ierr);
136*a12d352dSMatthew G. Knepley 
137*a12d352dSMatthew G. Knepley     ierr = DMCreateInterpolation(dm, adm, &Interp, NULL);CHKERRQ(ierr);
138*a12d352dSMatthew G. Knepley     ierr = DMGetGlobalVector(adm, &au);CHKERRQ(ierr);
139*a12d352dSMatthew G. Knepley     ierr = MatInterpolate(Interp, u, au);CHKERRQ(ierr);
140*a12d352dSMatthew G. Knepley     ierr = CheckError(adm, au, user.funcs);CHKERRQ(ierr);
141*a12d352dSMatthew G. Knepley     ierr = MatDestroy(&Interp);CHKERRQ(ierr);
142*a12d352dSMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
143*a12d352dSMatthew G. Knepley     ierr = DMDestroy(&dm);CHKERRQ(ierr);
144*a12d352dSMatthew G. Knepley     dm   = adm;
145*a12d352dSMatthew G. Knepley     u    = au;
146*a12d352dSMatthew G. Knepley   }
147*a12d352dSMatthew G. Knepley   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
148*a12d352dSMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
149*a12d352dSMatthew G. Knepley   ierr = PetscFinalize();
150*a12d352dSMatthew G. Knepley   return ierr;
151*a12d352dSMatthew G. Knepley }
152*a12d352dSMatthew G. Knepley 
153*a12d352dSMatthew G. Knepley /*TEST
154*a12d352dSMatthew G. Knepley 
155*a12d352dSMatthew G. Knepley   test:
156*a12d352dSMatthew G. Knepley     suffix: 0
157*a12d352dSMatthew G. Knepley     args: -num_refine 4 -petscspace_degree 3 \
158*a12d352dSMatthew G. Knepley           -dm_plex_dim 1 -dm_plex_box_faces 5 -dm_plex_transform_type refine_1d -dm_plex_hash_location -dm_view
159*a12d352dSMatthew G. Knepley 
160*a12d352dSMatthew G. Knepley TEST*/
161