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: mpiexec -n 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;
56: PetscTruth flg;
58: PetscInitialize(&argc,&argv,(char *)0,help);
59: time(&start);
60: PetscRandomCreate(PETSC_COMM_WORLD,&ran);
61: #if defined(PETSC_HAVE_SPRNG)
62: PetscRandomSetType(ran,PETSCSPRNG);
63: #elif defined(PETSC_HAVE_RAND)
64: PetscRandomSetType(ran,PETSCRAND);
65: #endif
66: PetscRandomSetFromOptions(ran);
69: MPI_Comm_size(PETSC_COMM_WORLD, &np); /* number of nodes */
70: MPI_Comm_rank(PETSC_COMM_WORLD, &myid); /* my ranking */
72: PetscOptionsHasName(PETSC_NULL, "-check_generators", &flg);
73: if (flg){
74: PetscRandomGetValue(ran,&r);
75: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] rval: %g\n",myid,r);
76: PetscSynchronizedFlush(PETSC_COMM_WORLD);
77: }
78:
79: hinfo.n = 31;
80: hinfo.r = 0.04;
81: hinfo.dt = 1.0/12; /* a month as a period */
82: hinfo.totalNumSim = 1000;
83: PetscOptionsGetInt(PETSC_NULL,"-num_of_stocks",&(hinfo.n),PETSC_NULL);
84: 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);
85: PetscOptionsGetReal(PETSC_NULL,"-interest_rate",&(hinfo.r),PETSC_NULL);
86: PetscOptionsGetReal(PETSC_NULL,"-time_interval",&(hinfo.dt),PETSC_NULL);
87: PetscOptionsGetInt(PETSC_NULL,"-num_of_simulations",&(hinfo.totalNumSim),PETSC_NULL);
89: n = hinfo.n;
90: r = hinfo.r;
91: dt = hinfo.dt;
92: totalNumSim = hinfo.totalNumSim;
93: vol = hinfo.vol = (double *)malloc(sizeof(double)*(2*n+1));
94: St0 = hinfo.St0 = hinfo.vol + n;
95: readData(PETSC_COMM_WORLD,&hinfo);
97: numdim = n*(n+1)/2;
98: if (numdim%2 == 1){
99: numdim++;
100: }
101: eps = (double *)malloc(sizeof(double)*numdim);
103: myNumSim = divWork(myid,totalNumSim,np);
105: x = 0;
106: for (i=0;i<myNumSim;i++){
107: stdNormalArray(eps,numdim,ran);
108: x += basketPayoff(vol,St0,n,r,dt,eps);
109: }
111: MPI_Reduce(&x, &totalx, 1, MPI_DOUBLE, MPI_SUM,0,PETSC_COMM_WORLD);
112: time(&stop);
113: if (myid == 0){
114: payoff = exp(-r*dt*n)*(totalx/totalNumSim);
115: PetscPrintf(PETSC_COMM_SELF,"Option price = $%.3f "
116: "using %ds of %s computation with %d %s "
117: "for %d stocks, %d trading period per year, "
118: "%.2f%% interest rate\n",
119: payoff,(int)(stop - start),"parallel",np,"processors",n,
120: (int)(1/dt),r);
121: }
122:
123: free(vol);
124: free(eps);
125: PetscRandomDestroy(ran);
126: PetscFinalize();
127: return 0;
128: }
130: void stdNormalArray(double *eps, unsigned long size, PetscRandom ran)
131: {
132: int i;
133: double u1,u2,t;
136: for (i=0;i<size;i+=2){
137: PetscRandomGetValue(ran,&u1);
138: PetscRandomGetValue(ran,&u2);
139:
140: t = sqrt(-2*log(u1));
141: eps[i] = t * cos(2*PI*u2);
142: eps[i+1] = t * sin(2*PI*u2);
143: }
144: }
147: double basketPayoff(double vol[], double St0[], int n, double r,double dt, double eps[])
148: {
149: double Stk[MAXBSIZE], temp;
150: double payoff;
151: int maxk,i,j;
152: int pointcount=0;
153:
154: for (i=0;i<n;i++) {
155: Stk[i] = St0[i];
156: }
158: for (i=0;i<n;i++){
159: maxk = 0;
160: for (j=0;j<(n-i);j++){
161: Stk[j] = mcVal(Stk[j],r,vol[j],dt,eps[pointcount++]);
162: if ((Stk[j]/St0[j]) > (Stk[maxk]/St0[maxk])){
163: maxk = j;
164: }
165: }
166: exchange(Stk+j-1,Stk+maxk);
167: exchange(St0+j-1,St0+maxk);
168: exchange(vol+j-1,vol+maxk);
169: }
170:
171: payoff = 0;
172: for (i=0;i<n;i++){
173: temp = (Stk[i]/St0[i]) - 1 ;
174: if (temp > 0) payoff += temp;
175: }
176: return payoff;
177: }
181: PetscErrorCode readData(MPI_Comm comm,himaInfo *hinfo)
182: {
183: int i;
184: FILE *fd;
185: char temp[50];
187: PetscMPIInt rank;
188: double *v = hinfo->vol, *t = hinfo->St0;
189: int num=hinfo->n;
190:
192: MPI_Comm_rank(comm,&rank);
193: if (!rank){
194: PetscFOpen(PETSC_COMM_SELF,DATAFILENAME,"r",&fd);
195: for (i=0;i<num;i++){
196: fscanf(fd,"%s%lf%lf",temp,v+i,t+i);
197: }
198: fclose(fd);
199: }
200: MPI_Bcast(v,2*num,MPI_DOUBLE,0,PETSC_COMM_WORLD);
201: //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]);
202: return(0);
203: }
205: void exchange(double *a, double *b)
206: {
207: double t;
208:
209: t = *a;
210: *a = *b;
211: *b = t;
212: }
214: double mcVal(double St, double r, double vol, double dt, double eps)
215: {
216: return (St * exp((r-0.5*vol*vol)*dt + vol*sqrt(dt)*eps));
217: }
219: unsigned long divWork(int id, unsigned long num, int np)
220: {
221: unsigned long numit;
223: numit = (unsigned long)(((double)num)/np);
224: numit++;
225: return numit;
226: }