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