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