xref: /petsc/src/ts/tutorials/ex15.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d. Modified from ex13.c for illustrating how to solve DAEs. \n";
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    u_t = uxx + uyy
5c4762a1bSJed Brown    0 < x < 1, 0 < y < 1;
6c4762a1bSJed Brown    At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125
7c4762a1bSJed Brown            u(x,y) = 0.0           if r >= .125
8c4762a1bSJed Brown 
9c4762a1bSJed Brown    Boundary conditions:
10c4762a1bSJed Brown    Drichlet BC:
11c4762a1bSJed Brown    At x=0, x=1, y=0, y=1: u = 0.0
12c4762a1bSJed Brown 
13c4762a1bSJed Brown    Neumann BC:
14c4762a1bSJed Brown    At x=0, x=1: du(x,y,t)/dx = 0
15c4762a1bSJed Brown    At y=0, y=1: du(x,y,t)/dy = 0
16c4762a1bSJed Brown 
17c4762a1bSJed Brown    mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
18c4762a1bSJed Brown          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -ts_monitor_draw_solution
19c4762a1bSJed Brown          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9
20c4762a1bSJed Brown 
21c4762a1bSJed Brown */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown #include <petscdm.h>
24c4762a1bSJed Brown #include <petscdmda.h>
25c4762a1bSJed Brown #include <petscts.h>
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /*
28c4762a1bSJed Brown    User-defined data structures and routines
29c4762a1bSJed Brown */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /* AppCtx: used by FormIFunction() and FormIJacobian() */
32c4762a1bSJed Brown typedef struct {
33c4762a1bSJed Brown   DM        da;
34c4762a1bSJed Brown   PetscInt  nstencilpts;         /* number of stencil points: 5 or 9 */
35c4762a1bSJed Brown   PetscReal c;
36c4762a1bSJed Brown   PetscInt  boundary;            /* Type of boundary condition */
37c4762a1bSJed Brown   PetscBool viewJacobian;
38c4762a1bSJed Brown } AppCtx;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown extern PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
41c4762a1bSJed Brown extern PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
42c4762a1bSJed Brown extern PetscErrorCode FormInitialSolution(Vec,void*);
43c4762a1bSJed Brown 
44c4762a1bSJed Brown int main(int argc,char **argv)
45c4762a1bSJed Brown {
46c4762a1bSJed Brown   TS             ts;                   /* nonlinear solver */
47c4762a1bSJed Brown   Vec            u,r;                  /* solution, residual vectors */
48c4762a1bSJed Brown   Mat            J,Jmf = NULL;   /* Jacobian matrices */
49c4762a1bSJed Brown   PetscErrorCode ierr;
50c4762a1bSJed Brown   DM             da;
51c4762a1bSJed Brown   PetscReal      dt;
52c4762a1bSJed Brown   AppCtx         user;              /* user-defined work context */
53c4762a1bSJed Brown   SNES           snes;
54c4762a1bSJed Brown   PetscInt       Jtype; /* Jacobian type
55c4762a1bSJed Brown                             0: user provide Jacobian;
56c4762a1bSJed Brown                             1: slow finite difference;
57c4762a1bSJed Brown                             2: fd with coloring; */
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
60c4762a1bSJed Brown   /* Initialize user application context */
61c4762a1bSJed Brown   user.da           = NULL;
62c4762a1bSJed Brown   user.nstencilpts  = 5;
63c4762a1bSJed Brown   user.c            = -30.0;
64c4762a1bSJed Brown   user.boundary     = 0;  /* 0: Drichlet BC; 1: Neumann BC */
65c4762a1bSJed Brown   user.viewJacobian = PETSC_FALSE;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nstencilpts",&user.nstencilpts,NULL);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian);CHKERRQ(ierr);
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
73c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74c4762a1bSJed Brown   if (user.nstencilpts == 5) {
75c4762a1bSJed Brown     ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
76c4762a1bSJed Brown   } else if (user.nstencilpts == 9) {
77c4762a1bSJed Brown     ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
78*98921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"nstencilpts %d is not supported",user.nstencilpts);
79c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
81c4762a1bSJed Brown   user.da = da;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84c4762a1bSJed Brown      Extract global vectors from DMDA;
85c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&u);CHKERRQ(ierr);
87c4762a1bSJed Brown   ierr = VecDuplicate(u,&r);CHKERRQ(ierr);
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90c4762a1bSJed Brown      Create timestepping solver context
91c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = TSSetIFunction(ts,r,FormIFunction,&user);CHKERRQ(ierr);
97c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown      Set initial conditions
102c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103c4762a1bSJed Brown   ierr = FormInitialSolution(u,&user);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = TSSetSolution(ts,u);CHKERRQ(ierr);
105c4762a1bSJed Brown   dt   = .01;
106c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109c4762a1bSJed Brown    Set Jacobian evaluation routine
110c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111c4762a1bSJed Brown   ierr = DMSetMatType(da,MATAIJ);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr  = DMCreateMatrix(da,&J);CHKERRQ(ierr);
113c4762a1bSJed Brown   Jtype = 0;
114c4762a1bSJed Brown   ierr  = PetscOptionsGetInt(NULL,NULL, "-Jtype",&Jtype,NULL);CHKERRQ(ierr);
115c4762a1bSJed Brown   if (Jtype == 0) { /* use user provided Jacobian evaluation routine */
116*98921bdaSJacob Faibussowitsch     if (user.nstencilpts != 5) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"user Jacobian routine FormIJacobian() does not support nstencilpts=%D",user.nstencilpts);
117c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,J,J,FormIJacobian,&user);CHKERRQ(ierr);
118c4762a1bSJed Brown   } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */
119c4762a1bSJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
120c4762a1bSJed Brown     ierr = MatCreateSNESMF(snes,&Jmf);CHKERRQ(ierr);
121c4762a1bSJed Brown     if (Jtype == 1) { /* slow finite difference J; */
122c4762a1bSJed Brown       ierr = SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefault,NULL);CHKERRQ(ierr);
123c4762a1bSJed Brown     } else if (Jtype == 2) { /* Use coloring to compute  finite difference J efficiently */
124c4762a1bSJed Brown       ierr = SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefaultColor,0);CHKERRQ(ierr);
125c4762a1bSJed Brown     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Jtype is not supported");
126c4762a1bSJed Brown   }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129c4762a1bSJed Brown    Sets various TS parameters from user options
130c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134c4762a1bSJed Brown      Solve nonlinear system
135c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136c4762a1bSJed Brown   ierr = TSSolve(ts,u);CHKERRQ(ierr);
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139c4762a1bSJed Brown      Free work space.
140c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = MatDestroy(&Jmf);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   ierr = PetscFinalize();
149c4762a1bSJed Brown   return ierr;
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /* --------------------------------------------------------------------- */
153c4762a1bSJed Brown /*
154c4762a1bSJed Brown   FormIFunction = Udot - RHSFunction
155c4762a1bSJed Brown */
156c4762a1bSJed Brown PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
157c4762a1bSJed Brown {
158c4762a1bSJed Brown   PetscErrorCode ierr;
159c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ctx;
160c4762a1bSJed Brown   DM             da   = (DM)user->da;
161c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
162c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
163c4762a1bSJed Brown   PetscScalar    u,uxx,uyy,**uarray,**f,**udot;
164c4762a1bSJed Brown   Vec            localU;
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   PetscFunctionBeginUser;
167c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
171c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
172c4762a1bSJed Brown   if (user->nstencilpts == 9 && hx != hy) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"hx must equal hy when nstencilpts = 9 for this example");
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /*
175c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
176c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
177c4762a1bSJed Brown      By placing code between these two statements, computations can be
178c4762a1bSJed Brown      done while messages are in transition.
179c4762a1bSJed Brown   */
180c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
181c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /* Get pointers to vector data */
184c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&uarray);CHKERRQ(ierr);
185c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,Udot,&udot);CHKERRQ(ierr);
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* Get local grid boundaries */
189c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
192c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
193c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
194c4762a1bSJed Brown       /* Boundary conditions */
195c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
196c4762a1bSJed Brown         if (user->boundary == 0) { /* Drichlet BC */
197c4762a1bSJed Brown           f[j][i] = uarray[j][i]; /* F = U */
198c4762a1bSJed Brown         } else {                  /* Neumann BC */
199c4762a1bSJed Brown           if (i == 0 && j == 0) {              /* SW corner */
200c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j+1][i+1];
201c4762a1bSJed Brown           } else if (i == Mx-1 && j == 0) {    /* SE corner */
202c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j+1][i-1];
203c4762a1bSJed Brown           } else if (i == 0 && j == My-1) {    /* NW corner */
204c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j-1][i+1];
205c4762a1bSJed Brown           } else if (i == Mx-1 && j == My-1) { /* NE corner */
206c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j-1][i-1];
207c4762a1bSJed Brown           } else if (i == 0) {                  /* Left */
208c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j][i+1];
209c4762a1bSJed Brown           } else if (i == Mx-1) {               /* Right */
210c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j][i-1];
211c4762a1bSJed Brown           } else if (j == 0) {                 /* Bottom */
212c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j+1][i];
213c4762a1bSJed Brown           } else if (j == My-1) {               /* Top */
214c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j-1][i];
215c4762a1bSJed Brown           }
216c4762a1bSJed Brown         }
217c4762a1bSJed Brown       } else { /* Interior */
218c4762a1bSJed Brown         u = uarray[j][i];
219c4762a1bSJed Brown         /* 5-point stencil */
220c4762a1bSJed Brown         uxx = (-2.0*u + uarray[j][i-1] + uarray[j][i+1]);
221c4762a1bSJed Brown         uyy = (-2.0*u + uarray[j-1][i] + uarray[j+1][i]);
222c4762a1bSJed Brown         if (user->nstencilpts == 9) {
223c4762a1bSJed Brown           /* 9-point stencil: assume hx=hy */
224c4762a1bSJed Brown           uxx = 2.0*uxx/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0;
225c4762a1bSJed Brown           uyy = 2.0*uyy/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0;
226c4762a1bSJed Brown         }
227c4762a1bSJed Brown         f[j][i] = udot[j][i] - (uxx*sx + uyy*sy);
228c4762a1bSJed Brown       }
229c4762a1bSJed Brown     }
230c4762a1bSJed Brown   }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* Restore vectors */
233c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&uarray);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,Udot,&udot);CHKERRQ(ierr);
236c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = PetscLogFlops(11.0*ym*xm);CHKERRQ(ierr);
238c4762a1bSJed Brown   PetscFunctionReturn(0);
239c4762a1bSJed Brown }
240c4762a1bSJed Brown 
241c4762a1bSJed Brown /* --------------------------------------------------------------------- */
242c4762a1bSJed Brown /*
243c4762a1bSJed Brown   FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot
244c4762a1bSJed Brown   This routine is not used with option '-use_coloring'
245c4762a1bSJed Brown */
246c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
247c4762a1bSJed Brown {
248c4762a1bSJed Brown   PetscErrorCode ierr;
249c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym,nc;
250c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
251c4762a1bSJed Brown   DM             da    = (DM)user->da;
252c4762a1bSJed Brown   MatStencil     col[5],row;
253c4762a1bSJed Brown   PetscScalar    vals[5],hx,hy,sx,sy;
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   PetscFunctionBeginUser;
256c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
257c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
260c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
263c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
264c4762a1bSJed Brown       nc    = 0;
265c4762a1bSJed Brown       row.j = j; row.i = i;
266c4762a1bSJed Brown       if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) {
267c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0;
268c4762a1bSJed Brown 
269c4762a1bSJed Brown       } else if (user->boundary > 0 && i == 0) {  /* Left Neumann */
270c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
271c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i+1; vals[nc++] = -1.0;
272c4762a1bSJed Brown       } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
273c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
274c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i-1; vals[nc++] = -1.0;
275c4762a1bSJed Brown       } else if (user->boundary > 0 && j == 0) {  /* Bottom Neumann */
276c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i; vals[nc++] = 1.0;
277c4762a1bSJed Brown         col[nc].j = j+1; col[nc].i = i; vals[nc++] = -1.0;
278c4762a1bSJed Brown       } else if (user->boundary > 0 && j == My-1) { /* Top Neumann */
279c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i;  vals[nc++] = 1.0;
280c4762a1bSJed Brown         col[nc].j = j-1; col[nc].i = i;  vals[nc++] = -1.0;
281c4762a1bSJed Brown       } else {   /* Interior */
282c4762a1bSJed Brown         col[nc].j = j-1; col[nc].i = i;   vals[nc++] = -sy;
283c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i-1; vals[nc++] = -sx;
284c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i;   vals[nc++] = 2.0*(sx + sy) + a;
285c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i+1; vals[nc++] = -sx;
286c4762a1bSJed Brown         col[nc].j = j+1; col[nc].i = i;   vals[nc++] = -sy;
287c4762a1bSJed Brown       }
288c4762a1bSJed Brown       ierr = MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES);CHKERRQ(ierr);
289c4762a1bSJed Brown     }
290c4762a1bSJed Brown   }
291c4762a1bSJed Brown   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
292c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293c4762a1bSJed Brown   if (J != Jpre) {
294c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
295c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
296c4762a1bSJed Brown   }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   if (user->viewJacobian) {
299c4762a1bSJed Brown     ierr = PetscPrintf(PetscObjectComm((PetscObject)Jpre),"Jpre:\n");CHKERRQ(ierr);
300c4762a1bSJed Brown     ierr = MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
301c4762a1bSJed Brown   }
302c4762a1bSJed Brown   PetscFunctionReturn(0);
303c4762a1bSJed Brown }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown /* ------------------------------------------------------------------- */
306c4762a1bSJed Brown PetscErrorCode FormInitialSolution(Vec U,void *ptr)
307c4762a1bSJed Brown {
308c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ptr;
309c4762a1bSJed Brown   DM             da   =user->da;
310c4762a1bSJed Brown   PetscReal      c    =user->c;
311c4762a1bSJed Brown   PetscErrorCode ierr;
312c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
313c4762a1bSJed Brown   PetscScalar    **u;
314c4762a1bSJed Brown   PetscReal      hx,hy,x,y,r;
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   PetscFunctionBeginUser;
317c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
320c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1);
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   /* Get pointers to vector data */
323c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   /* Get local grid boundaries */
326c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
329c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
330c4762a1bSJed Brown     y = j*hy;
331c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
332c4762a1bSJed Brown       x = i*hx;
333c4762a1bSJed Brown       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
334c4762a1bSJed Brown       if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
335c4762a1bSJed Brown       else u[j][i] = 0.0;
336c4762a1bSJed Brown     }
337c4762a1bSJed Brown   }
338c4762a1bSJed Brown 
339c4762a1bSJed Brown   /* Restore vectors */
340c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
341c4762a1bSJed Brown   PetscFunctionReturn(0);
342c4762a1bSJed Brown }
343c4762a1bSJed Brown 
344c4762a1bSJed Brown /*TEST
345c4762a1bSJed Brown 
346c4762a1bSJed Brown     test:
347c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor
348c4762a1bSJed Brown 
349c4762a1bSJed Brown     test:
350c4762a1bSJed Brown       suffix: 2
351c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor
352c4762a1bSJed Brown 
353c4762a1bSJed Brown     test:
354c4762a1bSJed Brown       suffix: 3
355c4762a1bSJed Brown       requires: !single
356c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
357c4762a1bSJed Brown 
358c4762a1bSJed Brown     test:
359c4762a1bSJed Brown       suffix: 4
360c4762a1bSJed Brown       requires: !single
361c4762a1bSJed Brown       nsize: 2
362c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
363c4762a1bSJed Brown 
364c4762a1bSJed Brown     test:
365c4762a1bSJed Brown       suffix: 5
366c4762a1bSJed Brown       nsize: 1
367c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor
368c4762a1bSJed Brown 
369c4762a1bSJed Brown TEST*/
370