hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_rand_bb_inc (g05xd)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_rand_bb_inc (g05xd) computes scaled increments of sample paths of a free or non-free Wiener process, where the sample paths are constructed by a Brownian bridge algorithm. The initialization function nag_rand_bb_inc_init (g05xc) must be called prior to the first call to nag_rand_bb_inc (g05xd).

Syntax

[z, b, ifail] = g05xd(npaths, dif, z, c, rcomm, 'rcord', rcord, 'd', d, 'a', a)
[z, b, ifail] = nag_rand_bb_inc(npaths, dif, z, c, rcomm, 'rcord', rcord, 'd', d, 'a', a)

Description

For details on the Brownian bridge algorithm and the bridge construction order see Brownian Bridge in the G05 Chapter Introduction and Description in nag_rand_bb_inc_init (g05xc). Recall that the terms Wiener process (or free Wiener process) and Brownian motion are often used interchangeably, while a non-free Wiener process (also known as a Brownian bridge process) refers to a process which is forced to terminate at a given point.
Fix two times t0<T, let ti 1iN  be any set of time points satisfying t0<t1<t2<<tN<T, and let Xt0, Xti 1iN , XT denote a d-dimensional Wiener sample path at these time points.
The Brownian bridge increments generator uses the Brownian bridge algorithm to construct sample paths for the (free or non-free) Wiener process X, and then uses this to compute the scaled Wiener increments
Xt1 - Xt0 t1 - t0 , Xt2 - Xt1 t2 - t1 ,, XtN - XtN-1 tN - tN-1 , XT - XtN T - tN  
The example program in Example shows how these increments can be used to compute a numerical solution to a stochastic differential equation (SDE) driven by a (free or non-free) Wiener process.

References

Glasserman P (2004) Monte Carlo Methods in Financial Engineering Springer

Parameters

Note: the following variable is used in the parameter descriptions: N=ntimes, the length of the array times passed to the initialization function nag_rand_bb_inc_init (g05xc).

Compulsory Input Parameters

1:     npaths int64int32nag_int scalar
The number of Wiener sample paths.
Constraint: npaths1.
2:     difd – double array
The difference between the terminal value and starting value of the Wiener process. If a=0, dif is ignored.
3:     zldz: – double array
The first dimension, ldz, of the array z must satisfy
  • if rcord=1, ldzd×N+1-a;
  • if rcord=2, ldznpaths.
The second dimension of the array z must be at least npaths if rcord=1 and at least d×n+1-a if rcord=2.
The Normal random numbers used to construct the sample paths.
If quasi-random numbers are used, the d×N+1-a-dimensional quasi-random points should be stored in successive rows of Z.
4:     cldc: – double array
The first dimension of the array c must be at least d.
The second dimension of the array c must be at least d.
The lower triangular Cholesky factorization C such that CCT gives the covariance matrix of the Wiener process. Elements of C above the diagonal are not referenced.
5:     rcomm* – double array
Note: the dimension of this array is dictated by the requirements of associated functions that must have been previously called. This array must be the same array passed as argument rcomm in the previous call to nag_rand_bb_inc_init (g05xc) or nag_rand_bb_inc (g05xd).
Communication array as returned by the last call to nag_rand_bb_inc_init (g05xc) or nag_rand_bb_inc (g05xd). This array must not be directly modified.

Optional Input Parameters

1:     rcord int64int32nag_int scalar
Default: 1
The order in which Normal random numbers are stored in z and in which the generated values are returned in b.
Constraint: rcord=1 or 2.
2:     d int64int32nag_int scalar
Default: the dimension of the array dif and the first dimension of the array c. (An error is raised if these dimensions are not equal.)
The dimension of each Wiener sample path.
Constraint: d1.
3:     a int64int32nag_int scalar
Default: 1
If a=0, a free Wiener process is created and dif is ignored.
If a=1, a non-free Wiener process is created where dif is the difference between the terminal value and the starting value of the process.
Constraint: a=0 or 1.

