Actual source code: ex2.c

petsc-master 2016-09-24
Report Typos and Errors
  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 exchange(PetscReal*, PetscReal*);
 21: PetscReal basketPayoff(PetscReal[], PetscReal[], PetscInt, PetscReal,PetscReal, PetscReal[]);
 22: void 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: */

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

 46:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 47:   PetscRandomCreate(PETSC_COMM_WORLD,&ran);
 48:   PetscRandomSetFromOptions(ran);

 50:   MPI_Comm_size(PETSC_COMM_WORLD, &size);       /* number of nodes */
 51:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);     /* my ranking */

 53:   hinfo.n           = 31;
 54:   hinfo.r           = 0.04;
 55:   hinfo.dt          = 1.0/12;   /* a month as a period */
 56:   hinfo.totalNumSim = 1000;

 58:   PetscOptionsGetInt(NULL,NULL,"-num_of_stocks",&(hinfo.n),NULL);
 59:   if (hinfo.n <1 || hinfo.n > 31) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only 31 stocks listed in stock.txt. num_of_stocks %D must between 1 and 31",hinfo.n);
 60:   PetscOptionsGetReal(NULL,NULL,"-interest_rate",&(hinfo.r),NULL);
 61:   PetscOptionsGetReal(NULL,NULL,"-time_interval",&(hinfo.dt),NULL);
 62:   PetscOptionsGetInt(NULL,NULL,"-num_of_simulations",&(hinfo.totalNumSim),NULL);

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

 73:   numdim = n*(n+1)/2;
 74:   if (numdim%2 == 1) numdim++;
 75:   PetscMalloc1(numdim,&eps);

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

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

 85:   MPI_Reduce(&x, &totalx, 1, MPIU_REAL, MPIU_SUM,0,PETSC_COMM_WORLD);
 86:   /* payoff = exp(-r*dt*n)*(totalx/totalNumSim);
 87:   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",
 88:    payoff,(int)(stop - start),"parallel",size,"processors",n,(int)(1/dt),r); */

 90:   PetscFree(vol);
 91:   PetscFree(eps);
 92:   PetscRandomDestroy(&ran);
 93:   PetscFinalize();
 94:   return ierr;
 95: }

 97: void stdNormalArray(PetscReal *eps, PetscInt numdim, PetscRandom ran)
 98: {
 99:   PetscInt       i;
100:   PetscScalar    u1,u2,t;

103:   for (i=0; i<numdim; i+=2) {
104:     PetscRandomGetValue(ran,&u1);CHKERRABORT(PETSC_COMM_WORLD,ierr);
105:     PetscRandomGetValue(ran,&u2);CHKERRABORT(PETSC_COMM_WORLD,ierr);

107:     t        = PetscSqrtReal(-2*PetscLogReal(PetscRealPart(u1)));
108:     eps[i]   = t * PetscCosReal(2*PETSC_PI*PetscRealPart(u2));
109:     eps[i+1] = t * PetscSinReal(2*PETSC_PI*PetscRealPart(u2));
110:   }
111: }


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

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

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

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

144: PetscErrorCode readData(MPI_Comm comm,himaInfo *hinfo)
145: {
146:   PetscInt       i;
147:   FILE           *fd;
148:   char           temp[50];
150:   PetscMPIInt    rank;
151:   PetscReal      *v = hinfo->vol, *t = hinfo->St0;
152:   PetscInt       num=hinfo->n;

155:   MPI_Comm_rank(comm,&rank);
156:   if (!rank) {
157:     PetscFOpen(PETSC_COMM_SELF,DATAFILENAME,"r",&fd);
158:     for (i=0;i<num;i++) {
159:       double vv,tt;
160:       if (fscanf(fd,"%s%lf%lf",temp,&vv,&tt) != 3) SETERRQ(PETSC_COMM_SELF,1,"Badly formatted input file\n");
161:       v[i] = vv;
162:       t[i] = tt;
163:     }
164:     fclose(fd);
165:   }
166:   MPI_Bcast(v,2*num,MPIU_REAL,0,PETSC_COMM_WORLD);
167:   /* 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]); */
168:   return(0);
169: }

171: void exchange(PetscReal *a, PetscReal *b)
172: {
173:   PetscReal t;

175:   t  = *a;
176:   *a = *b;
177:   *b = t;
178: }

180: PetscReal mcVal(PetscReal St, PetscReal r, PetscReal vol, PetscReal dt, PetscReal eps)
181: {
182:   return (St * PetscExpReal((r-0.5*vol*vol)*dt + vol*PetscSqrtReal(dt)*eps));
183: }

185: PetscInt divWork(PetscMPIInt id, PetscInt num, PetscMPIInt size)
186: {
187:   PetscInt numit;

189:   numit = (PetscInt)(((PetscReal)num)/size);
190:   numit++;
191:   return numit;
192: }