Actual source code: ex2.c

  1: static char help[] = "Tests PetscRandom functions.\n\n";

  3: #include <petscsys.h>

  5: #define PETSC_MAXBSIZE     40
  6: #define DATAFILENAME "ex2_stock.txt"

  8: struct himaInfoTag {
  9:   PetscInt    n;
 10:   PetscReal   r;
 11:   PetscReal   dt;
 12:   PetscInt    totalNumSim;
 13:   PetscReal   *St0;
 14:   PetscReal   *vol;
 15: };
 16: typedef struct himaInfoTag himaInfo;

 18: PetscErrorCode readData(MPI_Comm,himaInfo *);
 19: PetscReal mcVal(PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
 20: void exchangeVal(PetscReal*, PetscReal*);
 21: PetscReal basketPayoff(PetscReal[], PetscReal[], PetscInt, PetscReal,PetscReal, PetscReal[]);
 22: PetscErrorCode stdNormalArray(PetscReal*, PetscInt,PetscRandom);
 23: PetscInt divWork(PetscMPIInt, PetscInt, PetscMPIInt);

 25: /*
 26:    Contributed by Xiaoyan Zeng <zengxia@iit.edu> and Liu, Kwong Ip" <kiliu@math.hkbu.edu.hk>

 28:    Example of usage:
 29:      mpiexec -n 4 ./ex2 -num_of_stocks 30 -interest_rate 0.4 -time_interval 0.01 -num_of_simulations 10000
 30: */

 32: int main(int argc, char *argv[])
 33: {
 34:   PetscReal      r,dt;
 35:   PetscInt       n;
 36:   unsigned long  i,myNumSim,totalNumSim,numdim;
 37:   PetscReal      *vol, *St0, x, totalx;
 38:   PetscMPIInt    size,rank;
 39:   PetscReal      *eps;
 40:   himaInfo       hinfo;
 41:   PetscRandom    ran;

 43:   PetscInitialize(&argc,&argv,(char*)0,help);
 44:   PetscRandomCreate(PETSC_COMM_WORLD,&ran);
 45:   PetscRandomSetFromOptions(ran);

 47:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 48:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 50:   hinfo.n           = 31;
 51:   hinfo.r           = 0.04;
 52:   hinfo.dt          = 1.0/12;   /* a month as a period */
 53:   hinfo.totalNumSim = 1000;

 55:   PetscOptionsGetInt(NULL,NULL,"-num_of_stocks",&(hinfo.n),NULL);
 57:   PetscOptionsGetReal(NULL,NULL,"-interest_rate",&(hinfo.r),NULL);
 58:   PetscOptionsGetReal(NULL,NULL,"-time_interval",&(hinfo.dt),NULL);
 59:   PetscOptionsGetInt(NULL,NULL,"-num_of_simulations",&(hinfo.totalNumSim),NULL);

 61:   n           = hinfo.n;
 62:   r           = hinfo.r;
 63:   dt          = hinfo.dt;
 64:   totalNumSim = hinfo.totalNumSim;
 65:   PetscMalloc1(2*n+1,&hinfo.vol);
 66:   vol         = hinfo.vol;
 67:   St0         = hinfo.St0 = hinfo.vol + n;
 68:   readData(PETSC_COMM_WORLD,&hinfo);

 70:   numdim = n*(n+1)/2;
 71:   if (numdim%2 == 1) numdim++;
 72:   PetscMalloc1(numdim,&eps);

 74:   myNumSim = divWork(rank,totalNumSim,size);

 76:   x = 0;
 77:   for (i=0; i<myNumSim; i++) {
 78:     stdNormalArray(eps,numdim,ran);
 79:     x += basketPayoff(vol,St0,n,r,dt,eps);
 80:   }

 82:   MPI_Reduce(&x, &totalx, 1, MPIU_REAL, MPIU_SUM,0,PETSC_COMM_WORLD);
 83:   /* payoff = exp(-r*dt*n)*(totalx/totalNumSim);
 84:   PetscPrintf(PETSC_COMM_WORLD,"Option price = $%.3f using %ds of %s computation with %d %s for %d stocks, %d trading period per year, %.2f%% interest rate\n",
 85:    payoff,(int)(stop - start),"parallel",size,"processors",n,(int)(1/dt),r); */

 87:   PetscFree(vol);
 88:   PetscFree(eps);
 89:   PetscRandomDestroy(&ran);
 90:   PetscFinalize();
 91:   return 0;
 92: }

 94: PetscErrorCode stdNormalArray(PetscReal *eps, PetscInt numdim, PetscRandom ran)
 95: {
 96:   PetscInt       i;
 97:   PetscScalar    u1,u2;
 98:   PetscReal      t;

100:   for (i=0; i<numdim; i+=2) {
101:     PetscRandomGetValue(ran,&u1);
102:     PetscRandomGetValue(ran,&u2);

104:     t        = PetscSqrtReal(-2*PetscLogReal(PetscRealPart(u1)));
105:     eps[i]   = t * PetscCosReal(2*PETSC_PI*PetscRealPart(u2));
106:     eps[i+1] = t * PetscSinReal(2*PETSC_PI*PetscRealPart(u2));
107:   }
108:   return 0;
109: }

111: PetscReal basketPayoff(PetscReal vol[], PetscReal St0[], PetscInt n, PetscReal r,PetscReal dt, PetscReal eps[])
112: {
113:   PetscReal Stk[PETSC_MAXBSIZE], temp;
114:   PetscReal payoff;
115:   PetscInt  maxk,i,j;
116:   PetscInt  pointcount=0;

118:   for (i=0;i<n;i++) Stk[i] = St0[i];

120:   for (i=0;i<n;i++) {
121:     maxk = 0;
122:     for (j=0;j<(n-i);j++) {
123:       Stk[j] = mcVal(Stk[j],r,vol[j],dt,eps[pointcount++]);
124:       if ((Stk[j]/St0[j]) > (Stk[maxk]/St0[maxk])) maxk = j;
125:     }
126:     exchangeVal(Stk+j-1,Stk+maxk);
127:     exchangeVal(St0+j-1,St0+maxk);
128:     exchangeVal(vol+j-1,vol+maxk);
129:   }

131:   payoff = 0;
132:   for (i=0; i<n; i++) {
133:     temp = (Stk[i]/St0[i]) - 1;
134:     if (temp > 0) payoff += temp;
135:   }
136:   return payoff;
137: }

139: PetscErrorCode readData(MPI_Comm comm,himaInfo *hinfo)
140: {
141:   PetscInt       i;
142:   FILE           *fd;
143:   char           temp[50];
144:   PetscMPIInt    rank;
145:   PetscReal      *v = hinfo->vol, *t = hinfo->St0;
146:   PetscInt       num=hinfo->n;

149:   MPI_Comm_rank(comm,&rank);
150:   if (rank == 0) {
151:     PetscFOpen(PETSC_COMM_SELF,DATAFILENAME,"r",&fd);
152:     for (i=0;i<num;i++) {
153:       double vv,tt;
155:       v[i] = vv;
156:       t[i] = tt;
157:     }
158:     fclose(fd);
159:   }
160:   MPI_Bcast(v,2*num,MPIU_REAL,0,PETSC_COMM_WORLD);
161:   /* PetscPrintf(PETSC_COMM_SELF,"[%d] vol %g, ... %g; St0 %g, ... %g\n",rank,hinfo->vol[0],hinfo->vol[num-1],hinfo->St0 [0],hinfo->St0[num-1]); */
162:   return 0;
163: }

165: void exchangeVal(PetscReal *a, PetscReal *b)
166: {
167:   PetscReal t;

169:   t  = *a;
170:   *a = *b;
171:   *b = t;
172: }

174: PetscReal mcVal(PetscReal St, PetscReal r, PetscReal vol, PetscReal dt, PetscReal eps)
175: {
176:   return (St * PetscExpReal((r-0.5*vol*vol)*dt + vol*PetscSqrtReal(dt)*eps));
177: }

179: PetscInt divWork(PetscMPIInt id, PetscInt num, PetscMPIInt size)
180: {
181:   PetscInt numit;

183:   numit = (PetscInt)(((PetscReal)num)/size);
184:   numit++;
185:   return numit;
186: }

188: /*TEST

190:    test:
191:       nsize: 2
192:       output_file: output/ex1_1.out
193:       localrunfiles: ex2_stock.txt

195: TEST*/