Output Parameters

1:     zldz: – double array
The first dimension, ldz, of the array z will be
  • if rcord=1, ldz=d×N+1-a;
  • if rcord=2, ldz=npaths.
The second dimension of the array z will be npaths if rcord=1 and at least d×n+1-a if rcord=2.
The Normal random numbers premultiplied by c.
2:     bldb: – double array
The first dimension, ldb, of the array b will be
  • if rcord=1, ldb=d×N+1;
  • if rcord=2, ldb=npaths.
The second dimension of the array b will be npaths if rcord=1 and at least d×n+1 if rcord=2.
The scaled Wiener increments.
Let Xp,ik denote the kth dimension of the ith point of the pth sample path where 1kd, 1in+1 and 1pnpaths. The increment Xp,ik-Xp,i-1kti-ti-1 is stored at Bp,k+i-1×d.
3:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:
   ifail=1
On entry, rcomm was not initialized or has been corrupted.
   ifail=2
Constraint: npaths1.
   ifail=3
On entry, rcord=_ was an illegal value.
   ifail=4
Constraint: d1.
   ifail=5
Constraint: a=0​ or ​1.
   ifail=6
Constraint: ldzd×ntimes+1-a.
Constraint: ldznpaths.
   ifail=7
ldc is too small.
   ifail=8
Constraint: ldbd×ntimes+1.
Constraint: ldbnpaths.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

Not applicable.

Further Comments

None.

Example

The scaled Wiener increments produced by this function can be used to compute numerical solutions to stochastic differential equations (SDEs) driven by (free or non-free) Wiener processes. Consider an SDE of the form
d Yt = ft,Yt dt + σt,Yt dXt  
on the interval t0,T  where Xt t0 tT  is a (free or non-free) Wiener process and f and σ are suitable functions. A numerical solution to this SDE can be obtained by the Euler–Maruyama method. For any discretization t0 < t1 < t2 << tN+1=T  of t0,T , set
Y ti+1 = Y ti + f ti,Yti ti+1 - ti + σ ti,Yti Xti+1 - Xti  
for i=1,,N so that YtN+1  is an approximation to YT. The scaled Wiener increments produced by nag_rand_bb_inc (g05xd) can be used in the Euler–Maruyama scheme outlined above by writing
Yti+1 = Yti + ti+1 - ti f ti,Yti + σ ti,Yti Xti+1 - Xti ti+1 - ti .  
The following example program uses this method to solve the SDE for geometric Brownian motion
d St = rSt dt + σSt dXt  
where X is a Wiener process, and compares the results against the analytic solution
ST = S0 exp r-σ2/2 T + σXT .  
Quasi-random variates are used to construct the Wiener increments.
function g05xd_example


fprintf('g05xd example results\n\n');

% We wish to solve the stochastic differential equation (SDE)
%   dSt = r*St*dt + sigma*St*dXt
% where X is a one dimensional Wiener process.
% This means we have
%     a = 0,   d = 1,  c = 1
% Set the other parameters of the SDE and the Euler-Maruyama scheme

% Initial value of the process
s0    = 1;
r     = 0.05;
sigma = 0.12;

% Number of paths to simulate
npaths = int64(3);

% The time interval [t0,T] on which to solve the SDE
t0   = 0;
tend = 1;

% The time steps in the discretization of [t0,T]
tn = 20;
t = t0 + (1:tn)*(tend-t0)/(tn+1);

% Make the bridge construction order
bgord = int64(3);
move  = zeros(0, 1, 'int64');
[times, ifail] = g05xe( ...
                        t0, tend, t, move, 'bgord', bgord);

% Generate the Z values
d = 1;
a = int64(0);
[z] = get_z(npaths, d, a, tn);

% Initialize the generator
[rcomm, ifail] = g05xc( ...
                        t0, tend, times);

