xref: /petsc/src/dm/tests/ex3.c (revision 9566063d113dddea24716c546802770db7481bc0)
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*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
21*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
22*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
23c4762a1bSJed Brown 
24*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
25*9566063dSJacob Faibussowitsch   PetscCall(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   */
29*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-distribute",&flg,NULL));
30c4762a1bSJed Brown   if (flg) {
31*9566063dSJacob Faibussowitsch     PetscCall(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 */
37*9566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,M,1,1,localnodes,&da));
38*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
39*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
40*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(localnodes));
41*9566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&global));
42*9566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da,&local));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* Set up display to show combined wave graph */
45*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"Entire Solution",20,480,800,200,&viewer));
46*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
47*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* determine starting point of each processor */
50*9566063dSJacob Faibussowitsch   PetscCall(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);
55*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Local Portion of Solution",xbase,200,width,200,&viewer_private));
56*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(viewer_private,0,&draw));
57*9566063dSJacob Faibussowitsch   PetscCall(PetscDrawSetDoubleBuffer(draw));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* Initialize the array */
60*9566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(local,&localsize));
61*9566063dSJacob Faibussowitsch   PetscCall(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 
68*9566063dSJacob Faibussowitsch   PetscCall(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 */
78*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
79*9566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     /*Extract local array */
82*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(local,&localptr));
83*9566063dSJacob Faibussowitsch     PetscCall(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     }
90*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(global,&globalptr));
91*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(local,&localptr));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     /* View my part of Wave */
94*9566063dSJacob Faibussowitsch     PetscCall(VecView(global,viewer_private));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown     /* View global Wave */
97*9566063dSJacob Faibussowitsch     PetscCall(VecView(global,viewer));
98c4762a1bSJed Brown   }
99c4762a1bSJed Brown 
100*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
101*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
102*9566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer_private));
103*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
104*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
105c4762a1bSJed Brown 
106*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
107b122ec5aSJacob 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