Actual source code: ex2.c

  2: static char help[] = "Tests various 1-dimensional DA routines.\n\n";

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

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

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

 23:   /* Create viewers */
 24:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",280,480,600,200,&viewer);
 25:   PetscViewerDrawGetDraw(viewer,0,&draw);
 26:   PetscDrawSetDoubleBuffer(draw);

 28:   /* Readoptions */
 29:   PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
 30:   PetscOptionsGetInt(PETSC_NULL,"-w",&w,PETSC_NULL);
 31:   PetscOptionsGetInt(PETSC_NULL,"-s",&s,PETSC_NULL);

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

 39:   /* Set global vector; send ghost points to local vectors */
 40:   value = 1;
 41:   VecSet(global,value);
 42:   DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 43:   DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

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

 51:   VecView(global,viewer);
 52:   PetscPrintf(PETSC_COMM_WORLD,"\nGlobal Vector:\n");
 53:   VecView(global,PETSC_VIEWER_STDOUT_WORLD);
 54:   PetscPrintf(PETSC_COMM_WORLD,"\n");

 56:   /* Send ghost points to local vectors */
 57:   DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 58:   DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 60:   PetscOptionsHasName(PETSC_NULL,"-local_print",&flg);
 61:   if (flg) {
 62:     PetscViewer sviewer;
 63:     PetscSynchronizedPrintf(PETSC_COMM_WORLD,"\nLocal Vector: processor %d\n",rank);
 64:     PetscViewerGetSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
 65:     VecView(local,sviewer);
 66:     PetscViewerRestoreSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
 67:     PetscSynchronizedFlush(PETSC_COMM_WORLD);
 68:   }

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