xref: /petsc/src/tao/complementarity/tutorials/blackscholes.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 /*
92c4762a1bSJed Brown   User-defined application context - contains data needed by the
93c4762a1bSJed Brown   application-provided call-back routines, FormFunction(), and FormJacobian().
94c4762a1bSJed Brown */
95c4762a1bSJed Brown 
96c4762a1bSJed Brown typedef struct {
97c4762a1bSJed Brown   PetscReal *Vt1;                /* Value of the option at time T + dt */
98c4762a1bSJed Brown   PetscReal *c;                  /* Constant -- (r - D)S */
99c4762a1bSJed Brown   PetscReal *d;                  /* Constant -- -0.5(sigma**2)(S**alpha) */
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   PetscReal rate;                /* Interest rate */
102c4762a1bSJed Brown   PetscReal sigma, alpha, delta; /* Underlying asset properties */
103c4762a1bSJed Brown   PetscReal strike, expiry;      /* Option contract properties */
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   PetscReal es;                  /* Finite value used for maximum asset value */
106c4762a1bSJed Brown   PetscReal ds, dt;              /* Discretization properties */
107c4762a1bSJed Brown   PetscInt  ms, mt;               /* Number of elements */
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   DM        dm;
110c4762a1bSJed Brown } AppCtx;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown /* -------- User-defined Routines --------- */
113c4762a1bSJed Brown 
114c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
115c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *);
116c4762a1bSJed Brown PetscErrorCode ComputeVariableBounds(Tao, Vec, Vec, void*);
117c4762a1bSJed Brown 
118c4762a1bSJed Brown int main(int argc, char **argv)
119c4762a1bSJed Brown {
120c4762a1bSJed Brown   Vec            x;             /* solution vector */
121c4762a1bSJed Brown   Vec            c;             /* Constraints function vector */
122c4762a1bSJed Brown   Mat            J;                  /* jacobian matrix */
123c4762a1bSJed Brown   PetscBool      flg;         /* A return variable when checking for user options */
124c4762a1bSJed Brown   Tao            tao;          /* Tao solver context */
125c4762a1bSJed Brown   AppCtx         user;            /* user-defined work context */
126c4762a1bSJed Brown   PetscInt       i, j;
127c4762a1bSJed Brown   PetscInt       xs,xm,gxs,gxm;
128c4762a1bSJed Brown   PetscReal      sval = 0;
129c4762a1bSJed Brown   PetscReal      *x_array;
130c4762a1bSJed Brown   Vec            localX;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* Initialize PETSc, TAO */
133*327415f7SBarry Smith   PetscFunctionBeginUser;
1349566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /*
137c4762a1bSJed Brown      Initialize the user-defined application context with reasonable
138c4762a1bSJed Brown      values for the American option to price
139c4762a1bSJed Brown   */
140c4762a1bSJed Brown   user.rate = 0.04;
141c4762a1bSJed Brown   user.sigma = 0.40;
142c4762a1bSJed Brown   user.alpha = 2.00;
143c4762a1bSJed Brown   user.delta = 0.01;
144c4762a1bSJed Brown   user.strike = 10.0;
145c4762a1bSJed Brown   user.expiry = 1.0;
146c4762a1bSJed Brown   user.mt = 10;
147c4762a1bSJed Brown   user.ms = 150;
148c4762a1bSJed Brown   user.es = 100.0;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /* Read in alternative values for the American option to price */
1519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL, "-alpha", &user.alpha, &flg));
1529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL, "-delta", &user.delta, &flg));
1539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL, "-es", &user.es, &flg));
1549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL, "-expiry", &user.expiry, &flg));
1559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-ms", &user.ms, &flg));
1569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-mt", &user.mt, &flg));
1579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL, "-rate", &user.rate, &flg));
1589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL, "-sigma", &user.sigma, &flg));
1599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL, "-strike", &user.strike, &flg));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* Check that the options set are allowable (needs to be done) */
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   user.ms++;
1649566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.ms,1,1,NULL,&user.dm));
1659566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dm));
1669566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dm));
167c4762a1bSJed Brown   /* Create appropriate vectors and matrices */
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user.dm,&xs,NULL,NULL,&xm,NULL,NULL));
1709566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user.dm,&gxs,NULL,NULL,&gxm,NULL,NULL));
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dm,&x));
173c4762a1bSJed Brown   /*
174c4762a1bSJed Brown      Finish filling in the user-defined context with the values for
175c4762a1bSJed Brown      dS, dt, and allocating space for the constants
176c4762a1bSJed Brown   */
177c4762a1bSJed Brown   user.ds = user.es / (user.ms-1);
178c4762a1bSJed Brown   user.dt = user.expiry / user.mt;
179c4762a1bSJed Brown 
1809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(gxm,&(user.Vt1)));
1819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(gxm,&(user.c)));
1829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(gxm,&(user.d)));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /*
185c4762a1bSJed Brown      Calculate the values for the constant.  Vt1 begins with the ending
186c4762a1bSJed Brown      boundary condition.
187c4762a1bSJed Brown   */
188c4762a1bSJed Brown   for (i=0; i<gxm; i++) {
189c4762a1bSJed Brown     sval = (gxs+i)*user.ds;
190c4762a1bSJed Brown     user.Vt1[i] = PetscMax(user.strike - sval, 0);
191c4762a1bSJed Brown     user.c[i] = (user.delta - user.rate)*sval;
192c4762a1bSJed Brown     user.d[i] = -0.5*user.sigma*user.sigma*PetscPowReal(sval, user.alpha);
193c4762a1bSJed Brown   }
194c4762a1bSJed Brown   if (gxs+gxm==user.ms) {
195c4762a1bSJed Brown     user.Vt1[gxm-1] = 0;
196c4762a1bSJed Brown   }
1979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &c));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /*
200c4762a1bSJed Brown      Allocate the matrix used by TAO for the Jacobian.  Each row of
201c4762a1bSJed Brown      the Jacobian matrix will have at most three elements.
202c4762a1bSJed Brown   */
2039566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(user.dm,&J));
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /* The TAO code begins here */
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /* Create TAO solver and set desired solution method  */
2089566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
2099566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOSSILS));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* Set routines for constraints function and Jacobian evaluation */
2129566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user));
2139566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user));
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* Set the variable bounds */
2169566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBoundsRoutine(tao,ComputeVariableBounds,(void*)&user));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* Set initial solution guess */
2199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&x_array));
220c4762a1bSJed Brown   for (i=0; i< xm; i++)
221c4762a1bSJed Brown     x_array[i] = user.Vt1[i-gxs+xs];
2229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&x_array));
223c4762a1bSJed Brown   /* Set data structure */
2249566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, x));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* Set routines for function and Jacobian evaluation */
2279566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* Iteratively solve the linear complementarity problems  */
230c4762a1bSJed Brown   for (i = 1; i < user.mt; i++) {
231c4762a1bSJed Brown 
232c4762a1bSJed Brown     /* Solve the current version */
2339566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown     /* Update Vt1 with the solution */
2369566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(user.dm,&localX));
2379566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(user.dm,x,INSERT_VALUES,localX));
2389566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(user.dm,x,INSERT_VALUES,localX));
2399566063dSJacob Faibussowitsch     PetscCall(VecGetArray(localX,&x_array));
240c4762a1bSJed Brown     for (j = 0; j < gxm; j++) {
241c4762a1bSJed Brown       user.Vt1[j] = x_array[j];
242c4762a1bSJed Brown     }
2439566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x,&x_array));
2449566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(user.dm,&localX));
245c4762a1bSJed Brown   }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /* Free TAO data structures */
2489566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   /* Free PETSc data structures */
2519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&c));
2539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
2549566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dm));
255c4762a1bSJed Brown   /* Free user-defined workspace */
2569566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.Vt1));
2579566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.c));
2589566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.d));
259c4762a1bSJed Brown 
2609566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
261b122ec5aSJacob Faibussowitsch   return 0;
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
264c4762a1bSJed Brown /* -------------------------------------------------------------------- */
265c4762a1bSJed Brown PetscErrorCode ComputeVariableBounds(Tao tao, Vec xl, Vec xu, void*ctx)
266c4762a1bSJed Brown {
267c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ctx;
268c4762a1bSJed Brown   PetscInt       i;
269c4762a1bSJed Brown   PetscInt       xs,xm;
270c4762a1bSJed Brown   PetscInt       ms = user->ms;
271c4762a1bSJed Brown   PetscReal      sval=0.0,*xl_array,ub= PETSC_INFINITY;
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   /* Set the variable bounds */
2749566063dSJacob Faibussowitsch   PetscCall(VecSet(xu, ub));
2759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,NULL,NULL,&xm,NULL,NULL));
276c4762a1bSJed Brown 
2779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xl,&xl_array));
278c4762a1bSJed Brown   for (i = 0; i < xm; i++) {
279c4762a1bSJed Brown     sval = (xs+i)*user->ds;
280c4762a1bSJed Brown     xl_array[i] = PetscMax(user->strike - sval, 0);
281c4762a1bSJed Brown   }
2829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xl,&xl_array));
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   if (xs==0) {
2859566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xu,&xl_array));
286c4762a1bSJed Brown     xl_array[0] = PetscMax(user->strike, 0);
2879566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xu,&xl_array));
288c4762a1bSJed Brown   }
289c4762a1bSJed Brown   if (xs+xm==ms) {
2909566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xu,&xl_array));
291c4762a1bSJed Brown     xl_array[xm-1] = 0;
2929566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xu,&xl_array));
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   return 0;
296c4762a1bSJed Brown }
297c4762a1bSJed Brown /* -------------------------------------------------------------------- */
298c4762a1bSJed Brown 
299c4762a1bSJed Brown /*
300c4762a1bSJed Brown     FormFunction - Evaluates gradient of f.
301c4762a1bSJed Brown 
302c4762a1bSJed Brown     Input Parameters:
303c4762a1bSJed Brown .   tao  - the Tao context
304c4762a1bSJed Brown .   X    - input vector
305c4762a1bSJed Brown .   ptr  - optional user-defined context, as set by TaoAppSetConstraintRoutine()
306c4762a1bSJed Brown 
307c4762a1bSJed Brown     Output Parameters:
308c4762a1bSJed Brown .   F - vector containing the newly evaluated gradient
309c4762a1bSJed Brown */
310c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao tao, Vec X, Vec F, void *ptr)
311c4762a1bSJed Brown {
312c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
313c4762a1bSJed Brown   PetscReal      *x, *f;
314c4762a1bSJed Brown   PetscReal      *Vt1 = user->Vt1, *c = user->c, *d = user->d;
315c4762a1bSJed Brown   PetscReal      rate = user->rate;
316c4762a1bSJed Brown   PetscReal      dt = user->dt, ds = user->ds;
317c4762a1bSJed Brown   PetscInt       ms = user->ms;
318c4762a1bSJed Brown   PetscInt       i, xs,xm,gxs,gxm;
319c4762a1bSJed Brown   Vec            localX,localF;
320c4762a1bSJed Brown   PetscReal      zero=0.0;
321c4762a1bSJed Brown 
3229566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->dm,&localX));
3239566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->dm,&localF));
3249566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->dm,X,INSERT_VALUES,localX));
3259566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->dm,X,INSERT_VALUES,localX));
3269566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(user->dm,&xs,NULL,NULL,&xm,NULL,NULL));
3279566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,NULL,NULL,&gxm,NULL,NULL));
3289566063dSJacob Faibussowitsch   PetscCall(VecSet(F, zero));
329c4762a1bSJed Brown   /*
330c4762a1bSJed Brown      The problem size is smaller than the discretization because of the
331c4762a1bSJed Brown      two fixed elements (V(0,T) = E and V(Send,T) = 0.
332c4762a1bSJed Brown   */
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   /* Get pointers to the vector data */
3359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &x));
3369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localF, &f));
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   /* Left Boundary */
339c4762a1bSJed Brown   if (gxs==0) {
340c4762a1bSJed Brown     f[0] = x[0]-user->strike;
341c4762a1bSJed Brown   } else {
342c4762a1bSJed Brown     f[0] = 0;
343c4762a1bSJed Brown   }
344c4762a1bSJed Brown 
345c4762a1bSJed Brown   /* All points in the interior */
346c4762a1bSJed Brown   /*  for (i=gxs+1;i<gxm-1;i++) { */
347c4762a1bSJed Brown   for (i=1;i<gxm-1;i++) {
348c4762a1bSJed 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]) +
349c4762a1bSJed Brown            (d[i]/(2*ds*ds))*(x[i+1] -2*x[i] + x[i-1] + Vt1[i+1] - 2*Vt1[i] + Vt1[i-1]);
350c4762a1bSJed Brown   }
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   /* Right boundary */
353c4762a1bSJed Brown   if (gxs+gxm==ms) {
354c4762a1bSJed Brown     f[xm-1]=x[gxm-1];
355c4762a1bSJed Brown   } else {
356c4762a1bSJed Brown     f[xm-1]=0;
357c4762a1bSJed Brown   }
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   /* Restore vectors */
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &x));
3619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localF, &f));
3629566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(user->dm,localF,INSERT_VALUES,F));
3639566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(user->dm,localF,INSERT_VALUES,F));
3649566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->dm,&localX));
3659566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->dm,&localF));
3669566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(24.0*(gxm-2)));
367c4762a1bSJed Brown   /*
368c4762a1bSJed Brown   info=VecView(F,PETSC_VIEWER_STDOUT_WORLD);
369c4762a1bSJed Brown   */
370c4762a1bSJed Brown   return 0;
371c4762a1bSJed Brown }
372c4762a1bSJed Brown 
373c4762a1bSJed Brown /* ------------------------------------------------------------------- */
374c4762a1bSJed Brown /*
375c4762a1bSJed Brown    FormJacobian - Evaluates Jacobian matrix.
376c4762a1bSJed Brown 
377c4762a1bSJed Brown    Input Parameters:
378c4762a1bSJed Brown .  tao  - the Tao context
379c4762a1bSJed Brown .  X    - input vector
380c4762a1bSJed Brown .  ptr  - optional user-defined context, as set by TaoSetJacobian()
381c4762a1bSJed Brown 
382c4762a1bSJed Brown    Output Parameters:
383c4762a1bSJed Brown .  J    - Jacobian matrix
384c4762a1bSJed Brown */
385c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao tao, Vec X, Mat J, Mat tJPre, void *ptr)
386c4762a1bSJed Brown {
387c4762a1bSJed Brown   AppCtx         *user = (AppCtx *) ptr;
388c4762a1bSJed Brown   PetscReal      *c = user->c, *d = user->d;
389c4762a1bSJed Brown   PetscReal      rate = user->rate;
390c4762a1bSJed Brown   PetscReal      dt = user->dt, ds = user->ds;
391c4762a1bSJed Brown   PetscInt       ms = user->ms;
392c4762a1bSJed Brown   PetscReal      val[3];
393c4762a1bSJed Brown   PetscInt       col[3];
394c4762a1bSJed Brown   PetscInt       i;
395c4762a1bSJed Brown   PetscInt       gxs,gxm;
396c4762a1bSJed Brown   PetscBool      assembled;
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   /* Set various matrix options */
3999566063dSJacob Faibussowitsch   PetscCall(MatSetOption(J,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
4009566063dSJacob Faibussowitsch   PetscCall(MatAssembled(J,&assembled));
4019566063dSJacob Faibussowitsch   if (assembled) PetscCall(MatZeroEntries(J));
402c4762a1bSJed Brown 
4039566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(user->dm,&gxs,NULL,NULL,&gxm,NULL,NULL));
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   if (gxs==0) {
406c4762a1bSJed Brown     i = 0;
407c4762a1bSJed Brown     col[0] = 0;
408c4762a1bSJed Brown     val[0]=1.0;
4099566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J,1,&i,1,col,val,INSERT_VALUES));
410c4762a1bSJed Brown   }
411c4762a1bSJed Brown   for (i=1; i < gxm-1; i++) {
412c4762a1bSJed Brown     col[0] = gxs + i - 1;
413c4762a1bSJed Brown     col[1] = gxs + i;
414c4762a1bSJed Brown     col[2] = gxs + i + 1;
415c4762a1bSJed Brown     val[0] = -c[i]/(4*ds) + d[i]/(2*ds*ds);
416c4762a1bSJed Brown     val[1] = 1.0/dt + rate - d[i]/(ds*ds);
417c4762a1bSJed Brown     val[2] =  c[i]/(4*ds) + d[i]/(2*ds*ds);
4189566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J,1,&col[1],3,col,val,INSERT_VALUES));
419c4762a1bSJed Brown   }
420c4762a1bSJed Brown   if (gxs+gxm==ms) {
421c4762a1bSJed Brown     i = ms-1;
422c4762a1bSJed Brown     col[0] = i;
423c4762a1bSJed Brown     val[0]=1.0;
4249566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J,1,&i,1,col,val,INSERT_VALUES));
425c4762a1bSJed Brown   }
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   /* Assemble the Jacobian matrix */
4289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
4299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
4309566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(18.0*(gxm)+5));
431c4762a1bSJed Brown   return 0;
432c4762a1bSJed Brown }
433c4762a1bSJed Brown 
434c4762a1bSJed Brown /*TEST
435c4762a1bSJed Brown 
436c4762a1bSJed Brown    build:
437c4762a1bSJed Brown       requires: !complex
438c4762a1bSJed Brown 
439c4762a1bSJed Brown    test:
440c4762a1bSJed Brown       args: -tao_monitor -tao_type ssils -tao_gttol 1.e-5
441c4762a1bSJed Brown       requires: !single
442c4762a1bSJed Brown 
443c4762a1bSJed Brown    test:
444c4762a1bSJed Brown       suffix: 2
445c4762a1bSJed Brown       args: -tao_monitor -tao_type ssfls -tao_max_it 10 -tao_gttol 1.e-5
446c4762a1bSJed Brown       requires: !single
447c4762a1bSJed Brown 
448c4762a1bSJed Brown    test:
449c4762a1bSJed Brown       suffix: 3
450c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type subvec -tao_gttol 1.e-5
451c4762a1bSJed Brown       requires: !single
452c4762a1bSJed Brown 
453c4762a1bSJed Brown    test:
454c4762a1bSJed Brown       suffix: 4
455c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type mask -tao_gttol 1.e-5
456c4762a1bSJed Brown       requires: !single
457c4762a1bSJed Brown 
458c4762a1bSJed Brown    test:
459c4762a1bSJed Brown       suffix: 5
460c4762a1bSJed Brown       args: -tao_monitor -tao_type asils -tao_subset_type matrixfree -pc_type jacobi -tao_max_it 6 -tao_gttol 1.e-5
461c4762a1bSJed Brown       requires: !single
462c4762a1bSJed Brown 
463c4762a1bSJed Brown    test:
464c4762a1bSJed Brown       suffix: 6
465c4762a1bSJed Brown       args: -tao_monitor -tao_type asfls -tao_subset_type subvec -tao_max_it 10 -tao_gttol 1.e-5
466c4762a1bSJed Brown       requires: !single
467c4762a1bSJed Brown 
468c4762a1bSJed Brown    test:
469c4762a1bSJed Brown       suffix: 7
470c4762a1bSJed Brown       args: -tao_monitor -tao_type asfls -tao_subset_type mask -tao_max_it 10 -tao_gttol 1.e-5
471c4762a1bSJed Brown       requires: !single
472c4762a1bSJed Brown 
473c4762a1bSJed Brown TEST*/
474