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