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