Actual source code: ex32.c
petsc-3.8.4 2018-03-24
1: static char help[] = "Tests DMDA ghost coordinates\n\n";
3: #include <petscdm.h>
4: #include <petscdmda.h>
6: static PetscErrorCode CompareGhostedCoords(Vec gc1,Vec gc2)
7: {
9: PetscReal nrm,gnrm;
10: Vec tmp;
13: VecDuplicate(gc1,&tmp);
14: VecWAXPY(tmp,-1.0,gc1,gc2);
15: VecNorm(tmp,NORM_INFINITY,&nrm);
16: MPI_Allreduce(&nrm,&gnrm,1,MPIU_REAL,MPIU_MAX,PETSC_COMM_WORLD);
17: PetscPrintf(PETSC_COMM_WORLD,"norm of difference of ghosted coordinates %8.2e\n",gnrm);
18: VecDestroy(&tmp);
19: return(0);
20: }
22: static PetscErrorCode TestQ2Q1DA(void)
23: {
24: DM Q2_da,Q1_da,cda;
25: PetscInt mx,my,mz;
26: Vec coords,gcoords,gcoords2;
29: mx = 7;
30: my = 11;
31: mz = 13;
32: 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);
33: DMSetFromOptions(Q2_da);
34: DMSetUp(Q2_da);
35: DMDASetUniformCoordinates(Q2_da,-1.0,1.0,-2.0,2.0,-3.0,3.0);
36: DMGetCoordinates(Q2_da,&coords);
37: 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);
38: DMSetFromOptions(Q1_da);
39: DMSetUp(Q1_da);
40: DMSetCoordinates(Q1_da,coords);
42: /* Get ghost coordinates one way */
43: DMGetCoordinatesLocal(Q1_da,&gcoords);
45: /* And another */
46: DMGetCoordinates(Q1_da,&coords);
47: DMGetCoordinateDM(Q1_da,&cda);
48: DMGetLocalVector(cda,&gcoords2);
49: DMGlobalToLocalBegin(cda,coords,INSERT_VALUES,gcoords2);
50: DMGlobalToLocalEnd(cda,coords,INSERT_VALUES,gcoords2);
52: CompareGhostedCoords(gcoords,gcoords2);
53: DMRestoreLocalVector(cda,&gcoords2);
55: VecScale(coords,10.0);
56: VecScale(gcoords,10.0);
57: DMGetCoordinatesLocal(Q1_da,&gcoords2);
58: CompareGhostedCoords(gcoords,gcoords2);
60: DMDestroy(&Q2_da);
61: DMDestroy(&Q1_da);
62: return(0);
63: }
65: int main(int argc,char **argv)
66: {
69: PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
70: TestQ2Q1DA();
71: PetscFinalize();
72: return ierr;
73: }