xref: /petsc/src/dm/impls/stag/tests/ex53.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1 static char help[] = "Test local-to-local for DMStag.\n\n";
2 
3 #include <petscdmstag.h>
4 
main(int argc,char ** argv)5 int main(int argc, char **argv)
6 {
7   DM          dm;
8   PetscInt    dim, start, end, i;
9   PetscBool   flg;
10   Vec         g, l1, l2;
11   PetscMPIInt rank;
12   PetscScalar value;
13   PetscReal   norm, work;
14 
15   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
17 
18   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, &flg));
19   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Supply -dim option");
20 
21   if (dim == 1) PetscCall(DMStagCreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 64, 1, 1, DMSTAG_STENCIL_BOX, 1, NULL, &dm));
22   else if (dim == 2) PetscCall(DMStagCreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, 8, 8, PETSC_DECIDE, PETSC_DECIDE, 1, 1, 1, DMSTAG_STENCIL_BOX, 1, NULL, NULL, &dm));
23   else if (dim == 3) PetscCall(DMStagCreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, 4, 4, 4, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, 1, 1, DMSTAG_STENCIL_BOX, 1, NULL, NULL, NULL, &dm));
24   else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "dim must be 1, 2, or 3");
25   PetscCall(DMSetFromOptions(dm));
26   PetscCall(DMSetUp(dm));
27 
28   PetscCall(DMCreateGlobalVector(dm, &g));
29   PetscCall(DMCreateLocalVector(dm, &l1));
30   PetscCall(VecDuplicate(l1, &l2));
31 
32   PetscCall(VecSet(l1, 0.0));
33   PetscCall(VecSet(l2, 0.0));
34 
35   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
36   PetscCall(VecGetOwnershipRange(g, &start, &end));
37   for (i = start; i < end; ++i) {
38     value = rank + i;
39     PetscCall(VecSetValues(g, 1, &i, &value, INSERT_VALUES));
40   }
41   PetscCall(VecAssemblyBegin(g));
42   PetscCall(VecAssemblyEnd(g));
43 
44   PetscCall(DMGlobalToLocalBegin(dm, g, INSERT_VALUES, l1));
45   PetscCall(DMGlobalToLocalEnd(dm, g, INSERT_VALUES, l1));
46 
47   PetscCall(DMLocalToLocalBegin(dm, l1, INSERT_VALUES, l2));
48   PetscCall(DMLocalToLocalEnd(dm, l1, INSERT_VALUES, l2));
49 
50   /* l1 and l2 must be same. */
51   PetscCall(VecAXPY(l2, -1.0, l1));
52   PetscCall(VecNorm(l2, NORM_MAX, &work));
53   PetscCallMPI(MPIU_Allreduce(&work, &norm, 1, MPIU_REAL, MPIU_MAX, PETSC_COMM_WORLD));
54   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "norm = %g\n", (double)norm));
55 
56   PetscCall(VecDestroy(&g));
57   PetscCall(VecDestroy(&l1));
58   PetscCall(VecDestroy(&l2));
59   PetscCall(DMDestroy(&dm));
60   PetscCall(PetscFinalize());
61   return 0;
62 }
63 
64 /*TEST
65 
66    test:
67       suffix: 1
68       nsize: 2
69       args: -dim 1
70 
71    test:
72       suffix: 2
73       nsize: 3
74       args: -dim 1 -stag_boundary_type_x none
75       output_file: output/ex53_1.out
76 
77    test:
78       suffix: 3
79       nsize: 4
80       args: -dim 1 -stag_boundary_type_x periodic
81       output_file: output/ex53_1.out
82 
83    test:
84       suffix: 4
85       nsize: 4
86       args: -dim 2
87       output_file: output/ex53_1.out
88 
89    test:
90       suffix: 5
91       nsize: 4
92       args: -dim 2 -stag_boundary_type_x none -stag_stencil_type star
93       output_file: output/ex53_1.out
94 
95    test:
96       suffix: 6
97       nsize: 6
98       args: -dim 2 -stag_boundary_type_y periodic -stag_stencil_width 2 -stag_dof_0 0 -stag_dof_1 1 -stag_dof_2 0
99       output_file: output/ex53_1.out
100 
101    test:
102       suffix: 7
103       nsize: 8
104       args: -dim 3
105       output_file: output/ex53_1.out
106 
107    test:
108       suffix: 8
109       nsize: 8
110       args: -dim 3 -stag_boundary_type_x none -stag_boundary_type_y periodic
111       output_file: output/ex53_1.out
112 
113    test:
114       suffix: 9
115       nsize: 12
116       args: -dim 3 -stag_boundary_type_x none -stag_boundary_type_y none -stag_boundary_type_z none -stag_stencil_type star -stag_dof_0 0 -stag_dof_1 0 -stag_dof_2 0 -stag_dof_3 1
117       output_file: output/ex53_1.out
118 
119 TEST*/
120