00001 #include <iostream>
00002 #include <fstream>
00003 #include <iomanip>
00004 #include <sstream>
00005 #include <string>
00006 #include <math.h>
00007 #include "random.h"
00008 #include "parametermap.h"
00009 #include "cheezyparser.h"
00010 #include "constants.h"
00011
00012 #define NTYPES 6
00013
00014 #undef USE_OLD_INPUT_SCHEME
00015
00016 using namespace std;
00017
00018 string flowprof("transverse");
00019
00020 void getCoreMom(double ans[4], const double &temp, const double &mass, const double vflow[4]);
00021 void getCorePos(double ans[4], const double &Rx, const double &Ry, const double &Rz, const double &tau);
00022 void getFlowVelocity(double ans[4], const double r[4], const double R, const double vmax=1./sqrt(3.));
00023 void getHaloPos(double ans[4], const double r[4], const double p[4], const double &mass_w, const double &tau_w);
00024 void getHaloMom(double ans[4], const double p_res[4], const double mass, const double mass_res);
00025 bool passesCut( const double p[4], const double etamin, const double etamax,
00026 const double ymin, const double ymax, const double pTmin, const double pTmax);
00027 double dabs(const double& x){if (x>=0) {return x;} else {return -x;}}
00028
00029
00031 void getHelp(void){
00032 cout<<"\nUsage: chum <options> [inputFile.dat]"<<endl;
00033 cout<<" -h, -help, --help : Prints this message"<<endl;
00034 exit(0);
00035 }
00036
00037 CRandom Randomizer(197);
00038
00050 int main(int argc, char* argv[]){
00051
00052 cout << "*** CHUM: the Core-Halo Ur-Model ***"<<endl;
00053 cout << endl;
00054
00055 MESSAGE << CMessage::warning;
00056
00057
00058 if (argc==1) getHelp();
00059 string paramFile("");
00060 vector<string> modeList;
00061 for (int iarg = 1; iarg<argc; ++iarg){
00062 string sarg(argv[iarg]);
00063 if (sarg=="-help") getHelp();
00064 if (sarg=="--help") getHelp();
00065 if (sarg=="-h") getHelp();
00066 if (sarg.substr(0,1)=="-") modeList.push_back(sarg);
00067 else paramFile = sarg;
00068 }
00069
00070
00071 double p[4],r[4],vflow[4];
00072 double mass[NTYPES]={938.3,939.6,139.58,139.58,134.9766,493.7};
00073 int ident[NTYPES]={2212,2112,211,-211,111,321};
00074 double T, Rx, Ry, Rz, tau_fo, tau_w, fraction;
00075 double etamin,etamax,ymax,ymin,pTmax,pTmin;
00076 double R_Au, mass_w;
00077 int id,iblankline,iseed;
00078 string filename;
00079 int NUMPARTICLES = 100000;
00080
00081
00082 tau_w = 23.0;
00083 R_Au = 1.2*pow(197.,1./3.);
00084 id = 2;
00085 mass_w= 782.6;
00086
00087
00088 etamin=0.0;
00089 ymin=0.0;
00090 pTmin=0.0;
00091 etamax=200.0;
00092 ymax=200.0;
00093 pTmax=2e6;
00094
00095
00096 #ifdef USE_OLD_INPUT_SCHEME
00097 cout << "\n";
00098 cout << "What is the temperature? (in MeV)\n";
00099 cin >> T;
00100 cout << "What are Rx, Ry, Rz of the core and tau_fo of the rxn? (in fm)\n";
00101 cin >> Rx >> Ry >> Rz >> tau_fo;
00102 cout << "What fraction of particles come from the core? \n";
00103 cin >> fraction;
00104 #else
00105 if (argc==1){
00106 cout<<"\nUsage: chum [input parameter file]"<<endl;
00107 exit(false);
00108 }
00109 ifstream pmapFile(argv[1]);
00110 if (!pmapFile.good()){
00111 cerr<< "Cannot open file "<<argv[1]<<endl;
00112 exit(-1);
00113 }
00114 parameterMap pmap;
00115 pmapFile>>pmap;
00116 T = parameter::getD(pmap,"T",175.0);
00117 Rx = parameter::getD(pmap,"Rx",4.0);
00118 Ry = parameter::getD(pmap,"Ry",4.0);
00119 Rz = parameter::getD(pmap,"Rz",4.0);
00120 tau_fo = parameter::getD(pmap,"tau_fo",10.0);
00121 fraction = parameter::getD(pmap,"f",1.0);
00122 id = parameter::getI(pmap,"id",id);
00123 flowprof = parameter::getS(pmap,"flowprofile",flowprof);
00124 etamin = parameter::getD(pmap,"etamin",etamin);
00125 etamax = parameter::getD(pmap,"etamax",etamax);
00126 ymin = parameter::getD(pmap,"ymin",ymin);
00127 ymax = parameter::getD(pmap,"ymax",ymax);
00128 pTmin = parameter::getD(pmap,"pTmin",pTmin);
00129 pTmax = parameter::getD(pmap,"pTmax",pTmax);
00130 filename = parameter::getS(pmap,"filename",(string)"thermal_gauss01.dat");
00131 iseed = parameter::getI(pmap,"initial_seed",197);
00132 NUMPARTICLES = parameter::getI(pmap,"num_particles",100000);
00133 #endif
00134
00135 Randomizer.reset(iseed);
00136
00137 cout << "Core radii are " << Rx << " fm, " << Ry << " fm, " << Rz << " fm\n";
00138 cout << "Core freezeout duration is " << tau_fo << " fm/c" <<endl;
00139 cout << "Core fraction is " << fraction << "\n";
00140 cout << "Particle type is " << ident[id] << " with mass " << mass[id] <<" MeVc^2"<<endl;
00141 cout << "Accepting particles w/ "<< etamin << " < eta < "<<etamax<<endl;
00142 cout << "Accepting particles w/ "<< ymin << " < y < "<<ymax<<endl;
00143 cout << "Accepting particles w/ "<< pTmin << " < pT < "<<pTmax<<" MeV"<<endl;
00144 cout << "Output filename "<< filename <<endl;
00145 cout << "Random seed "<< iseed <<endl;
00146 cout << "Number of particles generated "<< NUMPARTICLES <<endl;
00147
00148 if (flowprof==(string)"noflow") cout << "Flow turned off" << endl;
00149 else if (flowprof==(string)"radial") cout << "Radial flow" << endl;
00150 else if (flowprof==(string)"transverse") cout << "Transverse flow only" << endl;
00151 else if (flowprof==(string)"bjorken") cout << "Longitudinal flow only" << endl;
00152 else {
00153 cerr << "Illegal flow profile:" << flowprof << endl;
00154 exit(false);
00155 }
00156
00157 ofstream fout(filename.c_str());
00158 for(iblankline=0;iblankline<3;iblankline++){
00159 fout << "line number " << iblankline << ", blah blah blah\n";
00160 }
00161 fout << 1 << " "
00162 << NUMPARTICLES << " "
00163 << setw(9) << setprecision(6) << 0.0 << " "
00164 << setw(9) << setprecision(6) << 0.0 << "\n";
00165
00166 int i=0;
00167 while (i<NUMPARTICLES) {
00168 getCorePos(r,Rx,Ry,Rz,tau_fo);
00169 getFlowVelocity(vflow,r,R_Au);
00170 if (Randomizer.ran()<=fraction){
00171 getCoreMom(p,T,mass[id],vflow);
00172 } else {
00173 getCoreMom(p,T,mass_w,vflow);
00174 getHaloPos(r,r,p,mass_w,tau_w);
00175 getHaloMom(p,p,mass[id],mass_w);
00176 }
00177 if (passesCut(p,etamin,etamax,ymin,ymax,pTmin,pTmax)){
00178
00179
00180
00181 fout << i+1 << " "<< ident[id] << " "
00182 << setw(9) << p[1]/1000.0 << " "
00183 << setw(9) << p[2]/1000.0 << " "
00184 << setw(9) << p[3]/1000.0 << " "
00185 << setw(9) << p[0]/1000.0 << " "
00186 << setw(9) << mass[id]/1000.0 << " "
00187 << setw(9) << r[1] << " "
00188 << setw(9) << r[2] << " "
00189 << setw(9) << r[3] << " "
00190 << setw(9) << r[0] << "\n";
00191 i++;
00192 }
00193 }
00194 return 0;
00195 }
00196
00197
00198
00199
00200
00201
00202
00203 void getCoreMom(double ans[4], const double &temp, const double &mass, const double vflow[4])
00204 {
00205 double Tx(0.), x(0.), mu(0.), gammainv2(1.), Txm(0.);
00206 double vmag = sqrt(vflow[1]*vflow[1]+vflow[2]*vflow[2]+vflow[3]*vflow[3]);
00207 x = Randomizer.ran_exp();
00208 mu=2.*Randomizer.ran()-1.;
00209 Tx = temp*x;
00210 Txm=Tx+mass;
00211 gammainv2 = 1.0-vmag*vmag*mu*mu;
00212 double pmag=(vmag*mu*Txm+sqrt(Txm*Txm-gammainv2*mass*mass))/gammainv2;
00213
00214 double sthet=sqrt(1.0-mu*mu);
00215 double phi=2.0*PI*Randomizer.ran();
00216
00217 ans[0]=sqrt(mass*mass+pmag*pmag);
00218 ans[1]=pmag*sthet*cos(phi);
00219 ans[2]=pmag*sthet*sin(phi);
00220 ans[3]=pmag*mu;
00221 }
00222
00223
00224
00225
00226
00227 void getCorePos(double ans[4], const double &Rx, const double &Ry, const double &Rz, const double &tau)
00228 {
00229 double tmp;
00230 ans[0] = Randomizer.ran_exp();
00231 Randomizer.gauss2(&(ans[1]),&tmp);
00232 Randomizer.gauss2(&(ans[2]),&(ans[3]));
00233 ans[0]=tau*ans[0];
00234 ans[1]=Rx*ans[1];
00235 ans[2]=Ry*ans[2];
00236 ans[3]=Rz*ans[3];
00237 }
00238
00239
00240
00241
00242
00243 void getFlowVelocity(double ans[4], const double r[4], const double R, const double vmax)
00244 {
00245 ans[0]=1.0; for (int i=1;i<4;++i)ans[i]=0.0;
00246 double alpha = vmax/R;
00247
00248 if (flowprof=="radial") {
00249
00250
00251 double rmag = sqrt(r[1]*r[1]+r[2]*r[2]+r[3]*r[3]);
00252 if (rmag>=R) { for (int i=1;i<4;++i){ans[i] = vmax;} }
00253 else{ for (int i=1;i<4;++i){ans[i] = alpha*r[i];} }
00254 }
00255 if (flowprof=="transverse") {
00256
00257
00258 double rT = sqrt(r[1]*r[1]+r[2]*r[2]);
00259 if (rT>=R) { for (int i=1;i<3;++i){ans[i] = vmax;} }
00260 else{ for (int i=1;i<3;++i){ans[i] = alpha*r[i];} }
00261 ans[3]=0.0;
00262 }
00263 if (flowprof=="bjorken") {
00264
00265
00266 if (r[3]>=R) { ans[3] = vmax;}
00267 else { ans[3] = alpha*r[3]; }
00268 ans[1]=0.0; ans[2]=0.0;
00269 }
00270
00271 }
00272
00273
00274
00275
00276
00277 void getHaloPos(double ans[4], const double r[4], const double p[4], const double &mass_w, const double &tau_w)
00278 {
00279 double t = Randomizer.ran_exp();
00280 for (int i=0;i<4;++i){ans[i]=r[i]+p[i]*tau_w*t/p[0];}
00281 }
00282
00283
00284
00285
00286
00287 void getHaloMom(double ans[4], const double p_res[4], const double mass, const double mass_res)
00288 {
00289
00290 int iparticle(static_cast<int>(10*Randomizer.ran()/3));
00291
00292
00293 double M = mass_res;
00294 double M1 = mass;
00295 double M2 = mass;
00296 double M3 = mass;
00297 double ES1, ES2, P1, P2, Cos12, Z;
00298
00299 do {
00300
00301 do {
00302 ES1 = Randomizer.ran() * (M - M2 - M3 - M1) + M1;
00303 ES2 = Randomizer.ran() * (M - M1 - M3 - M2) + M2;
00304 } while (ES1+ES2 > M);
00305
00306 P1 = sqrt(ES1*ES1 - M1*M1);
00307 P2 = sqrt(ES2*ES2 - M2*M2);
00308 Z = M - ES1 - ES2;
00309 Z *= Z;
00310 Cos12 = (Z - P1*P1 - P2*P2 - M3*M3)/(2*P1*P2);
00311 } while ((Cos12 < -1.0) || (Cos12 > 1.0));
00312
00313 double Pxr2 = P2 * sqrt(1-Cos12*Cos12);
00314 double Pzr2 = P2*Cos12;
00315 double Pxr3 = - Pxr2;
00316 double Pzr3 = - (P1 + Pzr2);
00317 double P3 = sqrt(Pxr3*Pxr3+Pzr3*Pzr3);
00318 double ES3 = sqrt(M3*M3+P3*P3);
00319
00320
00321 double Phi = Randomizer.ran() * 2 * PI;
00322 double Ksi = Randomizer.ran() * 2 * PI;
00323 double CosTh = Randomizer.ran() * 2.0 - 1.0;
00324
00325 double sp = sin(Phi);
00326 double cp = cos(Phi);
00327 double sk = sin(Ksi);
00328 double ck = cos(Ksi);
00329 double st = sqrt(1.0-CosTh*CosTh);
00330 double ct = CosTh;
00331
00332
00333 double Pxp1 = - st*ck * P1;
00334 double Pyp1 = st*sk * P1;
00335 double Pzp1 = ct * P1;
00336
00337 double Pxp2 = (cp*ct*ck - sp*sk) * Pxr2 + (-st*ck) * Pzr2;
00338 double Pyp2 = (-cp*ct*sk - sp*ck) * Pxr2 + (st*sk) * Pzr2;
00339 double Pzp2 = cp*st * Pxr2 + ct * Pzr2;
00340
00341 double Pxp3 = (cp*ct*ck - sp*sk) * Pxr3 + (-st*ck) * Pzr3;
00342 double Pyp3 = (-cp*ct*sk - sp*ck) * Pxr3 + (st*sk) * Pzr3;
00343 double Pzp3 = cp*st * Pxr3 + ct * Pzr3;
00344
00345 double Vx = p_res[1]/p_res[0];
00346 double Vy = p_res[2]/p_res[0];
00347 double Vz = p_res[3]/p_res[0];
00348
00349 ES1 = sqrt(M1*M1+Pxp1*Pxp1+Pyp1*Pyp1+Pzp1*Pzp1);
00350 ES2 = sqrt(M2*M2+Pxp2*Pxp2+Pyp2*Pyp2+Pzp2*Pzp2);
00351 ES3 = sqrt(M3*M3+Pxp3*Pxp3+Pyp3*Pyp3+Pzp3*Pzp3);
00352
00353 double V2 = Vx*Vx + Vy*Vy + Vz*Vz;
00354 double Gamma = 1.0/sqrt(1-V2);
00355
00356
00357 double VP = Vx*Pxp1 + Vy*Pyp1 + Vz*Pzp1;
00358 double gvp = (Gamma - 1.0) * (1.0/V2) * VP;
00359
00360 double Px1 = Pxp1 + (gvp + Gamma * ES1) * Vx;
00361 double Py1 = Pyp1 + (gvp + Gamma * ES1) * Vy;
00362 double Pz1 = Pzp1 + (gvp + Gamma * ES1) * Vz;
00363
00364 VP = Vx*Pxp2 + Vy*Pyp2 + Vz*Pzp2;
00365 gvp = (Gamma - 1.0) * (1.0/V2) * VP;
00366
00367 double Px2 = Pxp2 + (gvp + Gamma * ES2) * Vx;
00368 double Py2 = Pyp2 + (gvp + Gamma * ES2) * Vy;
00369 double Pz2 = Pzp2 + (gvp + Gamma * ES2) * Vz;
00370
00371 VP = Vx*Pxp3 + Vy*Pyp3 + Vz*Pzp3;
00372 gvp = (Gamma - 1.0) * (1.0/V2) * VP;
00373
00374 double Px3 = Pxp3 + (gvp + Gamma * ES3) * Vx;
00375 double Py3 = Pyp3 + (gvp + Gamma * ES3) * Vy;
00376 double Pz3 = Pzp3 + (gvp + Gamma * ES3) * Vz;
00377
00378 ES1 = sqrt(M1*M1+Px1*Px1+Py1*Py1+Pz1*Pz1);
00379 ES2 = sqrt(M2*M2+Px2*Px2+Py2*Py2+Pz2*Pz2);
00380 ES3 = sqrt(M3*M3+Px3*Px3+Py3*Py3+Pz3*Pz3);
00381
00382 if (iparticle==1) {
00383 ans[0]=ES1;ans[1]=Px1;ans[2]=Py1;ans[3]=Pz1;
00384 } else if (iparticle==2) {
00385 ans[0]=ES2;ans[1]=Px2;ans[2]=Py2;ans[3]=Pz2;
00386 } else {
00387 ans[0]=ES3;ans[1]=Px3;ans[2]=Py3;ans[3]=Pz3;
00388 }
00389
00390 }
00391
00392
00393
00394
00395
00396 bool passesCut( const double p[4], const double etamin, const double etamax, const double ymin, const double ymax, const double pTmin, const double pTmax)
00397 {
00398 double y=0.5*log((p[0]+p[3])/(p[0]-p[3]));
00399 double pT=sqrt(p[1]*p[1]+p[2]*p[2]);
00400
00401
00402
00403 bool goody = (dabs(y)>=ymin)&&(dabs(y)<=ymax);
00404 bool goodpT = (pT>=pTmin)&&(pT<=pTmax);
00405 bool goodeta= (dabs(y)>=etamin)&&(dabs(y)<=etamax);
00406
00407
00408
00409 return goody&&goodpT&&goodeta;
00410 }
00411