Actual source code: ex32.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: static char help[] = "Tests DMDA ghost coordinates\n\n";

  3: #include <petscdm.h>
  4: #include <petscdmda.h>

  8: static PetscErrorCode CompareGhostedCoords(Vec gc1,Vec gc2)
  9: {
 11:   PetscReal      nrm,gnrm;
 12:   Vec            tmp;

 15:   VecDuplicate(gc1,&tmp);
 16:   VecWAXPY(tmp,-1.0,gc1,gc2);
 17:   VecNorm(tmp,NORM_INFINITY,&nrm);
 18:   MPI_Allreduce(&nrm,&gnrm,1,MPIU_REAL,MPIU_MAX,PETSC_COMM_WORLD);
 19:   PetscPrintf(PETSC_COMM_WORLD,"norm of difference of ghosted coordinates %8.2e\n",gnrm);
 20:   VecDestroy(&tmp);
 21:   return(0);
 22: }

 26: static PetscErrorCode TestQ2Q1DA(void)
 27: {
 28:   DM             Q2_da,Q1_da,cda;
 29:   PetscInt       mx,my,mz;
 30:   Vec            coords,gcoords,gcoords2;

 33:   mx   = 7;
 34:   my   = 11;
 35:   mz   = 13;
 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,2,0,0,0,&Q2_da);
 37:   DMDASetUniformCoordinates(Q2_da,-1.0,1.0,-2.0,2.0,-3.0,3.0);
 38:   DMGetCoordinates(Q2_da,&coords);
 39:   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);
 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: }

 67: int main(int argc,char **argv)
 68: {

 71:   PetscInitialize(&argc,&argv,0,help);
 72:   TestQ2Q1DA();
 73:   PetscFinalize();
 74:   return 0;
 75: }