Actual source code: ex39.c

petsc-3.4.5 2014-06-29
  2: static char help[] = "Tests mirror boundary conditions in 1-d.\n\n";

  4: #include <petscdmda.h>

  8: int main(int argc,char **argv)
  9: {
 11:   PetscInt       M = 6,stencil_width = 1, dof = 1,m,xstart,i,j;
 12:   DM             da;
 13:   Vec            global,local;
 14:   PetscScalar    **vglobal;

 16:   PetscInitialize(&argc,&argv,(char*)0,help);

 19:   PetscOptionsGetInt(0,"-stencil_width",&stencil_width,0);
 20:   PetscOptionsGetInt(0,"-dof",&dof,0);

 22:   DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_MIRROR,M,dof,stencil_width,NULL,&da);
 23:   DMDAGetCorners(da,&xstart,0,0,&m,0,0);

 25:   DMCreateGlobalVector(da,&global);
 26:   DMDAVecGetArrayDOF(da,global,&vglobal);
 27:   for (i=xstart; i<xstart+m; i++) {
 28:     for (j=0; j<dof; j++) {
 29:       vglobal[i][j] = 100*(i+1) + j;
 30:     }
 31:   }
 32:   DMDAVecRestoreArrayDOF(da,global,&vglobal);

 34:   DMCreateLocalVector(da,&local);
 35:   DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 36:   DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 38:   PetscSequentialPhaseBegin(PETSC_COMM_WORLD,1);
 39:   VecView(local,PETSC_VIEWER_STDOUT_SELF);
 40:   PetscSequentialPhaseEnd(PETSC_COMM_WORLD,1);
 41:   VecView(global,PETSC_VIEWER_STDOUT_WORLD);

 43:   DMDestroy(&da);
 44:   VecDestroy(&local);
 45:   VecDestroy(&global);

 47:   PetscFinalize();
 48:   return 0;
 49: }