Actual source code: ex2.c
1: static char help[] = "Tests PetscRandom functions.\n\n";
3: #include <stdlib.h>
4: #include <stdio.h>
5: #include <string.h>
6: #include <time.h>
7: #include <sys/types.h>
9: #include "petsc.h"
10: #include "petscsys.h"
12: #define MAXBSIZE 40
13: #define PI 3.1415926535897
14: #define DATAFILENAME "ex2_stock.txt"
16: struct himaInfoTag {
17: int n;
18: double r;
19: double dt;
20: int totalNumSim;
21: double *St0;
22: double *vol;
23: };
24: typedef struct himaInfoTag himaInfo;
26: /* function protype */
27: PetscErrorCode readData(MPI_Comm comm,himaInfo *hinfo);
28: double mcVal(double St, double r, double vol, double dt, double eps);
29: void exchange(double *a, double *b);
30: double basketPayoff(double vol[], double St0[], int n, double r,double dt, double eps[]);
31: void stdNormalArray(double *eps, unsigned long size,PetscRandom ran);
32: unsigned long divWork(int id, unsigned long num, int np);
34: /*
35: Contributed by Xiaoyan Zeng <zengxia@iit.edu> and Liu, Kwong Ip" <kiliu@math.hkbu.edu.hk>
37: Example of usage:
38: mpirun -np 4 ./ex2 -num_of_stocks 30 -interest_rate 0.4 -time_interval 0.01 -num_of_simulations 10000
39: */
43: int main(int argc, char *argv[])
44: {
45: double r,dt;
46: int n;
47: unsigned long i,myNumSim,totalNumSim,numdim;
48: double payoff;
49: double *vol, *St0, x, totalx;
50: int np,myid;
51: time_t start,stop;
52: double *eps;
53: himaInfo hinfo;
54: PetscRandom ran;
57: PetscInitialize(&argc,&argv,(char *)0,help);
58: time(&start);
59: PetscRandomCreate(PETSC_COMM_WORLD,&ran);
60: #if defined(PETSC_HAVE_SPRNG)
61: PetscRandomSetType(ran,SPRNG);
62: #elif defined(PETSC_HAVE_RAND)
63: PetscRandomSetType(ran,PETSCRAND);
64: #endif
65: PetscRandomSetFromOptions(ran);
68: MPI_Comm_size(PETSC_COMM_WORLD, &np); /* number of nodes */
69: MPI_Comm_rank(PETSC_COMM_WORLD, &myid); /* my ranking */
70:
71: hinfo.n = 31;
72: hinfo.r = 0.04;
73: hinfo.dt = 1.0/12; /* a month as a period */
74: hinfo.totalNumSim = 1000;
75: PetscOptionsGetInt(PETSC_NULL,"-num_of_stocks",&(hinfo.n),PETSC_NULL);
76: if (hinfo.n <1 || hinfo.n > 31) SETERRQ1(PETSC_ERR_SUP,"Only 31 stocks listed in stock.txt. num_of_stocks %D must between 1 and 31",hinfo.n);
77: PetscOptionsGetReal(PETSC_NULL,"-interest_rate",&(hinfo.r),PETSC_NULL);
78: PetscOptionsGetReal(PETSC_NULL,"-time_interval",&(hinfo.dt),PETSC_NULL);
79: PetscOptionsGetInt(PETSC_NULL,"-num_of_simulations",&(hinfo.totalNumSim),PETSC_NULL);
81: n = hinfo.n;
82: r = hinfo.r;
83: dt = hinfo.dt;
84: totalNumSim = hinfo.totalNumSim;
85: vol = hinfo.vol = (double *)malloc(sizeof(double)*(2*n+1));
86: St0 = hinfo.St0 = hinfo.vol + n;
87: readData(PETSC_COMM_WORLD,&hinfo);
89: numdim = n*(n+1)/2;
90: if (numdim%2 == 1){
91: numdim++;
92: }
93: eps = (double *)malloc(sizeof(double)*numdim);
95: myNumSim = divWork(myid,totalNumSim,np);
97: x = 0;
98: for (i=0;i<myNumSim;i++){
99: stdNormalArray(eps,numdim,ran);
100: x += basketPayoff(vol,St0,n,r,dt,eps);
101: }
103: MPI_Reduce(&x, &totalx, 1, MPI_DOUBLE, MPI_SUM,0,PETSC_COMM_WORLD);
104: time(&stop);
105: if (myid == 0){
106: payoff = exp(-r*dt*n)*(totalx/totalNumSim);
107: PetscPrintf(PETSC_COMM_SELF,"Option price = $%.3f "
108: "using %ds of %s computation with %d %s "
109: "for %d stocks, %d trading period per year, "
110: "%.2f%% interest rate\n",
111: payoff,(int)(stop - start),"parallel",np,"processors",n,
112: (int)(1/dt),r);
113: }
114:
115: free(vol);
116: free(eps);
117: PetscRandomDestroy(ran);
118: PetscFinalize();
119: return 0;
120: }
122: void stdNormalArray(double *eps, unsigned long size, PetscRandom ran)
123: {
124: int i;
125: double u1,u2,t;
128: for (i=0;i<size;i+=2){
129: PetscRandomGetValue(ran,&u1);
130: PetscRandomGetValue(ran,&u2);
131:
132: t = sqrt(-2*log(u1));
133: eps[i] = t * cos(2*PI*u2);
134: eps[i+1] = t * sin(2*PI*u2);
135: }
136: }
139: double basketPayoff(double vol[], double St0[], int n, double r,double dt, double eps[])
140: {
141: double Stk[MAXBSIZE], temp;
142: double payoff;
143: int maxk,i,j;
144: int pointcount=0;
145:
146: for (i=0;i<n;i++) {
147: Stk[i] = St0[i];
148: }
150: for (i=0;i<n;i++){
151: maxk = 0;
152: for (j=0;j<(n-i);j++){
153: Stk[j] = mcVal(Stk[j],r,vol[j],dt,eps[pointcount++]);
154: if ((Stk[j]/St0[j]) > (Stk[maxk]/St0[maxk])){
155: maxk = j;
156: }
157: }
158: exchange(Stk+j-1,Stk+maxk);
159: exchange(St0+j-1,St0+maxk);
160: exchange(vol+j-1,vol+maxk);
161: }
162:
163: payoff = 0;
164: for (i=0;i<n;i++){
165: temp = (Stk[i]/St0[i]) - 1 ;
166: if (temp > 0) payoff += temp;
167: }
168: return payoff;
169: }
173: PetscErrorCode readData(MPI_Comm comm,himaInfo *hinfo)
174: {
175: int i;
176: FILE *fd;
177: char temp[50];
179: PetscMPIInt rank;
180: double *v = hinfo->vol, *t = hinfo->St0;
181: int num=hinfo->n;
182:
184: MPI_Comm_rank(comm,&rank);
185: if (!rank){
186: PetscFOpen(PETSC_COMM_SELF,DATAFILENAME,"r",&fd);
187: for (i=0;i<num;i++){
188: fscanf(fd,"%s%lf%lf",temp,v+i,t+i);
189: }
190: fclose(fd);
191: }
192: MPI_Bcast(v,2*num,MPI_DOUBLE,0,PETSC_COMM_WORLD);
193: //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]);
194: return(0);
195: }
197: void exchange(double *a, double *b)
198: {
199: double t;
200:
201: t = *a;
202: *a = *b;
203: *b = t;
204: }
206: double mcVal(double St, double r, double vol, double dt, double eps)
207: {
208: return (St * exp((r-0.5*vol*vol)*dt + vol*sqrt(dt)*eps));
209: }
211: unsigned long divWork(int id, unsigned long num, int np)
212: {
213: unsigned long numit;
215: numit = (unsigned long)(((double)num)/np);
216: numit++;
217: return numit;
218: }