Actual source code: ex2.c
1: const char help[] = "Simple example to get equally space points in high-order elements (and XGC mirror)";
3: #include <petscfe.h>
4: #include <petscdmplex.h>
5: static PetscErrorCode x(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
6: {
7: PetscFunctionBegin;
8: u[0] = x[0];
9: PetscFunctionReturn(PETSC_SUCCESS);
10: }
12: int main(int argc, char **argv)
13: {
14: PetscInt dim = 2, cells[] = {1, 1, 1};
15: DM K;
16: PetscReal radius = 2, lo[] = {-radius, -radius, -radius}, hi[] = {radius, radius, radius};
17: DMBoundaryType periodicity[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
18: PetscFE fe;
19: PetscErrorCode (*initu[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
20: Vec X;
22: PetscFunctionBeginUser;
23: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
24: PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for PETSCDUALSPACELAGRANGE test", "none");
25: PetscCall(PetscOptionsRangeInt("-dim", "The spatial dimension", "ex1.c", dim, &dim, NULL, 0, 3));
26: PetscOptionsEnd();
28: PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_SELF, dim, PETSC_FALSE, cells, lo, hi, periodicity, PETSC_TRUE, 0, PETSC_TRUE, &K));
29: PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, PETSC_DECIDE, &fe));
30: PetscCall(DMSetField(K, 0, NULL, (PetscObject)fe));
31: PetscCall(PetscFEDestroy(&fe));
32: PetscCall(DMCreateDS(K));
34: initu[0] = x;
35: PetscCall(DMCreateGlobalVector(K, &X));
36: PetscCall(DMProjectFunction(K, 0.0, initu, NULL, INSERT_ALL_VALUES, X));
37: PetscCall(DMViewFromOptions(K, NULL, "-dual_dm_view"));
38: PetscCall(VecViewFromOptions(X, NULL, "-dual_vec_view"));
39: PetscCall(DMDestroy(&K));
40: PetscCall(VecDestroy(&X));
41: PetscCall(PetscFinalize());
42: return 0;
43: }
45: /*TEST
47: testset:
48: filter: grep -v DM_
49: diff_args: -j
50: args: -petscspace_degree 4 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1 -dual_dm_view -dual_vec_view
51: test:
52: requires: !single
53: suffix: 0
54: test:
55: requires: !single
56: suffix: 3d
57: args: -dim 3
59: TEST*/