xref: /petsc/src/dm/tests/ex32.c (revision 9566063d113dddea24716c546802770db7481bc0)
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*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(gc1,&tmp));
13*9566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(tmp,-1.0,gc1,gc2));
14*9566063dSJacob Faibussowitsch   PetscCall(VecNorm(tmp,NORM_INFINITY,&nrm));
15*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&nrm,&gnrm,1,MPIU_REAL,MPIU_MAX,PETSC_COMM_WORLD));
16*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"norm of difference of ghosted coordinates %8.2e\n",gnrm));
17*9566063dSJacob Faibussowitsch   PetscCall(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*9566063dSJacob Faibussowitsch   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*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(Q2_da));
33*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(Q2_da));
34*9566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(Q2_da,-1.0,1.0,-2.0,2.0,-3.0,3.0));
35*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(Q2_da,&coords));
36*9566063dSJacob Faibussowitsch   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*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(Q1_da));
38*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(Q1_da));
39*9566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinates(Q1_da,coords));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Get ghost coordinates one way */
42*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(Q1_da,&gcoords));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* And another */
45*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(Q1_da,&coords));
46*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(Q1_da,&cda));
47*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(cda,&gcoords2));
48*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cda,coords,INSERT_VALUES,gcoords2));
49*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cda,coords,INSERT_VALUES,gcoords2));
50c4762a1bSJed Brown 
51*9566063dSJacob Faibussowitsch   PetscCall(CompareGhostedCoords(gcoords,gcoords2));
52*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(cda,&gcoords2));
53c4762a1bSJed Brown 
54*9566063dSJacob Faibussowitsch   PetscCall(VecScale(coords,10.0));
55*9566063dSJacob Faibussowitsch   PetscCall(VecScale(gcoords,10.0));
56*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(Q1_da,&gcoords2));
57*9566063dSJacob Faibussowitsch   PetscCall(CompareGhostedCoords(gcoords,gcoords2));
58c4762a1bSJed Brown 
59*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&Q2_da));
60*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&Q1_da));
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown int main(int argc,char **argv)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown 
67*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,0,help));
68*9566063dSJacob Faibussowitsch   PetscCall(TestQ2Q1DA());
69*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
70b122ec5aSJacob 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