xref: /petsc/src/dm/tests/ex12.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
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   DM             da;
17   Vec            global,local;
18   PetscScalar    *globalptr,*localptr;
19   PetscReal      h,k;
20 
21   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
22   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
23   PetscCall(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL));
24 
25   /* Set up the array */
26   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,M,w,s,NULL,&da));
27   PetscCall(DMSetFromOptions(da));
28   PetscCall(DMSetUp(da));
29   PetscCall(DMCreateGlobalVector(da,&global));
30   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
31   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
32 
33   /* Make copy of local array for doing updates */
34   PetscCall(DMCreateLocalVector(da,&local));
35 
36   /* determine starting point of each processor */
37   PetscCall(VecGetOwnershipRange(global,&mybase,&myend));
38 
39   /* Initialize the Array */
40   PetscCall(VecGetLocalSize (global,&globalsize));
41   PetscCall(VecGetArray (global,&globalptr));
42 
43   for (i=0; i<globalsize; i++) {
44     j = i + mybase;
45 
46     globalptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
47   }
48 
49   PetscCall(VecRestoreArray(global,&localptr));
50 
51   /* Assign Parameters */
52   h= 1.0/M;
53   k= h*h/2.2;
54   PetscCall(VecGetLocalSize(local,&localsize));
55 
56   for (j=0; j<time_steps; j++) {
57 
58     /* Global to Local */
59     PetscCall(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
60     PetscCall(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
61 
62     /*Extract local array */
63     PetscCall(VecGetArray(local,&localptr));
64     PetscCall(VecGetArray (global,&globalptr));
65 
66     /* Update Locally - Make array of new values */
67     /* Note: I don't do anything for the first and last entry */
68     for (i=1; i< localsize-1; i++) {
69       globalptr[i-1] = localptr[i] + (k/(h*h)) * (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
70     }
71 
72     PetscCall(VecRestoreArray (global,&globalptr));
73     PetscCall(VecRestoreArray(local,&localptr));
74 
75     /* View Wave */
76     /* Set Up Display to Show Heat Graph */
77 #if defined(PETSC_USE_SOCKET_VIEWER)
78     PetscCall(VecView(global,PETSC_VIEWER_SOCKET_WORLD));
79 #endif
80   }
81 
82   PetscCall(VecDestroy(&local));
83   PetscCall(VecDestroy(&global));
84   PetscCall(DMDestroy(&da));
85   PetscCall(PetscFinalize());
86   return 0;
87 }
88