1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves the 1-dimensional wave equation.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscdraw.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown int main(int argc,char **argv) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown PetscMPIInt rank,size; 11c4762a1bSJed Brown PetscInt M = 60,time_steps = 100, localsize,j,i,mybase,myend,width,xbase,*localnodes = NULL; 12c4762a1bSJed Brown DM da; 13c4762a1bSJed Brown PetscViewer viewer,viewer_private; 14c4762a1bSJed Brown PetscDraw draw; 15c4762a1bSJed Brown Vec local,global; 16c4762a1bSJed Brown PetscScalar *localptr,*globalptr; 17c4762a1bSJed Brown PetscReal a,h,k; 18c4762a1bSJed Brown PetscBool flg = PETSC_FALSE; 19c4762a1bSJed Brown 20*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 225f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 23c4762a1bSJed Brown 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL)); 26c4762a1bSJed Brown /* 27c4762a1bSJed Brown Test putting two nodes on each processor, exact last processor gets the rest 28c4762a1bSJed Brown */ 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-distribute",&flg,NULL)); 30c4762a1bSJed Brown if (flg) { 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size,&localnodes)); 32c4762a1bSJed Brown for (i=0; i<size-1; i++) localnodes[i] = 2; 33c4762a1bSJed Brown localnodes[size-1] = M - 2*(size-1); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* Set up the array */ 375f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,M,1,1,localnodes,&da)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(localnodes)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&global)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(da,&local)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Set up display to show combined wave graph */ 455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"Entire Solution",20,480,800,200,&viewer)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(viewer,0,&draw)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* determine starting point of each processor */ 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(global,&mybase,&myend)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* set up display to show my portion of the wave */ 53c4762a1bSJed Brown xbase = (int)((mybase)*((800.0 - 4.0*size)/M) + 4.0*rank); 54c4762a1bSJed Brown width = (int)((myend-mybase)*800./M); 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Local Portion of Solution",xbase,200,width,200,&viewer_private)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDraw(viewer_private,0,&draw)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetDoubleBuffer(draw)); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Initialize the array */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(local,&localsize)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(global,&globalptr)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown for (i=1; i<localsize-1; i++) { 64c4762a1bSJed Brown j = (i-1)+mybase; 65c4762a1bSJed Brown globalptr[i-1] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 2; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(global,&globalptr)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Assign Parameters */ 71c4762a1bSJed Brown a= 1.0; 72c4762a1bSJed Brown h= 1.0/M; 73c4762a1bSJed Brown k= h; 74c4762a1bSJed Brown 75c4762a1bSJed Brown for (j=0; j<time_steps; j++) { 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Global to Local */ 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /*Extract local array */ 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(local,&localptr)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(global,&globalptr)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Update Locally - Make array of new values */ 86c4762a1bSJed Brown /* Note: I don't do anything for the first and last entry */ 87c4762a1bSJed Brown for (i=1; i< localsize-1; i++) { 88c4762a1bSJed Brown globalptr[i-1] = .5*(localptr[i+1]+localptr[i-1]) - (k / (2.0*a*h)) * (localptr[i+1] - localptr[i-1]); 89c4762a1bSJed Brown } 905f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(global,&globalptr)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(local,&localptr)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* View my part of Wave */ 945f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(global,viewer_private)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* View global Wave */ 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(global,viewer)); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer_private)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&local)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&global)); 105c4762a1bSJed Brown 106*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 107*b122ec5aSJacob Faibussowitsch return 0; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown /*TEST 111c4762a1bSJed Brown 112c4762a1bSJed Brown test: 113c4762a1bSJed Brown nsize: 3 114c4762a1bSJed Brown args: -time 50 -nox 115c4762a1bSJed Brown requires: x 116c4762a1bSJed Brown 117c4762a1bSJed Brown TEST*/ 118