xref: /petsc/src/dm/tests/ex12.c (revision 8ebe3e4e9e00d86ece2e9fcd0cc84910b0ad437c)
1 
2 /*
3    Simple example to show how PETSc programs can be run from MATLAB.
4   See ex12.m.
5 */
6 
7 static char help[] = "Solves the one dimensional heat equation.\n\n";
8 
9 #include <petscdm.h>
10 #include <petscdmda.h>
11 
12 int main(int argc,char **argv)
13 {
14   PetscMPIInt    rank,size;
15   PetscInt       M = 14,time_steps = 20,w=1,s=1,localsize,j,i,mybase,myend,globalsize;
16   PetscErrorCode ierr;
17   DM             da;
18   Vec            global,local;
19   PetscScalar    *globalptr,*localptr;
20   PetscReal      h,k;
21 
22   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
23   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
24   ierr = PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);CHKERRQ(ierr);
25 
26   /* Set up the array */
27   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,M,w,s,NULL,&da);CHKERRQ(ierr);
28   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
29   ierr = DMSetUp(da);CHKERRQ(ierr);
30   ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr);
31   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
32   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
33 
34   /* Make copy of local array for doing updates */
35   ierr = DMCreateLocalVector(da,&local);CHKERRQ(ierr);
36 
37   /* determine starting point of each processor */
38   ierr = VecGetOwnershipRange(global,&mybase,&myend);CHKERRQ(ierr);
39 
40   /* Initialize the Array */
41   ierr = VecGetLocalSize (global,&globalsize);CHKERRQ(ierr);
42   ierr = VecGetArray (global,&globalptr);CHKERRQ(ierr);
43 
44   for (i=0; i<globalsize; i++) {
45     j = i + mybase;
46 
47     globalptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
48   }
49 
50   ierr = VecRestoreArray(global,&localptr);CHKERRQ(ierr);
51 
52   /* Assign Parameters */
53   h= 1.0/M;
54   k= h*h/2.2;
55   ierr = VecGetLocalSize(local,&localsize);CHKERRQ(ierr);
56 
57   for (j=0; j<time_steps; j++) {
58 
59     /* Global to Local */
60     ierr = DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
61     ierr = DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
62 
63     /*Extract local array */
64     ierr = VecGetArray(local,&localptr);CHKERRQ(ierr);
65     ierr = VecGetArray (global,&globalptr);CHKERRQ(ierr);
66 
67     /* Update Locally - Make array of new values */
68     /* Note: I don't do anything for the first and last entry */
69     for (i=1; i< localsize-1; i++) {
70       globalptr[i-1] = localptr[i] + (k/(h*h)) * (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
71     }
72 
73     ierr = VecRestoreArray (global,&globalptr);CHKERRQ(ierr);
74     ierr = VecRestoreArray(local,&localptr);CHKERRQ(ierr);
75 
76     /* View Wave */
77     /* Set Up Display to Show Heat Graph */
78 #if defined(PETSC_USE_SOCKET_VIEWER)
79     ierr = VecView(global,PETSC_VIEWER_SOCKET_WORLD);CHKERRQ(ierr);
80 #endif
81   }
82 
83   ierr = VecDestroy(&local);CHKERRQ(ierr);
84   ierr = VecDestroy(&global);CHKERRQ(ierr);
85   ierr = DMDestroy(&da);CHKERRQ(ierr);
86   ierr = PetscFinalize();
87   return ierr;
88 }
89 
90