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