xref: /petsc/src/ts/tests/ex4.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown /*
2*c4762a1bSJed Brown        The Problem:
3*c4762a1bSJed Brown            Solve the convection-diffusion equation:
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown              u_t+a*(u_x+u_y)=epsilon*(u_xx+u_yy)
6*c4762a1bSJed Brown              u=0   at x=0, y=0
7*c4762a1bSJed Brown              u_x=0 at x=1
8*c4762a1bSJed Brown              u_y=0 at y=1
9*c4762a1bSJed Brown              u = exp(-20.0*(pow(x-0.5,2.0)+pow(y-0.5,2.0))) at t=0
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown        This program tests the routine of computing the Jacobian by the
12*c4762a1bSJed Brown        finite difference method as well as PETSc.
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown */
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown static char help[] = "Solve the convection-diffusion equation. \n\n";
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown #include <petscts.h>
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown typedef struct
21*c4762a1bSJed Brown {
22*c4762a1bSJed Brown   PetscInt  m;          /* the number of mesh points in x-direction */
23*c4762a1bSJed Brown   PetscInt  n;          /* the number of mesh points in y-direction */
24*c4762a1bSJed Brown   PetscReal dx;         /* the grid space in x-direction */
25*c4762a1bSJed Brown   PetscReal dy;         /* the grid space in y-direction */
26*c4762a1bSJed Brown   PetscReal a;          /* the convection coefficient    */
27*c4762a1bSJed Brown   PetscReal epsilon;    /* the diffusion coefficient     */
28*c4762a1bSJed Brown   PetscReal tfinal;
29*c4762a1bSJed Brown } Data;
30*c4762a1bSJed Brown 
31*c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
32*c4762a1bSJed Brown extern PetscErrorCode Initial(Vec,void*);
33*c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
34*c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
35*c4762a1bSJed Brown extern PetscErrorCode PostStep(TS);
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown int main(int argc,char **argv)
38*c4762a1bSJed Brown {
39*c4762a1bSJed Brown   PetscErrorCode ierr;
40*c4762a1bSJed Brown   PetscInt       time_steps=100,iout,NOUT=1;
41*c4762a1bSJed Brown   PetscMPIInt    size;
42*c4762a1bSJed Brown   Vec            global;
43*c4762a1bSJed Brown   PetscReal      dt,ftime,ftime_original;
44*c4762a1bSJed Brown   TS             ts;
45*c4762a1bSJed Brown   PetscViewer    viewfile;
46*c4762a1bSJed Brown   Mat            J = 0;
47*c4762a1bSJed Brown   Vec            x;
48*c4762a1bSJed Brown   Data           data;
49*c4762a1bSJed Brown   PetscInt       mn;
50*c4762a1bSJed Brown   PetscBool      flg;
51*c4762a1bSJed Brown   MatColoring    mc;
52*c4762a1bSJed Brown   ISColoring     iscoloring;
53*c4762a1bSJed Brown   MatFDColoring  matfdcoloring        = 0;
54*c4762a1bSJed Brown   PetscBool      fd_jacobian_coloring = PETSC_FALSE;
55*c4762a1bSJed Brown   SNES           snes;
56*c4762a1bSJed Brown   KSP            ksp;
57*c4762a1bSJed Brown   PC             pc;
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
60*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   /* set data */
63*c4762a1bSJed Brown   data.m       = 9;
64*c4762a1bSJed Brown   data.n       = 9;
65*c4762a1bSJed Brown   data.a       = 1.0;
66*c4762a1bSJed Brown   data.epsilon = 0.1;
67*c4762a1bSJed Brown   data.dx      = 1.0/(data.m+1.0);
68*c4762a1bSJed Brown   data.dy      = 1.0/(data.n+1.0);
69*c4762a1bSJed Brown   mn           = (data.m)*(data.n);
70*c4762a1bSJed Brown   ierr         = PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);CHKERRQ(ierr);
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown   /* set initial conditions */
73*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&global);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = VecSetSizes(global,PETSC_DECIDE,mn);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = VecSetFromOptions(global);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = Initial(global,&data);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = VecDuplicate(global,&x);CHKERRQ(ierr);
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   /* create timestep context */
80*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = TSMonitorSet(ts,Monitor,&data,NULL);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr);
83*c4762a1bSJed Brown   dt             = 0.1;
84*c4762a1bSJed Brown   ftime_original = data.tfinal = 1.0;
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,time_steps);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,ftime_original);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = TSSetSolution(ts,global);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   /* set user provided RHSFunction and RHSJacobian */
93*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&data);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,mn,mn);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(J,5,NULL);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(J,5,NULL,5,NULL);CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-ts_fd",&flg);CHKERRQ(ierr);
101*c4762a1bSJed Brown   if (!flg) {
102*c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ts,J,J,RHSJacobian,&data);CHKERRQ(ierr);
103*c4762a1bSJed Brown   } else {
104*c4762a1bSJed Brown     ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
105*c4762a1bSJed Brown     ierr = PetscOptionsHasName(NULL,NULL,"-fd_color",&fd_jacobian_coloring);CHKERRQ(ierr);
106*c4762a1bSJed Brown     if (fd_jacobian_coloring) { /* Use finite differences with coloring */
107*c4762a1bSJed Brown       /* Get data structure of J */
108*c4762a1bSJed Brown       PetscBool pc_diagonal;
109*c4762a1bSJed Brown       ierr = PetscOptionsHasName(NULL,NULL,"-pc_diagonal",&pc_diagonal);CHKERRQ(ierr);
110*c4762a1bSJed Brown       if (pc_diagonal) { /* the preconditioner of J is a diagonal matrix */
111*c4762a1bSJed Brown         PetscInt    rstart,rend,i;
112*c4762a1bSJed Brown         PetscScalar zero=0.0;
113*c4762a1bSJed Brown         ierr = MatGetOwnershipRange(J,&rstart,&rend);CHKERRQ(ierr);
114*c4762a1bSJed Brown         for (i=rstart; i<rend; i++) {
115*c4762a1bSJed Brown           ierr = MatSetValues(J,1,&i,1,&i,&zero,INSERT_VALUES);CHKERRQ(ierr);
116*c4762a1bSJed Brown         }
117*c4762a1bSJed Brown         ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
118*c4762a1bSJed Brown         ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
119*c4762a1bSJed Brown       } else {
120*c4762a1bSJed Brown         /* Fill the structure using the expensive SNESComputeJacobianDefault. Temporarily set up the TS so we can call this function */
121*c4762a1bSJed Brown         ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
122*c4762a1bSJed Brown         ierr = TSSetUp(ts);CHKERRQ(ierr);
123*c4762a1bSJed Brown         ierr = SNESComputeJacobianDefault(snes,x,J,J,ts);CHKERRQ(ierr);
124*c4762a1bSJed Brown       }
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown       /* create coloring context */
127*c4762a1bSJed Brown       ierr = MatColoringCreate(J,&mc);CHKERRQ(ierr);
128*c4762a1bSJed Brown       ierr = MatColoringSetType(mc,MATCOLORINGSL);CHKERRQ(ierr);
129*c4762a1bSJed Brown       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
130*c4762a1bSJed Brown       ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr);
131*c4762a1bSJed Brown       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
132*c4762a1bSJed Brown       ierr = MatFDColoringCreate(J,iscoloring,&matfdcoloring);CHKERRQ(ierr);
133*c4762a1bSJed Brown       ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);CHKERRQ(ierr);
134*c4762a1bSJed Brown       ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr);
135*c4762a1bSJed Brown       ierr = MatFDColoringSetUp(J,iscoloring,matfdcoloring);CHKERRQ(ierr);
136*c4762a1bSJed Brown       ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);CHKERRQ(ierr);
137*c4762a1bSJed Brown       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
138*c4762a1bSJed Brown     } else { /* Use finite differences (slow) */
139*c4762a1bSJed Brown       ierr = SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);CHKERRQ(ierr);
140*c4762a1bSJed Brown     }
141*c4762a1bSJed Brown   }
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown   /* Pick up a Petsc preconditioner */
144*c4762a1bSJed Brown   /* one can always set method or preconditioner during the run time */
145*c4762a1bSJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = PCSetType(pc,PCJACOBI);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = TSSetUp(ts);CHKERRQ(ierr);
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   /* Test TSSetPostStep() */
155*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-test_PostStep",&flg);CHKERRQ(ierr);
156*c4762a1bSJed Brown   if (flg) {
157*c4762a1bSJed Brown     ierr = TSSetPostStep(ts,PostStep);CHKERRQ(ierr);
158*c4762a1bSJed Brown   }
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-NOUT",&NOUT,NULL);CHKERRQ(ierr);
161*c4762a1bSJed Brown   for (iout=1; iout<=NOUT; iout++) {
162*c4762a1bSJed Brown     ierr = TSSetMaxSteps(ts,time_steps);CHKERRQ(ierr);
163*c4762a1bSJed Brown     ierr = TSSetMaxTime(ts,iout*ftime_original/NOUT);CHKERRQ(ierr);
164*c4762a1bSJed Brown     ierr = TSSolve(ts,global);CHKERRQ(ierr);
165*c4762a1bSJed Brown     ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
166*c4762a1bSJed Brown     ierr = TSSetTime(ts,ftime);CHKERRQ(ierr);
167*c4762a1bSJed Brown     ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
168*c4762a1bSJed Brown   }
169*c4762a1bSJed Brown   /* Interpolate solution at tfinal */
170*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&global);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = TSInterpolate(ts,ftime_original,global);CHKERRQ(ierr);
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-matlab_view",&flg);CHKERRQ(ierr);
174*c4762a1bSJed Brown   if (flg) { /* print solution into a MATLAB file */
175*c4762a1bSJed Brown     ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,"out.m",&viewfile);CHKERRQ(ierr);
176*c4762a1bSJed Brown     ierr = PetscViewerPushFormat(viewfile,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
177*c4762a1bSJed Brown     ierr = VecView(global,viewfile);CHKERRQ(ierr);
178*c4762a1bSJed Brown     ierr = PetscViewerPopFormat(viewfile);CHKERRQ(ierr);
179*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&viewfile);CHKERRQ(ierr);
180*c4762a1bSJed Brown   }
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   /* free the memories */
183*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = VecDestroy(&global);CHKERRQ(ierr);
185*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
187*c4762a1bSJed Brown   if (fd_jacobian_coloring) {ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr);}
188*c4762a1bSJed Brown   ierr = PetscFinalize();
189*c4762a1bSJed Brown   return ierr;
190*c4762a1bSJed Brown }
191*c4762a1bSJed Brown 
192*c4762a1bSJed Brown /* -------------------------------------------------------------------*/
193*c4762a1bSJed Brown /* the initial function */
194*c4762a1bSJed Brown PetscReal f_ini(PetscReal x,PetscReal y)
195*c4762a1bSJed Brown {
196*c4762a1bSJed Brown   PetscReal f;
197*c4762a1bSJed Brown 
198*c4762a1bSJed Brown   f=PetscExpReal(-20.0*(PetscPowRealInt(x-0.5,2)+PetscPowRealInt(y-0.5,2)));
199*c4762a1bSJed Brown   return f;
200*c4762a1bSJed Brown }
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown PetscErrorCode Initial(Vec global,void *ctx)
203*c4762a1bSJed Brown {
204*c4762a1bSJed Brown   Data           *data = (Data*)ctx;
205*c4762a1bSJed Brown   PetscInt       m,row,col;
206*c4762a1bSJed Brown   PetscReal      x,y,dx,dy;
207*c4762a1bSJed Brown   PetscScalar    *localptr;
208*c4762a1bSJed Brown   PetscInt       i,mybase,myend,locsize;
209*c4762a1bSJed Brown   PetscErrorCode ierr;
210*c4762a1bSJed Brown 
211*c4762a1bSJed Brown   PetscFunctionBeginUser;
212*c4762a1bSJed Brown   /* make the local  copies of parameters */
213*c4762a1bSJed Brown   m  = data->m;
214*c4762a1bSJed Brown   dx = data->dx;
215*c4762a1bSJed Brown   dy = data->dy;
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown   /* determine starting point of each processor */
218*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(global,&mybase,&myend);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = VecGetLocalSize(global,&locsize);CHKERRQ(ierr);
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown   /* Initialize the array */
222*c4762a1bSJed Brown   ierr = VecGetArray(global,&localptr);CHKERRQ(ierr);
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   for (i=0; i<locsize; i++) {
225*c4762a1bSJed Brown     row         = 1+(mybase+i)-((mybase+i)/m)*m;
226*c4762a1bSJed Brown     col         = (mybase+i)/m+1;
227*c4762a1bSJed Brown     x           = dx*row;
228*c4762a1bSJed Brown     y           = dy*col;
229*c4762a1bSJed Brown     localptr[i] = f_ini(x,y);
230*c4762a1bSJed Brown   }
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown   ierr = VecRestoreArray(global,&localptr);CHKERRQ(ierr);
233*c4762a1bSJed Brown   PetscFunctionReturn(0);
234*c4762a1bSJed Brown }
235*c4762a1bSJed Brown 
236*c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
237*c4762a1bSJed Brown {
238*c4762a1bSJed Brown   VecScatter        scatter;
239*c4762a1bSJed Brown   IS                from,to;
240*c4762a1bSJed Brown   PetscInt          i,n,*idx,nsteps,maxsteps;
241*c4762a1bSJed Brown   Vec               tmp_vec;
242*c4762a1bSJed Brown   PetscErrorCode    ierr;
243*c4762a1bSJed Brown   const PetscScalar *tmp;
244*c4762a1bSJed Brown 
245*c4762a1bSJed Brown   PetscFunctionBeginUser;
246*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&nsteps);CHKERRQ(ierr);
247*c4762a1bSJed Brown   /* display output at selected time steps */
248*c4762a1bSJed Brown   ierr = TSGetMaxSteps(ts, &maxsteps);CHKERRQ(ierr);
249*c4762a1bSJed Brown   if (nsteps % 10 != 0) PetscFunctionReturn(0);
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown   /* Get the size of the vector */
252*c4762a1bSJed Brown   ierr = VecGetSize(global,&n);CHKERRQ(ierr);
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown   /* Set the index sets */
255*c4762a1bSJed Brown   ierr = PetscMalloc1(n,&idx);CHKERRQ(ierr);
256*c4762a1bSJed Brown   for (i=0; i<n; i++) idx[i]=i;
257*c4762a1bSJed Brown 
258*c4762a1bSJed Brown   /* Create local sequential vectors */
259*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);CHKERRQ(ierr);
260*c4762a1bSJed Brown 
261*c4762a1bSJed Brown   /* Create scatter context */
262*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
263*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
264*c4762a1bSJed Brown   ierr = VecScatterCreate(global,from,tmp_vec,to,&scatter);CHKERRQ(ierr);
265*c4762a1bSJed Brown   ierr = VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
266*c4762a1bSJed Brown   ierr = VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
267*c4762a1bSJed Brown 
268*c4762a1bSJed Brown   ierr = VecGetArrayRead(tmp_vec,&tmp);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"At t[%D] =%14.2e u= %14.2e at the center \n",nsteps,(double)time,(double)PetscRealPart(tmp[n/2]));CHKERRQ(ierr);
270*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(tmp_vec,&tmp);CHKERRQ(ierr);
271*c4762a1bSJed Brown 
272*c4762a1bSJed Brown   ierr = PetscFree(idx);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = ISDestroy(&from);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = ISDestroy(&to);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr);
276*c4762a1bSJed Brown   ierr = VecDestroy(&tmp_vec);CHKERRQ(ierr);
277*c4762a1bSJed Brown   PetscFunctionReturn(0);
278*c4762a1bSJed Brown }
279*c4762a1bSJed Brown 
280*c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ptr)
281*c4762a1bSJed Brown {
282*c4762a1bSJed Brown   Data           *data = (Data*)ptr;
283*c4762a1bSJed Brown   PetscScalar    v[5];
284*c4762a1bSJed Brown   PetscInt       idx[5],i,j,row;
285*c4762a1bSJed Brown   PetscErrorCode ierr;
286*c4762a1bSJed Brown   PetscInt       m,n,mn;
287*c4762a1bSJed Brown   PetscReal      dx,dy,a,epsilon,xc,xl,xr,yl,yr;
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown   PetscFunctionBeginUser;
290*c4762a1bSJed Brown   m       = data->m;
291*c4762a1bSJed Brown   n       = data->n;
292*c4762a1bSJed Brown   mn      = m*n;
293*c4762a1bSJed Brown   dx      = data->dx;
294*c4762a1bSJed Brown   dy      = data->dy;
295*c4762a1bSJed Brown   a       = data->a;
296*c4762a1bSJed Brown   epsilon = data->epsilon;
297*c4762a1bSJed Brown 
298*c4762a1bSJed Brown   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
299*c4762a1bSJed Brown   xl = 0.5*a/dx+epsilon/(dx*dx);
300*c4762a1bSJed Brown   xr = -0.5*a/dx+epsilon/(dx*dx);
301*c4762a1bSJed Brown   yl = 0.5*a/dy+epsilon/(dy*dy);
302*c4762a1bSJed Brown   yr = -0.5*a/dy+epsilon/(dy*dy);
303*c4762a1bSJed Brown 
304*c4762a1bSJed Brown   row    = 0;
305*c4762a1bSJed Brown   v[0]   = xc;  v[1] = xr;  v[2] = yr;
306*c4762a1bSJed Brown   idx[0] = 0; idx[1] = 2; idx[2] = m;
307*c4762a1bSJed Brown   ierr   = MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
308*c4762a1bSJed Brown 
309*c4762a1bSJed Brown   row    = m-1;
310*c4762a1bSJed Brown   v[0]   = 2.0*xl; v[1] = xc;    v[2] = yr;
311*c4762a1bSJed Brown   idx[0] = m-2;  idx[1] = m-1; idx[2] = m-1+m;
312*c4762a1bSJed Brown   ierr   = MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
313*c4762a1bSJed Brown 
314*c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
315*c4762a1bSJed Brown     row    = i;
316*c4762a1bSJed Brown     v[0]   = xl;    v[1] = xc;  v[2] = xr;    v[3] = yr;
317*c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1; idx[3] = i+m;
318*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);CHKERRQ(ierr);
319*c4762a1bSJed Brown   }
320*c4762a1bSJed Brown 
321*c4762a1bSJed Brown   for (j=1; j<n-1; j++) {
322*c4762a1bSJed Brown     row    = j*m;
323*c4762a1bSJed Brown     v[0]   = xc;    v[1] = xr;    v[2] = yl;      v[3] = yr;
324*c4762a1bSJed Brown     idx[0] = j*m; idx[1] = j*m; idx[2] = j*m-m; idx[3] = j*m+m;
325*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);CHKERRQ(ierr);
326*c4762a1bSJed Brown 
327*c4762a1bSJed Brown     row    = j*m+m-1;
328*c4762a1bSJed Brown     v[0]   = xc;        v[1] = 2.0*xl;      v[2] = yl;          v[3] = yr;
329*c4762a1bSJed 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;
330*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);CHKERRQ(ierr);
331*c4762a1bSJed Brown 
332*c4762a1bSJed Brown     for (i=1; i<m-1; i++) {
333*c4762a1bSJed Brown       row    = j*m+i;
334*c4762a1bSJed Brown       v[0]   = xc;      v[1] = xl;        v[2] = xr;        v[3] = yl; v[4]=yr;
335*c4762a1bSJed Brown       idx[0] = j*m+i; idx[1] = j*m+i-1; idx[2] = j*m+i+1; idx[3] = j*m+i-m;
336*c4762a1bSJed Brown       idx[4] = j*m+i+m;
337*c4762a1bSJed Brown       ierr   = MatSetValues(A,1,&row,5,idx,v,INSERT_VALUES);CHKERRQ(ierr);
338*c4762a1bSJed Brown     }
339*c4762a1bSJed Brown   }
340*c4762a1bSJed Brown 
341*c4762a1bSJed Brown   row    = mn-m;
342*c4762a1bSJed Brown   v[0]   = xc;     v[1] = xr;       v[2] = 2.0*yl;
343*c4762a1bSJed Brown   idx[0] = mn-m; idx[1] = mn-m+1; idx[2] = mn-m-m;
344*c4762a1bSJed Brown   ierr   = MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
345*c4762a1bSJed Brown 
346*c4762a1bSJed Brown   row    = mn-1;
347*c4762a1bSJed Brown   v[0]   = xc;     v[1] = 2.0*xl; v[2] = 2.0*yl;
348*c4762a1bSJed Brown   idx[0] = mn-1; idx[1] = mn-2; idx[2] = mn-1-m;
349*c4762a1bSJed Brown   ierr   = MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
350*c4762a1bSJed Brown 
351*c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
352*c4762a1bSJed Brown     row    = mn-m+i;
353*c4762a1bSJed Brown     v[0]   = xl;         v[1] = xc;       v[2] = xr;         v[3] = 2.0*yl;
354*c4762a1bSJed Brown     idx[0] = mn-m+i-1; idx[1] = mn-m+i; idx[2] = mn-m+i+1; idx[3] = mn-m+i-m;
355*c4762a1bSJed Brown     ierr   = MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);CHKERRQ(ierr);
356*c4762a1bSJed Brown   }
357*c4762a1bSJed Brown 
358*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
359*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
360*c4762a1bSJed Brown 
361*c4762a1bSJed Brown   PetscFunctionReturn(0);
362*c4762a1bSJed Brown }
363*c4762a1bSJed Brown 
364*c4762a1bSJed Brown /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
365*c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
366*c4762a1bSJed Brown {
367*c4762a1bSJed Brown   Data              *data = (Data*)ctx;
368*c4762a1bSJed Brown   PetscInt          m,n,mn;
369*c4762a1bSJed Brown   PetscReal         dx,dy;
370*c4762a1bSJed Brown   PetscReal         xc,xl,xr,yl,yr;
371*c4762a1bSJed Brown   PetscReal         a,epsilon;
372*c4762a1bSJed Brown   PetscScalar       *outptr;
373*c4762a1bSJed Brown   const PetscScalar *inptr;
374*c4762a1bSJed Brown   PetscInt          i,j,len;
375*c4762a1bSJed Brown   PetscErrorCode    ierr;
376*c4762a1bSJed Brown   IS                from,to;
377*c4762a1bSJed Brown   PetscInt          *idx;
378*c4762a1bSJed Brown   VecScatter        scatter;
379*c4762a1bSJed Brown   Vec               tmp_in,tmp_out;
380*c4762a1bSJed Brown 
381*c4762a1bSJed Brown   PetscFunctionBeginUser;
382*c4762a1bSJed Brown   m       = data->m;
383*c4762a1bSJed Brown   n       = data->n;
384*c4762a1bSJed Brown   mn      = m*n;
385*c4762a1bSJed Brown   dx      = data->dx;
386*c4762a1bSJed Brown   dy      = data->dy;
387*c4762a1bSJed Brown   a       = data->a;
388*c4762a1bSJed Brown   epsilon = data->epsilon;
389*c4762a1bSJed Brown 
390*c4762a1bSJed Brown   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
391*c4762a1bSJed Brown   xl = 0.5*a/dx+epsilon/(dx*dx);
392*c4762a1bSJed Brown   xr = -0.5*a/dx+epsilon/(dx*dx);
393*c4762a1bSJed Brown   yl = 0.5*a/dy+epsilon/(dy*dy);
394*c4762a1bSJed Brown   yr = -0.5*a/dy+epsilon/(dy*dy);
395*c4762a1bSJed Brown 
396*c4762a1bSJed Brown   /* Get the length of parallel vector */
397*c4762a1bSJed Brown   ierr = VecGetSize(globalin,&len);CHKERRQ(ierr);
398*c4762a1bSJed Brown 
399*c4762a1bSJed Brown   /* Set the index sets */
400*c4762a1bSJed Brown   ierr = PetscMalloc1(len,&idx);CHKERRQ(ierr);
401*c4762a1bSJed Brown   for (i=0; i<len; i++) idx[i]=i;
402*c4762a1bSJed Brown 
403*c4762a1bSJed Brown   /* Create local sequential vectors */
404*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,len,&tmp_in);CHKERRQ(ierr);
405*c4762a1bSJed Brown   ierr = VecDuplicate(tmp_in,&tmp_out);CHKERRQ(ierr);
406*c4762a1bSJed Brown 
407*c4762a1bSJed Brown   /* Create scatter context */
408*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
409*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
410*c4762a1bSJed Brown   ierr = VecScatterCreate(globalin,from,tmp_in,to,&scatter);CHKERRQ(ierr);
411*c4762a1bSJed Brown   ierr = VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
412*c4762a1bSJed Brown   ierr = VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
413*c4762a1bSJed Brown   ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr);
414*c4762a1bSJed Brown 
415*c4762a1bSJed Brown   /*Extract income array - include ghost points */
416*c4762a1bSJed Brown   ierr = VecGetArrayRead(tmp_in,&inptr);CHKERRQ(ierr);
417*c4762a1bSJed Brown 
418*c4762a1bSJed Brown   /* Extract outcome array*/
419*c4762a1bSJed Brown   ierr = VecGetArray(tmp_out,&outptr);CHKERRQ(ierr);
420*c4762a1bSJed Brown 
421*c4762a1bSJed Brown   outptr[0]   = xc*inptr[0]+xr*inptr[1]+yr*inptr[m];
422*c4762a1bSJed Brown   outptr[m-1] = 2.0*xl*inptr[m-2]+xc*inptr[m-1]+yr*inptr[m-1+m];
423*c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
424*c4762a1bSJed Brown     outptr[i] = xc*inptr[i]+xl*inptr[i-1]+xr*inptr[i+1]+yr*inptr[i+m];
425*c4762a1bSJed Brown   }
426*c4762a1bSJed Brown 
427*c4762a1bSJed Brown   for (j=1; j<n-1; j++) {
428*c4762a1bSJed Brown     outptr[j*m] = xc*inptr[j*m]+xr*inptr[j*m+1]+ yl*inptr[j*m-m]+yr*inptr[j*m+m];
429*c4762a1bSJed 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];
430*c4762a1bSJed Brown     for (i=1; i<m-1; i++) {
431*c4762a1bSJed 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];
432*c4762a1bSJed Brown     }
433*c4762a1bSJed Brown   }
434*c4762a1bSJed Brown 
435*c4762a1bSJed Brown   outptr[mn-m] = xc*inptr[mn-m]+xr*inptr[mn-m+1]+2.0*yl*inptr[mn-m-m];
436*c4762a1bSJed Brown   outptr[mn-1] = 2.0*xl*inptr[mn-2]+xc*inptr[mn-1]+2.0*yl*inptr[mn-1-m];
437*c4762a1bSJed Brown   for (i=1; i<m-1; i++) {
438*c4762a1bSJed 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];
439*c4762a1bSJed Brown   }
440*c4762a1bSJed Brown 
441*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(tmp_in,&inptr);CHKERRQ(ierr);
442*c4762a1bSJed Brown   ierr = VecRestoreArray(tmp_out,&outptr);CHKERRQ(ierr);
443*c4762a1bSJed Brown 
444*c4762a1bSJed Brown   ierr = VecScatterCreate(tmp_out,from,globalout,to,&scatter);CHKERRQ(ierr);
445*c4762a1bSJed Brown   ierr = VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
446*c4762a1bSJed Brown   ierr = VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
447*c4762a1bSJed Brown 
448*c4762a1bSJed Brown   /* Destroy idx aand scatter */
449*c4762a1bSJed Brown   ierr = VecDestroy(&tmp_in);CHKERRQ(ierr);
450*c4762a1bSJed Brown   ierr = VecDestroy(&tmp_out);CHKERRQ(ierr);
451*c4762a1bSJed Brown   ierr = ISDestroy(&from);CHKERRQ(ierr);
452*c4762a1bSJed Brown   ierr = ISDestroy(&to);CHKERRQ(ierr);
453*c4762a1bSJed Brown   ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr);
454*c4762a1bSJed Brown 
455*c4762a1bSJed Brown   ierr = PetscFree(idx);CHKERRQ(ierr);
456*c4762a1bSJed Brown   PetscFunctionReturn(0);
457*c4762a1bSJed Brown }
458*c4762a1bSJed Brown 
459*c4762a1bSJed Brown PetscErrorCode PostStep(TS ts)
460*c4762a1bSJed Brown {
461*c4762a1bSJed Brown   PetscErrorCode ierr;
462*c4762a1bSJed Brown   PetscReal      t;
463*c4762a1bSJed Brown 
464*c4762a1bSJed Brown   PetscFunctionBeginUser;
465*c4762a1bSJed Brown   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
466*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"  PostStep, t: %g\n",(double)t);CHKERRQ(ierr);
467*c4762a1bSJed Brown   PetscFunctionReturn(0);
468*c4762a1bSJed Brown }
469*c4762a1bSJed Brown 
470*c4762a1bSJed Brown /*TEST
471*c4762a1bSJed Brown 
472*c4762a1bSJed Brown     test:
473*c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
474*c4762a1bSJed Brown       output_file: output/ex4.out
475*c4762a1bSJed Brown 
476*c4762a1bSJed Brown     test:
477*c4762a1bSJed Brown       suffix: 2
478*c4762a1bSJed Brown       args: -ts_fd -ts_type beuler
479*c4762a1bSJed Brown       nsize: 2
480*c4762a1bSJed Brown       output_file: output/ex4.out
481*c4762a1bSJed Brown 
482*c4762a1bSJed Brown     test:
483*c4762a1bSJed Brown       suffix: 3
484*c4762a1bSJed Brown       args: -ts_fd -ts_type cn
485*c4762a1bSJed Brown 
486*c4762a1bSJed Brown     test:
487*c4762a1bSJed Brown       suffix: 4
488*c4762a1bSJed Brown       args: -ts_fd -ts_type cn
489*c4762a1bSJed Brown       output_file: output/ex4_3.out
490*c4762a1bSJed Brown       nsize: 2
491*c4762a1bSJed Brown 
492*c4762a1bSJed Brown     test:
493*c4762a1bSJed Brown       suffix: 5
494*c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
495*c4762a1bSJed Brown       output_file: output/ex4.out
496*c4762a1bSJed Brown 
497*c4762a1bSJed Brown     test:
498*c4762a1bSJed Brown       suffix: 6
499*c4762a1bSJed Brown       args: -ts_type beuler -ts_fd -fd_color -mat_coloring_type sl
500*c4762a1bSJed Brown       output_file: output/ex4.out
501*c4762a1bSJed Brown       nsize: 2
502*c4762a1bSJed Brown 
503*c4762a1bSJed Brown     test:
504*c4762a1bSJed Brown       suffix: 7
505*c4762a1bSJed Brown       requires: !single
506*c4762a1bSJed Brown       args: -ts_fd -ts_type beuler -test_PostStep -ts_dt .1
507*c4762a1bSJed Brown 
508*c4762a1bSJed Brown     test:
509*c4762a1bSJed Brown       suffix: 8
510*c4762a1bSJed Brown       requires: !single
511*c4762a1bSJed Brown       args: -ts_type rk -ts_rk_type 5dp -ts_dt .01 -ts_adapt_type none -ts_view
512*c4762a1bSJed Brown 
513*c4762a1bSJed Brown TEST*/
514*c4762a1bSJed Brown 
515