Actual source code: ex2.c

  1: /*$Id: ex2.c,v 1.43 2001/08/07 03:04:42 balay Exp $*/

  3: static char help[] = "Tests various 1-dimensional DA routines.nn";

 5:  #include petscda.h
 6:  #include petscsys.h

  8: int main(int argc,char **argv)
  9: {
 10:   int          rank,M = 13,ierr,w=1,s=1,wrap=1;
 11:   DA           da;
 12:   PetscViewer  viewer;
 13:   Vec          local,global;
 14:   PetscScalar  value;
 15:   PetscDraw    draw;
 16:   PetscTruth   flg;

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

 20:   /* Create viewers */
 21:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",280,480,600,200,&viewer);
 22:   PetscViewerDrawGetDraw(viewer,0,&draw);
 23:   PetscDrawSetDoubleBuffer(draw);

 25:   /* Readoptions */
 26:   PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
 27:   PetscOptionsGetInt(PETSC_NULL,"-w",&w,PETSC_NULL);
 28:   PetscOptionsGetInt(PETSC_NULL,"-s",&s,PETSC_NULL);

 30:   /* Create distributed array and get vectors */
 31:   DACreate1d(PETSC_COMM_WORLD,(DAPeriodicType)wrap,M,w,s,PETSC_NULL,&da);
 32:   DAView(da,viewer);
 33:   DACreateGlobalVector(da,&global);
 34:   DACreateLocalVector(da,&local);

 36:   /* Set global vector; send ghost points to local vectors */
 37:   value = 1;
 38:   VecSet(&value,global);
 39:   DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 40:   DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 42:   /* Scale local vectors according to processor rank; pass to global vector */
 43:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 44:   value = rank+1;
 45:   VecScale(&value,local);
 46:   DALocalToGlobal(da,local,INSERT_VALUES,global);

 48:   VecView(global,viewer);
 49:   PetscPrintf(PETSC_COMM_WORLD,"nGlobal Vector:n");
 50:   VecView(global,PETSC_VIEWER_STDOUT_WORLD);
 51:   PetscPrintf(PETSC_COMM_WORLD,"n");

 53:   /* Send ghost points to local vectors */
 54:   DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 55:   DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 57:   PetscOptionsHasName(PETSC_NULL,"-local_print",&flg);
 58:   if (flg) {
 59:     PetscViewer sviewer;
 60:     PetscSynchronizedPrintf(PETSC_COMM_WORLD,"nLocal Vector: processor %dn",rank);
 61:     PetscViewerGetSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
 62:     VecView(local,sviewer);
 63:     PetscViewerRestoreSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
 64:     PetscSynchronizedFlush(PETSC_COMM_WORLD);
 65:   }

 67:   /* Free memory */
 68:   PetscViewerDestroy(viewer);
 69:   VecDestroy(global);
 70:   VecDestroy(local);
 71:   DADestroy(da);
 72:   PetscFinalize();
 73:   return 0;
 74: }
 75: