xref: /petsc/src/tao/complementarity/tutorials/blackscholes.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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:
61606c0280SSatish Balay + * - Huang and Pang, "Options Pricing and Linear Complementarity,"
62c4762a1bSJed Brown        Journal of Computational Finance, volume 2, number 3, 1998.
63606c0280SSatish 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 */
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL, "-alpha", &user.alpha, &flg));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL, "-delta", &user.delta, &flg));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL, "-es", &user.es, &flg));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL, "-expiry", &user.expiry, &flg));
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-ms", &user.ms, &flg));
166*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-mt", &user.mt, &flg));
167*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL, "-rate", &user.rate, &flg));
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL, "-sigma", &user.sigma, &flg));
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL, "-strike", &user.strike, &flg));
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /* Check that the options set are allowable (needs to be done) */
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   user.ms++;
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.ms,1,1,NULL,&user.dm));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(user.dm));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(user.dm));
177c4762a1bSJed Brown   /* Create appropriate vectors and matrices */
178c4762a1bSJed Brown 
179*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL));
180*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user.dm,&gxs,NULL,NULL,&gxm,NULL,NULL));
181c4762a1bSJed Brown 
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(user.dm,&x));
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 
190*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(gxm,&(user.Vt1)));
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(gxm,&(user.c)));
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(gxm,&(user.d)));
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   }
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x, &c));
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   */
213*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(user.dm,&J));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* The TAO code begins here */
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* Create TAO solver and set desired solution method  */
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD, &tao));
219*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOSSILS));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Set routines for constraints function and Jacobian evaluation */
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user));
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* Set the variable bounds */
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetVariableBoundsRoutine(tao,ComputeVariableBounds,(void*)&user));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* Set initial solution guess */
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&x_array));
230c4762a1bSJed Brown   for (i=0; i< xm; i++)
231c4762a1bSJed Brown     x_array[i] = user.Vt1[i-gxs+xs];
232*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&x_array));
233c4762a1bSJed Brown   /* Set data structure */
234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao, x));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* Set routines for function and Jacobian evaluation */
237*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
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 */
243*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSolve(tao));
244c4762a1bSJed Brown 
245c4762a1bSJed Brown     /* Update Vt1 with the solution */
246*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(user.dm,&localX));
247*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(user.dm,x,INSERT_VALUES,localX));
248*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(user.dm,x,INSERT_VALUES,localX));
249*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(localX,&x_array));
250c4762a1bSJed Brown     for (j = 0; j < gxm; j++) {
251c4762a1bSJed Brown       user.Vt1[j] = x_array[j];
252c4762a1bSJed Brown     }
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(x,&x_array));
254*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(user.dm,&localX));
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* Free TAO data structures */
258*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /* Free PETSc data structures */
261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
262*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&c));
263*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
264*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.dm));
265c4762a1bSJed Brown   /* Free user-defined workspace */
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.Vt1));
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.c));
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user.d));
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   PetscInt       i;
279c4762a1bSJed Brown   PetscInt       xs,xm;
280c4762a1bSJed Brown   PetscInt       ms = user->ms;
281c4762a1bSJed Brown   PetscReal      sval=0.0,*xl_array,ub= PETSC_INFINITY;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /* Set the variable bounds */
284*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(xu, ub));
285*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,NULL,NULL,&xm,NULL,NULL));
286c4762a1bSJed Brown 
287*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(xl,&xl_array));
288c4762a1bSJed Brown   for (i = 0; i < xm; i++) {
289c4762a1bSJed Brown     sval = (xs+i)*user->ds;
290c4762a1bSJed Brown     xl_array[i] = PetscMax(user->strike - sval, 0);
291c4762a1bSJed Brown   }
292*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(xl,&xl_array));
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   if (xs==0) {
295*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(xu,&xl_array));
296c4762a1bSJed Brown     xl_array[0] = PetscMax(user->strike, 0);
297*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(xu,&xl_array));
298c4762a1bSJed Brown   }
299c4762a1bSJed Brown   if (xs+xm==ms) {
300*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(xu,&xl_array));
301c4762a1bSJed Brown     xl_array[xm-1] = 0;
302*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(xu,&xl_array));
303c4762a1bSJed Brown   }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   return 0;
306c4762a1bSJed Brown }
307c4762a1bSJed Brown /* -------------------------------------------------------------------- */
308c4762a1bSJed Brown 
309c4762a1bSJed Brown /*
310c4762a1bSJed Brown     FormFunction - Evaluates gradient of f.
311c4762a1bSJed Brown 
312c4762a1bSJed Brown     Input Parameters:
313c4762a1bSJed Brown .   tao  - the Tao context
314c4762a1bSJed Brown .   X    - input vector
315c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoAppSetConstraintRoutine()
316c4762a1bSJed Brown 
317c4762a1bSJed Brown     Output Parameters:
318c4762a1bSJed Brown .   F - vector containing the newly evaluated gradient
319c4762a1bSJed Brown */
320c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec F, void *ptr)
321c4762a1bSJed Brown {
322c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
323c4762a1bSJed Brown   PetscReal      *x, *f;
324c4762a1bSJed Brown   PetscReal      *Vt1 = user->Vt1, *c = user->c, *d = user->d;
325c4762a1bSJed Brown   PetscReal      rate = user->rate;
326c4762a1bSJed Brown   PetscReal      dt = user->dt, ds = user->ds;
327c4762a1bSJed Brown   PetscInt       ms = user->ms;
328c4762a1bSJed Brown   PetscInt       i, xs,xm,gxs,gxm;
329c4762a1bSJed Brown   Vec            localX,localF;
330c4762a1bSJed Brown   PetscReal      zero=0.0;
331c4762a1bSJed Brown 
332*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(user->dm,&localX));
333*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(user->dm,&localF));
334*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
336*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(user->dm,&xs,NULL,NULL,&xm,NULL,NULL));
337*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,NULL,NULL,&gxm,NULL,NULL));
338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(F, zero));
339c4762a1bSJed Brown   /*
340c4762a1bSJed Brown      The problem size is smaller than the discretization because of the
341c4762a1bSJed Brown      two fixed elements (V(0,T) = E and V(Send,T) = 0.
342c4762a1bSJed Brown   */
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   /* Get pointers to the vector data */
345*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localX, &x));
346*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localF, &f));
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   /* Left Boundary */
349c4762a1bSJed Brown   if (gxs==0) {
350c4762a1bSJed Brown     f[0] = x[0]-user->strike;
351c4762a1bSJed Brown   } else {
352c4762a1bSJed Brown     f[0] = 0;
353c4762a1bSJed Brown   }
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   /* All points in the interior */
356c4762a1bSJed Brown   /*  for (i=gxs+1;i<gxm-1;i++) { */
357c4762a1bSJed Brown   for (i=1;i<gxm-1;i++) {
358c4762a1bSJed 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]) +
359c4762a1bSJed Brown            (d[i]/(2*ds*ds))*(x[i+1] -2*x[i] + x[i-1] + Vt1[i+1] - 2*Vt1[i] + Vt1[i-1]);
360c4762a1bSJed Brown   }
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   /* Right boundary */
363c4762a1bSJed Brown   if (gxs+gxm==ms) {
364c4762a1bSJed Brown     f[xm-1]=x[gxm-1];
365c4762a1bSJed Brown   } else {
366c4762a1bSJed Brown     f[xm-1]=0;
367c4762a1bSJed Brown   }
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   /* Restore vectors */
370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localX, &x));
371*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localF, &f));
372*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(user->dm,localF,INSERT_VALUES,F));
373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(user->dm,localF,INSERT_VALUES,F));
374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(user->dm,&localX));
375*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(user->dm,&localF));
376*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(24.0*(gxm-2)));
377c4762a1bSJed Brown   /*
378c4762a1bSJed Brown   info=VecView(F,PETSC_VIEWER_STDOUT_WORLD);
379c4762a1bSJed Brown   */
380c4762a1bSJed Brown   return 0;
381c4762a1bSJed Brown }
382c4762a1bSJed Brown 
383c4762a1bSJed Brown /* ------------------------------------------------------------------- */
384c4762a1bSJed Brown /*
385c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
386c4762a1bSJed Brown 
387c4762a1bSJed Brown    Input Parameters:
388c4762a1bSJed Brown .  tao  - the Tao context
389c4762a1bSJed Brown .  X    - input vector
390c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetJacobian()
391c4762a1bSJed Brown 
392c4762a1bSJed Brown    Output Parameters:
393c4762a1bSJed Brown .  J    - Jacobian matrix
394c4762a1bSJed Brown */
395c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao tao, Vec X, Mat J, Mat tJPre, void *ptr)
396c4762a1bSJed Brown {
397c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
398c4762a1bSJed Brown   PetscReal      *c = user->c, *d = user->d;
399c4762a1bSJed Brown   PetscReal      rate = user->rate;
400c4762a1bSJed Brown   PetscReal      dt = user->dt, ds = user->ds;
401c4762a1bSJed Brown   PetscInt       ms = user->ms;
402c4762a1bSJed Brown   PetscReal      val[3];
403c4762a1bSJed Brown   PetscInt       col[3];
404c4762a1bSJed Brown   PetscInt       i;
405c4762a1bSJed Brown   PetscInt       gxs,gxm;
406c4762a1bSJed Brown   PetscBool      assembled;
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   /* Set various matrix options */
409*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(J,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
410*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssembled(J,&assembled));
411*5f80ce2aSJacob Faibussowitsch   if (assembled) CHKERRQ(MatZeroEntries(J));
412c4762a1bSJed Brown 
413*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(user->dm,&gxs,NULL,NULL,&gxm,NULL,NULL));
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   if (gxs==0) {
416c4762a1bSJed Brown     i = 0;
417c4762a1bSJed Brown     col[0] = 0;
418c4762a1bSJed Brown     val[0]=1.0;
419*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(J,1,&i,1,col,val,INSERT_VALUES));
420c4762a1bSJed Brown   }
421c4762a1bSJed Brown   for (i=1; i < gxm-1; i++) {
422c4762a1bSJed Brown     col[0] = gxs + i - 1;
423c4762a1bSJed Brown     col[1] = gxs + i;
424c4762a1bSJed Brown     col[2] = gxs + i + 1;
425c4762a1bSJed Brown     val[0] = -c[i]/(4*ds) + d[i]/(2*ds*ds);
426c4762a1bSJed Brown     val[1] = 1.0/dt + rate - d[i]/(ds*ds);
427c4762a1bSJed Brown     val[2] =  c[i]/(4*ds) + d[i]/(2*ds*ds);
428*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(J,1,&col[1],3,col,val,INSERT_VALUES));
429c4762a1bSJed Brown   }
430c4762a1bSJed Brown   if (gxs+gxm==ms) {
431c4762a1bSJed Brown     i = ms-1;
432c4762a1bSJed Brown     col[0] = i;
433c4762a1bSJed Brown     val[0]=1.0;
434*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(J,1,&i,1,col,val,INSERT_VALUES));
435c4762a1bSJed Brown   }
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   /* Assemble the Jacobian matrix */
438*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
439*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
440*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(18.0*(gxm)+5));
441c4762a1bSJed Brown   return 0;
442c4762a1bSJed Brown }
443c4762a1bSJed Brown 
444c4762a1bSJed Brown /*TEST
445c4762a1bSJed Brown 
446c4762a1bSJed Brown    build:
447c4762a1bSJed Brown       requires: !complex
448c4762a1bSJed Brown 
449c4762a1bSJed Brown    test:
450c4762a1bSJed Brown       args: -tao_monitor -tao_type ssils -tao_gttol 1.e-5
451c4762a1bSJed Brown       requires: !single
452c4762a1bSJed Brown 
453c4762a1bSJed Brown    test:
454c4762a1bSJed Brown       suffix: 2
455c4762a1bSJed Brown       args: -tao_monitor -tao_type ssfls -tao_max_it 10 -tao_gttol 1.e-5
456c4762a1bSJed Brown       requires: !single
457c4762a1bSJed Brown 
458c4762a1bSJed Brown    test:
459c4762a1bSJed Brown       suffix: 3
460c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type subvec -tao_gttol 1.e-5
461c4762a1bSJed Brown       requires: !single
462c4762a1bSJed Brown 
463c4762a1bSJed Brown    test:
464c4762a1bSJed Brown       suffix: 4
465c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type mask -tao_gttol 1.e-5
466c4762a1bSJed Brown       requires: !single
467c4762a1bSJed Brown 
468c4762a1bSJed Brown    test:
469c4762a1bSJed Brown       suffix: 5
470c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type matrixfree -pc_type jacobi -tao_max_it 6 -tao_gttol 1.e-5
471c4762a1bSJed Brown       requires: !single
472c4762a1bSJed Brown 
473c4762a1bSJed Brown    test:
474c4762a1bSJed Brown       suffix: 6
475c4762a1bSJed Brown       args: -tao_monitor -tao_type asfls -tao_subset_type subvec -tao_max_it 10 -tao_gttol 1.e-5
476c4762a1bSJed Brown       requires: !single
477c4762a1bSJed Brown 
478c4762a1bSJed Brown    test:
479c4762a1bSJed Brown       suffix: 7
480c4762a1bSJed Brown       args: -tao_monitor -tao_type asfls -tao_subset_type mask -tao_max_it 10 -tao_gttol 1.e-5
481c4762a1bSJed Brown       requires: !single
482c4762a1bSJed Brown 
483c4762a1bSJed Brown TEST*/
484