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