xref: /petsc/src/ts/tests/ex4.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown /*
2c4762a1bSJed Brown        The Problem:
3c4762a1bSJed Brown            Solve the convection-diffusion equation:
4c4762a1bSJed Brown 
5c4762a1bSJed Brown              u_t+a*(u_x+u_y)=epsilon*(u_xx+u_yy)
6c4762a1bSJed Brown              u=0   at x=0, y=0
7c4762a1bSJed Brown              u_x=0 at x=1
8c4762a1bSJed Brown              u_y=0 at y=1
9c4762a1bSJed Brown              u = exp(-20.0*(pow(x-0.5,2.0)+pow(y-0.5,2.0))) at t=0
10c4762a1bSJed Brown 
11c4762a1bSJed Brown        This program tests the routine of computing the Jacobian by the
12c4762a1bSJed Brown        finite difference method as well as PETSc.
13c4762a1bSJed Brown 
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16c4762a1bSJed Brown static char help[] = "Solve the convection-diffusion equation. \n\n";
17c4762a1bSJed Brown 
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct
21c4762a1bSJed Brown {
22c4762a1bSJed Brown   PetscInt  m;          /* the number of mesh points in x-direction */
23c4762a1bSJed Brown   PetscInt  n;          /* the number of mesh points in y-direction */
24c4762a1bSJed Brown   PetscReal dx;         /* the grid space in x-direction */
25c4762a1bSJed Brown   PetscReal dy;         /* the grid space in y-direction */
26c4762a1bSJed Brown   PetscReal a;          /* the convection coefficient    */
27c4762a1bSJed Brown   PetscReal epsilon;    /* the diffusion coefficient     */
28c4762a1bSJed Brown   PetscReal tfinal;
29c4762a1bSJed Brown } Data;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
32c4762a1bSJed Brown extern PetscErrorCode Initial(Vec,void*);
33c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
34c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
35c4762a1bSJed Brown extern PetscErrorCode PostStep(TS);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown int main(int argc,char **argv)
38c4762a1bSJed Brown {
39c4762a1bSJed Brown   PetscInt       time_steps=100,iout,NOUT=1;
40c4762a1bSJed Brown   Vec            global;
41c4762a1bSJed Brown   PetscReal      dt,ftime,ftime_original;
42c4762a1bSJed Brown   TS             ts;
43c4762a1bSJed Brown   PetscViewer    viewfile;
44c4762a1bSJed Brown   Mat            J = 0;
45c4762a1bSJed Brown   Vec            x;
46c4762a1bSJed Brown   Data           data;
47c4762a1bSJed Brown   PetscInt       mn;
48c4762a1bSJed Brown   PetscBool      flg;
49c4762a1bSJed Brown   MatColoring    mc;
50c4762a1bSJed Brown   ISColoring     iscoloring;
51c4762a1bSJed Brown   MatFDColoring  matfdcoloring        = 0;
52c4762a1bSJed Brown   PetscBool      fd_jacobian_coloring = PETSC_FALSE;
53c4762a1bSJed Brown   SNES           snes;
54c4762a1bSJed Brown   KSP            ksp;
55c4762a1bSJed Brown   PC             pc;
56c4762a1bSJed Brown 
57*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* set data */
60c4762a1bSJed Brown   data.m       = 9;
61c4762a1bSJed Brown   data.n       = 9;
62c4762a1bSJed Brown   data.a       = 1.0;
63c4762a1bSJed Brown   data.epsilon = 0.1;
64c4762a1bSJed Brown   data.dx      = 1.0/(data.m+1.0);
65c4762a1bSJed Brown   data.dy      = 1.0/(data.n+1.0);
66c4762a1bSJed Brown   mn           = (data.m)*(data.n);
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* set initial conditions */
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&global));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(global,PETSC_DECIDE,mn));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(global));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(Initial(global,&data));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(global,&x));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* create timestep context */
775f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(TSMonitorSet(ts,Monitor,&data,NULL));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSEULER));
80c4762a1bSJed Brown   dt   = 0.1;
81c4762a1bSJed Brown   ftime_original = data.tfinal = 1.0;
82c4762a1bSJed Brown 
835f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,time_steps));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,ftime_original));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,global));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* set user provided RHSFunction and RHSJacobian */
905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&data));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,mn,mn));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqAIJSetPreallocation(J,5,NULL));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIAIJSetPreallocation(J,5,NULL,5,NULL));
96c4762a1bSJed Brown 
975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-ts_fd",&flg));
98c4762a1bSJed Brown   if (!flg) {
995f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobian(ts,J,J,RHSJacobian,&data));
100c4762a1bSJed Brown   } else {
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSNES(ts,&snes));
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsHasName(NULL,NULL,"-fd_color",&fd_jacobian_coloring));
103c4762a1bSJed Brown     if (fd_jacobian_coloring) { /* Use finite differences with coloring */
104c4762a1bSJed Brown       /* Get data structure of J */
105c4762a1bSJed Brown       PetscBool pc_diagonal;
1065f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscOptionsHasName(NULL,NULL,"-pc_diagonal",&pc_diagonal));
107c4762a1bSJed Brown       if (pc_diagonal) { /* the preconditioner of J is a diagonal matrix */
108c4762a1bSJed Brown         PetscInt    rstart,rend,i;
109c4762a1bSJed Brown         PetscScalar zero=0.0;
1105f80ce2aSJacob Faibussowitsch         CHKERRQ(MatGetOwnershipRange(J,&rstart,&rend));
111c4762a1bSJed Brown         for (i=rstart; i<rend; i++) {
1125f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValues(J,1,&i,1,&i,&zero,INSERT_VALUES));
113c4762a1bSJed Brown         }
1145f80ce2aSJacob Faibussowitsch         CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1155f80ce2aSJacob Faibussowitsch         CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
116c4762a1bSJed Brown       } else {
117c4762a1bSJed Brown         /* Fill the structure using the expensive SNESComputeJacobianDefault. Temporarily set up the TS so we can call this function */
1185f80ce2aSJacob Faibussowitsch         CHKERRQ(TSSetType(ts,TSBEULER));
1195f80ce2aSJacob Faibussowitsch         CHKERRQ(TSSetUp(ts));
1205f80ce2aSJacob Faibussowitsch         CHKERRQ(SNESComputeJacobianDefault(snes,x,J,J,ts));
121c4762a1bSJed Brown       }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown       /* create coloring context */
1245f80ce2aSJacob Faibussowitsch       CHKERRQ(MatColoringCreate(J,&mc));
1255f80ce2aSJacob Faibussowitsch       CHKERRQ(MatColoringSetType(mc,MATCOLORINGSL));
1265f80ce2aSJacob Faibussowitsch       CHKERRQ(MatColoringSetFromOptions(mc));
1275f80ce2aSJacob Faibussowitsch       CHKERRQ(MatColoringApply(mc,&iscoloring));
1285f80ce2aSJacob Faibussowitsch       CHKERRQ(MatColoringDestroy(&mc));
1295f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFDColoringCreate(J,iscoloring,&matfdcoloring));
1305f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts));
1315f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFDColoringSetFromOptions(matfdcoloring));
1325f80ce2aSJacob Faibussowitsch       CHKERRQ(MatFDColoringSetUp(J,iscoloring,matfdcoloring));
1335f80ce2aSJacob Faibussowitsch       CHKERRQ(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring));
1345f80ce2aSJacob Faibussowitsch       CHKERRQ(ISColoringDestroy(&iscoloring));
135c4762a1bSJed Brown     } else { /* Use finite differences (slow) */
1365f80ce2aSJacob Faibussowitsch       CHKERRQ(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL));
137c4762a1bSJed Brown     }
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* Pick up a Petsc preconditioner */
141c4762a1bSJed Brown   /* one can always set method or preconditioner during the run time */
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSNES(ts,&snes));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(pc,PCJACOBI));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
147c4762a1bSJed Brown 
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetUp(ts));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Test TSSetPostStep() */
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_PostStep",&flg));
153c4762a1bSJed Brown   if (flg) {
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetPostStep(ts,PostStep));
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown 
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-NOUT",&NOUT,NULL));
158c4762a1bSJed Brown   for (iout=1; iout<=NOUT; iout++) {
1595f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetMaxSteps(ts,time_steps));
1605f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetMaxTime(ts,iout*ftime_original/NOUT));
1615f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSolve(ts,global));
1625f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSolveTime(ts,&ftime));
1635f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetTime(ts,ftime));
1645f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetTimeStep(ts,dt));
165c4762a1bSJed Brown   }
166c4762a1bSJed Brown   /* Interpolate solution at tfinal */
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts,&global));
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(TSInterpolate(ts,ftime_original,global));
169c4762a1bSJed Brown 
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-matlab_view",&flg));
171c4762a1bSJed Brown   if (flg) { /* print solution into a MATLAB file */
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIOpen(PETSC_COMM_WORLD,"out.m",&viewfile));
1735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPushFormat(viewfile,PETSC_VIEWER_ASCII_MATLAB));
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(global,viewfile));
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerPopFormat(viewfile));
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewfile));
177c4762a1bSJed Brown   }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* free the memories */
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&global));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
1845f80ce2aSJacob Faibussowitsch   if (fd_jacobian_coloring) CHKERRQ(MatFDColoringDestroy(&matfdcoloring));
185*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
186*b122ec5aSJacob Faibussowitsch   return 0;
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown /* -------------------------------------------------------------------*/
190c4762a1bSJed Brown /* the initial function */
191c4762a1bSJed Brown PetscReal f_ini(PetscReal x,PetscReal y)
192c4762a1bSJed Brown {
193c4762a1bSJed Brown   PetscReal f;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   f=PetscExpReal(-20.0*(PetscPowRealInt(x-0.5,2)+PetscPowRealInt(y-0.5,2)));
196c4762a1bSJed Brown   return f;
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown PetscErrorCode Initial(Vec global,void *ctx)
200c4762a1bSJed Brown {
201c4762a1bSJed Brown   Data           *data = (Data*)ctx;
202c4762a1bSJed Brown   PetscInt       m,row,col;
203c4762a1bSJed Brown   PetscReal      x,y,dx,dy;
204c4762a1bSJed Brown   PetscScalar    *localptr;
205c4762a1bSJed Brown   PetscInt       i,mybase,myend,locsize;
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   PetscFunctionBeginUser;
208c4762a1bSJed Brown   /* make the local  copies of parameters */
209c4762a1bSJed Brown   m  = data->m;
210c4762a1bSJed Brown   dx = data->dx;
211c4762a1bSJed Brown   dy = data->dy;
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* determine starting point of each processor */
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(global,&mybase,&myend));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(global,&locsize));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* Initialize the array */
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(global,&localptr));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   for (i=0; i<locsize; i++) {
221c4762a1bSJed Brown     row         = 1+(mybase+i)-((mybase+i)/m)*m;
222c4762a1bSJed Brown     col         = (mybase+i)/m+1;
223c4762a1bSJed Brown     x           = dx*row;
224c4762a1bSJed Brown     y           = dy*col;
225c4762a1bSJed Brown     localptr[i] = f_ini(x,y);
226c4762a1bSJed Brown   }
227c4762a1bSJed Brown 
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(global,&localptr));
229c4762a1bSJed Brown   PetscFunctionReturn(0);
230c4762a1bSJed Brown }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
233c4762a1bSJed Brown {
234c4762a1bSJed Brown   VecScatter        scatter;
235c4762a1bSJed Brown   IS                from,to;
236c4762a1bSJed Brown   PetscInt          i,n,*idx,nsteps,maxsteps;
237c4762a1bSJed Brown   Vec               tmp_vec;
238c4762a1bSJed Brown   const PetscScalar *tmp;
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   PetscFunctionBeginUser;
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&nsteps));
242c4762a1bSJed Brown   /* display output at selected time steps */
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetMaxSteps(ts, &maxsteps));
244c4762a1bSJed Brown   if (nsteps % 10 != 0) PetscFunctionReturn(0);
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /* Get the size of the vector */
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(global,&n));
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* Set the index sets */
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&idx));
251c4762a1bSJed Brown   for (i=0; i<n; i++) idx[i]=i;
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   /* Create local sequential vectors */
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Create scatter context */
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from));
2585f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to));
2595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(global,from,tmp_vec,to,&scatter));
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD));
2615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD));
262c4762a1bSJed Brown 
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(tmp_vec,&tmp));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"At t[%D] =%14.2e u= %14.2e at the center \n",nsteps,(double)time,(double)PetscRealPart(tmp[n/2])));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(tmp_vec,&tmp));
266c4762a1bSJed Brown 
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&from));
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&to));
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&scatter));
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&tmp_vec));
272c4762a1bSJed Brown   PetscFunctionReturn(0);
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
275c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ptr)
276c4762a1bSJed Brown {
277c4762a1bSJed Brown   Data           *data = (Data*)ptr;
278c4762a1bSJed Brown   PetscScalar    v[5];
279c4762a1bSJed Brown   PetscInt       idx[5],i,j,row;
280c4762a1bSJed Brown   PetscInt       m,n,mn;
281c4762a1bSJed Brown   PetscReal      dx,dy,a,epsilon,xc,xl,xr,yl,yr;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   PetscFunctionBeginUser;
284c4762a1bSJed Brown   m       = data->m;
285c4762a1bSJed Brown   n       = data->n;
286c4762a1bSJed Brown   mn      = m*n;
287c4762a1bSJed Brown   dx      = data->dx;
288c4762a1bSJed Brown   dy      = data->dy;
289c4762a1bSJed Brown   a       = data->a;
290c4762a1bSJed Brown   epsilon = data->epsilon;
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
293c4762a1bSJed Brown   xl = 0.5*a/dx+epsilon/(dx*dx);
294c4762a1bSJed Brown   xr = -0.5*a/dx+epsilon/(dx*dx);
295c4762a1bSJed Brown   yl = 0.5*a/dy+epsilon/(dy*dy);
296c4762a1bSJed Brown   yr = -0.5*a/dy+epsilon/(dy*dy);
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   row    = 0;
299c4762a1bSJed Brown   v[0]   = xc;  v[1] = xr;  v[2] = yr;
300c4762a1bSJed Brown   idx[0] = 0; idx[1] = 2; idx[2] = m;
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   row    = m-1;
304c4762a1bSJed Brown   v[0]   = 2.0*xl; v[1] = xc;    v[2] = yr;
305c4762a1bSJed Brown   idx[0] = m-2;  idx[1] = m-1; idx[2] = m-1+m;
3065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
309c4762a1bSJed Brown     row    = i;
310c4762a1bSJed Brown     v[0]   = xl;    v[1] = xc;  v[2] = xr;    v[3] = yr;
311c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1; idx[3] = i+m;
3125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
313c4762a1bSJed Brown   }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   for (j=1; j<n-1; j++) {
316c4762a1bSJed Brown     row    = j*m;
317c4762a1bSJed Brown     v[0]   = xc;    v[1] = xr;    v[2] = yl;      v[3] = yr;
318c4762a1bSJed Brown     idx[0] = j*m; idx[1] = j*m; idx[2] = j*m-m; idx[3] = j*m+m;
3195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
320c4762a1bSJed Brown 
321c4762a1bSJed Brown     row    = j*m+m-1;
322c4762a1bSJed Brown     v[0]   = xc;        v[1] = 2.0*xl;      v[2] = yl;          v[3] = yr;
323c4762a1bSJed Brown     idx[0] = j*m+m-1; idx[1] = j*m+m-1-1; idx[2] = j*m+m-1-m; idx[3] = j*m+m-1+m;
3245f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
325c4762a1bSJed Brown 
326c4762a1bSJed Brown     for (i=1; i<m-1; i++) {
327c4762a1bSJed Brown       row    = j*m+i;
328c4762a1bSJed Brown       v[0]   = xc;      v[1] = xl;        v[2] = xr;        v[3] = yl; v[4]=yr;
329c4762a1bSJed Brown       idx[0] = j*m+i; idx[1] = j*m+i-1; idx[2] = j*m+i+1; idx[3] = j*m+i-m;
330c4762a1bSJed Brown       idx[4] = j*m+i+m;
3315f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A,1,&row,5,idx,v,INSERT_VALUES));
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown   }
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   row    = mn-m;
336c4762a1bSJed Brown   v[0]   = xc;     v[1] = xr;       v[2] = 2.0*yl;
337c4762a1bSJed Brown   idx[0] = mn-m; idx[1] = mn-m+1; idx[2] = mn-m-m;
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES));
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   row    = mn-1;
341c4762a1bSJed Brown   v[0]   = xc;     v[1] = 2.0*xl; v[2] = 2.0*yl;
342c4762a1bSJed Brown   idx[0] = mn-1; idx[1] = mn-2; idx[2] = mn-1-m;
3435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES));
344c4762a1bSJed Brown 
345c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
346c4762a1bSJed Brown     row    = mn-m+i;
347c4762a1bSJed Brown     v[0]   = xl;         v[1] = xc;       v[2] = xr;         v[3] = 2.0*yl;
348c4762a1bSJed Brown     idx[0] = mn-m+i-1; idx[1] = mn-m+i; idx[2] = mn-m+i+1; idx[3] = mn-m+i-m;
3495f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
350c4762a1bSJed Brown   }
351c4762a1bSJed Brown 
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   PetscFunctionReturn(0);
356c4762a1bSJed Brown }
357c4762a1bSJed Brown 
358c4762a1bSJed Brown /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
359c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
360c4762a1bSJed Brown {
361c4762a1bSJed Brown   Data              *data = (Data*)ctx;
362c4762a1bSJed Brown   PetscInt          m,n,mn;
363c4762a1bSJed Brown   PetscReal         dx,dy;
364c4762a1bSJed Brown   PetscReal         xc,xl,xr,yl,yr;
365c4762a1bSJed Brown   PetscReal         a,epsilon;
366c4762a1bSJed Brown   PetscScalar       *outptr;
367c4762a1bSJed Brown   const PetscScalar *inptr;
368c4762a1bSJed Brown   PetscInt          i,j,len;
369c4762a1bSJed Brown   IS                from,to;
370c4762a1bSJed Brown   PetscInt          *idx;
371c4762a1bSJed Brown   VecScatter        scatter;
372c4762a1bSJed Brown   Vec               tmp_in,tmp_out;
373c4762a1bSJed Brown 
374c4762a1bSJed Brown   PetscFunctionBeginUser;
375c4762a1bSJed Brown   m       = data->m;
376c4762a1bSJed Brown   n       = data->n;
377c4762a1bSJed Brown   mn      = m*n;
378c4762a1bSJed Brown   dx      = data->dx;
379c4762a1bSJed Brown   dy      = data->dy;
380c4762a1bSJed Brown   a       = data->a;
381c4762a1bSJed Brown   epsilon = data->epsilon;
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
384c4762a1bSJed Brown   xl = 0.5*a/dx+epsilon/(dx*dx);
385c4762a1bSJed Brown   xr = -0.5*a/dx+epsilon/(dx*dx);
386c4762a1bSJed Brown   yl = 0.5*a/dy+epsilon/(dy*dy);
387c4762a1bSJed Brown   yr = -0.5*a/dy+epsilon/(dy*dy);
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   /* Get the length of parallel vector */
3905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetSize(globalin,&len));
391c4762a1bSJed Brown 
392c4762a1bSJed Brown   /* Set the index sets */
3935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(len,&idx));
394c4762a1bSJed Brown   for (i=0; i<len; i++) idx[i]=i;
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   /* Create local sequential vectors */
3975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,len,&tmp_in));
3985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(tmp_in,&tmp_out));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   /* Create scatter context */
4015f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&from));
4025f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&to));
4035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(globalin,from,tmp_in,to,&scatter));
4045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD));
4055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD));
4065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&scatter));
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   /*Extract income array - include ghost points */
4095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(tmp_in,&inptr));
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   /* Extract outcome array*/
4125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(tmp_out,&outptr));
413c4762a1bSJed Brown 
414c4762a1bSJed Brown   outptr[0]   = xc*inptr[0]+xr*inptr[1]+yr*inptr[m];
415c4762a1bSJed Brown   outptr[m-1] = 2.0*xl*inptr[m-2]+xc*inptr[m-1]+yr*inptr[m-1+m];
416c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
417c4762a1bSJed Brown     outptr[i] = xc*inptr[i]+xl*inptr[i-1]+xr*inptr[i+1]+yr*inptr[i+m];
418c4762a1bSJed Brown   }
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   for (j=1; j<n-1; j++) {
421c4762a1bSJed Brown     outptr[j*m] = xc*inptr[j*m]+xr*inptr[j*m+1]+ yl*inptr[j*m-m]+yr*inptr[j*m+m];
422c4762a1bSJed Brown     outptr[j*m+m-1] = xc*inptr[j*m+m-1]+2.0*xl*inptr[j*m+m-1-1]+ yl*inptr[j*m+m-1-m]+yr*inptr[j*m+m-1+m];
423c4762a1bSJed Brown     for (i=1; i<m-1; i++) {
424c4762a1bSJed Brown       outptr[j*m+i] = xc*inptr[j*m+i]+xl*inptr[j*m+i-1]+xr*inptr[j*m+i+1]+yl*inptr[j*m+i-m]+yr*inptr[j*m+i+m];
425c4762a1bSJed Brown     }
426c4762a1bSJed Brown   }
427c4762a1bSJed Brown 
428c4762a1bSJed Brown   outptr[mn-m] = xc*inptr[mn-m]+xr*inptr[mn-m+1]+2.0*yl*inptr[mn-m-m];
429c4762a1bSJed Brown   outptr[mn-1] = 2.0*xl*inptr[mn-2]+xc*inptr[mn-1]+2.0*yl*inptr[mn-1-m];
430c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
431c4762a1bSJed Brown     outptr[mn-m+i] = xc*inptr[mn-m+i]+xl*inptr[mn-m+i-1]+xr*inptr[mn-m+i+1]+2*yl*inptr[mn-m+i-m];
432c4762a1bSJed Brown   }
433c4762a1bSJed Brown 
4345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(tmp_in,&inptr));
4355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(tmp_out,&outptr));
436c4762a1bSJed Brown 
4375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterCreate(tmp_out,from,globalout,to,&scatter));
4385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD));
4395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD));
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /* Destroy idx aand scatter */
4425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&tmp_in));
4435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&tmp_out));
4445f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&from));
4455f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&to));
4465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScatterDestroy(&scatter));
447c4762a1bSJed Brown 
4485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
449c4762a1bSJed Brown   PetscFunctionReturn(0);
450c4762a1bSJed Brown }
451c4762a1bSJed Brown 
452c4762a1bSJed Brown PetscErrorCode PostStep(TS ts)
453c4762a1bSJed Brown {
454c4762a1bSJed Brown   PetscReal      t;
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   PetscFunctionBeginUser;
4575f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts,&t));
4585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"  PostStep, t: %g\n",(double)t));
459c4762a1bSJed Brown   PetscFunctionReturn(0);
460c4762a1bSJed Brown }
461c4762a1bSJed Brown 
462c4762a1bSJed Brown /*TEST
463c4762a1bSJed Brown 
464c4762a1bSJed Brown     test:
465c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
466c4762a1bSJed Brown       output_file: output/ex4.out
467c4762a1bSJed Brown 
468c4762a1bSJed Brown     test:
469c4762a1bSJed Brown       suffix: 2
470c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
471c4762a1bSJed Brown       nsize: 2
472c4762a1bSJed Brown       output_file: output/ex4.out
473c4762a1bSJed Brown 
474c4762a1bSJed Brown     test:
475c4762a1bSJed Brown       suffix: 3
476c4762a1bSJed Brown       args: -ts_fd -ts_type cn
477c4762a1bSJed Brown 
478c4762a1bSJed Brown     test:
479c4762a1bSJed Brown       suffix: 4
480c4762a1bSJed Brown       args: -ts_fd -ts_type cn
481c4762a1bSJed Brown       output_file: output/ex4_3.out
482c4762a1bSJed Brown       nsize: 2
483c4762a1bSJed Brown 
484c4762a1bSJed Brown     test:
485c4762a1bSJed Brown       suffix: 5
486c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
487c4762a1bSJed Brown       output_file: output/ex4.out
488c4762a1bSJed Brown 
489c4762a1bSJed Brown     test:
490c4762a1bSJed Brown       suffix: 6
491c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
492c4762a1bSJed Brown       output_file: output/ex4.out
493c4762a1bSJed Brown       nsize: 2
494c4762a1bSJed Brown 
495c4762a1bSJed Brown     test:
496c4762a1bSJed Brown       suffix: 7
497c4762a1bSJed Brown       requires: !single
498c4762a1bSJed Brown       args: -ts_fd -ts_type beuler -test_PostStep -ts_dt .1
499c4762a1bSJed Brown 
500c4762a1bSJed Brown     test:
501c4762a1bSJed Brown       suffix: 8
502c4762a1bSJed Brown       requires: !single
503c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 5dp -ts_dt .01 -ts_adapt_type none -ts_view
504c4762a1bSJed Brown 
505c4762a1bSJed Brown TEST*/
506