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