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 */
88434afd1SBarry Smith PetscSimplePointFn *funcs[2]; /* Functions to test */
9a12d352dSMatthew G. Knepley } AppCtx;
10a12d352dSMatthew G. Knepley
constant(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)11*2a8381b2SBarry Smith static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
12d71ae5a4SJacob Faibussowitsch {
13a12d352dSMatthew G. Knepley u[0] = 1.;
143ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
15a12d352dSMatthew G. Knepley }
16a12d352dSMatthew G. Knepley
linear(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)17*2a8381b2SBarry Smith static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
18d71ae5a4SJacob Faibussowitsch {
19a12d352dSMatthew G. Knepley u[0] = 2. * x[0] + 1.;
203ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
21a12d352dSMatthew G. Knepley }
22a12d352dSMatthew G. Knepley
quadratic(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)23*2a8381b2SBarry Smith static PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
24d71ae5a4SJacob Faibussowitsch {
25a12d352dSMatthew G. Knepley u[0] = 3. * x[0] * x[0] + 2. * x[0] + 1.;
263ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
27a12d352dSMatthew G. Knepley }
28a12d352dSMatthew G. Knepley
cubic(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)29*2a8381b2SBarry Smith static PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
30d71ae5a4SJacob Faibussowitsch {
31a12d352dSMatthew G. Knepley u[0] = 4. * x[0] * x[0] * x[0] + 3. * x[0] * x[0] + 2. * x[0] + 1.;
323ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
33a12d352dSMatthew G. Knepley }
34a12d352dSMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * options)35d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
36d71ae5a4SJacob Faibussowitsch {
37a12d352dSMatthew G. Knepley PetscFunctionBeginUser;
38a12d352dSMatthew G. Knepley options->Nr = 1;
39d0609cedSBarry Smith PetscOptionsBegin(comm, "", "1D Refinement Options", "DMPLEX");
409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-num_refine", "Refine cycles", "ex46.c", options->Nr, &options->Nr, NULL));
41d0609cedSBarry Smith PetscOptionsEnd();
423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
43a12d352dSMatthew G. Knepley }
44a12d352dSMatthew G. Knepley
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)45d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
46d71ae5a4SJacob Faibussowitsch {
47a12d352dSMatthew G. Knepley PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
499566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
519566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user));
529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
54a12d352dSMatthew G. Knepley }
55a12d352dSMatthew G. Knepley
SetupDiscretization(DM dm,AppCtx * user)56d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
57d71ae5a4SJacob Faibussowitsch {
58a12d352dSMatthew G. Knepley DM cdm = dm;
59a12d352dSMatthew G. Knepley PetscFE fe;
60a12d352dSMatthew G. Knepley PetscSpace sp;
61a12d352dSMatthew G. Knepley PetscInt dim, deg;
62a12d352dSMatthew G. Knepley
63a12d352dSMatthew G. Knepley PetscFunctionBeginUser;
649566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
659566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, -1, &fe));
669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "scalar"));
679566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
68ddd1d095SMatthew G. Knepley PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe));
699566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
70a12d352dSMatthew G. Knepley while (cdm) {
719566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
729566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
73a12d352dSMatthew G. Knepley }
749566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &sp));
759566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(sp, °, NULL));
76a12d352dSMatthew G. Knepley switch (deg) {
77d71ae5a4SJacob Faibussowitsch case 0:
78d71ae5a4SJacob Faibussowitsch user->funcs[0] = constant;
79d71ae5a4SJacob Faibussowitsch break;
80d71ae5a4SJacob Faibussowitsch case 1:
81d71ae5a4SJacob Faibussowitsch user->funcs[0] = linear;
82d71ae5a4SJacob Faibussowitsch break;
83d71ae5a4SJacob Faibussowitsch case 2:
84d71ae5a4SJacob Faibussowitsch user->funcs[0] = quadratic;
85d71ae5a4SJacob Faibussowitsch break;
86d71ae5a4SJacob Faibussowitsch case 3:
87d71ae5a4SJacob Faibussowitsch user->funcs[0] = cubic;
88d71ae5a4SJacob Faibussowitsch break;
89d71ae5a4SJacob Faibussowitsch default:
90d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Could not determine function to test for degree %" PetscInt_FMT, deg);
91a12d352dSMatthew G. Knepley }
92ddd1d095SMatthew G. Knepley user->funcs[1] = user->funcs[0];
939566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
95a12d352dSMatthew G. Knepley }
96a12d352dSMatthew G. Knepley
CheckError(DM dm,Vec u,PetscSimplePointFn * funcs[])978434afd1SBarry Smith static PetscErrorCode CheckError(DM dm, Vec u, PetscSimplePointFn *funcs[])
98d71ae5a4SJacob Faibussowitsch {
99a12d352dSMatthew G. Knepley PetscReal error, tol = PETSC_SMALL;
100a12d352dSMatthew G. Knepley MPI_Comm comm;
101a12d352dSMatthew G. Knepley
102a12d352dSMatthew G. Knepley PetscFunctionBeginUser;
1039566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, funcs, NULL, u, &error));
1049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1059566063dSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL at tolerance %g error %g\n", (double)tol, (double)error));
1069566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass at tolerance %g\n", (double)tol));
1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
108a12d352dSMatthew G. Knepley }
109a12d352dSMatthew G. Knepley
main(int argc,char ** argv)110d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
111d71ae5a4SJacob Faibussowitsch {
112a12d352dSMatthew G. Knepley DM dm;
113a12d352dSMatthew G. Knepley Vec u;
114a12d352dSMatthew G. Knepley AppCtx user;
115a12d352dSMatthew G. Knepley PetscInt cStart, cEnd, c, r;
116a12d352dSMatthew G. Knepley
117327415f7SBarry Smith PetscFunctionBeginUser;
1189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1199566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1209566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1219566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user));
1229566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u));
1239566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, user.funcs, NULL, INSERT_ALL_VALUES, u));
1249566063dSJacob Faibussowitsch PetscCall(CheckError(dm, u, user.funcs));
125a12d352dSMatthew G. Knepley for (r = 0; r < user.Nr; ++r) {
126a12d352dSMatthew G. Knepley DM adm;
127a12d352dSMatthew G. Knepley DMLabel adapt;
128a12d352dSMatthew G. Knepley Vec au;
129a12d352dSMatthew G. Knepley Mat Interp;
130a12d352dSMatthew G. Knepley
1319566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adapt));
1329566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(adapt, DM_ADAPT_COARSEN));
1339566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
134a12d352dSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) {
1359566063dSJacob Faibussowitsch if (c % 2) PetscCall(DMLabelSetValue(adapt, c, DM_ADAPT_REFINE));
136a12d352dSMatthew G. Knepley }
1379566063dSJacob Faibussowitsch PetscCall(DMAdaptLabel(dm, adapt, &adm));
1389566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adapt));
1399566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)adm, "Adapted Mesh"));
1409566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(adm, NULL, "-dm_view"));
141a12d352dSMatthew G. Knepley
1429566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dm, adm, &Interp, NULL));
1439566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(adm, &au));
1449566063dSJacob Faibussowitsch PetscCall(MatInterpolate(Interp, u, au));
1459566063dSJacob Faibussowitsch PetscCall(CheckError(adm, au, user.funcs));
1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp));
1479566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u));
1489566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
149a12d352dSMatthew G. Knepley dm = adm;
150a12d352dSMatthew G. Knepley u = au;
151a12d352dSMatthew G. Knepley }
1529566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u));
1539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
1549566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
155b122ec5aSJacob Faibussowitsch return 0;
156a12d352dSMatthew G. Knepley }
157a12d352dSMatthew G. Knepley
158a12d352dSMatthew G. Knepley /*TEST
159a12d352dSMatthew G. Knepley
160a12d352dSMatthew G. Knepley test:
161a12d352dSMatthew G. Knepley suffix: 0
162a12d352dSMatthew G. Knepley args: -num_refine 4 -petscspace_degree 3 \
163a12d352dSMatthew G. Knepley -dm_plex_dim 1 -dm_plex_box_faces 5 -dm_plex_transform_type refine_1d -dm_plex_hash_location -dm_view
164a12d352dSMatthew G. Knepley
165a12d352dSMatthew G. Knepley TEST*/
166