xref: /petsc/src/dm/tests/ex32.c (revision d8e47b638cf8f604a99e9678e1df24f82d959cd7)
1 static char help[] = "Tests DMDA ghost coordinates\n\n";
2 
3 #include <petscdm.h>
4 #include <petscdmda.h>
5 
CompareGhostedCoords(Vec gc1,Vec gc2)6 static PetscErrorCode CompareGhostedCoords(Vec gc1, Vec gc2)
7 {
8   PetscReal nrm, gnrm;
9   Vec       tmp;
10 
11   PetscFunctionBeginUser;
12   PetscCall(VecDuplicate(gc1, &tmp));
13   PetscCall(VecWAXPY(tmp, -1.0, gc1, gc2));
14   PetscCall(VecNorm(tmp, NORM_INFINITY, &nrm));
15   PetscCallMPI(MPIU_Allreduce(&nrm, &gnrm, 1, MPIU_REAL, MPIU_MAX, PETSC_COMM_WORLD));
16   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "norm of difference of ghosted coordinates %8.2e\n", (double)gnrm));
17   PetscCall(VecDestroy(&tmp));
18   PetscFunctionReturn(PETSC_SUCCESS);
19 }
20 
TestQ2Q1DA(void)21 static PetscErrorCode TestQ2Q1DA(void)
22 {
23   DM       Q2_da, Q1_da, cda;
24   PetscInt mx, my, mz;
25   Vec      coords, gcoords, gcoords2;
26 
27   PetscFunctionBeginUser;
28   mx = 7;
29   my = 11;
30   mz = 13;
31   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, mx, my, mz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 3, 2, 0, 0, 0, &Q2_da));
32   PetscCall(DMSetFromOptions(Q2_da));
33   PetscCall(DMSetUp(Q2_da));
34   PetscCall(DMDASetUniformCoordinates(Q2_da, -1.0, 1.0, -2.0, 2.0, -3.0, 3.0));
35   PetscCall(DMGetCoordinates(Q2_da, &coords));
36   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, mx, my, mz, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 3, 1, 0, 0, 0, &Q1_da));
37   PetscCall(DMSetFromOptions(Q1_da));
38   PetscCall(DMSetUp(Q1_da));
39   PetscCall(DMSetCoordinates(Q1_da, coords));
40 
41   /* Get ghost coordinates one way */
42   PetscCall(DMGetCoordinatesLocal(Q1_da, &gcoords));
43 
44   /* And another */
45   PetscCall(DMGetCoordinates(Q1_da, &coords));
46   PetscCall(DMGetCoordinateDM(Q1_da, &cda));
47   PetscCall(DMGetLocalVector(cda, &gcoords2));
48   PetscCall(DMGlobalToLocalBegin(cda, coords, INSERT_VALUES, gcoords2));
49   PetscCall(DMGlobalToLocalEnd(cda, coords, INSERT_VALUES, gcoords2));
50 
51   PetscCall(CompareGhostedCoords(gcoords, gcoords2));
52   PetscCall(DMRestoreLocalVector(cda, &gcoords2));
53 
54   PetscCall(VecScale(coords, 10.0));
55   PetscCall(VecScale(gcoords, 10.0));
56   PetscCall(DMGetCoordinatesLocal(Q1_da, &gcoords2));
57   PetscCall(CompareGhostedCoords(gcoords, gcoords2));
58 
59   PetscCall(DMDestroy(&Q2_da));
60   PetscCall(DMDestroy(&Q1_da));
61   PetscFunctionReturn(PETSC_SUCCESS);
62 }
63 
main(int argc,char ** argv)64 int main(int argc, char **argv)
65 {
66   PetscFunctionBeginUser;
67   PetscCall(PetscInitialize(&argc, &argv, 0, help));
68   PetscCall(TestQ2Q1DA());
69   PetscCall(PetscFinalize());
70   return 0;
71 }
72 
73 /*TEST
74 
75    test:
76       nsize: 2
77 
78 TEST*/
79