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