xref: /petsc/src/ts/tutorials/ex15.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   DM             da;
50c4762a1bSJed Brown   PetscReal      dt;
51c4762a1bSJed Brown   AppCtx         user;              /* user-defined work context */
52c4762a1bSJed Brown   SNES           snes;
53c4762a1bSJed Brown   PetscInt       Jtype; /* Jacobian type
54c4762a1bSJed Brown                             0: user provide Jacobian;
55c4762a1bSJed Brown                             1: slow finite difference;
56c4762a1bSJed Brown                             2: fd with coloring; */
57c4762a1bSJed Brown 
58*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
59c4762a1bSJed Brown   /* Initialize user application context */
60c4762a1bSJed Brown   user.da           = NULL;
61c4762a1bSJed Brown   user.nstencilpts  = 5;
62c4762a1bSJed Brown   user.c            = -30.0;
63c4762a1bSJed Brown   user.boundary     = 0;  /* 0: Drichlet BC; 1: Neumann BC */
64c4762a1bSJed Brown   user.viewJacobian = PETSC_FALSE;
65c4762a1bSJed Brown 
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nstencilpts",&user.nstencilpts,NULL));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
72c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73c4762a1bSJed Brown   if (user.nstencilpts == 5) {
745f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
75c4762a1bSJed Brown   } else if (user.nstencilpts == 9) {
765f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
7798921bdaSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"nstencilpts %d is not supported",user.nstencilpts);
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
80c4762a1bSJed Brown   user.da = da;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown      Extract global vectors from DMDA;
84c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&u));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&r));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89c4762a1bSJed Brown      Create timestepping solver context
90c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,r,FormIFunction,&user));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,1.0));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown      Set initial conditions
101c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialSolution(u,&user));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,u));
104c4762a1bSJed Brown   dt   = .01;
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown    Set Jacobian evaluation routine
109c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetMatType(da,MATAIJ));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(da,&J));
112c4762a1bSJed Brown   Jtype = 0;
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-Jtype",&Jtype,NULL));
114c4762a1bSJed Brown   if (Jtype == 0) { /* use user provided Jacobian evaluation routine */
1153c633725SBarry Smith     PetscCheck(user.nstencilpts == 5,PETSC_COMM_WORLD,PETSC_ERR_SUP,"user Jacobian routine FormIJacobian() does not support nstencilpts=%D",user.nstencilpts);
1165f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIJacobian(ts,J,J,FormIJacobian,&user));
117c4762a1bSJed Brown   } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSNES(ts,&snes));
1195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSNESMF(snes,&Jmf));
120c4762a1bSJed Brown     if (Jtype == 1) { /* slow finite difference J; */
1215f80ce2aSJacob Faibussowitsch       CHKERRQ(SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefault,NULL));
122c4762a1bSJed Brown     } else if (Jtype == 2) { /* Use coloring to compute  finite difference J efficiently */
1235f80ce2aSJacob Faibussowitsch       CHKERRQ(SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefaultColor,0));
124c4762a1bSJed Brown     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Jtype is not supported");
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128c4762a1bSJed Brown    Sets various TS parameters from user options
129c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133c4762a1bSJed Brown      Solve nonlinear system
134c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,u));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Free work space.
139c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Jmf));
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
146c4762a1bSJed Brown 
147*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
148*b122ec5aSJacob Faibussowitsch   return 0;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown /* --------------------------------------------------------------------- */
152c4762a1bSJed Brown /*
153c4762a1bSJed Brown   FormIFunction = Udot - RHSFunction
154c4762a1bSJed Brown */
155c4762a1bSJed Brown PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
156c4762a1bSJed Brown {
157c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ctx;
158c4762a1bSJed Brown   DM             da   = (DM)user->da;
159c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
160c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
161c4762a1bSJed Brown   PetscScalar    u,uxx,uyy,**uarray,**f,**udot;
162c4762a1bSJed Brown   Vec            localU;
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   PetscFunctionBeginUser;
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
169c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
1703c633725SBarry Smith   PetscCheck(user->nstencilpts != 9 || hx == hy,PETSC_COMM_WORLD,PETSC_ERR_SUP,"hx must equal hy when nstencilpts = 9 for this example");
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /*
173c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
174c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
175c4762a1bSJed Brown      By placing code between these two statements, computations can be
176c4762a1bSJed Brown      done while messages are in transition.
177c4762a1bSJed Brown   */
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Get pointers to vector data */
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,localU,&uarray));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Udot,&udot));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* Get local grid boundaries */
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
190c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
191c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
192c4762a1bSJed Brown       /* Boundary conditions */
193c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
194c4762a1bSJed Brown         if (user->boundary == 0) { /* Drichlet BC */
195c4762a1bSJed Brown           f[j][i] = uarray[j][i]; /* F = U */
196c4762a1bSJed Brown         } else {                  /* Neumann BC */
197c4762a1bSJed Brown           if (i == 0 && j == 0) {              /* SW corner */
198c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j+1][i+1];
199c4762a1bSJed Brown           } else if (i == Mx-1 && j == 0) {    /* SE corner */
200c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j+1][i-1];
201c4762a1bSJed Brown           } else if (i == 0 && j == My-1) {    /* NW corner */
202c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j-1][i+1];
203c4762a1bSJed Brown           } else if (i == Mx-1 && j == My-1) { /* NE corner */
204c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j-1][i-1];
205c4762a1bSJed Brown           } else if (i == 0) {                  /* Left */
206c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j][i+1];
207c4762a1bSJed Brown           } else if (i == Mx-1) {               /* Right */
208c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j][i-1];
209c4762a1bSJed Brown           } else if (j == 0) {                 /* Bottom */
210c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j+1][i];
211c4762a1bSJed Brown           } else if (j == My-1) {               /* Top */
212c4762a1bSJed Brown             f[j][i] = uarray[j][i] - uarray[j-1][i];
213c4762a1bSJed Brown           }
214c4762a1bSJed Brown         }
215c4762a1bSJed Brown       } else { /* Interior */
216c4762a1bSJed Brown         u = uarray[j][i];
217c4762a1bSJed Brown         /* 5-point stencil */
218c4762a1bSJed Brown         uxx = (-2.0*u + uarray[j][i-1] + uarray[j][i+1]);
219c4762a1bSJed Brown         uyy = (-2.0*u + uarray[j-1][i] + uarray[j+1][i]);
220c4762a1bSJed Brown         if (user->nstencilpts == 9) {
221c4762a1bSJed Brown           /* 9-point stencil: assume hx=hy */
222c4762a1bSJed 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;
223c4762a1bSJed 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;
224c4762a1bSJed Brown         }
225c4762a1bSJed Brown         f[j][i] = udot[j][i] - (uxx*sx + uyy*sy);
226c4762a1bSJed Brown       }
227c4762a1bSJed Brown     }
228c4762a1bSJed Brown   }
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   /* Restore vectors */
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&uarray));
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Udot,&udot));
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(11.0*ym*xm));
236c4762a1bSJed Brown   PetscFunctionReturn(0);
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown /* --------------------------------------------------------------------- */
240c4762a1bSJed Brown /*
241c4762a1bSJed Brown   FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot
242c4762a1bSJed Brown   This routine is not used with option '-use_coloring'
243c4762a1bSJed Brown */
244c4762a1bSJed Brown PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
245c4762a1bSJed Brown {
246c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym,nc;
247c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ctx;
248c4762a1bSJed Brown   DM             da    = (DM)user->da;
249c4762a1bSJed Brown   MatStencil     col[5],row;
250c4762a1bSJed Brown   PetscScalar    vals[5],hx,hy,sx,sy;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   PetscFunctionBeginUser;
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
257c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
260c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
261c4762a1bSJed Brown       nc    = 0;
262c4762a1bSJed Brown       row.j = j; row.i = i;
263c4762a1bSJed Brown       if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) {
264c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0;
265c4762a1bSJed Brown 
266c4762a1bSJed Brown       } else if (user->boundary > 0 && i == 0) {  /* Left Neumann */
267c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
268c4762a1bSJed Brown         col[nc].j = j; col[nc].i = i+1; vals[nc++] = -1.0;
269c4762a1bSJed Brown       } else if (user->boundary > 0 && i == Mx-1) { /* Right 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 && j == 0) {  /* Bottom Neumann */
273c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i; vals[nc++] = 1.0;
274c4762a1bSJed Brown         col[nc].j = j+1; col[nc].i = i; vals[nc++] = -1.0;
275c4762a1bSJed Brown       } else if (user->boundary > 0 && j == My-1) { /* Top 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 {   /* Interior */
279c4762a1bSJed Brown         col[nc].j = j-1; col[nc].i = i;   vals[nc++] = -sy;
280c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i-1; vals[nc++] = -sx;
281c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i;   vals[nc++] = 2.0*(sx + sy) + a;
282c4762a1bSJed Brown         col[nc].j = j;   col[nc].i = i+1; vals[nc++] = -sx;
283c4762a1bSJed Brown         col[nc].j = j+1; col[nc].i = i;   vals[nc++] = -sy;
284c4762a1bSJed Brown       }
2855f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES));
286c4762a1bSJed Brown     }
287c4762a1bSJed Brown   }
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
2895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
290c4762a1bSJed Brown   if (J != Jpre) {
2915f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
2925f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   if (user->viewJacobian) {
2965f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)Jpre),"Jpre:\n"));
2975f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD));
298c4762a1bSJed Brown   }
299c4762a1bSJed Brown   PetscFunctionReturn(0);
300c4762a1bSJed Brown }
301c4762a1bSJed Brown 
302c4762a1bSJed Brown /* ------------------------------------------------------------------- */
303c4762a1bSJed Brown PetscErrorCode FormInitialSolution(Vec U,void *ptr)
304c4762a1bSJed Brown {
305c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ptr;
306c4762a1bSJed Brown   DM             da   =user->da;
307c4762a1bSJed Brown   PetscReal      c    =user->c;
308c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
309c4762a1bSJed Brown   PetscScalar    **u;
310c4762a1bSJed Brown   PetscReal      hx,hy,x,y,r;
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   PetscFunctionBeginUser;
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
316c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1);
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /* Get pointers to vector data */
3195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,U,&u));
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   /* Get local grid boundaries */
3225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
325c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
326c4762a1bSJed Brown     y = j*hy;
327c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
328c4762a1bSJed Brown       x = i*hx;
329c4762a1bSJed Brown       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
330c4762a1bSJed Brown       if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
331c4762a1bSJed Brown       else u[j][i] = 0.0;
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown   }
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   /* Restore vectors */
3365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,U,&u));
337c4762a1bSJed Brown   PetscFunctionReturn(0);
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown /*TEST
341c4762a1bSJed Brown 
342c4762a1bSJed Brown     test:
343c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor
344c4762a1bSJed Brown 
345c4762a1bSJed Brown     test:
346c4762a1bSJed Brown       suffix: 2
347c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor
348c4762a1bSJed Brown 
349c4762a1bSJed Brown     test:
350c4762a1bSJed Brown       suffix: 3
351c4762a1bSJed Brown       requires: !single
352c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
353c4762a1bSJed Brown 
354c4762a1bSJed Brown     test:
355c4762a1bSJed Brown       suffix: 4
356c4762a1bSJed Brown       requires: !single
357c4762a1bSJed Brown       nsize: 2
358c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
359c4762a1bSJed Brown 
360c4762a1bSJed Brown     test:
361c4762a1bSJed Brown       suffix: 5
362c4762a1bSJed Brown       nsize: 1
363c4762a1bSJed Brown       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor
364c4762a1bSJed Brown 
365c4762a1bSJed Brown TEST*/
366