xref: /petsc/src/ts/tests/ex4.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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 
579566063dSJacob Faibussowitsch   PetscCall(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);
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* set initial conditions */
709566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&global));
719566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(global,PETSC_DECIDE,mn));
729566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(global));
739566063dSJacob Faibussowitsch   PetscCall(Initial(global,&data));
749566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(global,&x));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* create timestep context */
779566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
789566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts,Monitor,&data,NULL));
799566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSEULER));
80c4762a1bSJed Brown   dt   = 0.1;
81c4762a1bSJed Brown   ftime_original = data.tfinal = 1.0;
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,dt));
849566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,time_steps));
859566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime_original));
869566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
879566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,global));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* set user provided RHSFunction and RHSJacobian */
909566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&data));
919566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,mn,mn));
939566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
949566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J,5,NULL));
959566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(J,5,NULL,5,NULL));
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-ts_fd",&flg));
98c4762a1bSJed Brown   if (!flg) {
999566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts,J,J,RHSJacobian,&data));
100c4762a1bSJed Brown   } else {
1019566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts,&snes));
1029566063dSJacob Faibussowitsch     PetscCall(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;
1069566063dSJacob Faibussowitsch       PetscCall(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;
1109566063dSJacob Faibussowitsch         PetscCall(MatGetOwnershipRange(J,&rstart,&rend));
111c4762a1bSJed Brown         for (i=rstart; i<rend; i++) {
1129566063dSJacob Faibussowitsch           PetscCall(MatSetValues(J,1,&i,1,&i,&zero,INSERT_VALUES));
113c4762a1bSJed Brown         }
1149566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1159566063dSJacob Faibussowitsch         PetscCall(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 */
1189566063dSJacob Faibussowitsch         PetscCall(TSSetType(ts,TSBEULER));
1199566063dSJacob Faibussowitsch         PetscCall(TSSetUp(ts));
1209566063dSJacob Faibussowitsch         PetscCall(SNESComputeJacobianDefault(snes,x,J,J,ts));
121c4762a1bSJed Brown       }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown       /* create coloring context */
1249566063dSJacob Faibussowitsch       PetscCall(MatColoringCreate(J,&mc));
1259566063dSJacob Faibussowitsch       PetscCall(MatColoringSetType(mc,MATCOLORINGSL));
1269566063dSJacob Faibussowitsch       PetscCall(MatColoringSetFromOptions(mc));
1279566063dSJacob Faibussowitsch       PetscCall(MatColoringApply(mc,&iscoloring));
1289566063dSJacob Faibussowitsch       PetscCall(MatColoringDestroy(&mc));
1299566063dSJacob Faibussowitsch       PetscCall(MatFDColoringCreate(J,iscoloring,&matfdcoloring));
1309566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts));
1319566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetFromOptions(matfdcoloring));
1329566063dSJacob Faibussowitsch       PetscCall(MatFDColoringSetUp(J,iscoloring,matfdcoloring));
1339566063dSJacob Faibussowitsch       PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring));
1349566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&iscoloring));
135c4762a1bSJed Brown     } else { /* Use finite differences (slow) */
1369566063dSJacob Faibussowitsch       PetscCall(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 */
1429566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&snes));
1439566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes,&ksp));
1449566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp,&pc));
1459566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCJACOBI));
1469566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1499566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Test TSSetPostStep() */
1529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_PostStep",&flg));
153*1baa6e33SBarry Smith   if (flg) PetscCall(TSSetPostStep(ts,PostStep));
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-NOUT",&NOUT,NULL));
156c4762a1bSJed Brown   for (iout=1; iout<=NOUT; iout++) {
1579566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts,time_steps));
1589566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts,iout*ftime_original/NOUT));
1599566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts,global));
1609566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts,&ftime));
1619566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts,ftime));
1629566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts,dt));
163c4762a1bSJed Brown   }
164c4762a1bSJed Brown   /* Interpolate solution at tfinal */
1659566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&global));
1669566063dSJacob Faibussowitsch   PetscCall(TSInterpolate(ts,ftime_original,global));
167c4762a1bSJed Brown 
1689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-matlab_view",&flg));
169c4762a1bSJed Brown   if (flg) { /* print solution into a MATLAB file */
1709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD,"out.m",&viewfile));
1719566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewfile,PETSC_VIEWER_ASCII_MATLAB));
1729566063dSJacob Faibussowitsch     PetscCall(VecView(global,viewfile));
1739566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewfile));
1749566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewfile));
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* free the memories */
1789566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
1809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1829566063dSJacob Faibussowitsch   if (fd_jacobian_coloring) PetscCall(MatFDColoringDestroy(&matfdcoloring));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
184b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /* -------------------------------------------------------------------*/
188c4762a1bSJed Brown /* the initial function */
189c4762a1bSJed Brown PetscReal f_ini(PetscReal x,PetscReal y)
190c4762a1bSJed Brown {
191c4762a1bSJed Brown   PetscReal f;
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   f=PetscExpReal(-20.0*(PetscPowRealInt(x-0.5,2)+PetscPowRealInt(y-0.5,2)));
194c4762a1bSJed Brown   return f;
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown PetscErrorCode Initial(Vec global,void *ctx)
198c4762a1bSJed Brown {
199c4762a1bSJed Brown   Data           *data = (Data*)ctx;
200c4762a1bSJed Brown   PetscInt       m,row,col;
201c4762a1bSJed Brown   PetscReal      x,y,dx,dy;
202c4762a1bSJed Brown   PetscScalar    *localptr;
203c4762a1bSJed Brown   PetscInt       i,mybase,myend,locsize;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   PetscFunctionBeginUser;
206c4762a1bSJed Brown   /* make the local  copies of parameters */
207c4762a1bSJed Brown   m  = data->m;
208c4762a1bSJed Brown   dx = data->dx;
209c4762a1bSJed Brown   dy = data->dy;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* determine starting point of each processor */
2129566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(global,&mybase,&myend));
2139566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(global,&locsize));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* Initialize the array */
2169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(global,&localptr));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   for (i=0; i<locsize; i++) {
219c4762a1bSJed Brown     row         = 1+(mybase+i)-((mybase+i)/m)*m;
220c4762a1bSJed Brown     col         = (mybase+i)/m+1;
221c4762a1bSJed Brown     x           = dx*row;
222c4762a1bSJed Brown     y           = dy*col;
223c4762a1bSJed Brown     localptr[i] = f_ini(x,y);
224c4762a1bSJed Brown   }
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(global,&localptr));
227c4762a1bSJed Brown   PetscFunctionReturn(0);
228c4762a1bSJed Brown }
229c4762a1bSJed Brown 
230c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
231c4762a1bSJed Brown {
232c4762a1bSJed Brown   VecScatter        scatter;
233c4762a1bSJed Brown   IS                from,to;
234c4762a1bSJed Brown   PetscInt          i,n,*idx,nsteps,maxsteps;
235c4762a1bSJed Brown   Vec               tmp_vec;
236c4762a1bSJed Brown   const PetscScalar *tmp;
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   PetscFunctionBeginUser;
2399566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&nsteps));
240c4762a1bSJed Brown   /* display output at selected time steps */
2419566063dSJacob Faibussowitsch   PetscCall(TSGetMaxSteps(ts, &maxsteps));
242c4762a1bSJed Brown   if (nsteps % 10 != 0) PetscFunctionReturn(0);
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /* Get the size of the vector */
2459566063dSJacob Faibussowitsch   PetscCall(VecGetSize(global,&n));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /* Set the index sets */
2489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&idx));
249c4762a1bSJed Brown   for (i=0; i<n; i++) idx[i]=i;
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* Create local sequential vectors */
2529566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /* Create scatter context */
2559566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from));
2569566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to));
2579566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(global,from,tmp_vec,to,&scatter));
2589566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD));
2599566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD));
260c4762a1bSJed Brown 
2619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_vec,&tmp));
26263a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"At t[%" PetscInt_FMT "] =%14.2e u= %14.2e at the center \n",nsteps,(double)time,(double)PetscRealPart(tmp[n/2])));
2639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_vec,&tmp));
264c4762a1bSJed Brown 
2659566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
2669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
2679566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
2689566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
2699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_vec));
270c4762a1bSJed Brown   PetscFunctionReturn(0);
271c4762a1bSJed Brown }
272c4762a1bSJed Brown 
273c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ptr)
274c4762a1bSJed Brown {
275c4762a1bSJed Brown   Data           *data = (Data*)ptr;
276c4762a1bSJed Brown   PetscScalar    v[5];
277c4762a1bSJed Brown   PetscInt       idx[5],i,j,row;
278c4762a1bSJed Brown   PetscInt       m,n,mn;
279c4762a1bSJed Brown   PetscReal      dx,dy,a,epsilon,xc,xl,xr,yl,yr;
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   PetscFunctionBeginUser;
282c4762a1bSJed Brown   m       = data->m;
283c4762a1bSJed Brown   n       = data->n;
284c4762a1bSJed Brown   mn      = m*n;
285c4762a1bSJed Brown   dx      = data->dx;
286c4762a1bSJed Brown   dy      = data->dy;
287c4762a1bSJed Brown   a       = data->a;
288c4762a1bSJed Brown   epsilon = data->epsilon;
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
291c4762a1bSJed Brown   xl = 0.5*a/dx+epsilon/(dx*dx);
292c4762a1bSJed Brown   xr = -0.5*a/dx+epsilon/(dx*dx);
293c4762a1bSJed Brown   yl = 0.5*a/dy+epsilon/(dy*dy);
294c4762a1bSJed Brown   yr = -0.5*a/dy+epsilon/(dy*dy);
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   row    = 0;
297c4762a1bSJed Brown   v[0]   = xc;  v[1] = xr;  v[2] = yr;
298c4762a1bSJed Brown   idx[0] = 0; idx[1] = 2; idx[2] = m;
2999566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES));
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   row    = m-1;
302c4762a1bSJed Brown   v[0]   = 2.0*xl; v[1] = xc;    v[2] = yr;
303c4762a1bSJed Brown   idx[0] = m-2;  idx[1] = m-1; idx[2] = m-1+m;
3049566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES));
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
307c4762a1bSJed Brown     row    = i;
308c4762a1bSJed Brown     v[0]   = xl;    v[1] = xc;  v[2] = xr;    v[3] = yr;
309c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1; idx[3] = i+m;
3109566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
311c4762a1bSJed Brown   }
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   for (j=1; j<n-1; j++) {
314c4762a1bSJed Brown     row    = j*m;
315c4762a1bSJed Brown     v[0]   = xc;    v[1] = xr;    v[2] = yl;      v[3] = yr;
316c4762a1bSJed Brown     idx[0] = j*m; idx[1] = j*m; idx[2] = j*m-m; idx[3] = j*m+m;
3179566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
318c4762a1bSJed Brown 
319c4762a1bSJed Brown     row    = j*m+m-1;
320c4762a1bSJed Brown     v[0]   = xc;        v[1] = 2.0*xl;      v[2] = yl;          v[3] = yr;
321c4762a1bSJed 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;
3229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
323c4762a1bSJed Brown 
324c4762a1bSJed Brown     for (i=1; i<m-1; i++) {
325c4762a1bSJed Brown       row    = j*m+i;
326c4762a1bSJed Brown       v[0]   = xc;      v[1] = xl;        v[2] = xr;        v[3] = yl; v[4]=yr;
327c4762a1bSJed Brown       idx[0] = j*m+i; idx[1] = j*m+i-1; idx[2] = j*m+i+1; idx[3] = j*m+i-m;
328c4762a1bSJed Brown       idx[4] = j*m+i+m;
3299566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A,1,&row,5,idx,v,INSERT_VALUES));
330c4762a1bSJed Brown     }
331c4762a1bSJed Brown   }
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   row    = mn-m;
334c4762a1bSJed Brown   v[0]   = xc;     v[1] = xr;       v[2] = 2.0*yl;
335c4762a1bSJed Brown   idx[0] = mn-m; idx[1] = mn-m+1; idx[2] = mn-m-m;
3369566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES));
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   row    = mn-1;
339c4762a1bSJed Brown   v[0]   = xc;     v[1] = 2.0*xl; v[2] = 2.0*yl;
340c4762a1bSJed Brown   idx[0] = mn-1; idx[1] = mn-2; idx[2] = mn-1-m;
3419566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES));
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
344c4762a1bSJed Brown     row    = mn-m+i;
345c4762a1bSJed Brown     v[0]   = xl;         v[1] = xc;       v[2] = xr;         v[3] = 2.0*yl;
346c4762a1bSJed Brown     idx[0] = mn-m+i-1; idx[1] = mn-m+i; idx[2] = mn-m+i+1; idx[3] = mn-m+i-m;
3479566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES));
348c4762a1bSJed Brown   }
349c4762a1bSJed Brown 
3509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
3519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   PetscFunctionReturn(0);
354c4762a1bSJed Brown }
355c4762a1bSJed Brown 
356c4762a1bSJed Brown /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
357c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
358c4762a1bSJed Brown {
359c4762a1bSJed Brown   Data              *data = (Data*)ctx;
360c4762a1bSJed Brown   PetscInt          m,n,mn;
361c4762a1bSJed Brown   PetscReal         dx,dy;
362c4762a1bSJed Brown   PetscReal         xc,xl,xr,yl,yr;
363c4762a1bSJed Brown   PetscReal         a,epsilon;
364c4762a1bSJed Brown   PetscScalar       *outptr;
365c4762a1bSJed Brown   const PetscScalar *inptr;
366c4762a1bSJed Brown   PetscInt          i,j,len;
367c4762a1bSJed Brown   IS                from,to;
368c4762a1bSJed Brown   PetscInt          *idx;
369c4762a1bSJed Brown   VecScatter        scatter;
370c4762a1bSJed Brown   Vec               tmp_in,tmp_out;
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   PetscFunctionBeginUser;
373c4762a1bSJed Brown   m       = data->m;
374c4762a1bSJed Brown   n       = data->n;
375c4762a1bSJed Brown   mn      = m*n;
376c4762a1bSJed Brown   dx      = data->dx;
377c4762a1bSJed Brown   dy      = data->dy;
378c4762a1bSJed Brown   a       = data->a;
379c4762a1bSJed Brown   epsilon = data->epsilon;
380c4762a1bSJed Brown 
381c4762a1bSJed Brown   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
382c4762a1bSJed Brown   xl = 0.5*a/dx+epsilon/(dx*dx);
383c4762a1bSJed Brown   xr = -0.5*a/dx+epsilon/(dx*dx);
384c4762a1bSJed Brown   yl = 0.5*a/dy+epsilon/(dy*dy);
385c4762a1bSJed Brown   yr = -0.5*a/dy+epsilon/(dy*dy);
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /* Get the length of parallel vector */
3889566063dSJacob Faibussowitsch   PetscCall(VecGetSize(globalin,&len));
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   /* Set the index sets */
3919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(len,&idx));
392c4762a1bSJed Brown   for (i=0; i<len; i++) idx[i]=i;
393c4762a1bSJed Brown 
394c4762a1bSJed Brown   /* Create local sequential vectors */
3959566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,len,&tmp_in));
3969566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tmp_in,&tmp_out));
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   /* Create scatter context */
3999566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&from));
4009566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&to));
4019566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(globalin,from,tmp_in,to,&scatter));
4029566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD));
4039566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD));
4049566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   /*Extract income array - include ghost points */
4079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(tmp_in,&inptr));
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   /* Extract outcome array*/
4109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(tmp_out,&outptr));
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   outptr[0]   = xc*inptr[0]+xr*inptr[1]+yr*inptr[m];
413c4762a1bSJed Brown   outptr[m-1] = 2.0*xl*inptr[m-2]+xc*inptr[m-1]+yr*inptr[m-1+m];
414c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
415c4762a1bSJed Brown     outptr[i] = xc*inptr[i]+xl*inptr[i-1]+xr*inptr[i+1]+yr*inptr[i+m];
416c4762a1bSJed Brown   }
417c4762a1bSJed Brown 
418c4762a1bSJed Brown   for (j=1; j<n-1; j++) {
419c4762a1bSJed Brown     outptr[j*m] = xc*inptr[j*m]+xr*inptr[j*m+1]+ yl*inptr[j*m-m]+yr*inptr[j*m+m];
420c4762a1bSJed 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];
421c4762a1bSJed Brown     for (i=1; i<m-1; i++) {
422c4762a1bSJed 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];
423c4762a1bSJed Brown     }
424c4762a1bSJed Brown   }
425c4762a1bSJed Brown 
426c4762a1bSJed Brown   outptr[mn-m] = xc*inptr[mn-m]+xr*inptr[mn-m+1]+2.0*yl*inptr[mn-m-m];
427c4762a1bSJed Brown   outptr[mn-1] = 2.0*xl*inptr[mn-2]+xc*inptr[mn-1]+2.0*yl*inptr[mn-1-m];
428c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
429c4762a1bSJed 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];
430c4762a1bSJed Brown   }
431c4762a1bSJed Brown 
4329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(tmp_in,&inptr));
4339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(tmp_out,&outptr));
434c4762a1bSJed Brown 
4359566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(tmp_out,from,globalout,to,&scatter));
4369566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD));
4379566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD));
438c4762a1bSJed Brown 
439c4762a1bSJed Brown   /* Destroy idx aand scatter */
4409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_in));
4419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tmp_out));
4429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
4439566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
4449566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&scatter));
445c4762a1bSJed Brown 
4469566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
447c4762a1bSJed Brown   PetscFunctionReturn(0);
448c4762a1bSJed Brown }
449c4762a1bSJed Brown 
450c4762a1bSJed Brown PetscErrorCode PostStep(TS ts)
451c4762a1bSJed Brown {
452c4762a1bSJed Brown   PetscReal      t;
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   PetscFunctionBeginUser;
4559566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&t));
4569566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"  PostStep, t: %g\n",(double)t));
457c4762a1bSJed Brown   PetscFunctionReturn(0);
458c4762a1bSJed Brown }
459c4762a1bSJed Brown 
460c4762a1bSJed Brown /*TEST
461c4762a1bSJed Brown 
462c4762a1bSJed Brown     test:
463c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
464c4762a1bSJed Brown       output_file: output/ex4.out
465c4762a1bSJed Brown 
466c4762a1bSJed Brown     test:
467c4762a1bSJed Brown       suffix: 2
468c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
469c4762a1bSJed Brown       nsize: 2
470c4762a1bSJed Brown       output_file: output/ex4.out
471c4762a1bSJed Brown 
472c4762a1bSJed Brown     test:
473c4762a1bSJed Brown       suffix: 3
474c4762a1bSJed Brown       args: -ts_fd -ts_type cn
475c4762a1bSJed Brown 
476c4762a1bSJed Brown     test:
477c4762a1bSJed Brown       suffix: 4
478c4762a1bSJed Brown       args: -ts_fd -ts_type cn
479c4762a1bSJed Brown       output_file: output/ex4_3.out
480c4762a1bSJed Brown       nsize: 2
481c4762a1bSJed Brown 
482c4762a1bSJed Brown     test:
483c4762a1bSJed Brown       suffix: 5
484c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
485c4762a1bSJed Brown       output_file: output/ex4.out
486c4762a1bSJed Brown 
487c4762a1bSJed Brown     test:
488c4762a1bSJed Brown       suffix: 6
489c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
490c4762a1bSJed Brown       output_file: output/ex4.out
491c4762a1bSJed Brown       nsize: 2
492c4762a1bSJed Brown 
493c4762a1bSJed Brown     test:
494c4762a1bSJed Brown       suffix: 7
495c4762a1bSJed Brown       requires: !single
496c4762a1bSJed Brown       args: -ts_fd -ts_type beuler -test_PostStep -ts_dt .1
497c4762a1bSJed Brown 
498c4762a1bSJed Brown     test:
499c4762a1bSJed Brown       suffix: 8
500c4762a1bSJed Brown       requires: !single
501c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 5dp -ts_dt .01 -ts_adapt_type none -ts_view
502c4762a1bSJed Brown 
503c4762a1bSJed Brown TEST*/
504