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 769371c9d4SSatish Balay static char help[] = "This example demonstrates use of the TAO package to\n\ 77c4762a1bSJed Brown solve a linear complementarity problem for pricing American put options.\n\ 78c4762a1bSJed Brown The code uses backward differences in time and central differences in\n\ 79c4762a1bSJed Brown space. The command line options are:\n\ 80c4762a1bSJed Brown -rate <r>, where <r> = interest rate\n\ 81c4762a1bSJed Brown -sigma <s>, where <s> = volatility of the underlying\n\ 82c4762a1bSJed Brown -alpha <a>, where <a> = elasticity of the underlying\n\ 83c4762a1bSJed Brown -delta <d>, where <d> = dividend rate\n\ 84c4762a1bSJed Brown -strike <e>, where <e> = strike price\n\ 85c4762a1bSJed Brown -expiry <t>, where <t> = the expiration date\n\ 86c4762a1bSJed Brown -mt <tg>, where <tg> = number of grid points in time\n\ 87c4762a1bSJed Brown -ms <sg>, where <sg> = number of grid points in space\n\ 88c4762a1bSJed Brown -es <se>, where <se> = ending point of the space discretization\n\n"; 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown User-defined application context - contains data needed by the 92c4762a1bSJed Brown application-provided call-back routines, FormFunction(), and FormJacobian(). 93c4762a1bSJed Brown */ 94c4762a1bSJed Brown 95c4762a1bSJed Brown typedef struct { 96c4762a1bSJed Brown PetscReal *Vt1; /* Value of the option at time T + dt */ 97c4762a1bSJed Brown PetscReal *c; /* Constant -- (r - D)S */ 98c4762a1bSJed Brown PetscReal *d; /* Constant -- -0.5(sigma**2)(S**alpha) */ 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscReal rate; /* Interest rate */ 101c4762a1bSJed Brown PetscReal sigma, alpha, delta; /* Underlying asset properties */ 102c4762a1bSJed Brown PetscReal strike, expiry; /* Option contract properties */ 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscReal es; /* Finite value used for maximum asset value */ 105c4762a1bSJed Brown PetscReal ds, dt; /* Discretization properties */ 106c4762a1bSJed Brown PetscInt ms, mt; /* Number of elements */ 107c4762a1bSJed Brown 108c4762a1bSJed Brown DM dm; 109c4762a1bSJed Brown } AppCtx; 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* -------- User-defined Routines --------- */ 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscErrorCode FormConstraints(Tao, Vec, Vec, void *); 114c4762a1bSJed Brown PetscErrorCode FormJacobian(Tao, Vec, Mat, Mat, void *); 115c4762a1bSJed Brown PetscErrorCode ComputeVariableBounds(Tao, Vec, Vec, void *); 116c4762a1bSJed Brown 1179371c9d4SSatish Balay int main(int argc, char **argv) { 118c4762a1bSJed Brown Vec x; /* solution vector */ 119c4762a1bSJed Brown Vec c; /* Constraints function vector */ 120c4762a1bSJed Brown Mat J; /* jacobian matrix */ 121c4762a1bSJed Brown PetscBool flg; /* A return variable when checking for user options */ 122c4762a1bSJed Brown Tao tao; /* Tao solver context */ 123c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 124c4762a1bSJed Brown PetscInt i, j; 125c4762a1bSJed Brown PetscInt xs, xm, gxs, gxm; 126c4762a1bSJed Brown PetscReal sval = 0; 127c4762a1bSJed Brown PetscReal *x_array; 128c4762a1bSJed Brown Vec localX; 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* Initialize PETSc, TAO */ 131327415f7SBarry Smith PetscFunctionBeginUser; 1329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* 135c4762a1bSJed Brown Initialize the user-defined application context with reasonable 136c4762a1bSJed Brown values for the American option to price 137c4762a1bSJed Brown */ 138c4762a1bSJed Brown user.rate = 0.04; 139c4762a1bSJed Brown user.sigma = 0.40; 140c4762a1bSJed Brown user.alpha = 2.00; 141c4762a1bSJed Brown user.delta = 0.01; 142c4762a1bSJed Brown user.strike = 10.0; 143c4762a1bSJed Brown user.expiry = 1.0; 144c4762a1bSJed Brown user.mt = 10; 145c4762a1bSJed Brown user.ms = 150; 146c4762a1bSJed Brown user.es = 100.0; 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* Read in alternative values for the American option to price */ 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-alpha", &user.alpha, &flg)); 1509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-delta", &user.delta, &flg)); 1519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-es", &user.es, &flg)); 1529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-expiry", &user.expiry, &flg)); 1539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ms", &user.ms, &flg)); 1549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mt", &user.mt, &flg)); 1559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-rate", &user.rate, &flg)); 1569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-sigma", &user.sigma, &flg)); 1579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-strike", &user.strike, &flg)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* Check that the options set are allowable (needs to be done) */ 160c4762a1bSJed Brown 161c4762a1bSJed Brown user.ms++; 1629566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.ms, 1, 1, NULL, &user.dm)); 1639566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dm)); 1649566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dm)); 165c4762a1bSJed Brown /* Create appropriate vectors and matrices */ 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user.dm, &xs, NULL, NULL, &xm, NULL, NULL)); 1689566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user.dm, &gxs, NULL, NULL, &gxm, NULL, NULL)); 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dm, &x)); 171c4762a1bSJed Brown /* 172c4762a1bSJed Brown Finish filling in the user-defined context with the values for 173c4762a1bSJed Brown dS, dt, and allocating space for the constants 174c4762a1bSJed Brown */ 175c4762a1bSJed Brown user.ds = user.es / (user.ms - 1); 176c4762a1bSJed Brown user.dt = user.expiry / user.mt; 177c4762a1bSJed Brown 1789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(gxm, &(user.Vt1))); 1799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(gxm, &(user.c))); 1809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(gxm, &(user.d))); 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* 183c4762a1bSJed Brown Calculate the values for the constant. Vt1 begins with the ending 184c4762a1bSJed Brown boundary condition. 185c4762a1bSJed Brown */ 186c4762a1bSJed Brown for (i = 0; i < gxm; i++) { 187c4762a1bSJed Brown sval = (gxs + i) * user.ds; 188c4762a1bSJed Brown user.Vt1[i] = PetscMax(user.strike - sval, 0); 189c4762a1bSJed Brown user.c[i] = (user.delta - user.rate) * sval; 190c4762a1bSJed Brown user.d[i] = -0.5 * user.sigma * user.sigma * PetscPowReal(sval, user.alpha); 191c4762a1bSJed Brown } 192*ad540459SPierre Jolivet if (gxs + gxm == user.ms) user.Vt1[gxm - 1] = 0; 1939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &c)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* 196c4762a1bSJed Brown Allocate the matrix used by TAO for the Jacobian. Each row of 197c4762a1bSJed Brown the Jacobian matrix will have at most three elements. 198c4762a1bSJed Brown */ 1999566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.dm, &J)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* The TAO code begins here */ 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2049566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2059566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOSSILS)); 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* Set routines for constraints function and Jacobian evaluation */ 2089566063dSJacob Faibussowitsch PetscCall(TaoSetConstraintsRoutine(tao, c, FormConstraints, (void *)&user)); 2099566063dSJacob Faibussowitsch PetscCall(TaoSetJacobianRoutine(tao, J, J, FormJacobian, (void *)&user)); 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* Set the variable bounds */ 2129566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBoundsRoutine(tao, ComputeVariableBounds, (void *)&user)); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* Set initial solution guess */ 2159566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_array)); 2169371c9d4SSatish Balay for (i = 0; i < xm; i++) x_array[i] = user.Vt1[i - gxs + xs]; 2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_array)); 218c4762a1bSJed Brown /* Set data structure */ 2199566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* Set routines for function and Jacobian evaluation */ 2229566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Iteratively solve the linear complementarity problems */ 225c4762a1bSJed Brown for (i = 1; i < user.mt; i++) { 226c4762a1bSJed Brown /* Solve the current version */ 2279566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* Update Vt1 with the solution */ 2309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user.dm, &localX)); 2319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user.dm, x, INSERT_VALUES, localX)); 2329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user.dm, x, INSERT_VALUES, localX)); 2339566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x_array)); 234*ad540459SPierre Jolivet for (j = 0; j < gxm; j++) user.Vt1[j] = x_array[j]; 2359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_array)); 2369566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user.dm, &localX)); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 239c4762a1bSJed Brown /* Free TAO data structures */ 2409566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* Free PETSc data structures */ 2439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 2459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 2469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dm)); 247c4762a1bSJed Brown /* Free user-defined workspace */ 2489566063dSJacob Faibussowitsch PetscCall(PetscFree(user.Vt1)); 2499566063dSJacob Faibussowitsch PetscCall(PetscFree(user.c)); 2509566063dSJacob Faibussowitsch PetscCall(PetscFree(user.d)); 251c4762a1bSJed Brown 2529566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 253b122ec5aSJacob Faibussowitsch return 0; 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 2579371c9d4SSatish Balay PetscErrorCode ComputeVariableBounds(Tao tao, Vec xl, Vec xu, void *ctx) { 258c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 259c4762a1bSJed Brown PetscInt i; 260c4762a1bSJed Brown PetscInt xs, xm; 261c4762a1bSJed Brown PetscInt ms = user->ms; 262c4762a1bSJed Brown PetscReal sval = 0.0, *xl_array, ub = PETSC_INFINITY; 263c4762a1bSJed Brown 264c4762a1bSJed Brown /* Set the variable bounds */ 2659566063dSJacob Faibussowitsch PetscCall(VecSet(xu, ub)); 2669566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, NULL, NULL, &xm, NULL, NULL)); 267c4762a1bSJed Brown 2689566063dSJacob Faibussowitsch PetscCall(VecGetArray(xl, &xl_array)); 269c4762a1bSJed Brown for (i = 0; i < xm; i++) { 270c4762a1bSJed Brown sval = (xs + i) * user->ds; 271c4762a1bSJed Brown xl_array[i] = PetscMax(user->strike - sval, 0); 272c4762a1bSJed Brown } 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xl, &xl_array)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown if (xs == 0) { 2769566063dSJacob Faibussowitsch PetscCall(VecGetArray(xu, &xl_array)); 277c4762a1bSJed Brown xl_array[0] = PetscMax(user->strike, 0); 2789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xu, &xl_array)); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown if (xs + xm == ms) { 2819566063dSJacob Faibussowitsch PetscCall(VecGetArray(xu, &xl_array)); 282c4762a1bSJed Brown xl_array[xm - 1] = 0; 2839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xu, &xl_array)); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown return 0; 287c4762a1bSJed Brown } 288c4762a1bSJed Brown /* -------------------------------------------------------------------- */ 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* 291c4762a1bSJed Brown FormFunction - Evaluates gradient of f. 292c4762a1bSJed Brown 293c4762a1bSJed Brown Input Parameters: 294c4762a1bSJed Brown . tao - the Tao context 295c4762a1bSJed Brown . X - input vector 296c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoAppSetConstraintRoutine() 297c4762a1bSJed Brown 298c4762a1bSJed Brown Output Parameters: 299c4762a1bSJed Brown . F - vector containing the newly evaluated gradient 300c4762a1bSJed Brown */ 3019371c9d4SSatish Balay PetscErrorCode FormConstraints(Tao tao, Vec X, Vec F, void *ptr) { 302c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 303c4762a1bSJed Brown PetscReal *x, *f; 304c4762a1bSJed Brown PetscReal *Vt1 = user->Vt1, *c = user->c, *d = user->d; 305c4762a1bSJed Brown PetscReal rate = user->rate; 306c4762a1bSJed Brown PetscReal dt = user->dt, ds = user->ds; 307c4762a1bSJed Brown PetscInt ms = user->ms; 308c4762a1bSJed Brown PetscInt i, xs, xm, gxs, gxm; 309c4762a1bSJed Brown Vec localX, localF; 310c4762a1bSJed Brown PetscReal zero = 0.0; 311c4762a1bSJed Brown 3129566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localX)); 3139566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->dm, &localF)); 3149566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->dm, X, INSERT_VALUES, localX)); 3159566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->dm, X, INSERT_VALUES, localX)); 3169566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(user->dm, &xs, NULL, NULL, &xm, NULL, NULL)); 3179566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, NULL, NULL, &gxm, NULL, NULL)); 3189566063dSJacob Faibussowitsch PetscCall(VecSet(F, zero)); 319c4762a1bSJed Brown /* 320c4762a1bSJed Brown The problem size is smaller than the discretization because of the 321c4762a1bSJed Brown two fixed elements (V(0,T) = E and V(Send,T) = 0. 322c4762a1bSJed Brown */ 323c4762a1bSJed Brown 324c4762a1bSJed Brown /* Get pointers to the vector data */ 3259566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &x)); 3269566063dSJacob Faibussowitsch PetscCall(VecGetArray(localF, &f)); 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* Left Boundary */ 329c4762a1bSJed Brown if (gxs == 0) { 330c4762a1bSJed Brown f[0] = x[0] - user->strike; 331c4762a1bSJed Brown } else { 332c4762a1bSJed Brown f[0] = 0; 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 335c4762a1bSJed Brown /* All points in the interior */ 336c4762a1bSJed Brown /* for (i=gxs+1;i<gxm-1;i++) { */ 337c4762a1bSJed Brown for (i = 1; i < gxm - 1; i++) { 3389371c9d4SSatish Balay 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]) + (d[i] / (2 * ds * ds)) * (x[i + 1] - 2 * x[i] + x[i - 1] + Vt1[i + 1] - 2 * Vt1[i] + Vt1[i - 1]); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* Right boundary */ 342c4762a1bSJed Brown if (gxs + gxm == ms) { 343c4762a1bSJed Brown f[xm - 1] = x[gxm - 1]; 344c4762a1bSJed Brown } else { 345c4762a1bSJed Brown f[xm - 1] = 0; 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* Restore vectors */ 3499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &x)); 3509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localF, &f)); 3519566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(user->dm, localF, INSERT_VALUES, F)); 3529566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(user->dm, localF, INSERT_VALUES, F)); 3539566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localX)); 3549566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->dm, &localF)); 3559566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(24.0 * (gxm - 2))); 356c4762a1bSJed Brown /* 357c4762a1bSJed Brown info=VecView(F,PETSC_VIEWER_STDOUT_WORLD); 358c4762a1bSJed Brown */ 359c4762a1bSJed Brown return 0; 360c4762a1bSJed Brown } 361c4762a1bSJed Brown 362c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 363c4762a1bSJed Brown /* 364c4762a1bSJed Brown FormJacobian - Evaluates Jacobian matrix. 365c4762a1bSJed Brown 366c4762a1bSJed Brown Input Parameters: 367c4762a1bSJed Brown . tao - the Tao context 368c4762a1bSJed Brown . X - input vector 369c4762a1bSJed Brown . ptr - optional user-defined context, as set by TaoSetJacobian() 370c4762a1bSJed Brown 371c4762a1bSJed Brown Output Parameters: 372c4762a1bSJed Brown . J - Jacobian matrix 373c4762a1bSJed Brown */ 3749371c9d4SSatish Balay PetscErrorCode FormJacobian(Tao tao, Vec X, Mat J, Mat tJPre, void *ptr) { 375c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr; 376c4762a1bSJed Brown PetscReal *c = user->c, *d = user->d; 377c4762a1bSJed Brown PetscReal rate = user->rate; 378c4762a1bSJed Brown PetscReal dt = user->dt, ds = user->ds; 379c4762a1bSJed Brown PetscInt ms = user->ms; 380c4762a1bSJed Brown PetscReal val[3]; 381c4762a1bSJed Brown PetscInt col[3]; 382c4762a1bSJed Brown PetscInt i; 383c4762a1bSJed Brown PetscInt gxs, gxm; 384c4762a1bSJed Brown PetscBool assembled; 385c4762a1bSJed Brown 386c4762a1bSJed Brown /* Set various matrix options */ 3879566063dSJacob Faibussowitsch PetscCall(MatSetOption(J, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE)); 3889566063dSJacob Faibussowitsch PetscCall(MatAssembled(J, &assembled)); 3899566063dSJacob Faibussowitsch if (assembled) PetscCall(MatZeroEntries(J)); 390c4762a1bSJed Brown 3919566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(user->dm, &gxs, NULL, NULL, &gxm, NULL, NULL)); 392c4762a1bSJed Brown 393c4762a1bSJed Brown if (gxs == 0) { 394c4762a1bSJed Brown i = 0; 395c4762a1bSJed Brown col[0] = 0; 396c4762a1bSJed Brown val[0] = 1.0; 3979566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, &i, 1, col, val, INSERT_VALUES)); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown for (i = 1; i < gxm - 1; i++) { 400c4762a1bSJed Brown col[0] = gxs + i - 1; 401c4762a1bSJed Brown col[1] = gxs + i; 402c4762a1bSJed Brown col[2] = gxs + i + 1; 403c4762a1bSJed Brown val[0] = -c[i] / (4 * ds) + d[i] / (2 * ds * ds); 404c4762a1bSJed Brown val[1] = 1.0 / dt + rate - d[i] / (ds * ds); 405c4762a1bSJed Brown val[2] = c[i] / (4 * ds) + d[i] / (2 * ds * ds); 4069566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, &col[1], 3, col, val, INSERT_VALUES)); 407c4762a1bSJed Brown } 408c4762a1bSJed Brown if (gxs + gxm == ms) { 409c4762a1bSJed Brown i = ms - 1; 410c4762a1bSJed Brown col[0] = i; 411c4762a1bSJed Brown val[0] = 1.0; 4129566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, &i, 1, col, val, INSERT_VALUES)); 413c4762a1bSJed Brown } 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* Assemble the Jacobian matrix */ 4169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 4189566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(18.0 * (gxm) + 5)); 419c4762a1bSJed Brown return 0; 420c4762a1bSJed Brown } 421c4762a1bSJed Brown 422c4762a1bSJed Brown /*TEST 423c4762a1bSJed Brown 424c4762a1bSJed Brown build: 425c4762a1bSJed Brown requires: !complex 426c4762a1bSJed Brown 427c4762a1bSJed Brown test: 428c4762a1bSJed Brown args: -tao_monitor -tao_type ssils -tao_gttol 1.e-5 429c4762a1bSJed Brown requires: !single 430c4762a1bSJed Brown 431c4762a1bSJed Brown test: 432c4762a1bSJed Brown suffix: 2 433c4762a1bSJed Brown args: -tao_monitor -tao_type ssfls -tao_max_it 10 -tao_gttol 1.e-5 434c4762a1bSJed Brown requires: !single 435c4762a1bSJed Brown 436c4762a1bSJed Brown test: 437c4762a1bSJed Brown suffix: 3 438c4762a1bSJed Brown args: -tao_monitor -tao_type asils -tao_subset_type subvec -tao_gttol 1.e-5 439c4762a1bSJed Brown requires: !single 440c4762a1bSJed Brown 441c4762a1bSJed Brown test: 442c4762a1bSJed Brown suffix: 4 443c4762a1bSJed Brown args: -tao_monitor -tao_type asils -tao_subset_type mask -tao_gttol 1.e-5 444c4762a1bSJed Brown requires: !single 445c4762a1bSJed Brown 446c4762a1bSJed Brown test: 447c4762a1bSJed Brown suffix: 5 448c4762a1bSJed Brown args: -tao_monitor -tao_type asils -tao_subset_type matrixfree -pc_type jacobi -tao_max_it 6 -tao_gttol 1.e-5 449c4762a1bSJed Brown requires: !single 450c4762a1bSJed Brown 451c4762a1bSJed Brown test: 452c4762a1bSJed Brown suffix: 6 453c4762a1bSJed Brown args: -tao_monitor -tao_type asfls -tao_subset_type subvec -tao_max_it 10 -tao_gttol 1.e-5 454c4762a1bSJed Brown requires: !single 455c4762a1bSJed Brown 456c4762a1bSJed Brown test: 457c4762a1bSJed Brown suffix: 7 458c4762a1bSJed Brown args: -tao_monitor -tao_type asfls -tao_subset_type mask -tao_max_it 10 -tao_gttol 1.e-5 459c4762a1bSJed Brown requires: !single 460c4762a1bSJed Brown 461c4762a1bSJed Brown TEST*/ 462