xref: /petsc/src/dm/impls/plex/tests/ex46.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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);
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-num_refine", "Refine cycles", "ex46.c", options->Nr, &options->Nr, NULL));
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   PetscFunctionBeginUser;
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(*dm, user));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
55a12d352dSMatthew G. Knepley   PetscFunctionReturn(0);
56a12d352dSMatthew G. Knepley }
57a12d352dSMatthew G. Knepley 
58a12d352dSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
59a12d352dSMatthew G. Knepley {
60a12d352dSMatthew G. Knepley   DM             cdm = dm;
61a12d352dSMatthew G. Knepley   PetscFE        fe;
62a12d352dSMatthew G. Knepley   PetscSpace     sp;
63a12d352dSMatthew G. Knepley   PetscInt       dim, deg;
64a12d352dSMatthew G. Knepley 
65a12d352dSMatthew G. Knepley   PetscFunctionBeginUser;
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, -1, &fe));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, "scalar"));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
71a12d352dSMatthew G. Knepley   while (cdm) {
72*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm,cdm));
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
74a12d352dSMatthew G. Knepley   }
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEGetBasisSpace(fe, &sp));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSpaceGetDegree(sp, &deg, NULL));
77a12d352dSMatthew G. Knepley   switch (deg) {
78a12d352dSMatthew G. Knepley   case 0: user->funcs[0] = constant;break;
79a12d352dSMatthew G. Knepley   case 1: user->funcs[0] = linear;break;
80a12d352dSMatthew G. Knepley   case 2: user->funcs[0] = quadratic;break;
81a12d352dSMatthew G. Knepley   case 3: user->funcs[0] = cubic;break;
8298921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Could not determine function to test for degree %D", deg);
83a12d352dSMatthew G. Knepley   }
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
85a12d352dSMatthew G. Knepley   PetscFunctionReturn(0);
86a12d352dSMatthew G. Knepley }
87a12d352dSMatthew G. Knepley 
88a12d352dSMatthew G. Knepley static PetscErrorCode CheckError(DM dm, Vec u, PetscSimplePointFunc funcs[])
89a12d352dSMatthew G. Knepley {
90a12d352dSMatthew G. Knepley   PetscReal      error, tol = PETSC_SMALL;
91a12d352dSMatthew G. Knepley   MPI_Comm       comm;
92a12d352dSMatthew G. Knepley 
93a12d352dSMatthew G. Knepley   PetscFunctionBeginUser;
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMComputeL2Diff(dm, 0.0, funcs, NULL, u, &error));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
96*5f80ce2aSJacob Faibussowitsch   if (error > tol) CHKERRQ(PetscPrintf(comm, "Function tests FAIL at tolerance %g error %g\n", (double)tol,(double) error));
97*5f80ce2aSJacob Faibussowitsch   else             CHKERRQ(PetscPrintf(comm, "Function tests pass at tolerance %g\n", (double)tol));
98a12d352dSMatthew G. Knepley   PetscFunctionReturn(0);
99a12d352dSMatthew G. Knepley }
100a12d352dSMatthew G. Knepley 
101a12d352dSMatthew G. Knepley int main(int argc, char **argv)
102a12d352dSMatthew G. Knepley {
103a12d352dSMatthew G. Knepley   DM             dm;
104a12d352dSMatthew G. Knepley   Vec            u;
105a12d352dSMatthew G. Knepley   AppCtx         user;
106a12d352dSMatthew G. Knepley   PetscInt       cStart, cEnd, c, r;
107a12d352dSMatthew G. Knepley   PetscErrorCode ierr;
108a12d352dSMatthew G. Knepley 
109a12d352dSMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, &user));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &u));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMProjectFunction(dm, 0.0, user.funcs, NULL, INSERT_ALL_VALUES, u));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckError(dm, u, user.funcs));
116a12d352dSMatthew G. Knepley   for (r = 0; r < user.Nr; ++r) {
117a12d352dSMatthew G. Knepley     DM      adm;
118a12d352dSMatthew G. Knepley     DMLabel adapt;
119a12d352dSMatthew G. Knepley     Vec     au;
120a12d352dSMatthew G. Knepley     Mat     Interp;
121a12d352dSMatthew G. Knepley 
122*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adapt));
123*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelSetDefaultValue(adapt, DM_ADAPT_COARSEN));
124*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
125a12d352dSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
126*5f80ce2aSJacob Faibussowitsch       if (c % 2) CHKERRQ(DMLabelSetValue(adapt, c, DM_ADAPT_REFINE));
127a12d352dSMatthew G. Knepley     }
128*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAdaptLabel(dm, adapt, &adm));
129*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelDestroy(&adapt));
130*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) adm, "Adapted Mesh"));
131*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMViewFromOptions(adm, NULL, "-dm_view"));
132a12d352dSMatthew G. Knepley 
133*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateInterpolation(dm, adm, &Interp, NULL));
134*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalVector(adm, &au));
135*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatInterpolate(Interp, u, au));
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(CheckError(adm, au, user.funcs));
137*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Interp));
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreGlobalVector(dm, &u));
139*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&dm));
140a12d352dSMatthew G. Knepley     dm   = adm;
141a12d352dSMatthew G. Knepley     u    = au;
142a12d352dSMatthew G. Knepley   }
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &u));
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
145a12d352dSMatthew G. Knepley   ierr = PetscFinalize();
146a12d352dSMatthew G. Knepley   return ierr;
147a12d352dSMatthew G. Knepley }
148a12d352dSMatthew G. Knepley 
149a12d352dSMatthew G. Knepley /*TEST
150a12d352dSMatthew G. Knepley 
151a12d352dSMatthew G. Knepley   test:
152a12d352dSMatthew G. Knepley     suffix: 0
153a12d352dSMatthew G. Knepley     args: -num_refine 4 -petscspace_degree 3 \
154a12d352dSMatthew G. Knepley           -dm_plex_dim 1 -dm_plex_box_faces 5 -dm_plex_transform_type refine_1d -dm_plex_hash_location -dm_view
155a12d352dSMatthew G. Knepley 
156a12d352dSMatthew G. Knepley TEST*/
157