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;

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

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

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

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

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

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

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

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