xref: /petsc/src/sys/classes/random/tutorials/ex2.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests PetscRandom functions.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscsys.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown #define PETSC_MAXBSIZE     40
6*c4762a1bSJed Brown #define DATAFILENAME "ex2_stock.txt"
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown struct himaInfoTag {
9*c4762a1bSJed Brown   PetscInt    n;
10*c4762a1bSJed Brown   PetscReal   r;
11*c4762a1bSJed Brown   PetscReal   dt;
12*c4762a1bSJed Brown   PetscInt    totalNumSim;
13*c4762a1bSJed Brown   PetscReal   *St0;
14*c4762a1bSJed Brown   PetscReal   *vol;
15*c4762a1bSJed Brown };
16*c4762a1bSJed Brown typedef struct himaInfoTag himaInfo;
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown PetscErrorCode readData(MPI_Comm,himaInfo *);
19*c4762a1bSJed Brown PetscReal mcVal(PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
20*c4762a1bSJed Brown void exchangeVal(PetscReal*, PetscReal*);
21*c4762a1bSJed Brown PetscReal basketPayoff(PetscReal[], PetscReal[], PetscInt, PetscReal,PetscReal, PetscReal[]);
22*c4762a1bSJed Brown PetscErrorCode stdNormalArray(PetscReal*, PetscInt,PetscRandom);
23*c4762a1bSJed Brown PetscInt divWork(PetscMPIInt, PetscInt, PetscMPIInt);
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown /*
26*c4762a1bSJed Brown    Contributed by Xiaoyan Zeng <zengxia@iit.edu> and Liu, Kwong Ip" <kiliu@math.hkbu.edu.hk>
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown    Example of usage:
29*c4762a1bSJed Brown      mpiexec -n 4 ./ex2 -num_of_stocks 30 -interest_rate 0.4 -time_interval 0.01 -num_of_simulations 10000
30*c4762a1bSJed Brown */
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown int main(int argc, char *argv[])
33*c4762a1bSJed Brown {
34*c4762a1bSJed Brown   PetscReal      r,dt;
35*c4762a1bSJed Brown   PetscInt       n;
36*c4762a1bSJed Brown   unsigned long  i,myNumSim,totalNumSim,numdim;
37*c4762a1bSJed Brown   PetscReal      *vol, *St0, x, totalx;
38*c4762a1bSJed Brown   PetscMPIInt    size,rank;
39*c4762a1bSJed Brown   PetscReal      *eps;
40*c4762a1bSJed Brown   himaInfo       hinfo;
41*c4762a1bSJed Brown   PetscRandom    ran;
42*c4762a1bSJed Brown   PetscErrorCode ierr;
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
45*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&ran);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(ran);CHKERRQ(ierr);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);       /* number of nodes */
49*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);     /* my ranking */
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   hinfo.n           = 31;
52*c4762a1bSJed Brown   hinfo.r           = 0.04;
53*c4762a1bSJed Brown   hinfo.dt          = 1.0/12;   /* a month as a period */
54*c4762a1bSJed Brown   hinfo.totalNumSim = 1000;
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-num_of_stocks",&(hinfo.n),NULL);CHKERRQ(ierr);
57*c4762a1bSJed Brown   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);
58*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-interest_rate",&(hinfo.r),NULL);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-time_interval",&(hinfo.dt),NULL);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-num_of_simulations",&(hinfo.totalNumSim),NULL);CHKERRQ(ierr);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown   n           = hinfo.n;
63*c4762a1bSJed Brown   r           = hinfo.r;
64*c4762a1bSJed Brown   dt          = hinfo.dt;
65*c4762a1bSJed Brown   totalNumSim = hinfo.totalNumSim;
66*c4762a1bSJed Brown   ierr        = PetscMalloc1(2*n+1,&hinfo.vol);CHKERRQ(ierr);
67*c4762a1bSJed Brown   vol         = hinfo.vol;
68*c4762a1bSJed Brown   St0         = hinfo.St0 = hinfo.vol + n;
69*c4762a1bSJed Brown   ierr        = readData(PETSC_COMM_WORLD,&hinfo);CHKERRQ(ierr);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   numdim = n*(n+1)/2;
72*c4762a1bSJed Brown   if (numdim%2 == 1) numdim++;
73*c4762a1bSJed Brown   ierr = PetscMalloc1(numdim,&eps);CHKERRQ(ierr);
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown   myNumSim = divWork(rank,totalNumSim,size);
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   x = 0;
78*c4762a1bSJed Brown   for (i=0; i<myNumSim; i++) {
79*c4762a1bSJed Brown     ierr = stdNormalArray(eps,numdim,ran);CHKERRQ(ierr);
80*c4762a1bSJed Brown     x += basketPayoff(vol,St0,n,r,dt,eps);
81*c4762a1bSJed Brown   }
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown   ierr = MPI_Reduce(&x, &totalx, 1, MPIU_REAL, MPIU_SUM,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
84*c4762a1bSJed Brown   /* payoff = exp(-r*dt*n)*(totalx/totalNumSim);
85*c4762a1bSJed Brown   ierr = 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",
86*c4762a1bSJed Brown    payoff,(int)(stop - start),"parallel",size,"processors",n,(int)(1/dt),r);CHKERRQ(ierr); */
87*c4762a1bSJed Brown 
88*c4762a1bSJed Brown   ierr = PetscFree(vol);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = PetscFree(eps);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&ran);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = PetscFinalize();
92*c4762a1bSJed Brown   return ierr;
93*c4762a1bSJed Brown }
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown PetscErrorCode stdNormalArray(PetscReal *eps, PetscInt numdim, PetscRandom ran)
96*c4762a1bSJed Brown {
97*c4762a1bSJed Brown   PetscInt       i;
98*c4762a1bSJed Brown   PetscScalar    u1,u2;
99*c4762a1bSJed Brown   PetscReal      t;
100*c4762a1bSJed Brown   PetscErrorCode ierr;
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown   PetscFunctionBegin;
103*c4762a1bSJed Brown   for (i=0; i<numdim; i+=2) {
104*c4762a1bSJed Brown     ierr = PetscRandomGetValue(ran,&u1);CHKERRQ(ierr);
105*c4762a1bSJed Brown     ierr = PetscRandomGetValue(ran,&u2);CHKERRQ(ierr);
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown     t        = PetscSqrtReal(-2*PetscLogReal(PetscRealPart(u1)));
108*c4762a1bSJed Brown     eps[i]   = t * PetscCosReal(2*PETSC_PI*PetscRealPart(u2));
109*c4762a1bSJed Brown     eps[i+1] = t * PetscSinReal(2*PETSC_PI*PetscRealPart(u2));
110*c4762a1bSJed Brown   }
111*c4762a1bSJed Brown   PetscFunctionReturn(0);
112*c4762a1bSJed Brown }
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown PetscReal basketPayoff(PetscReal vol[], PetscReal St0[], PetscInt n, PetscReal r,PetscReal dt, PetscReal eps[])
116*c4762a1bSJed Brown {
117*c4762a1bSJed Brown   PetscReal Stk[PETSC_MAXBSIZE], temp;
118*c4762a1bSJed Brown   PetscReal payoff;
119*c4762a1bSJed Brown   PetscInt  maxk,i,j;
120*c4762a1bSJed Brown   PetscInt  pointcount=0;
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown   for (i=0;i<n;i++) Stk[i] = St0[i];
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown   for (i=0;i<n;i++) {
125*c4762a1bSJed Brown     maxk = 0;
126*c4762a1bSJed Brown     for (j=0;j<(n-i);j++) {
127*c4762a1bSJed Brown       Stk[j] = mcVal(Stk[j],r,vol[j],dt,eps[pointcount++]);
128*c4762a1bSJed Brown       if ((Stk[j]/St0[j]) > (Stk[maxk]/St0[maxk])) maxk = j;
129*c4762a1bSJed Brown     }
130*c4762a1bSJed Brown     exchangeVal(Stk+j-1,Stk+maxk);
131*c4762a1bSJed Brown     exchangeVal(St0+j-1,St0+maxk);
132*c4762a1bSJed Brown     exchangeVal(vol+j-1,vol+maxk);
133*c4762a1bSJed Brown   }
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown   payoff = 0;
136*c4762a1bSJed Brown   for (i=0; i<n; i++) {
137*c4762a1bSJed Brown     temp = (Stk[i]/St0[i]) - 1;
138*c4762a1bSJed Brown     if (temp > 0) payoff += temp;
139*c4762a1bSJed Brown   }
140*c4762a1bSJed Brown   return payoff;
141*c4762a1bSJed Brown }
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown PetscErrorCode readData(MPI_Comm comm,himaInfo *hinfo)
144*c4762a1bSJed Brown {
145*c4762a1bSJed Brown   PetscInt       i;
146*c4762a1bSJed Brown   FILE           *fd;
147*c4762a1bSJed Brown   char           temp[50];
148*c4762a1bSJed Brown   PetscErrorCode ierr;
149*c4762a1bSJed Brown   PetscMPIInt    rank;
150*c4762a1bSJed Brown   PetscReal      *v = hinfo->vol, *t = hinfo->St0;
151*c4762a1bSJed Brown   PetscInt       num=hinfo->n;
152*c4762a1bSJed Brown 
153*c4762a1bSJed Brown   PetscFunctionBeginUser;
154*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
155*c4762a1bSJed Brown   if (!rank) {
156*c4762a1bSJed Brown     ierr = PetscFOpen(PETSC_COMM_SELF,DATAFILENAME,"r",&fd);CHKERRQ(ierr);
157*c4762a1bSJed Brown     for (i=0;i<num;i++) {
158*c4762a1bSJed Brown       double vv,tt;
159*c4762a1bSJed Brown       if (fscanf(fd,"%s%lf%lf",temp,&vv,&tt) != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file\n");
160*c4762a1bSJed Brown       v[i] = vv;
161*c4762a1bSJed Brown       t[i] = tt;
162*c4762a1bSJed Brown     }
163*c4762a1bSJed Brown     fclose(fd);
164*c4762a1bSJed Brown   }
165*c4762a1bSJed Brown   ierr = MPI_Bcast(v,2*num,MPIU_REAL,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
166*c4762a1bSJed Brown   /* ierr = 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]); */
167*c4762a1bSJed Brown   PetscFunctionReturn(0);
168*c4762a1bSJed Brown }
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown void exchangeVal(PetscReal *a, PetscReal *b)
171*c4762a1bSJed Brown {
172*c4762a1bSJed Brown   PetscReal t;
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown   t  = *a;
175*c4762a1bSJed Brown   *a = *b;
176*c4762a1bSJed Brown   *b = t;
177*c4762a1bSJed Brown }
178*c4762a1bSJed Brown 
179*c4762a1bSJed Brown PetscReal mcVal(PetscReal St, PetscReal r, PetscReal vol, PetscReal dt, PetscReal eps)
180*c4762a1bSJed Brown {
181*c4762a1bSJed Brown   return (St * PetscExpReal((r-0.5*vol*vol)*dt + vol*PetscSqrtReal(dt)*eps));
182*c4762a1bSJed Brown }
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown PetscInt divWork(PetscMPIInt id, PetscInt num, PetscMPIInt size)
185*c4762a1bSJed Brown {
186*c4762a1bSJed Brown   PetscInt numit;
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   numit = (PetscInt)(((PetscReal)num)/size);
189*c4762a1bSJed Brown   numit++;
190*c4762a1bSJed Brown   return numit;
191*c4762a1bSJed Brown }
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown /*TEST
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown    build:
196*c4762a1bSJed Brown       requires: !comple
197*c4762a1bSJed Brown       output_file: output/ex1_1.out
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown    test:
200*c4762a1bSJed Brown       nsize: 2
201*c4762a1bSJed Brown       output_file: output/ex1_1.out
202*c4762a1bSJed Brown       localrunfiles: ex2_stock.txt
203*c4762a1bSJed Brown 
204*c4762a1bSJed Brown TEST*/
205