xref: /petsc/src/ts/tests/ex12.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] ="Tests PetscObjectSetOptions() for TS object\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*
5*c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
6*c4762a1bSJed Brown    Processors: n
7*c4762a1bSJed Brown */
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown /* ------------------------------------------------------------------------
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown    This program solves the PDE
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown                u * u_xx
14*c4762a1bSJed Brown          u_t = ---------
15*c4762a1bSJed Brown                2*(t+1)^2
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown     on the domain 0 <= x <= 1, with boundary conditions
18*c4762a1bSJed Brown          u(t,0) = t + 1,  u(t,1) = 2*t + 2,
19*c4762a1bSJed Brown     and initial condition
20*c4762a1bSJed Brown          u(0,x) = 1 + x*x.
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown     The exact solution is:
23*c4762a1bSJed Brown          u(t,x) = (1 + x*x) * (1 + t)
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown     Note that since the solution is linear in time and quadratic in x,
26*c4762a1bSJed Brown     the finite difference scheme actually computes the "exact" solution.
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown     We use by default the backward Euler method.
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ------------------------------------------------------------------------- */
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown /*
33*c4762a1bSJed Brown    Include "petscts.h" to use the PETSc timestepping routines. Note that
34*c4762a1bSJed Brown    this file automatically includes "petscsys.h" and other lower-level
35*c4762a1bSJed Brown    PETSc include files.
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown    Include the "petscdmda.h" to allow us to use the distributed array data
38*c4762a1bSJed Brown    structures to manage the parallel grid.
39*c4762a1bSJed Brown */
40*c4762a1bSJed Brown #include <petscts.h>
41*c4762a1bSJed Brown #include <petscdm.h>
42*c4762a1bSJed Brown #include <petscdmda.h>
43*c4762a1bSJed Brown #include <petscdraw.h>
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown /*
46*c4762a1bSJed Brown    User-defined application context - contains data needed by the
47*c4762a1bSJed Brown    application-provided callback routines.
48*c4762a1bSJed Brown */
49*c4762a1bSJed Brown typedef struct {
50*c4762a1bSJed Brown   MPI_Comm  comm;           /* communicator */
51*c4762a1bSJed Brown   DM        da;             /* distributed array data structure */
52*c4762a1bSJed Brown   Vec       localwork;      /* local ghosted work vector */
53*c4762a1bSJed Brown   Vec       u_local;        /* local ghosted approximate solution vector */
54*c4762a1bSJed Brown   Vec       solution;       /* global exact solution vector */
55*c4762a1bSJed Brown   PetscInt  m;              /* total number of grid points */
56*c4762a1bSJed Brown   PetscReal h;              /* mesh width: h = 1/(m-1) */
57*c4762a1bSJed Brown } AppCtx;
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown /*
60*c4762a1bSJed Brown    User-defined routines, provided below.
61*c4762a1bSJed Brown */
62*c4762a1bSJed Brown extern PetscErrorCode InitialConditions(Vec,AppCtx*);
63*c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
64*c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
65*c4762a1bSJed Brown extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown int main(int argc,char **argv)
68*c4762a1bSJed Brown {
69*c4762a1bSJed Brown   AppCtx         appctx;                 /* user-defined application context */
70*c4762a1bSJed Brown   TS             ts;                     /* timestepping context */
71*c4762a1bSJed Brown   Mat            A;                      /* Jacobian matrix data structure */
72*c4762a1bSJed Brown   Vec            u;                      /* approximate solution vector */
73*c4762a1bSJed Brown   PetscInt       time_steps_max = 100;  /* default max timesteps */
74*c4762a1bSJed Brown   PetscErrorCode ierr;
75*c4762a1bSJed Brown   PetscReal      dt;
76*c4762a1bSJed Brown   PetscReal      time_total_max = 100.0; /* default max total time */
77*c4762a1bSJed Brown   PetscOptions   options,optionscopy;
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80*c4762a1bSJed Brown      Initialize program and set problem parameters
81*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown   ierr = PetscOptionsCreate(&options);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = PetscOptionsSetValue(options,"-ts_monitor","ascii");CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = PetscOptionsSetValue(options,"-snes_monitor","ascii");CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = PetscOptionsSetValue(options,"-ksp_monitor","ascii");CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   appctx.comm = PETSC_COMM_WORLD;
91*c4762a1bSJed Brown   appctx.m    = 60;
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(options,NULL,"-M",&appctx.m,NULL);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   appctx.h    = 1.0/(appctx.m-1.0);
96*c4762a1bSJed Brown 
97*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98*c4762a1bSJed Brown      Create vector data structures
99*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   /*
102*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
103*c4762a1bSJed Brown      and to set up the ghost point communication pattern.  There are M
104*c4762a1bSJed Brown      total grid values spread equally among all the processors.
105*c4762a1bSJed Brown   */
106*c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,appctx.m,1,1,NULL,&appctx.da);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = PetscObjectSetOptions((PetscObject)appctx.da,options);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = DMSetFromOptions(appctx.da);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = DMSetUp(appctx.da);CHKERRQ(ierr);
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown   /*
112*c4762a1bSJed Brown      Extract global and local vectors from DMDA; we use these to store the
113*c4762a1bSJed Brown      approximate solution.  Then duplicate these for remaining vectors that
114*c4762a1bSJed Brown      have the same types.
115*c4762a1bSJed Brown   */
116*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(appctx.da,&u);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = DMCreateLocalVector(appctx.da,&appctx.u_local);CHKERRQ(ierr);
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown   /*
120*c4762a1bSJed Brown      Create local work vector for use in evaluating right-hand-side function;
121*c4762a1bSJed Brown      create global work vector for storing exact solution.
122*c4762a1bSJed Brown   */
123*c4762a1bSJed Brown   ierr = VecDuplicate(appctx.u_local,&appctx.localwork);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = VecDuplicate(u,&appctx.solution);CHKERRQ(ierr);
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127*c4762a1bSJed Brown      Create timestepping solver context; set callback routine for
128*c4762a1bSJed Brown      right-hand-side function evaluation.
129*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = PetscObjectSetOptions((PetscObject)ts,options);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr);
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137*c4762a1bSJed Brown      For nonlinear problems, the user can provide a Jacobian evaluation
138*c4762a1bSJed Brown      routine (or use a finite differencing approximation).
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine.
141*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = PetscObjectSetOptions((PetscObject)A,options);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,appctx.m,appctx.m);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,&appctx);CHKERRQ(ierr);
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151*c4762a1bSJed Brown      Set solution vector and initial timestep
152*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   dt   = appctx.h/2.0;
155*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158*c4762a1bSJed Brown      Customize timestepping solver:
159*c4762a1bSJed Brown        - Set the solution method to be the Backward Euler method.
160*c4762a1bSJed Brown        - Set timestepping duration info
161*c4762a1bSJed Brown      Then set runtime options, which can override these defaults.
162*c4762a1bSJed Brown      For example,
163*c4762a1bSJed Brown           -ts_max_steps <maxsteps> -ts_max_time <maxtime>
164*c4762a1bSJed Brown      to override the defaults set by TSSetMaxSteps()/TSSetMaxTime().
165*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,time_steps_max);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,time_total_max);CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174*c4762a1bSJed Brown      Solve the problem
175*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176*c4762a1bSJed Brown 
177*c4762a1bSJed Brown   /*
178*c4762a1bSJed Brown      Evaluate initial conditions
179*c4762a1bSJed Brown   */
180*c4762a1bSJed Brown   ierr = InitialConditions(u,&appctx);CHKERRQ(ierr);
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   /*
183*c4762a1bSJed Brown      Run the timestepping solver
184*c4762a1bSJed Brown   */
185*c4762a1bSJed Brown   ierr = TSSolve(ts,u);CHKERRQ(ierr);
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
189*c4762a1bSJed Brown      are no longer needed.
190*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191*c4762a1bSJed Brown 
192*c4762a1bSJed Brown   ierr = PetscObjectGetOptions((PetscObject)ts,&optionscopy);CHKERRQ(ierr);
193*c4762a1bSJed Brown   if (options != optionscopy) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"PetscObjectGetOptions() failed");
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
196*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
197*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
198*c4762a1bSJed Brown   ierr = DMDestroy(&appctx.da);CHKERRQ(ierr);
199*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.localwork);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.solution);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr = VecDestroy(&appctx.u_local);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = PetscOptionsDestroy(&options);CHKERRQ(ierr);
203*c4762a1bSJed Brown 
204*c4762a1bSJed Brown   /*
205*c4762a1bSJed Brown      Always call PetscFinalize() before exiting a program.  This routine
206*c4762a1bSJed Brown        - finalizes the PETSc libraries as well as MPI
207*c4762a1bSJed Brown        - provides summary and diagnostic information if certain runtime
208*c4762a1bSJed Brown          options are chosen (e.g., -log_view).
209*c4762a1bSJed Brown   */
210*c4762a1bSJed Brown   ierr = PetscFinalize();
211*c4762a1bSJed Brown   return ierr;
212*c4762a1bSJed Brown }
213*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
214*c4762a1bSJed Brown /*
215*c4762a1bSJed Brown    InitialConditions - Computes the solution at the initial time.
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown    Input Parameters:
218*c4762a1bSJed Brown    u - uninitialized solution vector (global)
219*c4762a1bSJed Brown    appctx - user-defined application context
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown    Output Parameter:
222*c4762a1bSJed Brown    u - vector with solution at initial time (global)
223*c4762a1bSJed Brown */
224*c4762a1bSJed Brown PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
225*c4762a1bSJed Brown {
226*c4762a1bSJed Brown   PetscScalar    *u_localptr,h = appctx->h,x;
227*c4762a1bSJed Brown   PetscInt       i,mybase,myend;
228*c4762a1bSJed Brown   PetscErrorCode ierr;
229*c4762a1bSJed Brown 
230*c4762a1bSJed Brown   /*
231*c4762a1bSJed Brown      Determine starting point of each processor's range of
232*c4762a1bSJed Brown      grid values.
233*c4762a1bSJed Brown   */
234*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(u,&mybase,&myend);CHKERRQ(ierr);
235*c4762a1bSJed Brown 
236*c4762a1bSJed Brown   /*
237*c4762a1bSJed Brown     Get a pointer to vector data.
238*c4762a1bSJed Brown     - For default PETSc vectors, VecGetArray() returns a pointer to
239*c4762a1bSJed Brown       the data array.  Otherwise, the routine is implementation dependent.
240*c4762a1bSJed Brown     - You MUST call VecRestoreArray() when you no longer need access to
241*c4762a1bSJed Brown       the array.
242*c4762a1bSJed Brown     - Note that the Fortran interface to VecGetArray() differs from the
243*c4762a1bSJed Brown       C version.  See the users manual for details.
244*c4762a1bSJed Brown   */
245*c4762a1bSJed Brown   ierr = VecGetArray(u,&u_localptr);CHKERRQ(ierr);
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown   /*
248*c4762a1bSJed Brown      We initialize the solution array by simply writing the solution
249*c4762a1bSJed Brown      directly into the array locations.  Alternatively, we could use
250*c4762a1bSJed Brown      VecSetValues() or VecSetValuesLocal().
251*c4762a1bSJed Brown   */
252*c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
253*c4762a1bSJed Brown     x = h*(PetscReal)i; /* current location in global grid */
254*c4762a1bSJed Brown     u_localptr[i-mybase] = 1.0 + x*x;
255*c4762a1bSJed Brown   }
256*c4762a1bSJed Brown 
257*c4762a1bSJed Brown   /*
258*c4762a1bSJed Brown      Restore vector
259*c4762a1bSJed Brown   */
260*c4762a1bSJed Brown   ierr = VecRestoreArray(u,&u_localptr);CHKERRQ(ierr);
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   return 0;
263*c4762a1bSJed Brown }
264*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
265*c4762a1bSJed Brown /*
266*c4762a1bSJed Brown    ExactSolution - Computes the exact solution at a given time.
267*c4762a1bSJed Brown 
268*c4762a1bSJed Brown    Input Parameters:
269*c4762a1bSJed Brown    t - current time
270*c4762a1bSJed Brown    solution - vector in which exact solution will be computed
271*c4762a1bSJed Brown    appctx - user-defined application context
272*c4762a1bSJed Brown 
273*c4762a1bSJed Brown    Output Parameter:
274*c4762a1bSJed Brown    solution - vector with the newly computed exact solution
275*c4762a1bSJed Brown */
276*c4762a1bSJed Brown PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
277*c4762a1bSJed Brown {
278*c4762a1bSJed Brown   PetscScalar    *s_localptr,h = appctx->h,x;
279*c4762a1bSJed Brown   PetscInt       i,mybase,myend;
280*c4762a1bSJed Brown   PetscErrorCode ierr;
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown   /*
283*c4762a1bSJed Brown      Determine starting and ending points of each processor's
284*c4762a1bSJed Brown      range of grid values
285*c4762a1bSJed Brown   */
286*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(solution,&mybase,&myend);CHKERRQ(ierr);
287*c4762a1bSJed Brown 
288*c4762a1bSJed Brown   /*
289*c4762a1bSJed Brown      Get a pointer to vector data.
290*c4762a1bSJed Brown   */
291*c4762a1bSJed Brown   ierr = VecGetArray(solution,&s_localptr);CHKERRQ(ierr);
292*c4762a1bSJed Brown 
293*c4762a1bSJed Brown   /*
294*c4762a1bSJed Brown      Simply write the solution directly into the array locations.
295*c4762a1bSJed Brown      Alternatively, we could use VecSetValues() or VecSetValuesLocal().
296*c4762a1bSJed Brown   */
297*c4762a1bSJed Brown   for (i=mybase; i<myend; i++) {
298*c4762a1bSJed Brown     x = h*(PetscReal)i;
299*c4762a1bSJed Brown     s_localptr[i-mybase] = (t + 1.0)*(1.0 + x*x);
300*c4762a1bSJed Brown   }
301*c4762a1bSJed Brown 
302*c4762a1bSJed Brown   /*
303*c4762a1bSJed Brown      Restore vector
304*c4762a1bSJed Brown   */
305*c4762a1bSJed Brown   ierr = VecRestoreArray(solution,&s_localptr);CHKERRQ(ierr);
306*c4762a1bSJed Brown   return 0;
307*c4762a1bSJed Brown }
308*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
309*c4762a1bSJed Brown /*
310*c4762a1bSJed Brown    RHSFunction - User-provided routine that evalues the right-hand-side
311*c4762a1bSJed Brown    function of the ODE.  This routine is set in the main program by
312*c4762a1bSJed Brown    calling TSSetRHSFunction().  We compute:
313*c4762a1bSJed Brown           global_out = F(global_in)
314*c4762a1bSJed Brown 
315*c4762a1bSJed Brown    Input Parameters:
316*c4762a1bSJed Brown    ts         - timesteping context
317*c4762a1bSJed Brown    t          - current time
318*c4762a1bSJed Brown    global_in  - vector containing the current iterate
319*c4762a1bSJed Brown    ctx        - (optional) user-provided context for function evaluation.
320*c4762a1bSJed Brown                 In this case we use the appctx defined above.
321*c4762a1bSJed Brown 
322*c4762a1bSJed Brown    Output Parameter:
323*c4762a1bSJed Brown    global_out - vector containing the newly evaluated function
324*c4762a1bSJed Brown */
325*c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec global_in,Vec global_out,void *ctx)
326*c4762a1bSJed Brown {
327*c4762a1bSJed Brown   AppCtx            *appctx   = (AppCtx*) ctx;     /* user-defined application context */
328*c4762a1bSJed Brown   DM                da        = appctx->da;        /* distributed array */
329*c4762a1bSJed Brown   Vec               local_in  = appctx->u_local;   /* local ghosted input vector */
330*c4762a1bSJed Brown   Vec               localwork = appctx->localwork; /* local ghosted work vector */
331*c4762a1bSJed Brown   PetscErrorCode    ierr;
332*c4762a1bSJed Brown   PetscInt          i,localsize;
333*c4762a1bSJed Brown   PetscMPIInt       rank,size;
334*c4762a1bSJed Brown   PetscScalar       *copyptr,sc;
335*c4762a1bSJed Brown   const PetscScalar *localptr;
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
338*c4762a1bSJed Brown      Get ready for local function computations
339*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
340*c4762a1bSJed Brown   /*
341*c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
342*c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
343*c4762a1bSJed Brown      By placing code between these two statements, computations can be
344*c4762a1bSJed Brown      done while messages are in transition.
345*c4762a1bSJed Brown   */
346*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
347*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
348*c4762a1bSJed Brown 
349*c4762a1bSJed Brown   /*
350*c4762a1bSJed Brown       Access directly the values in our local INPUT work array
351*c4762a1bSJed Brown   */
352*c4762a1bSJed Brown   ierr = VecGetArrayRead(local_in,&localptr);CHKERRQ(ierr);
353*c4762a1bSJed Brown 
354*c4762a1bSJed Brown   /*
355*c4762a1bSJed Brown       Access directly the values in our local OUTPUT work array
356*c4762a1bSJed Brown   */
357*c4762a1bSJed Brown   ierr = VecGetArray(localwork,&copyptr);CHKERRQ(ierr);
358*c4762a1bSJed Brown 
359*c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
360*c4762a1bSJed Brown 
361*c4762a1bSJed Brown   /*
362*c4762a1bSJed Brown       Evaluate our function on the nodes owned by this processor
363*c4762a1bSJed Brown   */
364*c4762a1bSJed Brown   ierr = VecGetLocalSize(local_in,&localsize);CHKERRQ(ierr);
365*c4762a1bSJed Brown 
366*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
367*c4762a1bSJed Brown      Compute entries for the locally owned part
368*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
369*c4762a1bSJed Brown 
370*c4762a1bSJed Brown   /*
371*c4762a1bSJed Brown      Handle boundary conditions: This is done by using the boundary condition
372*c4762a1bSJed Brown         u(t,boundary) = g(t,boundary)
373*c4762a1bSJed Brown      for some function g. Now take the derivative with respect to t to obtain
374*c4762a1bSJed Brown         u_{t}(t,boundary) = g_{t}(t,boundary)
375*c4762a1bSJed Brown 
376*c4762a1bSJed Brown      In our case, u(t,0) = t + 1, so that u_{t}(t,0) = 1
377*c4762a1bSJed Brown              and  u(t,1) = 2t+ 2, so that u_{t}(t,1) = 2
378*c4762a1bSJed Brown   */
379*c4762a1bSJed Brown   ierr = MPI_Comm_rank(appctx->comm,&rank);CHKERRQ(ierr);
380*c4762a1bSJed Brown   ierr = MPI_Comm_size(appctx->comm,&size);CHKERRQ(ierr);
381*c4762a1bSJed Brown   if (!rank)          copyptr[0]           = 1.0;
382*c4762a1bSJed Brown   if (rank == size-1) copyptr[localsize-1] = 2.0;
383*c4762a1bSJed Brown 
384*c4762a1bSJed Brown   /*
385*c4762a1bSJed Brown      Handle the interior nodes where the PDE is replace by finite
386*c4762a1bSJed Brown      difference operators.
387*c4762a1bSJed Brown   */
388*c4762a1bSJed Brown   for (i=1; i<localsize-1; i++) copyptr[i] =  localptr[i] * sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);
389*c4762a1bSJed Brown 
390*c4762a1bSJed Brown   /*
391*c4762a1bSJed Brown      Restore vectors
392*c4762a1bSJed Brown   */
393*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(local_in,&localptr);CHKERRQ(ierr);
394*c4762a1bSJed Brown   ierr = VecRestoreArray(localwork,&copyptr);CHKERRQ(ierr);
395*c4762a1bSJed Brown 
396*c4762a1bSJed Brown   /*
397*c4762a1bSJed Brown      Insert values from the local OUTPUT vector into the global
398*c4762a1bSJed Brown      output vector
399*c4762a1bSJed Brown   */
400*c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,global_out);CHKERRQ(ierr);
401*c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,global_out);CHKERRQ(ierr);
402*c4762a1bSJed Brown 
403*c4762a1bSJed Brown   return 0;
404*c4762a1bSJed Brown }
405*c4762a1bSJed Brown /* --------------------------------------------------------------------- */
406*c4762a1bSJed Brown /*
407*c4762a1bSJed Brown    RHSJacobian - User-provided routine to compute the Jacobian of
408*c4762a1bSJed Brown    the nonlinear right-hand-side function of the ODE.
409*c4762a1bSJed Brown 
410*c4762a1bSJed Brown    Input Parameters:
411*c4762a1bSJed Brown    ts - the TS context
412*c4762a1bSJed Brown    t - current time
413*c4762a1bSJed Brown    global_in - global input vector
414*c4762a1bSJed Brown    dummy - optional user-defined context, as set by TSetRHSJacobian()
415*c4762a1bSJed Brown 
416*c4762a1bSJed Brown    Output Parameters:
417*c4762a1bSJed Brown    AA - Jacobian matrix
418*c4762a1bSJed Brown    BB - optionally different preconditioning matrix
419*c4762a1bSJed Brown    str - flag indicating matrix structure
420*c4762a1bSJed Brown 
421*c4762a1bSJed Brown   Notes:
422*c4762a1bSJed Brown   RHSJacobian computes entries for the locally owned part of the Jacobian.
423*c4762a1bSJed Brown    - Currently, all PETSc parallel matrix formats are partitioned by
424*c4762a1bSJed Brown      contiguous chunks of rows across the processors.
425*c4762a1bSJed Brown    - Each processor needs to insert only elements that it owns
426*c4762a1bSJed Brown      locally (but any non-local elements will be sent to the
427*c4762a1bSJed Brown      appropriate processor during matrix assembly).
428*c4762a1bSJed Brown    - Always specify global row and columns of matrix entries when
429*c4762a1bSJed Brown      using MatSetValues().
430*c4762a1bSJed Brown    - Here, we set all entries for a particular row at once.
431*c4762a1bSJed Brown    - Note that MatSetValues() uses 0-based row and column numbers
432*c4762a1bSJed Brown      in Fortran as well as in C.
433*c4762a1bSJed Brown */
434*c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec global_in,Mat AA,Mat BB,void *ctx)
435*c4762a1bSJed Brown {
436*c4762a1bSJed Brown   AppCtx            *appctx  = (AppCtx*)ctx;    /* user-defined application context */
437*c4762a1bSJed Brown   Vec               local_in = appctx->u_local;   /* local ghosted input vector */
438*c4762a1bSJed Brown   DM                da       = appctx->da;        /* distributed array */
439*c4762a1bSJed Brown   PetscScalar       v[3],sc;
440*c4762a1bSJed Brown   const PetscScalar *localptr;
441*c4762a1bSJed Brown   PetscErrorCode    ierr;
442*c4762a1bSJed Brown   PetscInt          i,mstart,mend,mstarts,mends,idx[3],is;
443*c4762a1bSJed Brown 
444*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
445*c4762a1bSJed Brown      Get ready for local Jacobian computations
446*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
447*c4762a1bSJed Brown   /*
448*c4762a1bSJed Brown      Scatter ghost points to local vector, using the 2-step process
449*c4762a1bSJed Brown         DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
450*c4762a1bSJed Brown      By placing code between these two statements, computations can be
451*c4762a1bSJed Brown      done while messages are in transition.
452*c4762a1bSJed Brown   */
453*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
454*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,global_in,INSERT_VALUES,local_in);CHKERRQ(ierr);
455*c4762a1bSJed Brown 
456*c4762a1bSJed Brown   /*
457*c4762a1bSJed Brown      Get pointer to vector data
458*c4762a1bSJed Brown   */
459*c4762a1bSJed Brown   ierr = VecGetArrayRead(local_in,&localptr);CHKERRQ(ierr);
460*c4762a1bSJed Brown 
461*c4762a1bSJed Brown   /*
462*c4762a1bSJed Brown      Get starting and ending locally owned rows of the matrix
463*c4762a1bSJed Brown   */
464*c4762a1bSJed Brown   ierr   = MatGetOwnershipRange(BB,&mstarts,&mends);CHKERRQ(ierr);
465*c4762a1bSJed Brown   mstart = mstarts; mend = mends;
466*c4762a1bSJed Brown 
467*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
468*c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
469*c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
470*c4762a1bSJed Brown         contiguous chunks of rows across the processors.
471*c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
472*c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
473*c4762a1bSJed Brown         appropriate processor during matrix assembly).
474*c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
475*c4762a1bSJed Brown       - We can set matrix entries either using either
476*c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues().
477*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
478*c4762a1bSJed Brown 
479*c4762a1bSJed Brown   /*
480*c4762a1bSJed Brown      Set matrix rows corresponding to boundary data
481*c4762a1bSJed Brown   */
482*c4762a1bSJed Brown   if (mstart == 0) {
483*c4762a1bSJed Brown     v[0] = 0.0;
484*c4762a1bSJed Brown     ierr = MatSetValues(BB,1,&mstart,1,&mstart,v,INSERT_VALUES);CHKERRQ(ierr);
485*c4762a1bSJed Brown     mstart++;
486*c4762a1bSJed Brown   }
487*c4762a1bSJed Brown   if (mend == appctx->m) {
488*c4762a1bSJed Brown     mend--;
489*c4762a1bSJed Brown     v[0] = 0.0;
490*c4762a1bSJed Brown     ierr = MatSetValues(BB,1,&mend,1,&mend,v,INSERT_VALUES);CHKERRQ(ierr);
491*c4762a1bSJed Brown   }
492*c4762a1bSJed Brown 
493*c4762a1bSJed Brown   /*
494*c4762a1bSJed Brown      Set matrix rows corresponding to interior data.  We construct the
495*c4762a1bSJed Brown      matrix one row at a time.
496*c4762a1bSJed Brown   */
497*c4762a1bSJed Brown   sc = 1.0/(appctx->h*appctx->h*2.0*(1.0+t)*(1.0+t));
498*c4762a1bSJed Brown   for (i=mstart; i<mend; i++) {
499*c4762a1bSJed Brown     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
500*c4762a1bSJed Brown     is     = i - mstart + 1;
501*c4762a1bSJed Brown     v[0]   = sc*localptr[is];
502*c4762a1bSJed Brown     v[1]   = sc*(localptr[is+1] + localptr[is-1] - 4.0*localptr[is]);
503*c4762a1bSJed Brown     v[2]   = sc*localptr[is];
504*c4762a1bSJed Brown     ierr   = MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES);CHKERRQ(ierr);
505*c4762a1bSJed Brown   }
506*c4762a1bSJed Brown 
507*c4762a1bSJed Brown   /*
508*c4762a1bSJed Brown      Restore vector
509*c4762a1bSJed Brown   */
510*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(local_in,&localptr);CHKERRQ(ierr);
511*c4762a1bSJed Brown 
512*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
513*c4762a1bSJed Brown      Complete the matrix assembly process and set some options
514*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
515*c4762a1bSJed Brown   /*
516*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
517*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
518*c4762a1bSJed Brown      Computations can be done while messages are in transition
519*c4762a1bSJed Brown      by placing code between these two statements.
520*c4762a1bSJed Brown   */
521*c4762a1bSJed Brown   ierr = MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
522*c4762a1bSJed Brown   ierr = MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
523*c4762a1bSJed Brown   if (BB != AA) {
524*c4762a1bSJed Brown     ierr = MatAssemblyBegin(AA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
525*c4762a1bSJed Brown     ierr = MatAssemblyEnd(AA,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
526*c4762a1bSJed Brown   }
527*c4762a1bSJed Brown 
528*c4762a1bSJed Brown   /*
529*c4762a1bSJed Brown      Set and option to indicate that we will never add a new nonzero location
530*c4762a1bSJed Brown      to the matrix. If we do, it will generate an error.
531*c4762a1bSJed Brown   */
532*c4762a1bSJed Brown   ierr = MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
533*c4762a1bSJed Brown 
534*c4762a1bSJed Brown   return 0;
535*c4762a1bSJed Brown }
536*c4762a1bSJed Brown 
537*c4762a1bSJed Brown 
538*c4762a1bSJed Brown /*TEST
539*c4762a1bSJed Brown 
540*c4762a1bSJed Brown     test:
541*c4762a1bSJed Brown       requires: !single
542*c4762a1bSJed Brown 
543*c4762a1bSJed Brown TEST*/
544*c4762a1bSJed Brown 
545*c4762a1bSJed Brown 
546