Actual source code: ex14.c


  2: static char help[] = "Tests saving DMDA vectors to files.\n\n";

  4: /*
  5:     ex13.c reads in the DMDA and vector written by this program.
  6: */

  8: #include <petscdm.h>
  9: #include <petscdmda.h>

 11: int main(int argc,char **argv)
 12: {
 13:   PetscMPIInt    rank;
 14:   PetscInt       M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE, dof = 1;
 15:   DM             da;
 16:   Vec            local,global,natural;
 17:   PetscScalar    value;
 18:   PetscViewer    bviewer;

 20:   PetscInitialize(&argc,&argv,(char*)0,help);
 21:   PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);
 22:   PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);
 23:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 24:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 25:   PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL);

 27:   /* Create distributed array and get vectors */
 28:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,NULL,NULL,&da);
 29:   DMSetFromOptions(da);
 30:   DMSetUp(da);
 31:   DMCreateGlobalVector(da,&global);
 32:   DMCreateLocalVector(da,&local);

 34:   value = -3.0;
 35:   VecSet(global,value);
 36:   DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 37:   DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 39:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 40:   value = rank+1;
 41:   VecScale(local,value);
 42:   DMLocalToGlobalBegin(da,local,ADD_VALUES,global);
 43:   DMLocalToGlobalEnd(da,local,ADD_VALUES,global);

 45:   DMDACreateNaturalVector(da,&natural);
 46:   DMDAGlobalToNaturalBegin(da,global,INSERT_VALUES,natural);
 47:   DMDAGlobalToNaturalEnd(da,global,INSERT_VALUES,natural);

 49:   DMDASetFieldName(da,0,"First field");
 50:   /*  VecView(global,PETSC_VIEWER_DRAW_WORLD); */

 52:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"daoutput",FILE_MODE_WRITE,&bviewer);
 53:   DMView(da,bviewer);
 54:   VecView(global,bviewer);
 55:   PetscViewerDestroy(&bviewer);

 57:   /* Free memory */
 58:   VecDestroy(&local);
 59:   VecDestroy(&global);
 60:   VecDestroy(&natural);
 61:   DMDestroy(&da);
 62:   PetscFinalize();
 63:   return 0;
 64: }