xref: /petsc/src/tao/complementarity/tutorials/blackscholes.c (revision 606c028001d6169ccbb8dc4a5aa61aa3bda31a09)
1c4762a1bSJed Brown /**********************************************************************
2c4762a1bSJed Brown     American Put Options Pricing using the Black-Scholes Equation
3c4762a1bSJed Brown 
4c4762a1bSJed Brown    Background (European Options):
5c4762a1bSJed Brown      The standard European option is a contract where the holder has the right
6c4762a1bSJed Brown      to either buy (call option) or sell (put option) an underlying asset at
7c4762a1bSJed Brown      a designated future time and price.
8c4762a1bSJed Brown 
9c4762a1bSJed Brown      The classic Black-Scholes model begins with an assumption that the
10c4762a1bSJed Brown      price of the underlying asset behaves as a lognormal random walk.
11c4762a1bSJed Brown      Using this assumption and a no-arbitrage argument, the following
12c4762a1bSJed Brown      linear parabolic partial differential equation for the value of the
13c4762a1bSJed Brown      option results:
14c4762a1bSJed Brown 
15c4762a1bSJed Brown        dV/dt + 0.5(sigma**2)(S**alpha)(d2V/dS2) + (r - D)S(dV/dS) - rV = 0.
16c4762a1bSJed Brown 
17c4762a1bSJed Brown      Here, sigma is the volatility of the underling asset, alpha is a
18c4762a1bSJed Brown      measure of elasticity (typically two), D measures the dividend payments
19c4762a1bSJed Brown      on the underling asset, and r is the interest rate.
20c4762a1bSJed Brown 
21c4762a1bSJed Brown      To completely specify the problem, we need to impose some boundary
22c4762a1bSJed Brown      conditions.  These are as follows:
23c4762a1bSJed Brown 
24c4762a1bSJed Brown        V(S, T) = max(E - S, 0)
25c4762a1bSJed Brown        V(0, t) = E for all 0 <= t <= T
26c4762a1bSJed Brown        V(s, t) = 0 for all 0 <= t <= T and s->infinity
27c4762a1bSJed Brown 
28c4762a1bSJed Brown      where T is the exercise time time and E the strike price (price paid
29c4762a1bSJed Brown      for the contract).
30c4762a1bSJed Brown 
31c4762a1bSJed Brown      An explicit formula for the value of an European option can be
32c4762a1bSJed Brown      found.  See the references for examples.
33c4762a1bSJed Brown 
34c4762a1bSJed Brown    Background (American Options):
35c4762a1bSJed Brown      The American option is similar to its European counterpart.  The
36a5b23f4aSJose E. Roman      difference is that the holder of the American option can exercise
37c4762a1bSJed Brown      their right to buy or sell the asset at any time prior to the
38c4762a1bSJed Brown      expiration.  This additional ability introduce a free boundary into
39c4762a1bSJed Brown      the Black-Scholes equation which can be modeled as a linear
40c4762a1bSJed Brown      complementarity problem.
41c4762a1bSJed Brown 
42c4762a1bSJed Brown        0 <= -(dV/dt + 0.5(sigma**2)(S**alpha)(d2V/dS2) + (r - D)S(dV/dS) - rV)
43c4762a1bSJed Brown          complements
44c4762a1bSJed Brown        V(S,T) >= max(E-S,0)
45c4762a1bSJed Brown 
46c4762a1bSJed Brown      where the variables are the same as before and we have the same boundary
47c4762a1bSJed Brown      conditions.
48c4762a1bSJed Brown 
49c4762a1bSJed Brown      There is not explicit formula for calculating the value of an American
50c4762a1bSJed Brown      option.  Therefore, we discretize the above problem and solve the
51c4762a1bSJed Brown      resulting linear complementarity problem.
52c4762a1bSJed Brown 
53c4762a1bSJed Brown      We will use backward differences for the time variables and central
54c4762a1bSJed Brown      differences for the space variables.  Crank-Nicholson averaging will
55c4762a1bSJed Brown      also be used in the discretization.  The algorithm used by the code
56c4762a1bSJed Brown      solves for V(S,t) for a fixed t and then uses this value in the
57c4762a1bSJed Brown      calculation of V(S,t - dt).  The method stops when V(S,0) has been
58c4762a1bSJed Brown      found.
59c4762a1bSJed Brown 
60c4762a1bSJed Brown    References:
61*606c0280SSatish Balay + * - Huang and Pang, "Options Pricing and Linear Complementarity,"
62c4762a1bSJed Brown        Journal of Computational Finance, volume 2, number 3, 1998.
63*606c0280SSatish Balay - * - Wilmott, "Derivatives: The Theory and Practice of Financial Engineering,"
64c4762a1bSJed Brown        John Wiley and Sons, New York, 1998.
65c4762a1bSJed Brown ***************************************************************************/
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /*
68c4762a1bSJed Brown   Include "petsctao.h" so we can use TAO solvers.
69c4762a1bSJed Brown   Include "petscdmda.h" so that we can use distributed meshes (DMs) for managing
70c4762a1bSJed Brown   the parallel mesh.
71c4762a1bSJed Brown */
72c4762a1bSJed Brown 
73c4762a1bSJed Brown #include <petscdmda.h>
74c4762a1bSJed Brown #include <petsctao.h>
75c4762a1bSJed Brown 
76c4762a1bSJed Brown static char  help[] =
77c4762a1bSJed Brown "This example demonstrates use of the TAO package to\n\
78c4762a1bSJed Brown solve a linear complementarity problem for pricing American put options.\n\
79c4762a1bSJed Brown The code uses backward differences in time and central differences in\n\
80c4762a1bSJed Brown space.  The command line options are:\n\
81c4762a1bSJed Brown   -rate <r>, where <r> = interest rate\n\
82c4762a1bSJed Brown   -sigma <s>, where <s> = volatility of the underlying\n\
83c4762a1bSJed Brown   -alpha <a>, where <a> = elasticity of the underlying\n\
84c4762a1bSJed Brown   -delta <d>, where <d> = dividend rate\n\
85c4762a1bSJed Brown   -strike <e>, where <e> = strike price\n\
86c4762a1bSJed Brown   -expiry <t>, where <t> = the expiration date\n\
87c4762a1bSJed Brown   -mt <tg>, where <tg> = number of grid points in time\n\
88c4762a1bSJed Brown   -ms <sg>, where <sg> = number of grid points in space\n\
89c4762a1bSJed Brown   -es <se>, where <se> = ending point of the space discretization\n\n";
90c4762a1bSJed Brown 
91c4762a1bSJed Brown /*T
92c4762a1bSJed Brown    Concepts: TAO^Solving a complementarity problem
93c4762a1bSJed Brown    Routines: TaoCreate(); TaoDestroy();
94c4762a1bSJed Brown    Routines: TaoSetJacobianRoutine(); TaoAppSetConstraintRoutine();
95c4762a1bSJed Brown    Routines: TaoSetFromOptions();
96c4762a1bSJed Brown    Routines: TaoSolveApplication();
97c4762a1bSJed Brown    Routines: TaoSetVariableBoundsRoutine(); TaoSetInitialSolutionVec();
98c4762a1bSJed Brown    Processors: 1
99c4762a1bSJed Brown T*/
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /*
102c4762a1bSJed Brown   User-defined application context - contains data needed by the
103c4762a1bSJed Brown   application-provided call-back routines, FormFunction(), and FormJacobian().
104c4762a1bSJed Brown */
105c4762a1bSJed Brown 
106c4762a1bSJed Brown typedef struct {
107c4762a1bSJed Brown   PetscReal *Vt1;                /* Value of the option at time T + dt */
108c4762a1bSJed Brown   PetscReal *c;                  /* Constant -- (r - D)S */
109c4762a1bSJed Brown   PetscReal *d;                  /* Constant -- -0.5(sigma**2)(S**alpha) */
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   PetscReal rate;                /* Interest rate */
112c4762a1bSJed Brown   PetscReal sigma, alpha, delta; /* Underlying asset properties */
113c4762a1bSJed Brown   PetscReal strike, expiry;      /* Option contract properties */
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   PetscReal es;                  /* Finite value used for maximum asset value */
116c4762a1bSJed Brown   PetscReal ds, dt;              /* Discretization properties */
117c4762a1bSJed Brown   PetscInt  ms, mt;               /* Number of elements */
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   DM        dm;
120c4762a1bSJed Brown } AppCtx;
121c4762a1bSJed Brown 
122c4762a1bSJed Brown /* -------- User-defined Routines --------- */
123c4762a1bSJed Brown 
124c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
125c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *);
126c4762a1bSJed Brown PetscErrorCode ComputeVariableBounds(Tao, Vec, Vec, void*);
127c4762a1bSJed Brown 
128c4762a1bSJed Brown int main(int argc, char **argv)
129c4762a1bSJed Brown {
130c4762a1bSJed Brown   PetscErrorCode ierr;    /* used to check for functions returning nonzeros */
131c4762a1bSJed Brown   Vec            x;             /* solution vector */
132c4762a1bSJed Brown   Vec            c;             /* Constraints function vector */
133c4762a1bSJed Brown   Mat            J;                  /* jacobian matrix */
134c4762a1bSJed Brown   PetscBool      flg;         /* A return variable when checking for user options */
135c4762a1bSJed Brown   Tao            tao;          /* Tao solver context */
136c4762a1bSJed Brown   AppCtx         user;            /* user-defined work context */
137c4762a1bSJed Brown   PetscInt       i, j;
138c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm;
139c4762a1bSJed Brown   PetscReal      sval = 0;
140c4762a1bSJed Brown   PetscReal      *x_array;
141c4762a1bSJed Brown   Vec            localX;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* Initialize PETSc, TAO */
144c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, (char *)0, help);if (ierr) return ierr;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /*
147c4762a1bSJed Brown      Initialize the user-defined application context with reasonable
148c4762a1bSJed Brown      values for the American option to price
149c4762a1bSJed Brown   */
150c4762a1bSJed Brown   user.rate = 0.04;
151c4762a1bSJed Brown   user.sigma = 0.40;
152c4762a1bSJed Brown   user.alpha = 2.00;
153c4762a1bSJed Brown   user.delta = 0.01;
154c4762a1bSJed Brown   user.strike = 10.0;
155c4762a1bSJed Brown   user.expiry = 1.0;
156c4762a1bSJed Brown   user.mt = 10;
157c4762a1bSJed Brown   user.ms = 150;
158c4762a1bSJed Brown   user.es = 100.0;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* Read in alternative values for the American option to price */
161c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL, "-alpha", &user.alpha, &flg);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL, "-delta", &user.delta, &flg);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL, "-es", &user.es, &flg);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL, "-expiry", &user.expiry, &flg);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-ms", &user.ms, &flg);CHKERRQ(ierr);
166c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-mt", &user.mt, &flg);CHKERRQ(ierr);
167c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL, "-rate", &user.rate, &flg);CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL, "-sigma", &user.sigma, &flg);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL, "-strike", &user.strike, &flg);CHKERRQ(ierr);
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /* Check that the options set are allowable (needs to be done) */
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   user.ms++;
174c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.ms,1,1,NULL,&user.dm);CHKERRQ(ierr);
175c4762a1bSJed Brown   ierr = DMSetFromOptions(user.dm);CHKERRQ(ierr);
176c4762a1bSJed Brown   ierr = DMSetUp(user.dm);CHKERRQ(ierr);
177c4762a1bSJed Brown   /* Create appropriate vectors and matrices */
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   ierr = DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
180c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user.dm,&gxs,NULL,NULL,&gxm,NULL,NULL);CHKERRQ(ierr);
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   ierr = DMCreateGlobalVector(user.dm,&x);CHKERRQ(ierr);
183c4762a1bSJed Brown   /*
184c4762a1bSJed Brown      Finish filling in the user-defined context with the values for
185c4762a1bSJed Brown      dS, dt, and allocating space for the constants
186c4762a1bSJed Brown   */
187c4762a1bSJed Brown   user.ds = user.es / (user.ms-1);
188c4762a1bSJed Brown   user.dt = user.expiry / user.mt;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   ierr = PetscMalloc1(gxm,&(user.Vt1));CHKERRQ(ierr);
191c4762a1bSJed Brown   ierr = PetscMalloc1(gxm,&(user.c));CHKERRQ(ierr);
192c4762a1bSJed Brown   ierr = PetscMalloc1(gxm,&(user.d));CHKERRQ(ierr);
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /*
195c4762a1bSJed Brown      Calculate the values for the constant.  Vt1 begins with the ending
196c4762a1bSJed Brown      boundary condition.
197c4762a1bSJed Brown   */
198c4762a1bSJed Brown   for (i=0; i<gxm; i++) {
199c4762a1bSJed Brown     sval = (gxs+i)*user.ds;
200c4762a1bSJed Brown     user.Vt1[i] = PetscMax(user.strike - sval, 0);
201c4762a1bSJed Brown     user.c[i] = (user.delta - user.rate)*sval;
202c4762a1bSJed Brown     user.d[i] = -0.5*user.sigma*user.sigma*PetscPowReal(sval, user.alpha);
203c4762a1bSJed Brown   }
204c4762a1bSJed Brown   if (gxs+gxm==user.ms) {
205c4762a1bSJed Brown     user.Vt1[gxm-1] = 0;
206c4762a1bSJed Brown   }
207c4762a1bSJed Brown   ierr = VecDuplicate(x, &c);CHKERRQ(ierr);
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /*
210c4762a1bSJed Brown      Allocate the matrix used by TAO for the Jacobian.  Each row of
211c4762a1bSJed Brown      the Jacobian matrix will have at most three elements.
212c4762a1bSJed Brown   */
213c4762a1bSJed Brown   ierr = DMCreateMatrix(user.dm,&J);CHKERRQ(ierr);
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* The TAO code begins here */
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* Create TAO solver and set desired solution method  */
218c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD, &tao);CHKERRQ(ierr);
219c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOSSILS);CHKERRQ(ierr);
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Set routines for constraints function and Jacobian evaluation */
222c4762a1bSJed Brown   ierr = TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user);CHKERRQ(ierr);
223c4762a1bSJed Brown   ierr = TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user);CHKERRQ(ierr);
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* Set the variable bounds */
226c4762a1bSJed Brown   ierr = TaoSetVariableBoundsRoutine(tao,ComputeVariableBounds,(void*)&user);CHKERRQ(ierr);
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* Set initial solution guess */
229c4762a1bSJed Brown   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
230c4762a1bSJed Brown   for (i=0; i< xm; i++)
231c4762a1bSJed Brown     x_array[i] = user.Vt1[i-gxs+xs];
232c4762a1bSJed Brown   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
233c4762a1bSJed Brown   /* Set data structure */
234a82e8c82SStefano Zampini   ierr = TaoSetSolution(tao, x);CHKERRQ(ierr);
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* Set routines for function and Jacobian evaluation */
237c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /* Iteratively solve the linear complementarity problems  */
240c4762a1bSJed Brown   for (i = 1; i < user.mt; i++) {
241c4762a1bSJed Brown 
242c4762a1bSJed Brown     /* Solve the current version */
243c4762a1bSJed Brown     ierr = TaoSolve(tao);CHKERRQ(ierr);
244c4762a1bSJed Brown 
245c4762a1bSJed Brown     /* Update Vt1 with the solution */
246c4762a1bSJed Brown     ierr = DMGetLocalVector(user.dm,&localX);CHKERRQ(ierr);
247c4762a1bSJed Brown     ierr = DMGlobalToLocalBegin(user.dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
248c4762a1bSJed Brown     ierr = DMGlobalToLocalEnd(user.dm,x,INSERT_VALUES,localX);CHKERRQ(ierr);
249c4762a1bSJed Brown     ierr = VecGetArray(localX,&x_array);CHKERRQ(ierr);
250c4762a1bSJed Brown     for (j = 0; j < gxm; j++) {
251c4762a1bSJed Brown       user.Vt1[j] = x_array[j];
252c4762a1bSJed Brown     }
253c4762a1bSJed Brown     ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
254c4762a1bSJed Brown     ierr = DMRestoreLocalVector(user.dm,&localX);CHKERRQ(ierr);
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* Free TAO data structures */
258c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /* Free PETSc data structures */
261c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr = VecDestroy(&c);CHKERRQ(ierr);
263c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
265c4762a1bSJed Brown   /* Free user-defined workspace */
266c4762a1bSJed Brown   ierr = PetscFree(user.Vt1);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = PetscFree(user.c);CHKERRQ(ierr);
268c4762a1bSJed Brown   ierr = PetscFree(user.d);CHKERRQ(ierr);
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   ierr = PetscFinalize();
271c4762a1bSJed Brown   return ierr;
272c4762a1bSJed Brown }
273c4762a1bSJed Brown 
274c4762a1bSJed Brown /* -------------------------------------------------------------------- */
275c4762a1bSJed Brown PetscErrorCode ComputeVariableBounds(Tao tao, Vec xl, Vec xu, void*ctx)
276c4762a1bSJed Brown {
277c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ctx;
278c4762a1bSJed Brown   PetscErrorCode ierr;
279c4762a1bSJed Brown   PetscInt       i;
280c4762a1bSJed Brown   PetscInt       xs,xm;
281c4762a1bSJed Brown   PetscInt       ms = user->ms;
282c4762a1bSJed Brown   PetscReal      sval=0.0,*xl_array,ub= PETSC_INFINITY;
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /* Set the variable bounds */
285c4762a1bSJed Brown   ierr = VecSet(xu, ub);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   ierr = VecGetArray(xl,&xl_array);CHKERRQ(ierr);
289c4762a1bSJed Brown   for (i = 0; i < xm; i++) {
290c4762a1bSJed Brown     sval = (xs+i)*user->ds;
291c4762a1bSJed Brown     xl_array[i] = PetscMax(user->strike - sval, 0);
292c4762a1bSJed Brown   }
293c4762a1bSJed Brown   ierr = VecRestoreArray(xl,&xl_array);CHKERRQ(ierr);
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   if (xs==0) {
296c4762a1bSJed Brown     ierr = VecGetArray(xu,&xl_array);CHKERRQ(ierr);
297c4762a1bSJed Brown     xl_array[0] = PetscMax(user->strike, 0);
298c4762a1bSJed Brown     ierr = VecRestoreArray(xu,&xl_array);CHKERRQ(ierr);
299c4762a1bSJed Brown   }
300c4762a1bSJed Brown   if (xs+xm==ms) {
301c4762a1bSJed Brown     ierr = VecGetArray(xu,&xl_array);CHKERRQ(ierr);
302c4762a1bSJed Brown     xl_array[xm-1] = 0;
303c4762a1bSJed Brown     ierr = VecRestoreArray(xu,&xl_array);CHKERRQ(ierr);
304c4762a1bSJed Brown   }
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   return 0;
307c4762a1bSJed Brown }
308c4762a1bSJed Brown /* -------------------------------------------------------------------- */
309c4762a1bSJed Brown 
310c4762a1bSJed Brown /*
311c4762a1bSJed Brown     FormFunction - Evaluates gradient of f.
312c4762a1bSJed Brown 
313c4762a1bSJed Brown     Input Parameters:
314c4762a1bSJed Brown .   tao  - the Tao context
315c4762a1bSJed Brown .   X    - input vector
316c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoAppSetConstraintRoutine()
317c4762a1bSJed Brown 
318c4762a1bSJed Brown     Output Parameters:
319c4762a1bSJed Brown .   F - vector containing the newly evaluated gradient
320c4762a1bSJed Brown */
321c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec F, void *ptr)
322c4762a1bSJed Brown {
323c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
324c4762a1bSJed Brown   PetscReal      *x, *f;
325c4762a1bSJed Brown   PetscReal      *Vt1 = user->Vt1, *c = user->c, *d = user->d;
326c4762a1bSJed Brown   PetscReal      rate = user->rate;
327c4762a1bSJed Brown   PetscReal      dt = user->dt, ds = user->ds;
328c4762a1bSJed Brown   PetscInt       ms = user->ms;
329c4762a1bSJed Brown   PetscErrorCode ierr;
330c4762a1bSJed Brown   PetscInt       i, xs,xm,gxs,gxm;
331c4762a1bSJed Brown   Vec            localX,localF;
332c4762a1bSJed Brown   PetscReal      zero=0.0;
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   ierr = DMGetLocalVector(user->dm,&localX);CHKERRQ(ierr);
335c4762a1bSJed Brown   ierr = DMGetLocalVector(user->dm,&localF);CHKERRQ(ierr);
336c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
337c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
338c4762a1bSJed Brown   ierr = DMDAGetCorners(user->dm,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
339c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,NULL,NULL,&gxm,NULL,NULL);CHKERRQ(ierr);
340c4762a1bSJed Brown   ierr = VecSet(F, zero);CHKERRQ(ierr);
341c4762a1bSJed Brown   /*
342c4762a1bSJed Brown      The problem size is smaller than the discretization because of the
343c4762a1bSJed Brown      two fixed elements (V(0,T) = E and V(Send,T) = 0.
344c4762a1bSJed Brown   */
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   /* Get pointers to the vector data */
347c4762a1bSJed Brown   ierr = VecGetArray(localX, &x);CHKERRQ(ierr);
348c4762a1bSJed Brown   ierr = VecGetArray(localF, &f);CHKERRQ(ierr);
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   /* Left Boundary */
351c4762a1bSJed Brown   if (gxs==0) {
352c4762a1bSJed Brown     f[0] = x[0]-user->strike;
353c4762a1bSJed Brown   } else {
354c4762a1bSJed Brown     f[0] = 0;
355c4762a1bSJed Brown   }
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /* All points in the interior */
358c4762a1bSJed Brown   /*  for (i=gxs+1;i<gxm-1;i++) { */
359c4762a1bSJed Brown   for (i=1;i<gxm-1;i++) {
360c4762a1bSJed Brown     f[i] = (1.0/dt + rate)*x[i] - Vt1[i]/dt + (c[i]/(4*ds))*(x[i+1] - x[i-1] + Vt1[i+1] - Vt1[i-1]) +
361c4762a1bSJed Brown            (d[i]/(2*ds*ds))*(x[i+1] -2*x[i] + x[i-1] + Vt1[i+1] - 2*Vt1[i] + Vt1[i-1]);
362c4762a1bSJed Brown   }
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /* Right boundary */
365c4762a1bSJed Brown   if (gxs+gxm==ms) {
366c4762a1bSJed Brown     f[xm-1]=x[gxm-1];
367c4762a1bSJed Brown   } else {
368c4762a1bSJed Brown     f[xm-1]=0;
369c4762a1bSJed Brown   }
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   /* Restore vectors */
372c4762a1bSJed Brown   ierr = VecRestoreArray(localX, &x);CHKERRQ(ierr);
373c4762a1bSJed Brown   ierr = VecRestoreArray(localF, &f);CHKERRQ(ierr);
374c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(user->dm,localF,INSERT_VALUES,F);CHKERRQ(ierr);
375c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(user->dm,localF,INSERT_VALUES,F);CHKERRQ(ierr);
376c4762a1bSJed Brown   ierr = DMRestoreLocalVector(user->dm,&localX);CHKERRQ(ierr);
377c4762a1bSJed Brown   ierr = DMRestoreLocalVector(user->dm,&localF);CHKERRQ(ierr);
378ca0c957dSBarry Smith   ierr = PetscLogFlops(24.0*(gxm-2));CHKERRQ(ierr);
379c4762a1bSJed Brown   /*
380c4762a1bSJed Brown   info=VecView(F,PETSC_VIEWER_STDOUT_WORLD);
381c4762a1bSJed Brown   */
382c4762a1bSJed Brown   return 0;
383c4762a1bSJed Brown }
384c4762a1bSJed Brown 
385c4762a1bSJed Brown /* ------------------------------------------------------------------- */
386c4762a1bSJed Brown /*
387c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
388c4762a1bSJed Brown 
389c4762a1bSJed Brown    Input Parameters:
390c4762a1bSJed Brown .  tao  - the Tao context
391c4762a1bSJed Brown .  X    - input vector
392c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetJacobian()
393c4762a1bSJed Brown 
394c4762a1bSJed Brown    Output Parameters:
395c4762a1bSJed Brown .  J    - Jacobian matrix
396c4762a1bSJed Brown */
397c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao tao, Vec X, Mat J, Mat tJPre, void *ptr)
398c4762a1bSJed Brown {
399c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
400c4762a1bSJed Brown   PetscReal      *c = user->c, *d = user->d;
401c4762a1bSJed Brown   PetscReal      rate = user->rate;
402c4762a1bSJed Brown   PetscReal      dt = user->dt, ds = user->ds;
403c4762a1bSJed Brown   PetscInt       ms = user->ms;
404c4762a1bSJed Brown   PetscReal      val[3];
405c4762a1bSJed Brown   PetscErrorCode ierr;
406c4762a1bSJed Brown   PetscInt       col[3];
407c4762a1bSJed Brown   PetscInt       i;
408c4762a1bSJed Brown   PetscInt       gxs,gxm;
409c4762a1bSJed Brown   PetscBool      assembled;
410c4762a1bSJed Brown 
411c4762a1bSJed Brown   /* Set various matrix options */
412c4762a1bSJed Brown   ierr = MatSetOption(J,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
413c4762a1bSJed Brown   ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
414c4762a1bSJed Brown   if (assembled) {ierr = MatZeroEntries(J);CHKERRQ(ierr);}
415c4762a1bSJed Brown 
416c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(user->dm,&gxs,NULL,NULL,&gxm,NULL,NULL);CHKERRQ(ierr);
417c4762a1bSJed Brown 
418c4762a1bSJed Brown   if (gxs==0) {
419c4762a1bSJed Brown     i = 0;
420c4762a1bSJed Brown     col[0] = 0;
421c4762a1bSJed Brown     val[0]=1.0;
422c4762a1bSJed Brown     ierr = MatSetValues(J,1,&i,1,col,val,INSERT_VALUES);CHKERRQ(ierr);
423c4762a1bSJed Brown   }
424c4762a1bSJed Brown   for (i=1; i < gxm-1; i++) {
425c4762a1bSJed Brown     col[0] = gxs + i - 1;
426c4762a1bSJed Brown     col[1] = gxs + i;
427c4762a1bSJed Brown     col[2] = gxs + i + 1;
428c4762a1bSJed Brown     val[0] = -c[i]/(4*ds) + d[i]/(2*ds*ds);
429c4762a1bSJed Brown     val[1] = 1.0/dt + rate - d[i]/(ds*ds);
430c4762a1bSJed Brown     val[2] =  c[i]/(4*ds) + d[i]/(2*ds*ds);
431c4762a1bSJed Brown     ierr = MatSetValues(J,1,&col[1],3,col,val,INSERT_VALUES);CHKERRQ(ierr);
432c4762a1bSJed Brown   }
433c4762a1bSJed Brown   if (gxs+gxm==ms) {
434c4762a1bSJed Brown     i = ms-1;
435c4762a1bSJed Brown     col[0] = i;
436c4762a1bSJed Brown     val[0]=1.0;
437c4762a1bSJed Brown     ierr = MatSetValues(J,1,&i,1,col,val,INSERT_VALUES);CHKERRQ(ierr);
438c4762a1bSJed Brown   }
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   /* Assemble the Jacobian matrix */
441c4762a1bSJed Brown   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
442c4762a1bSJed Brown   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
443ca0c957dSBarry Smith   ierr = PetscLogFlops(18.0*(gxm)+5);CHKERRQ(ierr);
444c4762a1bSJed Brown   return 0;
445c4762a1bSJed Brown }
446c4762a1bSJed Brown 
447c4762a1bSJed Brown /*TEST
448c4762a1bSJed Brown 
449c4762a1bSJed Brown    build:
450c4762a1bSJed Brown       requires: !complex
451c4762a1bSJed Brown 
452c4762a1bSJed Brown    test:
453c4762a1bSJed Brown       args: -tao_monitor -tao_type ssils -tao_gttol 1.e-5
454c4762a1bSJed Brown       requires: !single
455c4762a1bSJed Brown 
456c4762a1bSJed Brown    test:
457c4762a1bSJed Brown       suffix: 2
458c4762a1bSJed Brown       args: -tao_monitor -tao_type ssfls -tao_max_it 10 -tao_gttol 1.e-5
459c4762a1bSJed Brown       requires: !single
460c4762a1bSJed Brown 
461c4762a1bSJed Brown    test:
462c4762a1bSJed Brown       suffix: 3
463c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type subvec -tao_gttol 1.e-5
464c4762a1bSJed Brown       requires: !single
465c4762a1bSJed Brown 
466c4762a1bSJed Brown    test:
467c4762a1bSJed Brown       suffix: 4
468c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type mask -tao_gttol 1.e-5
469c4762a1bSJed Brown       requires: !single
470c4762a1bSJed Brown 
471c4762a1bSJed Brown    test:
472c4762a1bSJed Brown       suffix: 5
473c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type matrixfree -pc_type jacobi -tao_max_it 6 -tao_gttol 1.e-5
474c4762a1bSJed Brown       requires: !single
475c4762a1bSJed Brown 
476c4762a1bSJed Brown    test:
477c4762a1bSJed Brown       suffix: 6
478c4762a1bSJed Brown       args: -tao_monitor -tao_type asfls -tao_subset_type subvec -tao_max_it 10 -tao_gttol 1.e-5
479c4762a1bSJed Brown       requires: !single
480c4762a1bSJed Brown 
481c4762a1bSJed Brown    test:
482c4762a1bSJed Brown       suffix: 7
483c4762a1bSJed Brown       args: -tao_monitor -tao_type asfls -tao_subset_type mask -tao_max_it 10 -tao_gttol 1.e-5
484c4762a1bSJed Brown       requires: !single
485c4762a1bSJed Brown 
486c4762a1bSJed Brown TEST*/
487