xref: /petsc/src/dm/tests/ex3.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
12c4762a1bSJed Brown   PetscInt       M = 60,time_steps = 100, localsize,j,i,mybase,myend,width,xbase,*localnodes = NULL;
13c4762a1bSJed Brown   DM             da;
14c4762a1bSJed Brown   PetscViewer    viewer,viewer_private;
15c4762a1bSJed Brown   PetscDraw      draw;
16c4762a1bSJed Brown   Vec            local,global;
17c4762a1bSJed Brown   PetscScalar    *localptr,*globalptr;
18c4762a1bSJed Brown   PetscReal      a,h,k;
19c4762a1bSJed Brown   PetscBool      flg = PETSC_FALSE;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
22*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
23*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
24c4762a1bSJed Brown 
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL));
27c4762a1bSJed Brown   /*
28c4762a1bSJed Brown       Test putting two nodes on each processor, exact last processor gets the rest
29c4762a1bSJed Brown   */
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-distribute",&flg,NULL));
31c4762a1bSJed Brown   if (flg) {
32*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(size,&localnodes));
33c4762a1bSJed Brown     for (i=0; i<size-1; i++) localnodes[i] = 2;
34c4762a1bSJed Brown     localnodes[size-1] = M - 2*(size-1);
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Set up the array */
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,M,1,1,localnodes,&da));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(localnodes));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&global));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateLocalVector(da,&local));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Set up display to show combined wave graph */
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"Entire Solution",20,480,800,200,&viewer));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawGetDraw(viewer,0,&draw));
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetDoubleBuffer(draw));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* determine starting point of each processor */
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(global,&mybase,&myend));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* set up display to show my portion of the wave */
54c4762a1bSJed Brown   xbase = (int)((mybase)*((800.0 - 4.0*size)/M) + 4.0*rank);
55c4762a1bSJed Brown   width = (int)((myend-mybase)*800./M);
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Local Portion of Solution",xbase,200,width,200,&viewer_private));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDrawGetDraw(viewer_private,0,&draw));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDrawSetDoubleBuffer(draw));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Initialize the array */
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(local,&localsize));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(global,&globalptr));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   for (i=1; i<localsize-1; i++) {
65c4762a1bSJed Brown     j           = (i-1)+mybase;
66c4762a1bSJed Brown     globalptr[i-1] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 2;
67c4762a1bSJed Brown   }
68c4762a1bSJed Brown 
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(global,&globalptr));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* Assign Parameters */
72c4762a1bSJed Brown   a= 1.0;
73c4762a1bSJed Brown   h= 1.0/M;
74c4762a1bSJed Brown   k= h;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   for (j=0; j<time_steps; j++) {
77c4762a1bSJed Brown 
78c4762a1bSJed Brown     /* Global to Local */
79*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
80*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown     /*Extract local array */
83*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(local,&localptr));
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(global,&globalptr));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown     /* Update Locally - Make array of new values */
87c4762a1bSJed Brown     /* Note: I don't do anything for the first and last entry */
88c4762a1bSJed Brown     for (i=1; i< localsize-1; i++) {
89c4762a1bSJed Brown       globalptr[i-1] = .5*(localptr[i+1]+localptr[i-1]) - (k / (2.0*a*h)) * (localptr[i+1] - localptr[i-1]);
90c4762a1bSJed Brown     }
91*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(global,&globalptr));
92*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(local,&localptr));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown     /* View my part of Wave */
95*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(global,viewer_private));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown     /* View global Wave */
98*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(global,viewer));
99c4762a1bSJed Brown   }
100c4762a1bSJed Brown 
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer_private));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&local));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&global));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   ierr = PetscFinalize();
108c4762a1bSJed Brown   return ierr;
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111c4762a1bSJed Brown /*TEST
112c4762a1bSJed Brown 
113c4762a1bSJed Brown     test:
114c4762a1bSJed Brown       nsize: 3
115c4762a1bSJed Brown       args: -time 50 -nox
116c4762a1bSJed Brown       requires: x
117c4762a1bSJed Brown 
118c4762a1bSJed Brown TEST*/
119