xref: /petsc/src/dm/tests/ex32.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(gc1,&tmp));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecWAXPY(tmp,-1.0,gc1,gc2));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(tmp,NORM_INFINITY,&nrm));
155f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(&nrm,&gnrm,1,MPIU_REAL,MPIU_MAX,PETSC_COMM_WORLD));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"norm of difference of ghosted coordinates %8.2e\n",gnrm));
175f80ce2aSJacob 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;
315f80ce2aSJacob 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));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(Q2_da));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(Q2_da));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(Q2_da,-1.0,1.0,-2.0,2.0,-3.0,3.0));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(Q2_da,&coords));
365f80ce2aSJacob 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));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(Q1_da));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(Q1_da));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetCoordinates(Q1_da,coords));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Get ghost coordinates one way */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(Q1_da,&gcoords));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* And another */
455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(Q1_da,&coords));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(Q1_da,&cda));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(cda,&gcoords2));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(cda,coords,INSERT_VALUES,gcoords2));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(cda,coords,INSERT_VALUES,gcoords2));
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(CompareGhostedCoords(gcoords,gcoords2));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(cda,&gcoords2));
53c4762a1bSJed Brown 
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(coords,10.0));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(gcoords,10.0));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(Q1_da,&gcoords2));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(CompareGhostedCoords(gcoords,gcoords2));
58c4762a1bSJed Brown 
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&Q2_da));
605f80ce2aSJacob 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 
67*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(TestQ2Q1DA());
69*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
70*b122ec5aSJacob Faibussowitsch   return 0;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*TEST
74c4762a1bSJed Brown 
75c4762a1bSJed Brown    test:
76c4762a1bSJed Brown       nsize: 2
77c4762a1bSJed Brown 
78c4762a1bSJed Brown TEST*/
79