% Get the scaled increments of the Wiener process
dif = zeros(d, 1);
c   = ones(d, 1);
[z, b, ifail] = g05xd( ...
                       npaths, dif, z, c, rcomm, 'a', a);

% Do the Euler-Maruyama time stepping
st = zeros(tn+1, npaths);

% Do first time step
st(1,:) = s0 + (t(1)-t0)*(r*s0+sigma*s0*b(1, :));
v = st(1,:);
for i=2:tn
  v = v + (t(i)-t(i-1))*(r*v+sigma*v.*b(i,:));
  st(i,:) = v;
end
% Do last time step
st(tn+1,:) = v + (tend-t(tn))*(r*v+sigma*v.*b(tn+1,:));

% Compute the analytic solution:
%    ST = S0*exp( (r-sigma^2/2)T + sigma WT )
analytic = s0*exp((r-0.5*sigma^2)*tend+sigma*sqrt(tend-t0)*z(1,:));

% Display the results
fprintf('\nEuler-Maruyama solution for Geometric Brownian motion\n');
fprintf('        Path 1     Path 2     Path 3\n');
for i = 1:tn+1
  fprintf('%2d  %10.4f %10.4f %10.4f\n', i, st(i, :));
end
fprintf('\nAnalytic solution at final time step\n');
fprintf('        Path 1     Path 2     Path 3\n');
fprintf('    %10.4f %10.4f %10.4f\n', analytic);



function [z] = get_z(npaths, d, a, ntimes)
  idim = d*(ntimes+1-a);

  % We now need to generate the input pseudorandom points

  % First initialize the base pseudorandom number generator
  state = initialize_prng(int64(6), int64(0), [int64(1023401)]);

  % Scrambled quasi-random sequences preserve the good discrepancy
  % properties of quasi-random sequences while counteracting the bias
  % some applications experience when using quasi-random sequences.
  % Initialize the scrambled quasi-random generator.
  [iref, state] = initialize_scrambled_qrng(int64(1), int64(2), ...
                                            idim, state);

  % Generate the quasi-random points from N(0,1)
  xmean = zeros(idim, 1);
  std   = ones(idim, 1);
  [z, iref, ifail] = g05yj( ...
                            xmean, std, npaths, iref);
  z = z';

function [state] = initialize_prng(genid, subid, seed)
  % Initialize the generator to a repeatable sequence
  [state, ifail] = g05kf( ...
                          genid, subid, seed);


function [iref, state] = initialize_scrambled_qrng(genid,stype,idim,state)
  iskip = int64(0);
  nsdigits = int64(32);
  [iref, state, ifail] = g05yn( ...
                                genid, stype, int64(idim), ...
                                iskip, nsdigits, state);
g05xd example results


Euler-Maruyama solution for Geometric Brownian motion
        Path 1     Path 2     Path 3
 1      0.9668     1.0367     0.9992
 2      0.9717     1.0254     1.0077
 3      0.9954     1.0333     1.0098
 4      0.9486     1.0226     0.9911
 5      0.9270     1.0113     1.0630
 6      0.8997     1.0127     1.0164
 7      0.8955     1.0138     1.0771
 8      0.8953     0.9953     1.0691
 9      0.8489     1.0462     1.0484
10      0.8449     1.0592     1.0429
11      0.8158     1.0233     1.0625
12      0.7997     1.0384     1.0729
13      0.8025     1.0138     1.0725
14      0.8187     1.0499     1.0554
15      0.8270     1.0459     1.0529
16      0.7914     1.0294     1.0783
17      0.8076     1.0224     1.0943
18      0.8208     1.0359     1.0773
19      0.8190     1.0326     1.0857
20      0.8217     1.0326     1.1095
21      0.8084     0.9695     1.1389

Analytic solution at final time step
        Path 1     Path 2     Path 3
        0.8079     0.9685     1.1389

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015