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