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