Class FeynmanKac
- All Implemented Interfaces:
Serializable,Cloneable
Class FeynmanKac solves the generalized Feynman-Kac PDE on a rectangular
grid with boundary conditions using a finite element Galerkin method. The generalized
Feynman-Kac differential equation has the form
$$
f_t +\mu(x,t) f_x +\frac{\sigma^2(x,t)}{2}f_{xx}-\kappa(x,t)f=\phi(f,x,t),
$$
where the initial data satisfies \(f(x,T)=p(x)\). The derivatives are
\(f_t=\frac{\partial{f}}{\partial{t}}\), etc. Method
computeCoefficients uses a finite element Galerkin method over the
rectangle
$$
[x_{\min},x_{\max}] \times [\bar{T},T]
$$
in \((x,t)\) to compute the approximate solution. The interval
\([x_{\min},x_{\max}]\) is decomposed with a grid
$$
(x_{\min}=)x_1 \lt x_2 \lt \ldots \lt x_m(=x_{\max})\,.
$$
On each subinterval the solution is represented by
$$
f(x,t) = f_ib_0(z)+f_{i+1}b_0(1-z)+h_i f_i' b_1(z)-h_i f_{i+1}' b_1(1-z)+
h_i^2f_i''b_2(z)+h_i^2f_{i+1}''b_2(1-z).
$$
The values
$$
f_i,f_i',f_i'',f_{i+1},f_{i+1}',f_{i+1}''
$$
are time-dependent coefficients associated with each interval. The basis functions
\(b_0,b_1,b_2\) are given for
$$
x \in [x_i, x_{i+1}], \, h_i=x_{i+1}-x_i, \, z=(x-x_i)/h_i, z \in [0,1],
$$
by
$$
b_0(z)=-6z^5+15z^4-10z^3+1=(1-z)^3(6z^2+3z+1),
$$
$$
b_1(z)=-3z^5+8z^4-6z^3+z=(1-z)^3z(3z+1),
$$
$$
b_2(z)=\frac{1}{2}(-z^5+3z^4-3z^3+z^2)=\frac{1}{2}(1-z)^3z^2.
$$
By adding the piece-wise definitions the unknown solution function may be arranged as the
series
$$
f(x,t) = \sum_{i=1}^{3m}y_i\beta_i(x), \; x \in [x_{\min},x_{\max}],
$$
where the time-dependent coefficients are defined by re-labeling:
$$
y_{3i-2} = f_i,\, y_{3i-1} = f_i',\, y_{3i} = f_i'',\, i=1,\ldots,m.
$$
The Galerkin principle is then applied. Using the provided initial and boundary conditions
leads to an Index 1 differential-algebraic equation for the coefficients
\(y_i \, i=1,\ldots,3m\).
This system is integrated using the variable order, variable step algorithm DDASLX
as noted in Hanson, R. and Krogh, F. (2008).
Solving Constrained Differential-Algebraic Systems Using Projections.
Solution values
and their time derivatives at a grid preceding time T,
expressed in units of time remaining, are given back by methods getSplineCoefficients
and getSplineCoefficientsPrime, respectively. To evaluate f or its
partials \(f_x, f_{xx}, f_{xxx}, f_t, f_{tx}, f_{txx}, f_{txxx}\) at any
time point in the grid, use method getSplineValue.
For further details of deriving and solving the system see Hanson, R. (2008) and Integrating Feynman-Kac Equations Using Hermite Quintic Finite Elements.
One useful application of the FeynmanKac class is financial analytics.
This is illustrated in Example 2, which solves a diffusion model for call options which, in the
special case \(\alpha = 2\) reduces to the Black-Scholes (BS) model.
Another useful application for the FeynmanKac class is the calculation of the Greeks,
i.e. various derivatives of Feynman-Kac solutions applicable to, e.g., the pricing of options and
related financial derivatives. In Example 5, the FeynmanKac class is used to
calculate eleven of the Greeks for the same diffusion model introduced in Example 2
in the special BS case. These Greeks are also calculated using the BS closed form Greek equations
(https://en.wikipedia.org/wiki/Greeks_(finance)).
The Feynman-Kac and BS solutions are output and compared.
Example 5 illustrates that the FeynmanKac class can be used to
explore the Greeks for a much wider class of financial models than can BS.
-
Nested Class Summary
Nested ClassesModifier and TypeClassDescriptionstatic interfacePublic interface for user supplied boundary coefficients and terminal condition the PDE must satisfy.static classThe boundary conditions are inconsistent.static classThe constraints are inconsistent.static classCorrector failed to converge.static classError test failure detected.static interfacePublic interface for non-zero forcing term in the Feynman-Kac equation.static classThe constraints at the initial point are inconsistent.static interfacePublic interface for adjustment of initial data or as an opportunity for output during the integration steps.static classIteration matrix is singular.static interfacePublic interface for user supplied PDE coefficients in the Feynman-Kac PDE.static classThe end value for the integration in time, tout, is not consistent with the current time value, t.static classThe current integration point in time and the end point are equal.static classDistance between starting time point and end point for the integration is too small.static classTolerance is too small.static classToo many iterations required by the DAE solver. -
Field Summary
FieldsModifier and TypeFieldDescriptionstatic final intUsed by methodsetStepControlMethodto indicate that the step control algorithm of the original Petzold code is used in the integration.static final intUsed by methodsetStepControlMethodto indicate that the step control method by Soederlind is used in the integration. -
Constructor Summary
ConstructorsConstructorDescriptionFeynmanKac(FeynmanKac.PdeCoefficients pdeCoeffs) Constructs a PDE solver to solve the Feynman-Kac PDE. -
Method Summary
Modifier and TypeMethodDescriptionvoidcomputeCoefficients(int numLeftBounds, int numRightBounds, FeynmanKac.Boundaries pdeBounds, double[] xGrid, double[] tGrid) Determines the coefficients of the Hermite quintic splines that represent an approximate solution for the Feynman-Kac PDE.double[]Returns absolute error tolerances.intReturns the number of quadrature points used in the Gauss-Legendre quadrature formula.doubleReturns the starting step size for the integration.intReturns the maximum order of the BDF formulas.doubleReturns the maximum internal step size used by the integrator.intReturns the maximum number of internal steps allowed.double[]Returns relative error tolerances.double[][]Returns the coefficients of the Hermite quintic splines that represent an approximate solution of the Feynman-Kac PDE.double[][]Returns the first derivatives of the Hermite quintic spline coefficients that represent an approximate solution of the Feynman-Kac PDE.double[]getSplineValue(double[] evaluationPoints, double[] coefficients, int ideriv) Evaluates for time value 0 or a time value in tGrid the derivative of the Hermite quintic spline interpolant at evaluation points within the range ofxGrid.intReturns the step control method used in the integration of the Feynman-Kac PDE.doubleReturns the barrier set for integration in the time direction.voidsetAbsoluteErrorTolerances(double atol) Sets the absolute error tolerances.voidsetAbsoluteErrorTolerances(double[] atol) Sets the absolute error tolerances.voidsetForcingTerm(FeynmanKac.ForcingTerm forceTerm) Sets the user-supplied method that computes approximations to the forcing term \(\phi(x)\) and its derivative \(\partial \phi/\partial y\) used in the FeynmanKac PDE.voidsetGaussLegendreDegree(int degree) Sets the number of quadrature points used in the Gauss-Legendre quadrature formula.voidsetInitialData(FeynmanKac.InitialData initData) Sets the user-supplied method for adjustment of initial data or as an opportunity for output during the integration steps.voidsetInitialStepsize(double initStepsize) Sets the starting stepsize for the integration.voidsetMaximumBDFOrder(int maxBDFOrder) Sets the maximum order of the BDF formulas.voidsetMaximumStepsize(double maximumStepsize) Sets the maximum internal step size used by the integrator.voidsetMaxSteps(int maxSteps) Sets the maximum number of internal steps allowed.voidsetRelativeErrorTolerances(double rtol) Sets the relative error tolerances.voidsetRelativeErrorTolerances(double[] rtol) Sets the relative error tolerances.voidsetStepControlMethod(int stepControlMethod) Sets the step control method used in the integration of the Feynman-Kac PDE.voidsetTimeBarrier(double timeBarrier) Sets a barrier for the integration in the time direction.voidsetTimeDependence(boolean[] timeFlag) Sets the time dependence of the coefficients, boundary conditions and function \(\phi\) in the Feynman Kac equation.
-
Field Details
-
METHOD_OF_SOEDERLIND
public static final int METHOD_OF_SOEDERLINDUsed by methodsetStepControlMethodto indicate that the step control method by Soederlind is used in the integration.- See Also:
-
METHOD_OF_PETZOLD
public static final int METHOD_OF_PETZOLDUsed by methodsetStepControlMethodto indicate that the step control algorithm of the original Petzold code is used in the integration.- See Also:
-
-
Constructor Details
-
FeynmanKac
Constructs a PDE solver to solve the Feynman-Kac PDE.- Parameters:
pdeCoeffs- Implementation of interface PdeCoefficients that computes the values of the Feynman-Kac coefficients at a given point(t,x).
-
-
Method Details
-
getSplineCoefficients
public double[][] getSplineCoefficients()Returns the coefficients of the Hermite quintic splines that represent an approximate solution of the Feynman-Kac PDE.- Returns:
- a
doublearray of dimension(tGrid.length+1) by (3*xGrid.length)containing the coefficients of the Hermite quintic spline representation of the approximate solution for the Feynman-Kac PDE at time points0, tGrid[0],...,tGrid[tGrid.length-1]. SettingntGrid = tGrid.lengthandnxGrid = xGrid.lengththe approximate solution is given by $$ f(x,t) = \sum_{j=0}^{3*\text{nxGrid}-1}y_{ij}\beta_j(x) \; \mbox{for} \; t=\text{tGrid}[i-1], i=1,\ldots,\text{ntGrid}. $$ The representation for the initial data att=0is $$ p(x) = \sum_{j=0}^{3*\text{nxGrid}-1} y_{0j}\beta_j(x)\,. $$ The(ntGrid+1) by (3*nxGrid)matrix $$ (y_{ij})_{i=0,\ldots,\text{ntGrid}}^{j=0,\ldots,3*\text{nxGrid}-1} $$ is stored row-wise in the returned array.
-
getSplineCoefficientsPrime
public double[][] getSplineCoefficientsPrime()Returns the first derivatives of the Hermite quintic spline coefficients that represent an approximate solution of the Feynman-Kac PDE.- Returns:
- a
doublearray of dimension(tGrid.length+1) by (3*xGrid.length)containing the first derivatives (in time) of the coefficients of the Hermite quintic spline representation of the approximate solution for the Feynman-Kac PDE at time points0, tGrid[0],...,tGrid[tGrid.length-1]. The approximate solution itself is given by $$ f_t(x,\bar{t}) = \sum_{j=0}^{3*\text{nxGrid}-1}y_{ij}' \beta_j(x) \quad \mbox{for} \; \bar{t} = \text{tGrid}[i-1], i=1,\ldots,\text{ntGrid}, $$ and $$ f_t(x,\bar{t})=\sum_{j=0}^{3*\text{nxGrid}-1}y_{0j}^\prime\beta_j(x) \quad \mbox{for} \; \bar{t} = 0 \, . $$ The(ntGrid+1) by (3*nxGrid)matrix $$ (y_{ij}')_{i=0,\ldots,\text{ntGrid}}^{j=0,\ldots,3*\text{nxGrid}-1} $$ is stored row-wise in the returned array.
-
setAbsoluteErrorTolerances
public void setAbsoluteErrorTolerances(double[] atol) Sets the absolute error tolerances.- Parameters:
atol- adoublearray of length3*xGrid.lengthspecifying the absolute error tolerances for the row-wise solutions returned by methodgetSplineCoefficients. All entries inatolmust be greater than or equal to zero. Also, not all entries inatolandrtolare allowed to be equal to 0 simultaneously. Default value is 1.0e-5 for each solution component.
-
setRelativeErrorTolerances
public void setRelativeErrorTolerances(double[] rtol) Sets the relative error tolerances.- Parameters:
rtol- adoublearray of length3*xGrid.lengthspecifying the relative error tolerances for the solution. All entries inrtolmust be greater than or equal to zero. Also, not all entries inatolandrtolare allowed to be equal to 0 simultaneously. Default value is 1.0e-5 for each solution component.
-
setAbsoluteErrorTolerances
public void setAbsoluteErrorTolerances(double atol) Sets the absolute error tolerances.- Parameters:
atol- adoublescalar specifying the absolute error tolerances for the row-wise solutions returned by methodgetSplineCoefficients. The tolerance valueatolis applied to all3*xGrid.lengthsolution components.atolmust be greater than or equal to zero. Also, not all entries inatolandrtolare allowed to be equal to 0 simultaneously. Default value is 1.0e-5 for each solution component.
-
setRelativeErrorTolerances
public void setRelativeErrorTolerances(double rtol) Sets the relative error tolerances.- Parameters:
rtol- adoublescalar specifying the relative error tolerances for the solution. The tolerance valuertolis applied to all3*xGrid.lengthsolution components.rtolmust be greater than or equal to zero. Also, not all entries inatolandrtolare allowed to be equal to 0 simultaneously. Default value is 1.0e-5 for each solution component.
-
getAbsoluteErrorTolerances
public double[] getAbsoluteErrorTolerances()Returns absolute error tolerances.- Returns:
- a
doublearray of length3*xGrid.lengthcontaining absolute error tolerances for the solutions.
-
getRelativeErrorTolerances
public double[] getRelativeErrorTolerances()Returns relative error tolerances.- Returns:
- a
doublearray of length3*xGrid.lengthcontaining relative error tolerances for the solutions.
-
setGaussLegendreDegree
public void setGaussLegendreDegree(int degree) Sets the number of quadrature points used in the Gauss-Legendre quadrature formula.- Parameters:
degree- anint, the degree of the polynomial used in the Gauss-Legendre quadrature. It is required thatdegreeis greater than or equal to 6. The default value is 6.
-
getGaussLegendreDegree
public int getGaussLegendreDegree()Returns the number of quadrature points used in the Gauss-Legendre quadrature formula.- Returns:
- an
int, the degree of the polynomial used in the Gauss-Legendre quadrature.
-
setMaximumStepsize
public void setMaximumStepsize(double maximumStepsize) Sets the maximum internal step size used by the integrator.- Parameters:
maximumStepsize- a positive scalar of typedouble, the maximum internal step size. Default value isDouble.MAX_VALUE, the largest possible machine number.
-
getMaximumStepsize
public double getMaximumStepsize()Returns the maximum internal step size used by the integrator.- Returns:
- a
double, the maximum internal step size.
-
setInitialStepsize
public void setInitialStepsize(double initStepsize) Sets the starting stepsize for the integration.- Parameters:
initStepsize- adouble, the starting stepsize used by the integrator. Must be less than zero since the integration is internally done fromt=0tot=tGrid[tGrid.length-1]in a negative direction. The default isinitStepsize = -1.1102230246252e-16.
-
getInitialStepsize
public double getInitialStepsize()Returns the starting step size for the integration.- Returns:
- a
double, the starting step size used in the integrator.
-
setMaxSteps
public void setMaxSteps(int maxSteps) Sets the maximum number of internal steps allowed.- Parameters:
maxSteps- anintspecifying the maximum number of internal steps allowed between each output point of the integration.maxStepsmust be positive. The default value is 500000.
-
getMaxSteps
public int getMaxSteps()Returns the maximum number of internal steps allowed.- Returns:
- an
intspecifying the maximum number of internal steps allowed between each output point of the integration.
-
setStepControlMethod
public void setStepControlMethod(int stepControlMethod) Sets the step control method used in the integration of the Feynman-Kac PDE.- Parameters:
stepControlMethod- anintscalar specifying the step control method to be used in the integration. If this member function is not calledstepControlMethodis set toMETHOD_OF_SOEDERLINDby default.stepControlMethodDescription METHOD_OF_SOEDERLINDUse method of Soederlind. METHOD_OF_PETZOLDUse method from the original Petzold code DASSL.
-
getStepControlMethod
public int getStepControlMethod()Returns the step control method used in the integration of the Feynman-Kac PDE.- Returns:
- an
intscalar specifying which step control method to be used.Return Value Description METHOD_OF_SOEDERLINDMethod of Soederlind METHOD_OF_PETZOLDMethod from the original Petzold code DASSL
-
setMaximumBDFOrder
public void setMaximumBDFOrder(int maxBDFOrder) Sets the maximum order of the BDF formulas.- Parameters:
maxBDFOrder- anintspecifying the maximum order of the backward differentiation formulas used in the integrator. It is required thatmaxBDFOrderis greater than zero and smaller than 6. The default ismaxBDFOrder = 5.
-
getMaximumBDFOrder
public int getMaximumBDFOrder()Returns the maximum order of the BDF formulas.- Returns:
- an
int, the maximum order of the backward differentiation formulas used in the integrator.
-
setTimeBarrier
public void setTimeBarrier(double timeBarrier) Sets a barrier for the integration in the time direction.- Parameters:
timeBarrier- adouble, controls whether the integrator should integrate in the time direction beyond a special point,timeBarrier, and then interpolate to get the Hermite quintic spline coefficients and its derivatives at the points intGrid. It is required thattimeBarrierbe greater than or equal totGrid[tGrid.length-1]. The default istimeBarrier = tGrid[tGrid.length-1].
-
getTimeBarrier
public double getTimeBarrier()Returns the barrier set for integration in the time direction.- Returns:
- a
double, the time point beyond which the integrator should not integrate.
-
setInitialData
Sets the user-supplied method for adjustment of initial data or as an opportunity for output during the integration steps.- Parameters:
initData- anInitialDataobject specifying the user-defined function for adjustment of initial data or the object can be used as an opportunity for output during the integration steps. If this member function is not called, no adjustment of initial data or output during the integration steps is done.
-
setForcingTerm
Sets the user-supplied method that computes approximations to the forcing term \(\phi(x)\) and its derivative \(\partial \phi/\partial y\) used in the FeynmanKac PDE.- Parameters:
forceTerm- aForcingTerminterface specifying the user defined function used for computation of the forcing term \(\phi(f,x,t)\) and its derivative \(\partial \phi/\partial y\). If this member function is not called it is assumed that \(\phi(f,x,t)\) is identically zero.
-
setTimeDependence
public void setTimeDependence(boolean[] timeFlag) Sets the time dependence of the coefficients, boundary conditions and function \(\phi\) in the Feynman Kac equation.- Parameters:
timeFlag- abooleanvector of length 7 indicating time dependencies in the Feynman-Kac PDE.Index Time dependency of 0\(\sigma'\) 1\(\sigma\) 2\(\mu\) 3\(\kappa\) 4Left boundary conditions 5Right boundary conditions 6\(\phi\) timeFlag[i] = trueindicates that the associated value is time-dependent,timeFlag[i] = falseindicates that the associated value is not time-dependent. By default,timeFlag[i] = falsefor \(i=0,\ldots,6\).
-
computeCoefficients
public void computeCoefficients(int numLeftBounds, int numRightBounds, FeynmanKac.Boundaries pdeBounds, double[] xGrid, double[] tGrid) throws FeynmanKac.ToleranceTooSmallException, FeynmanKac.TooManyIterationsException, FeynmanKac.ErrorTestException, FeynmanKac.CorrectorConvergenceException, FeynmanKac.IterationMatrixSingularException, FeynmanKac.TimeIntervalTooSmallException, FeynmanKac.TcurrentTstopInconsistentException, FeynmanKac.TEqualsToutException, FeynmanKac.InitialConstraintsException, FeynmanKac.ConstraintsInconsistentException, SingularMatrixException, FeynmanKac.BoundaryInconsistentException Determines the coefficients of the Hermite quintic splines that represent an approximate solution for the Feynman-Kac PDE.- Parameters:
numLeftBounds- anint, the number of left boundary conditions. It is required that
1\(\le\)numLeftBounds\(\le\)3.numRightBounds- anint, the number of right boundary conditions. It is required that
1\(\le\)numRightBounds\(\le\)3.pdeBounds- Implementation of interface Boundaries that computes the boundary coefficients and terminal condition for given(t,x).xGrid- adoublearray containing the breakpoints for the Hermite quintic splines used in thexdiscretization. The length ofxGridmust be at least 2,xGrid.length >= 2, and the elements inxGridmust be in strictly increasing order.tGrid- adoublearray containing the set of time points (in time-remaining units) where an approximate solution is returned. The elements in arraytGridmust be positive and in strictly increasing order.- Throws:
FeynmanKac.ToleranceTooSmallException- is thrown if the absolute or relative error tolerances used in the integrator are too small.FeynmanKac.TooManyIterationsException- is thrown if the integrator needs too many iteration steps.FeynmanKac.ErrorTestException- is thrown if the error test used in the integrator failed repeatedly.FeynmanKac.CorrectorConvergenceException- is thrown if the corrector failed to converge repeatedly.FeynmanKac.IterationMatrixSingularException- thrown if one of the iteration matrices used in the integrator is singular.FeynmanKac.TimeIntervalTooSmallException- is thrown if the distance between an intermediate starting and end point for the integration is too small.FeynmanKac.TcurrentTstopInconsistentException- is thrown if during the integration the current integration time and given stepsize is inconsistent with the endpoint of the integration.FeynmanKac.TEqualsToutException- is thrown if during the integration process the actual integration time and the end time of the integration are identical.FeynmanKac.InitialConstraintsException- is thrown if at the initial integration point some of the constraints are inconsistent.FeynmanKac.ConstraintsInconsistentException- is thrown if during the integration process the constraints for the actual time point and given stepsize are inconsistent.SingularMatrixException- is thrown if one of the matrices used outside the integrator is singular.FeynmanKac.BoundaryInconsistentException- is thrown if the boundary conditions are inconsistent.
-
getSplineValue
public double[] getSplineValue(double[] evaluationPoints, double[] coefficients, int ideriv) Evaluates for time value 0 or a time value in tGrid the derivative of the Hermite quintic spline interpolant at evaluation points within the range ofxGrid.- Parameters:
evaluationPoints- adoublearray containing the points inx-direction at which the Hermite quintic spline representing the approximate solution to the Feynman-Kac PDE or one of its derivatives is to be evaluated. It is required that all elements in arrayevaluationPointsare greater than or equal toxGrid[0]and less than or equal toxGrid[xGrid.length-1].coefficients- adoublearray of length3*xGrid.lengthcontaining the coefficients of the Hermite quintic spline representing the approximate solution f or \(f_t\) to the Feynman-Kac PDE. These coefficients are the rows of the arrayssplineCoeffsandsplineCoeffsPrimereturned by methodsgetSplineCoefficientsandgetSplineCoefficientsPrime. If the user wants to compute approximate solutions \(f\) or \(f_x,f_{xx},f_{xxx}\) to the Feynman-Kac PDE at time point0, one must assign rowsplineCoeffs[0]to arraycoefficients. If the user wants to compute these approximate solutions for time pointst=tGrid[i], i=0,...,tGrid.length-1, one must assign rowsplineCoeffs[i+1]to arraycoefficients. The same reasoning applies to the computation of approximate solutions \(f_t\) and \(f_{tx},f_{txx},f_{txxx}\) and assignment of rows of arraygetSplineCoefficientsPrimeto arraycoefficients.ideriv- anintspecifying the derivative to be computed. It must be 0, 1, 2 or 3.- Returns:
- a
doublearray containing the derivative of orderiderivof the Hermite quintic spline representing the approximate solution f or \(f_t\) to the Feynman Kac PDE atevaluationPoints. Ifideriv=0, then the spline values are returned. Ifideriv=1, then the first derivative is returned, etc.
-