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