CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/GeneratorInterface/CosmicMuonGenerator/src/CosmicMuonGenerator.cc

Go to the documentation of this file.
00001 
00002 // modified by P. Biallass 29.03.2006 to implement new cosmic generator (CMSCGEN.cc) and new normalization of flux (CMSCGENnorm.cc)
00004 // 04.12.2008 sonne: replaced Min/MaxE by Min/MaxP to get cos_sf/ug scripts working again
00005 // 20.04.2009 sonne: Implemented mechanism to read in multi muon events and propagate each muon
00006 #define sim_cxx
00007 
00008 
00009 #include "GeneratorInterface/CosmicMuonGenerator/interface/CosmicMuonGenerator.h"
00010 
00011 
00012 void CosmicMuonGenerator::runCMG(){
00013   initialize();
00014   for (unsigned int iGen=0; iGen<NumberOfEvents; ++iGen){ nextEvent(); }
00015   terminate();
00016 }
00017 
00018 void CosmicMuonGenerator::initialize(CLHEP::HepRandomEngine *rng){
00019   if (delRanGen)
00020     delete RanGen;
00021   if (!rng) {
00022     RanGen = new CLHEP::HepJamesRandom;
00023     RanGen->setSeed(RanSeed, 0); //set seed for Random Generator (seed can be controled by config-file)
00024     delRanGen = true;
00025   } else {
00026     RanGen = rng;
00027     delRanGen = false;
00028   }
00029   checkIn();
00030   if (NumberOfEvents > 0){
00031     // set up "surface geometry" dimensions
00032     double RadiusTargetEff = RadiusOfTarget; //get this from cfg-file
00033     double Z_DistTargetEff = ZDistOfTarget;  //get this from cfg-file
00034     //double Z_CentrTargetEff = ZCentrOfTarget;  //get this from cfg-file
00035     if(TrackerOnly==true){
00036     RadiusTargetEff = RadiusTracker;
00037     Z_DistTargetEff = Z_DistTracker;
00038     }
00039     Target3dRadius = sqrt(RadiusTargetEff*RadiusTargetEff + Z_DistTargetEff*Z_DistTargetEff) + MinStepSize;
00040     if (Debug) std::cout << "  radius of sphere  around  target = " << Target3dRadius << " mm" << std::endl;
00041 
00042     if (MinTheta > 90.*Deg2Rad) //upgoing muons from neutrinos
00043       SurfaceRadius = (RadiusCMS)*(-tan(MinTheta)) + MinStepSize;
00044     else
00045       SurfaceRadius = (SurfaceOfEarth+PlugWidth+RadiusTargetEff)*tan(MaxTheta) + Target3dRadius;  
00046     if (Debug) std::cout << "  starting point radius at Surface + PlugWidth = " << SurfaceRadius << " mm" << std::endl;
00047 
00048     OneMuoEvt.PlugVx = PlugVx;
00049     OneMuoEvt.PlugVz = PlugVz;
00050     OneMuoEvt.RhoAir = RhoAir;
00051     OneMuoEvt.RhoWall = RhoWall;
00052     OneMuoEvt.RhoRock = RhoRock;
00053     OneMuoEvt.RhoClay = RhoClay;
00054     OneMuoEvt.RhoPlug = RhoPlug;
00055     //set energy and angle limits for CMSCGEN, give same seed as above 
00056     if (MinTheta >= 90.*Deg2Rad) //upgoing muons from neutrinos
00057       Cosmics->initializeNuMu(MinP, MaxP, MinTheta, MaxTheta, MinEnu, MaxEnu, 
00058                               MinPhi, MaxPhi, NuProdAlt, RanGen);
00059     else 
00060       Cosmics->initialize(MinP, MaxP, MinTheta, MaxTheta, RanGen, TIFOnly_constant, TIFOnly_linear);
00061    
00062 #if ROOT_INTERACTIVE
00063   // book histos
00064   TH1D* ene = new TH1D("ene","generated energy",210,0.,1050.);
00065   TH1D* the = new TH1D("the","generated theta",90,0.,90.);
00066   TH1D* phi = new TH1D("phi","generated phi",120,0.,360.);
00067   TH3F* ver = new TH3F("ver","Z-X-Y coordinates",100,-25.,25.,20,-10.,10.,40,-10.,10.);
00068 #endif
00069     if (EventDisplay) initEvDis();
00070     std::cout << std::endl;
00071     
00072     if (MultiMuon) {
00073       MultiIn = 0;
00074       
00075       std::cout << "MultiMuonFileName.c_str()=" << MultiMuonFileName.c_str() << std::endl;
00076       MultiIn = new TFile( MultiMuonFileName.c_str() );
00077       
00078       if (!MultiIn) std::cout << "MultiMuon=True: MultiMuonFileName='"
00079                               << MultiMuonFileName.c_str() << "' does not exist" << std::endl;
00080       else std::cout << "MultiMuonFile: " << MultiMuonFileName.c_str() << " opened!" << std::endl;
00081       //MultiTree = (TTree*) gDirectory->Get("sim");
00082       MultiTree = (TTree*) MultiIn->Get("sim");
00083       SimTree = new sim(MultiTree);
00084       SimTree->Init(MultiTree);
00085       SimTreeEntries = SimTree->fChain->GetEntriesFast();
00086       std::cout << "SimTreeEntries=" << SimTreeEntries << std::endl;
00087 
00088       if (MultiMuonFileFirstEvent <= 0)
00089         SimTree_jentry = 0;
00090       else 
00091         SimTree_jentry = MultiMuonFileFirstEvent - 1; //1=1st evt (SimTree_jentry=0)
00092       
00093       NcloseMultiMuonEvents = 0;
00094       NskippedMultiMuonEvents = 0;
00095     }
00096     
00097     if (!MultiMuon || (MultiMuon && MultiIn)) NotInitialized = false;
00098      
00099   }
00100 }
00101 
00102 void CosmicMuonGenerator::nextEvent(){
00103 
00104   double E = 0.; double Theta = 0.; double Phi = 0.; double RxzV = 0.; double PhiV = 0.;
00105   if (int(Nsel)%100 == 0) std::cout << "    generated " << int(Nsel) << " events" << std::endl;
00106   // generate cosmic (E,theta,phi)
00107   bool   notSelected = true;
00108   while (notSelected){
00109     bool   badMomentumGenerated = true;
00110     while (badMomentumGenerated){
00111       
00112       if (MinTheta > 90.*Deg2Rad) //upgoing muons from neutrinos
00113         Cosmics->generateNuMu();
00114       else
00115         Cosmics->generate(); //dice one event now
00116       
00117       E = sqrt(Cosmics->momentum_times_charge()*Cosmics->momentum_times_charge() + MuonMass*MuonMass);
00118       Theta = TMath::ACos( Cosmics->cos_theta() ) ; //angle has to be in RAD here
00119           Ngen+=1.;   //count number of initial cosmic events (in surface area), vertices will be added later
00120           badMomentumGenerated = false;
00121           Phi = RanGen->flat()*(MaxPhi-MinPhi) + MinPhi;
00122     }
00123     Norm->events_n100cos(E, Theta); //test if this muon is in normalization range
00124     Ndiced += 1; //one more cosmic is diced
00125 
00126     // generate vertex
00127     double Nver = 0.;
00128     bool   badVertexGenerated = true;
00129     while (badVertexGenerated){
00130       RxzV = sqrt(RanGen->flat())*SurfaceRadius;
00131       PhiV = RanGen->flat()*TwoPi;
00132       // check phi range (for a sphere with Target3dRadius around the target)
00133       double dPhi = Pi; if (RxzV > Target3dRadius) dPhi = asin(Target3dRadius/RxzV);
00134       double rotPhi = PhiV + Pi; if (rotPhi > TwoPi) rotPhi -= TwoPi;
00135       double disPhi = std::fabs(rotPhi - Phi); if (disPhi > Pi) disPhi = TwoPi - disPhi;
00136       if (disPhi < dPhi || AcptAllMu) badVertexGenerated = false;
00137       Nver+=1.;
00138     }
00139     Ngen += (Nver-1.); //add number of generated vertices to initial cosmic events
00140     
00141     // complete event at surface
00142     int id =  13; // mu-
00143     if (Cosmics->momentum_times_charge() >0.) id = -13; // mu+
00144     double absMom = sqrt(E*E - MuonMass*MuonMass);
00145     double verMom = absMom*cos(Theta);
00146     double horMom = absMom*sin(Theta);
00147     double Px = horMom*sin(Phi); // [GeV/c]
00148     double Py = -verMom;         // [GeV/c]
00149     double Pz = horMom*cos(Phi); // [GeV/c]
00150     double Vx = RxzV*sin(PhiV);  // [mm]
00151 
00152     double Vy;
00153     if (MinTheta > 90.*Deg2Rad) //upgoing muons from neutrinos
00154       Vy = -RadiusCMS;
00155     else
00156       Vy = SurfaceOfEarth + PlugWidth;  // [mm]
00157 
00158     double Vz = RxzV*cos(PhiV);  // [mm]
00159     double T0 = (RanGen->flat()*(MaxT0-MinT0) + MinT0)*SpeedOfLight; // [mm/c];
00160 
00161     Id_at = id;
00162     Px_at = Px; Py_at = Py; Pz_at = Pz; E_at = E; //M_at = MuonMass;
00163     Vx_at = Vx; Vy_at = Vy; Vz_at = Vz; T0_at = T0;
00164 
00165     OneMuoEvt.create(id, Px, Py, Pz, E, MuonMass, Vx, Vy, Vz, T0); 
00166     // if angles are ok, propagate to target
00167     if (goodOrientation()) { 
00168       if (MinTheta > 90.*Deg2Rad) //upgoing muons from neutrinos
00169         OneMuoEvt.propagate(0., RadiusOfTarget, ZDistOfTarget, ZCentrOfTarget, TrackerOnly, MTCCHalf);
00170       else
00171         OneMuoEvt.propagate(ElossScaleFactor, RadiusOfTarget, ZDistOfTarget, ZCentrOfTarget, TrackerOnly, MTCCHalf);
00172     }
00173 
00174     if ( (OneMuoEvt.hitTarget() && sqrt(OneMuoEvt.e()*OneMuoEvt.e() - MuonMass*MuonMass) > MinP_CMS)
00175         || AcptAllMu==true){
00176       Nsel+=1.; //count number of generated and accepted events  
00177       notSelected = false;
00178     }
00179   }
00180 
00181   EventWeight = 1.;
00182 
00183   //just one outgoing particle at SurFace
00184   Id_sf.resize(1); 
00185   Px_sf.resize(1); 
00186   Py_sf.resize(1); 
00187   Pz_sf.resize(1);
00188   E_sf.resize(1); 
00189   //M_sf.resize(1);
00190   Vx_sf.resize(1);
00191   Vy_sf.resize(1);
00192   Vz_sf.resize(1);
00193   T0_sf.resize(1);
00194 
00195   Id_sf[0] = Id_at;
00196   Px_sf[0] = Px_at; Py_sf[0] = Py_at; Pz_sf[0] = Pz_at; E_sf[0] = E_at; //M_fs[0] = MuonMass;
00197   Vx_sf[0] = Vx_at; Vy_sf[0] = Vy_at; Vz_sf[0] = Vz_at; T0_sf[0] = T0_at;
00198   
00199 
00200   //just one particle at UnderGround  
00201   Id_ug.resize(1); 
00202   Px_ug.resize(1); 
00203   Py_ug.resize(1); 
00204   Pz_ug.resize(1);
00205   E_ug.resize(1); 
00206   //M_ug.resize(1);
00207   Vx_ug.resize(1);
00208   Vy_ug.resize(1);
00209   Vz_ug.resize(1);
00210   T0_ug.resize(1);
00211 
00212   Id_ug[0] = OneMuoEvt.id();
00213   Px_ug[0] = OneMuoEvt.px(); 
00214   Py_ug[0] = OneMuoEvt.py(); 
00215   Pz_ug[0] = OneMuoEvt.pz();
00216   E_ug[0] = OneMuoEvt.e(); 
00217   //M_ug[0] = OneMuoEvt.m();
00218   Vx_ug[0] = OneMuoEvt.vx(); 
00219   Vy_ug[0] = OneMuoEvt.vy(); 
00220   Vz_ug[0] = OneMuoEvt.vz();
00221   T0_ug[0] = OneMuoEvt.t0();
00222 
00223   // plot variables of selected events
00224 #if ROOT_INTERACTIVE
00225   ene->Fill(OneMuoEvt.e());
00226   the->Fill((OneMuoEvt.theta()*Rad2Deg));
00227   phi->Fill((OneMuoEvt.phi()*Rad2Deg));
00228   ver->Fill((OneMuoEvt.vz()/1000.),(OneMuoEvt.vx()/1000.),(OneMuoEvt.vy()/1000.));
00229 #endif
00230   if (Debug){
00231     std::cout << "new event" << std::endl;
00232     std::cout << "  Px,Py,Pz,E,m = " << OneMuoEvt.px() << ", " << OneMuoEvt.py() << ", "
00233          << OneMuoEvt.pz() << ", " << OneMuoEvt.e() << ", " << OneMuoEvt.m() << " GeV" << std::endl;
00234     std::cout << "  Vx,Vy,Vz,t0  = " << OneMuoEvt.vx() << ", " << OneMuoEvt.vy() << ", " 
00235          << OneMuoEvt.vz() << ", " << OneMuoEvt.t0() << " mm" << std::endl;
00236   }
00237   if (EventDisplay) displayEv();
00238   
00239 }
00240 
00241 
00242 
00243 bool CosmicMuonGenerator::nextMultiEvent() {
00244 
00245   if (Debug) std::cout << "\nEntered CosmicMuonGenerator::nextMultiEvent()" << std::endl;
00246   bool EvtRejected = true;
00247   bool MuInMaxDist = false;
00248   double MinDist; //[mm] 
00249 
00250   while (EvtRejected) {
00251 
00252     //read in event from SimTree
00253     //ULong64_t ientry = SimTree->LoadTree(SimTree_jentry);
00254     Long64_t ientry = SimTree->GetEntry(SimTree_jentry);
00255     std::cout << "CosmicMuonGenerator::nextMultiEvent(): SimTree_jentry=" << SimTree_jentry
00256       //<< " ientry=" << ientry 
00257               << " SimTreeEntries=" << SimTreeEntries << std::endl;
00258     if (ientry < 0) return false; //stop run
00259     if (SimTree_jentry < SimTreeEntries) {
00260       SimTree_jentry++;
00261     }
00262     else {
00263       std::cout << "CosmicMuonGenerator.cc::nextMultiEvent: No more events in file!" << std::endl;
00264       return false; //stop run
00265     }
00266     
00267 
00268 
00269     int nmuons = SimTree->shower_nParticlesWritten;
00270     if (nmuons<MultiMuonNmin) {
00271       std::cout << "CosmicMuonGenerator.cc: Warning!  Less than " << MultiMuonNmin << 
00272         " muons in event!" << std::endl;
00273       std::cout << "trying next event from file" << std::endl;
00274       NskippedMultiMuonEvents++;
00275       continue; //EvtRejected while loop: get next event from file
00276     }
00277 
00278 
00279 
00280     Px_mu.resize(nmuons); Py_mu.resize(nmuons); Pz_mu.resize(nmuons);
00281     P_mu.resize(nmuons);
00282 
00283     MinDist = 99999.e9; //[mm] 
00284     double MuMuDist;
00285     MuInMaxDist = false;
00286     //check if at least one muon pair closer than 30m at surface
00287     int NmuPmin = 0;
00288     for (int imu=0; imu<nmuons; ++imu) {
00289 
00290       Px_mu[imu] = -SimTree->particle__Px[imu]*sin(NorthCMSzDeltaPhi)
00291         + SimTree->particle__Py[imu]*cos(NorthCMSzDeltaPhi);
00292       Pz_mu[imu] = SimTree->particle__Px[imu]*cos(NorthCMSzDeltaPhi)
00293         + SimTree->particle__Py[imu]*sin(NorthCMSzDeltaPhi);
00294       Py_mu[imu] = -SimTree->particle__Pz[imu]; //Corsika down going particles defined in -z direction!
00295       P_mu[imu] = sqrt(Px_mu[imu]*Px_mu[imu] + Py_mu[imu]*Py_mu[imu] + Pz_mu[imu]*Pz_mu[imu]);
00296       
00297       if (P_mu[imu] < MinP_CMS  && AcptAllMu==false) continue;
00298       else if (SimTree->particle__ParticleID[imu] != 5 &&
00299                SimTree->particle__ParticleID[imu] != 6) continue;
00300       else NmuPmin++;
00301       
00302       for (int jmu=0; jmu<imu; ++jmu) {
00303         if (P_mu[jmu] < MinP_CMS && AcptAllMu==false) continue;
00304         if (SimTree->particle__ParticleID[imu] != 5 &&
00305             SimTree->particle__ParticleID[imu] != 6) continue;
00306         MuMuDist = sqrt( (SimTree->particle__x[imu]-SimTree->particle__x[jmu])*
00307                          (SimTree->particle__x[imu]-SimTree->particle__x[jmu]) 
00308                          +(SimTree->particle__y[imu]-SimTree->particle__y[jmu])*
00309                          (SimTree->particle__y[imu]-SimTree->particle__y[jmu])
00310                          )*10.; //CORSIKA [cm] to CMSCGEN [mm]
00311         if (MuMuDist < MinDist) MinDist = MuMuDist;
00312         if (MuMuDist < 2.*Target3dRadius) MuInMaxDist = true;
00313       }
00314     }
00315     if (MultiMuonNmin>=2) {
00316       if (MuInMaxDist) { 
00317         NcloseMultiMuonEvents++;
00318       }
00319       else {
00320         std::cout << "CosmicMuonGenerator.cc: Warning! No muon pair closer than " 
00321                   << 2.*Target3dRadius/1000. << "m   MinDist=" << MinDist/1000. << "m at surface" << std::endl;
00322         std::cout << "Fraction of too wide opening angle multi muon events: "
00323                   << 1 - double(NcloseMultiMuonEvents)/SimTree_jentry << std::endl;
00324         std::cout << "NcloseMultiMuonEvents=" << NcloseMultiMuonEvents << std::endl;
00325         std::cout << "trying next event from file" << std::endl;
00326         NskippedMultiMuonEvents++;
00327         continue; //EvtRejected while loop: get next event from file
00328       }
00329     }
00330   
00331     if (NmuPmin < MultiMuonNmin && AcptAllMu==false) { //take single muon events consistently into account
00332       NskippedMultiMuonEvents++;
00333       continue; //EvtRejected while loop: get next event from file
00334     }
00335 
00336     if (Debug) 
00337       if (MultiMuonNmin>=2)
00338         std::cout << "start trial do loop: MuMuDist=" << MinDist/1000. << "[m]   Nmuons=" 
00339                   << nmuons << "  NcloseMultiMuonEvents=" << NcloseMultiMuonEvents 
00340                   << "  NskippedMultiMuonEvents=" << NskippedMultiMuonEvents << std::endl;  
00341     
00342     
00343     //int primary_id = SimTree->run_ParticleID;
00344     Id_at = SimTree->shower_EventID;
00345     
00346     double M_at = 0.;
00347     //if (Id_at == 13) {
00348     Id_at = 2212; //convert from Corsika to HepPDT
00349     M_at = 938.272e-3; //[GeV] mass
00350     //}
00351     
00352     E_at = SimTree->shower_Energy;
00353     Theta_at = SimTree->shower_Theta;
00354     double phi_at = SimTree->shower_Phi - NorthCMSzDeltaPhi; //rotate by almost 90 degrees
00355     if (phi_at < -Pi) phi_at +=TwoPi; //bring into interval (-Pi,Pi]
00356     else if (phi_at > Pi) phi_at -= TwoPi; 
00357     double P_at = sqrt(E_at*E_at - M_at*M_at);
00358     //need to rotate about 90degrees around x->N axis => y<=>z, 
00359     //then rotate new x-z-plane from x->North to x->LHC centre
00360     Px_at = P_at*sin(Theta_at)*sin(phi_at);
00361     Py_at = -P_at*cos(Theta_at);
00362     Pz_at = P_at*sin(Theta_at)*cos(phi_at);
00363 
00364     //compute maximal theta of secondary muons
00365     double theta_mu_max = Theta_at;
00366     double theta_mu_min = Theta_at;
00367 
00368     double phi_rel_min = 0.; //phi_mu_min - phi_at
00369     double phi_rel_max = 0.; //phi_mu_max - phi_at
00370 
00371     //std::cout << "SimTree->shower_Energy=" << SimTree->shower_Energy <<std::endl;
00372 
00373     Theta_mu.resize(nmuons);
00374     for (int imu=0; imu<nmuons; ++imu) {
00375       Theta_mu[imu] = acos(-Py_mu[imu]/P_mu[imu]);
00376       if (Theta_mu[imu]>theta_mu_max) theta_mu_max = Theta_mu[imu];
00377       if (Theta_mu[imu]<theta_mu_min) theta_mu_min = Theta_mu[imu];
00378 
00379       double phi_mu = atan2(Px_mu[imu],Pz_mu[imu]); // in (-Pi,Pi]
00380       double phi_rel = phi_mu - phi_at;
00381       if (phi_rel < -Pi) phi_rel += TwoPi; //bring into interval (-Pi,Pi] 
00382       else if (phi_rel > Pi) phi_rel -= TwoPi;
00383       if (phi_rel < phi_rel_min) phi_rel_min = phi_rel;
00384       else if (phi_rel > phi_rel_max) phi_rel_max =phi_rel;
00385 
00386 
00387     }
00388     
00389     
00390     double h_sf = SurfaceOfEarth + PlugWidth; //[mm]
00391 
00392     double R_at = h_sf*tan(Theta_at);
00393     
00394     double JdRxzV_dR_trans = 1.; 
00395     double JdPhiV_dPhi_trans = 1.; 
00396     double JdR_trans_sqrt = 1.;
00397     
00398     //chose random vertex Phi and Rxz weighted to speed up and smoothen
00399     double R_mu_max = (h_sf+Target3dRadius)*tan(theta_mu_max);
00400     double R_max = std::min(SurfaceRadius, R_mu_max);
00401     double R_mu_min = (h_sf-Target3dRadius)*tan(theta_mu_min);
00402     double R_min = std::max(0., R_mu_min);
00403     
00404     if (R_at>SurfaceRadius) {
00405       std::cout << "CosmicMuonGenerator.cc: Warning! R_at=" << R_at 
00406            << " > SurfaceRadius=" << SurfaceRadius <<std::endl;
00407     }
00408     
00409     //do phase space transformation for horizontal radius R
00410     
00411     //determine first phase space limits
00412     
00413     double psR1min = R_min + 0.25*(R_max-R_min); 
00414     double psR1max = std::min(SurfaceRadius,R_max-0.25*(R_max-R_min)); //no R's beyond R_max
00415     double psR1 = psR1max - psR1min;
00416 
00417     double psR2min = R_min;
00418     double psR2max = R_max;
00419     double psR2 = psR2max - psR2min;
00420 
00421     double psR3min = 0.;
00422     double psR3max = SurfaceRadius;
00423     double psR3 = psR3max - psR3min;
00424 
00425     //double psall = psR1+psR2+psR3;
00426     double psRall = psR3;
00427 
00428     double fR1=psR1/psRall, fR2=psR2/psRall, fR3=psR3/psRall; //f1+f2+f3=130%
00429     double pR1=0.25, pR2=0.7, pR3=0.05;
00430 
00431 
00432     //do phase space transformation for azimuthal angle phi
00433     double psPh1 = 0.5*(phi_rel_max - phi_rel_min);
00434     double psPh2 = phi_rel_max - phi_rel_min;
00435     double psPh3 = TwoPi;
00436     double psPhall = psPh3;
00437     
00438     double fPh1=psPh1/psPhall, fPh2=psPh2/psPhall, fPh3=psPh3/psPhall; //(f1+f2+f3=TwoPi+f1+f2)/(TwoPi+f1+f2) 
00439 
00440     double pPh1=0.25, pPh2=0.7, pPh3=0.05;
00441 
00442     int temp = 0;
00443     bool btemp = false;
00444   
00445     Trials = 0; //global int trials
00446     double trials = 0.; //local weighted trials
00447     Vx_mu.resize(nmuons); Vy_mu.resize(nmuons); Vz_mu.resize(nmuons);
00448     int NmuHitTarget = 0;
00449     while (NmuHitTarget < MultiMuonNmin) { 
00450    
00451       NmuHitTarget = 0; //re-initialize every loop iteration
00452       double Nver = 0.;
00453 
00454       
00455       //chose phase space class
00456       double RxzV;
00457       double which_R_class = RanGen->flat();
00458       if (which_R_class < pR1) { //pR1% in psR1
00459         RxzV = psR1min + psR1 * RanGen->flat();
00460         JdRxzV_dR_trans = fR1/pR1 * SurfaceRadius/psR1;
00461       }
00462       else if (which_R_class < pR1+pR2) { //further pR2% in psR2 
00463         RxzV = psR2min + psR2 * RanGen->flat();   
00464         JdRxzV_dR_trans = fR2/pR2 * SurfaceRadius/psR2;
00465       }
00466       else { //remaining pR3% in psR3=[0., R_max]
00467         RxzV = psR3min + psR3 * RanGen->flat();
00468         JdRxzV_dR_trans = fR3/pR3 * SurfaceRadius/psR3;
00469       }
00470       
00471       JdR_trans_sqrt = 2.*RxzV/SurfaceRadius; //flat in sqrt(r) space
00472 
00473       //chose phase space class
00474       double PhiV;
00475       double which_phi_class = RanGen->flat();
00476       if (which_phi_class < pPh1) { //pPh1% in psPh1 
00477         PhiV = phi_at + phi_rel_min + psPh1 * RanGen->flat();
00478         JdPhiV_dPhi_trans = fPh1/pPh1 * TwoPi/psPh1;
00479       }
00480       else if (which_phi_class < pPh1+pPh2) { //further pPh2% in psPh2
00481         PhiV = phi_at + phi_rel_min + psPh2 * RanGen->flat();     
00482         JdPhiV_dPhi_trans = fPh2/pPh2 * TwoPi/psPh2;
00483       }
00484       else { //remaining pPh3% in psPh3=[-Pi,Pi]
00485         PhiV = phi_at + phi_rel_min + psPh3 * RanGen->flat();
00486         JdPhiV_dPhi_trans = fPh3/pPh3 * TwoPi/psPh3;
00487       }
00488 
00489       //shuffle PhiV into [-Pi,+Pi] interval
00490       if (PhiV < -Pi) PhiV+=TwoPi;
00491       else if (PhiV > Pi) PhiV-=TwoPi;
00492       
00493       
00494       Nver++;
00495       trials += JdR_trans_sqrt * JdRxzV_dR_trans * JdPhiV_dPhi_trans;
00496       Trials++;
00497       if (trials > max_Trials) break; //while (Id_sf.size() < 2) loop 
00498       Ngen += (Nver-1.); //add number of generated vertices to initial cosmic events
00499          
00500    
00501       Vx_at = RxzV*sin(PhiV); // [mm]
00502 
00503       Vy_at = h_sf; // [mm] (SurfaceOfEarth + PlugWidth; Determine primary particle height below)
00504       //Vy_at = SimTree->shower_StartingAltitude*10. + h_sf; // [mm]
00505       //std::cout << "SimTree->shower_StartingAltitude*10=" << SimTree->shower_StartingAltitude*10 <<std::endl;
00506       Vz_at = RxzV*cos(PhiV); // [mm]
00507       
00508       int NmuHitTargetSphere = 0;
00509       for (int imu=0; imu<nmuons; ++imu) {
00510         
00511         Vx_mu[imu] = Vx_at + (-SimTree->particle__x[imu]*sin(NorthCMSzDeltaPhi)
00512                               +SimTree->particle__y[imu]*cos(NorthCMSzDeltaPhi) )*10; //[mm] (Corsika cm to CMSCGEN mm)
00513         Vy_mu[imu] = h_sf; //[mm] fixed at surface + PlugWidth
00514         Vz_mu[imu] = Vz_at + ( SimTree->particle__x[imu]*cos(NorthCMSzDeltaPhi)
00515                                +SimTree->particle__y[imu]*sin(NorthCMSzDeltaPhi) )*10; //[mm] (Corsika cm to CMSCGEN mm)
00516 
00517         
00518         //add atmospheric height to primary particle (default SimTree->shower_StartingAltitude = 0.)
00519         double pt_sec = sqrt(Px_mu[imu]*Px_mu[imu]+Pz_mu[imu]*Pz_mu[imu]);
00520         double theta_sec = atan2(std::fabs(Py_mu[imu]),pt_sec);
00521         double r_sec = sqrt((Vx_mu[imu]-Vx_at)*(Vx_mu[imu]-Vx_at)
00522                             +(Vz_mu[imu]-Vz_at)*(Vz_mu[imu]-Vz_at));
00523         double h_prod =  r_sec * tan(theta_sec);
00524         if (h_prod + h_sf > Vy_at) Vy_at = h_prod + h_sf;
00525 
00526         //only muons
00527         if (SimTree->particle__ParticleID[imu] != 5 &&
00528             SimTree->particle__ParticleID[imu] != 6) continue;
00529 
00530         if (P_mu[imu] < MinP_CMS && AcptAllMu==false) continue; 
00531 
00532         //check here if at least 2 muons make it to the target sphere
00533         double Vxz_mu = sqrt(Vx_mu[imu]*Vx_mu[imu] + Vz_mu[imu]*Vz_mu[imu]);
00534         theta_mu_max = atan((Vxz_mu+Target3dRadius)/(h_sf-Target3dRadius));
00535         theta_mu_min = atan((Vxz_mu-Target3dRadius)/(h_sf+Target3dRadius));
00536         if (Theta_mu[imu] > theta_mu_min && Theta_mu[imu] < theta_mu_max) {
00537           
00538           // check phi range (for a sphere with Target3dRadius around the target)
00539           double dPhi = Pi; if (Vxz_mu > Target3dRadius) dPhi = asin(Target3dRadius/Vxz_mu);
00540           double PhiPmu = atan2(Px_mu[imu],Pz_mu[imu]); //muon phi
00541           double PhiVmu = atan2(Vx_mu[imu],Vz_mu[imu]); //muon phi
00542           double rotPhi = PhiVmu + Pi; if (rotPhi > Pi) rotPhi -= TwoPi;
00543           double disPhi = std::fabs(rotPhi - PhiPmu); if (disPhi > Pi) disPhi = TwoPi - disPhi;
00544           if (disPhi < dPhi) {
00545             NmuHitTargetSphere++;
00546           }
00547 
00548         }
00549 
00550       } //end imu for loop
00551 
00552       
00553       if (temp > 0) btemp = true;
00554 
00555 
00556       if (NmuHitTargetSphere < MultiMuonNmin) continue; //while (Id_sf.size() < 2) loop
00557     
00558       //nmuons outgoing particle at SurFace
00559       Id_sf.clear(); 
00560       Px_sf.clear(); 
00561       Py_sf.clear(); 
00562       Pz_sf.clear();
00563       E_sf.clear(); 
00564       //M_sf_out.clear();
00565       Vx_sf.clear(); 
00566       Vy_sf.clear(); 
00567       Vz_sf.clear();
00568       T0_sf.clear();
00569       
00570       //nmuons particles at UnderGround  
00571       Id_ug.clear(); 
00572       Px_ug.clear(); 
00573       Py_ug.clear(); 
00574       Pz_ug.clear();
00575       E_ug.clear(); 
00576       //M_ug.clear();
00577       Vx_ug.clear();
00578       Vy_ug.clear();
00579       Vz_ug.clear();
00580       T0_ug.clear();
00581       
00582       int Id_sf_this =0;
00583       double Px_sf_this =0., Py_sf_this=0., Pz_sf_this=0.;
00584       double E_sf_this=0.;
00585       //double M_sf_this=0.;
00586       double Vx_sf_this=0., Vy_sf_this=0., Vz_sf_this=0.;
00587       double T0_sf_this=0.;
00588       
00589       T0_at = SimTree->shower_GH_t0 * SpeedOfLight; // [mm]
00590 
00591       for (int imu=0; imu<nmuons; ++imu) {
00592 
00593         if (P_mu[imu] < MinP_CMS && AcptAllMu==false) continue;
00594         //for the time being only muons
00595         if (SimTree->particle__ParticleID[imu] != 5 &&
00596             SimTree->particle__ParticleID[imu] != 6) continue;
00597         
00598         Id_sf_this = SimTree->particle__ParticleID[imu];
00599         if (Id_sf_this == 5) Id_sf_this = -13; //mu+
00600         else if (Id_sf_this == 6) Id_sf_this = 13; //mu-
00601 
00602         else if (Id_sf_this == 1) Id_sf_this = 22; //gamma
00603         else if (Id_sf_this == 2) Id_sf_this = -11; //e+
00604         else if (Id_sf_this == 3) Id_sf_this = 11; //e-
00605         else if (Id_sf_this == 7) Id_sf_this = 111; //pi0
00606         else if (Id_sf_this == 8) Id_sf_this = 211; //pi+
00607         else if (Id_sf_this == 9) Id_sf_this = -211; //pi-
00608         else if (Id_sf_this == 10) Id_sf_this = 130; //KL0
00609         else if (Id_sf_this == 11) Id_sf_this = 321; //K+
00610         else if (Id_sf_this == 12) Id_sf_this = -321; //K-
00611         else if (Id_sf_this == 13) Id_sf_this = 2112; //n
00612         else if (Id_sf_this == 14) Id_sf_this = 2212; //p
00613         else if (Id_sf_this == 15) Id_sf_this = -2212; //pbar
00614         else if (Id_sf_this == 16) Id_sf_this = 310; //Ks0
00615         else if (Id_sf_this == 17) Id_sf_this = 221; //eta
00616         else if (Id_sf_this == 18) Id_sf_this = 3122; //Lambda
00617         
00618         else {
00619           std::cout << "CosmicMuonGenerator.cc: Warning! Muon Id=" << Id_sf_this 
00620                << " from file read in" <<std::endl;
00621           Id_sf_this = 99999; //trouble
00622         }
00623         
00624         T0_sf_this = SimTree->particle__Time[imu] * SpeedOfLight; //in [mm]
00625 
00626         Px_sf_this = Px_mu[imu];
00627         Py_sf_this = Py_mu[imu]; //Corsika down going particles defined in -z direction!
00628         Pz_sf_this = Pz_mu[imu];
00629         E_sf_this = sqrt(P_mu[imu]*P_mu[imu] + MuonMass*MuonMass);
00630         Vx_sf_this = Vx_mu[imu];
00631         Vy_sf_this = Vy_mu[imu]; //[mm] fixed at surface + PlugWidth
00632         Vz_sf_this = Vz_mu[imu];
00633 
00634 
00635         OneMuoEvt.create(Id_sf_this, Px_sf_this, Py_sf_this, Pz_sf_this, E_sf_this, MuonMass, Vx_sf_this, Vy_sf_this, Vz_sf_this, T0_sf_this); 
00636         // if angles are ok, propagate to target
00637         if (goodOrientation()) {
00638           OneMuoEvt.propagate(ElossScaleFactor, RadiusOfTarget, ZDistOfTarget, ZCentrOfTarget, TrackerOnly, MTCCHalf);
00639         }
00640         
00641         if ( (OneMuoEvt.hitTarget() 
00642               && sqrt(OneMuoEvt.e()*OneMuoEvt.e() - MuonMass*MuonMass) > MinP_CMS)
00643              || AcptAllMu==true ) {
00644           
00645           Id_sf.push_back(Id_sf_this);
00646           Px_sf.push_back(Px_sf_this);
00647           Py_sf.push_back(Py_sf_this);
00648           Pz_sf.push_back(Pz_sf_this);
00649           E_sf.push_back(E_sf_this);
00650           //M_sf.push_back(M_sf_this);
00651           Vx_sf.push_back(Vx_sf_this);
00652           Vy_sf.push_back(Vy_sf_this);
00653           Vz_sf.push_back(Vz_sf_this);
00654           T0_sf.push_back(T0_sf_this);
00655           //T0_sf.push_back(0.); //synchronised arrival for 100% efficient full simulation tests
00656           
00657           Id_ug.push_back(OneMuoEvt.id());
00658           Px_ug.push_back(OneMuoEvt.px());
00659           Py_ug.push_back(OneMuoEvt.py());
00660           Pz_ug.push_back(OneMuoEvt.pz());
00661           E_ug.push_back(OneMuoEvt.e());
00662           //M_sf.push_back(OneMuoEvt.m());
00663           Vx_ug.push_back(OneMuoEvt.vx());
00664           Vy_ug.push_back(OneMuoEvt.vy());
00665           Vz_ug.push_back(OneMuoEvt.vz());
00666           T0_ug.push_back(OneMuoEvt.t0());
00667 
00668           NmuHitTarget++;
00669         }
00670       }
00671       
00672 
00673     } // while (Id_sf.size() < 2); //end of do loop
00674     
00675 
00676     if (trials > max_Trials) {
00677       std::cout << "CosmicMuonGenerator.cc: Warning! trials reach max_trials=" << max_Trials
00678                 << " without accepting event!" << std::endl;
00679       if (Debug) {
00680         std::cout << " N(mu)=" << Id_ug.size();
00681         if (Id_ug.size()>=1) 
00682           std::cout << " E[0]=" << E_ug[0] << " theta=" 
00683                     << acos(-Py_ug[0]/sqrt(Px_ug[0]*Px_ug[0]+Py_ug[0]*Py_ug[0]+Pz_ug[0]*Pz_ug[0]))
00684                     << " R_xz=" << sqrt(Vx_sf[0]*Vx_sf[0]+Vy_sf[0]*Vy_sf[0]);
00685         std::cout << " Theta_at=" << Theta_at << std::endl; 
00686       }
00687       std::cout << "Unweighted int num of Trials = " << Trials << std::endl;
00688       std::cout << "trying next event (" << SimTree_jentry << ") from file" << std::endl;
00689       NskippedMultiMuonEvents++;
00690       continue; //EvtRejected while loop: get next event from file
00691     }
00692     else {
00693       if (NmuHitTarget < MultiMuonNmin) {
00694         std::cout << "CosmicMuonGenerator.cc: Warning! less than " << MultiMuonNmin << 
00695           " muons hit target: N(mu=)" << NmuHitTarget << std::endl;
00696         std::cout << "trying next event (" << SimTree_jentry << ") from file" << std::endl;
00697         NskippedMultiMuonEvents++;
00698         continue; //EvtRejected while loop: get next event from file
00699       }
00700       else { //if (MuInMaxDist) {
00701 
00702         //re-adjust T0's of surviving muons shifted to trigger time box
00703         //(possible T0 increase due to propagation (loss of energy/speed + travelled distance))
00704         double T0_ug_min, T0_ug_max;
00705         T0_ug_min = T0_ug_max = T0_ug[0];
00706         double Tbox = (MaxT0 - MinT0) * SpeedOfLight; // [mm]
00707         double minDeltaT0 = 2*Tbox;
00708         for (unsigned int imu=0; imu<Id_ug.size(); ++imu) {
00709           double T0_this = T0_ug[imu];
00710           if (T0_this < T0_ug_min) T0_ug_min = T0_this;
00711           if (T0_this > T0_ug_max) T0_ug_max = T0_this;
00712           if (Debug) std::cout << "imu=" << imu << " T0_this=" << T0_this 
00713                                << " P=" << sqrt(pow(Px_ug[imu],2) + pow(Py_ug[imu],2) + pow(Pz_ug[imu],2)) 
00714                                << std::endl;
00715           for (unsigned int jmu=0; jmu<imu; ++jmu) {
00716             if (std::fabs(T0_ug[imu]-T0_ug[jmu]) < minDeltaT0) minDeltaT0 = std::fabs(T0_ug[imu]-T0_ug[jmu]);
00717           }
00718         }
00719         
00720         if (int(Id_ug.size()) >= MultiMuonNmin && MultiMuonNmin>=2 && minDeltaT0 > Tbox)
00721           continue; //EvtRejected while loop: get next event from file
00722 
00723         double T0_min = T0_ug_min +MinT0*SpeedOfLight; //-12.5ns * c [mm]
00724         double T0_max = T0_ug_max +MaxT0*SpeedOfLight; //+12.5ns * c [mm]
00725 
00726         //ckeck if >= NmuMin in time box, else throw new random number + augment evt weight
00727         int TboxTrials = 0;
00728         int NmuInTbox;
00729         double T0_offset, T0diff;
00730         for (int tboxtrial=0; tboxtrial<1000; ++tboxtrial) { //max 1000 trials
00731           T0_offset = RanGen->flat()*(T0_max -T0_min); // [mm]
00732           TboxTrials++;
00733           T0diff = T0_offset - T0_max; // [mm] 
00734           NmuInTbox = 0;
00735           for (unsigned int imu=0; imu<Id_ug.size(); ++imu) {
00736             if (T0_ug[imu]+T0diff > MinT0*SpeedOfLight && T0_ug[imu]+T0diff < MaxT0*SpeedOfLight)
00737               NmuInTbox++;
00738           }
00739           if (NmuInTbox >= MultiMuonNmin) break;
00740 
00741         }
00742         if (NmuInTbox < MultiMuonNmin) continue; //EvtRejected while loop: get next event from file
00743 
00744 
00745         if (Debug) std::cout << "initial T0_at=" << T0_at << " T0_min=" << T0_min << " T0_max=" << T0_max 
00746                         << " T0_offset=" << T0_offset;
00747         T0_at += T0diff; //[mm]
00748         if (Debug) std::cout << " T0diff=" << T0diff << std::endl;
00749         for (unsigned int imu=0; imu<Id_ug.size(); ++imu) { //adjust @ surface + underground
00750           if (Debug) std::cout << "before: T0_sf[" << imu << "]=" << T0_sf[imu] << "  T0_ug=" <<  T0_ug[imu]; 
00751           T0_sf[imu] += T0diff;
00752           T0_ug[imu] += T0diff;
00753           if (Debug) 
00754             std::cout << "   after: T0_sf[" << imu << "]=" << T0_sf[imu] << "  T0_ug=" <<  T0_ug[imu] << std::endl;
00755         }       
00756         if (Debug) std::cout << "T0diff=" << T0diff << " T0_at=" << T0_at << std::endl;
00757 
00758 
00759 
00760         Nsel += 1;
00761         EventWeight = JdR_trans_sqrt * JdRxzV_dR_trans * JdPhiV_dPhi_trans 
00762           / (trials * TboxTrials);
00763         EvtRejected = false;
00764         if (Debug) std::cout << "CosmicMuonGenerator.cc: Theta_at=" << Theta_at << " phi_at=" << phi_at 
00765                              << " Px_at=" << Px_at << " Py_at=" << Py_at << " Pz_at=" << Pz_at 
00766                              << " T0_at=" << T0_at
00767                              << " Vx_at=" << Vx_at << " Vy_at=" << Vy_at << " Vz_at=" << Vz_at 
00768                              << " EventWeight=" << EventWeight << " Nmuons=" << Id_sf.size() << std::endl;
00769       }
00770     }
00771     
00772     
00773   } //while loop EvtRejected
00774   
00775 
00776   return true; //write event to HepMC;
00777 
00778 }
00779 
00780 
00781 
00782 void CosmicMuonGenerator::terminate(){
00783   if (NumberOfEvents > 0){
00784     std::cout << std::endl;
00785     std::cout << "*********************************************************" << std::endl;
00786     std::cout << "*********************************************************" << std::endl;
00787     std::cout << "***                                                   ***" << std::endl;
00788     std::cout << "***    C O S M I C   M U O N   S T A T I S T I C S    ***" << std::endl;
00789     std::cout << "***                                                   ***" << std::endl;
00790     std::cout << "*********************************************************" << std::endl;
00791     std::cout << "*********************************************************" << std::endl;
00792     std::cout << std::endl;  
00793     std::cout << "       number of initial cosmic events:  " << int(Ngen) << std::endl;
00794     std::cout << "       number of actually diced events:  " << int(Ndiced) << std::endl;
00795     std::cout << "       number of generated and accepted events:  " << int(Nsel) << std::endl;
00796     double selEff = Nsel/Ngen; // selection efficiency
00797     std::cout << "       event selection efficiency:  " << selEff*100. << "%" << std::endl;
00798     int n100cos =  Norm->events_n100cos(0., 0.); //get final amount of cosmics in defined range for normalisation of flux
00799     std::cout << "       events with ~100 GeV and 1 - cos(theta) < 1/2pi: " << n100cos << std::endl;
00800     std::cout << std::endl;
00801     std::cout << "       momentum range: " << MinP             << " ... " << MaxP << " GeV" << std::endl;
00802     std::cout << "       theta  range:   " << MinTheta*Rad2Deg << " ... " << MaxTheta*Rad2Deg << " deg" << std::endl; 
00803     std::cout << "       phi    range:   " << MinPhi*Rad2Deg   << " ... " << MaxPhi*Rad2Deg << " deg" << std::endl;
00804     std::cout << "       time   range:   " << MinT0            << " ... " << MaxT0 << " ns" << std::endl;
00805     std::cout << "       energy  loss:   " << ElossScaleFactor*100. << "%" << std::endl;
00806     std::cout << std::endl;
00807     double area = 1.e-6*Pi*SurfaceRadius*SurfaceRadius; // area on surface [m^2] 
00808     if (MinTheta > 90.*Deg2Rad) //upgoing muons from neutrinos)
00809       std::cout << "       area of initial cosmics at CMS detector bottom surface:   " << area << " m^2" << std::endl;
00810     else
00811       std::cout << "       area of initial cosmics on Surface + PlugWidth:   " << area << " m^2" << std::endl;
00812     std::cout << "       depth of CMS detector (from Surface):   " << SurfaceOfEarth/1000 << " m" << std::endl;
00813        
00814     //at least 100 evts., and
00815     //downgoing inside theta parametersisation range
00816     //or upgoing neutrino muons 
00817     if( (n100cos>0 && MaxTheta<84.26*Deg2Rad) 
00818        || MinTheta>90.*Deg2Rad) {
00819       // rate: corrected for area and selection-Eff. and normalized to known flux, integration over solid angle (dOmega) is implicit
00820       // flux is normalised with respect to known flux of vertical 100GeV muons in area at suface level 
00821       // rate seen by detector is lower than rate at surface area, so has to be corrected for selection-Eff.
00822       // normalisation factor has unit [1/s/m^2] 
00823       // rate = N/time --> normalization factor gives 1/runtime/area 
00824       // normalization with respect to number of actually diced events (Ndiced)
00825 
00826       if (MinTheta > 90.*Deg2Rad) {//upgoing muons from neutrinos) 
00827         double Omega = (cos(MinTheta)-cos(MaxTheta)) * (MaxPhi-MinPhi);
00828         //EventRate = (Ndiced * 3.e-13) * Omega * area*1.e4 * selEff;//area in cm, flux=3.e-13cm^-2s^-1sr^-1
00829         EventRate = (Ndiced * 3.e-13) * Omega * 4.*RadiusOfTarget*ZDistOfTarget*1.e-2 * selEff;//area in cm, flux=3.e-13cm^-2s^-1sr^-1
00830         rateErr_stat = EventRate/sqrt( (double) Ndiced);  // stat. rate error 
00831         rateErr_syst = EventRate/3.e-13 * 1.0e-13;  // syst. rate error, from error of known flux 
00832       }
00833       else {
00834         EventRate= (Ndiced * Norm->norm(n100cos)) * area * selEff;
00835         rateErr_stat = EventRate/sqrt( (double) n100cos);  // stat. rate error 
00836         rateErr_syst = EventRate/2.63e-3 * 0.06e-3;  // syst. rate error, from error of known flux 
00837       }
00838 
00839       // normalisation in region 1.-cos(theta) < 1./(2.*Pi), if MaxTheta even lower correct for this
00840       if(MaxTheta<0.572){
00841         double spacean = 2.*Pi*(1.-cos(MaxTheta));
00842         EventRate= (Ndiced * Norm->norm(n100cos)) * area * selEff * spacean;
00843         rateErr_stat = EventRate/sqrt( (double) n100cos);  // rate error 
00844         rateErr_syst = EventRate/2.63e-3 * 0.06e-3;  // syst. rate error, from error of known flux 
00845       }
00846 
00847     }else{
00848       EventRate=Nsel; //no info as no muons at 100 GeV
00849       rateErr_stat =Nsel;
00850       rateErr_syst =Nsel;
00851       std::cout << std::endl;
00852       if (MinP > 100.)
00853         std::cout << " !!! MinP > 100 GeV. Cannot apply normalisation!" << std::endl;
00854       else if (MaxTheta > 84.26*Deg2Rad)
00855         std::cout << " !!! Note: generated cosmics exceed parameterisation. No flux calculated!" << std::endl;
00856 
00857       else
00858         std::cout << " !!! Not enough statistics to apply normalisation (rate=1 +- 1) !!!" << std::endl;
00859     } 
00860     
00861     std::cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << std::endl;
00862     std::cout << "       rate is " << EventRate << " +-" << rateErr_stat <<" (stat) " << "+-" << 
00863       rateErr_syst << " (syst) " <<" muons per second" << std::endl;
00864     if(EventRate!=0) std::cout << "       number of events corresponds to " << Nsel/EventRate << " s" << std::endl;  //runtime at CMS = Nsel/rate
00865     std::cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << std::endl;
00866     std::cout << std::endl;
00867     std::cout << "*********************************************************" << std::endl;
00868     std::cout << "*********************************************************" << std::endl;
00869   }
00870 }
00871 
00872 void CosmicMuonGenerator::checkIn(){
00873   if (MinP < 0.){ NumberOfEvents = 0;
00874     std::cout << "  CMG-ERR: min.energy is out of range (0 GeV ... inf]" << std::endl << std::endl; }
00875   if (MaxP < 0.){ NumberOfEvents = 0;
00876     std::cout << "  CMG-ERR: max.energy is out of range (0 GeV ... inf]" << std::endl << std::endl; }
00877   if (MaxP <= MinP){ NumberOfEvents = 0;
00878     std::cout << "  CMG-ERR: max.energy is not greater than min.energy" << std::endl << std::endl; }
00879   if (MinTheta < 0.){ NumberOfEvents = 0;
00880     std::cout << "  CMG-ERR: min.theta is out of range [0 deg ... 90 deg)" << std::endl << std::endl; }
00881   if (MaxTheta < 0.){ NumberOfEvents = 0;
00882     std::cout << "  CMG-ERR: max.theta is out of range [0 deg ... 90 deg)" << std::endl << std::endl; }
00883   if (MaxTheta <= MinTheta){ NumberOfEvents = 0;
00884     std::cout << "  CMG-ERR: max.theta is not greater than min.theta" << std::endl << std::endl; }
00885   if (MinPhi < 0.){ NumberOfEvents = 0;
00886     std::cout << "  CMG-ERR: min.phi is out of range [0 deg ... 360 deg]" << std::endl << std::endl; }
00887   if (MaxPhi < 0.){ NumberOfEvents = 0;
00888     std::cout << "  CMG-ERR: max.phi is out of range [0 deg ... 360 deg]" << std::endl << std::endl; }
00889   if (MaxPhi <= MinPhi){ NumberOfEvents = 0;
00890     std::cout << "  CMG-ERR: max.phi is not greater than min.phi" << std::endl << std::endl; }
00891   if (MaxT0 <= MinT0){ NumberOfEvents = 0;
00892     std::cout << "  CMG-ERR: max.t0 is not greater than min.t0" << std::endl << std::endl; }
00893   if (ElossScaleFactor < 0.){ NumberOfEvents = 0;
00894     std::cout << "  CMG-ERR: E-loss scale factor is out of range [0 ... inf)" << std::endl << std::endl; }
00895   if (MinEnu < 0.){ NumberOfEvents = 0;
00896     std::cout << "  CMG-ERR: min.Enu is out of range [0 GeV ... inf]" << std::endl << std::endl; }
00897   if (MaxEnu < 0.){ NumberOfEvents = 0;
00898     std::cout << "  CMG-ERR: max.Enu is out of range [0 GeV ... inf]" << std::endl << std::endl; }
00899   if (MaxEnu <= MinEnu){ NumberOfEvents = 0;
00900     std::cout << "  CMG-ERR: max.Enu is not greater than min.Enu" << std::endl << std::endl; }
00901 
00902 }
00903 
00904 bool CosmicMuonGenerator::goodOrientation(){
00905   // check angular range (for a sphere with Target3dRadius around the target)
00906   bool goodAngles = false;
00907   bool phiaccepted = false;
00908   bool thetaaccepted = false;
00909   double RxzV = sqrt(OneMuoEvt.vx()*OneMuoEvt.vx() + OneMuoEvt.vz()*OneMuoEvt.vz());
00910 
00911   double rVY;
00912   if (MinTheta > 90.*Deg2Rad) //upgoing muons from neutrinos
00913     rVY = -sqrt(RxzV*RxzV + RadiusCMS*RadiusCMS);
00914   else 
00915     rVY = sqrt(RxzV*RxzV + (SurfaceOfEarth+PlugWidth)*(SurfaceOfEarth+PlugWidth));
00916 
00917   double Phi = OneMuoEvt.phi();
00918   double PhiV = atan2(OneMuoEvt.vx(),OneMuoEvt.vz()) + Pi; if (PhiV > TwoPi) PhiV -= TwoPi;
00919   double disPhi = std::fabs(PhiV - Phi); if (disPhi > Pi) disPhi = TwoPi - disPhi;
00920   double dPhi = Pi; if (RxzV > Target3dRadius) dPhi = asin(Target3dRadius/RxzV);
00921   if (disPhi < dPhi) phiaccepted = true;
00922   double Theta = OneMuoEvt.theta();
00923   double ThetaV = asin(RxzV/rVY);
00924   double dTheta = Pi; if (std::fabs(rVY) > Target3dRadius) dTheta = asin(Target3dRadius/std::fabs(rVY));
00925   //std::cout << "    dPhi = " <<   dPhi << "  (" <<   Phi << " <p|V> " <<   PhiV << ")" << std::endl;
00926   //std::cout << "  dTheta = " << dTheta << "  (" << Theta << " <p|V> " << ThetaV << ")" << std::endl;
00927 
00928   if (!phiaccepted && RxzV < Target3dRadius)
00929   //if (RxzV < Target3dRadius)
00930     std::cout << "Rejected phi=" << Phi << "  PhiV="  << PhiV 
00931          << "  dPhi=" << dPhi << "  disPhi=" << disPhi
00932          << "  RxzV=" << RxzV << "  Target3dRadius=" << Target3dRadius 
00933          << "  Theta=" << Theta << std::endl;
00934 
00935   if (std::fabs(Theta-ThetaV) < dTheta) thetaaccepted = true;
00936   if (phiaccepted && thetaaccepted) goodAngles = true;
00937   return goodAngles;
00938 }
00939 
00940 void CosmicMuonGenerator::initEvDis(){
00941 #if ROOT_INTERACTIVE
00942   float rCMS = RadiusCMS/1000.;
00943   float zCMS = Z_DistCMS/1000.;
00944   if(TrackerOnly==true){
00945     rCMS = RadiusTracker/1000.;
00946     zCMS = Z_DistTracker/1000.;
00947 }
00948   TH2F* disXY = new TH2F("disXY","X-Y view",160,-rCMS,rCMS,160,-rCMS,rCMS);
00949   TH2F* disZY = new TH2F("disZY","Z-Y view",150,-zCMS,zCMS,160,-rCMS,rCMS);
00950   gStyle->SetPalette(1,0);
00951   gStyle->SetMarkerColor(1);
00952   gStyle->SetMarkerSize(1.5);
00953   TCanvas *disC = new TCanvas("disC","Cosmic Muon Event Display",0,0,800,410);
00954   disC->Divide(2,1);
00955   disC->cd(1);
00956   gPad->SetTicks(1,1);
00957   disXY->SetMinimum(log10(MinP));
00958   disXY->SetMaximum(log10(MaxP));
00959   disXY->GetXaxis()->SetLabelSize(0.05);
00960   disXY->GetXaxis()->SetTitleSize(0.05);
00961   disXY->GetXaxis()->SetTitleOffset(1.0);
00962   disXY->GetXaxis()->SetTitle("X [m]");
00963   disXY->GetYaxis()->SetLabelSize(0.05);
00964   disXY->GetYaxis()->SetTitleSize(0.05);
00965   disXY->GetYaxis()->SetTitleOffset(0.8);
00966   disXY->GetYaxis()->SetTitle("Y [m]");
00967   disC->cd(2);
00968   gPad->SetGrid(1,1);
00969   gPad->SetTicks(1,1);
00970   disZY->SetMinimum(log10(MinP));
00971   disZY->SetMaximum(log10(MaxP));
00972   disZY->GetXaxis()->SetLabelSize(0.05);
00973   disZY->GetXaxis()->SetTitleSize(0.05);
00974   disZY->GetXaxis()->SetTitleOffset(1.0);
00975   disZY->GetXaxis()->SetTitle("Z [m]");
00976   disZY->GetYaxis()->SetLabelSize(0.05);
00977   disZY->GetYaxis()->SetTitleSize(0.05);
00978   disZY->GetYaxis()->SetTitleOffset(0.8);
00979   disZY->GetYaxis()->SetTitle("Y [m]");
00980 #endif
00981 }
00982 
00983 void CosmicMuonGenerator::displayEv(){
00984 #if ROOT_INTERACTIVE
00985   double RadiusDet=RadiusCMS;
00986   double Z_DistDet=Z_DistCMS;
00987   if(TrackerOnly==true){
00988     RadiusDet = RadiusTracker;
00989     Z_DistDet = Z_DistTracker;
00990   }
00991   disXY->Reset();
00992   disZY->Reset();
00993   TMarker* InteractionPoint = new TMarker(0.,0.,2);
00994   TArc* r8m = new TArc(0.,0.,(RadiusDet/1000.));
00995   TLatex* logEaxis = new TLatex(); logEaxis->SetTextSize(0.05);
00996   float energy = float(OneMuoEvt.e());
00997   float verX = float(OneMuoEvt.vx()/1000.); // [m]
00998   float verY = float(OneMuoEvt.vy()/1000.); // [m]
00999   float verZ = float(OneMuoEvt.vz()/1000.); // [m]
01000   float dirX = float(OneMuoEvt.px())/std::fabs(OneMuoEvt.py());
01001   float dirY = float(OneMuoEvt.py())/std::fabs(OneMuoEvt.py());
01002   float dirZ = float(OneMuoEvt.pz())/std::fabs(OneMuoEvt.py());
01003   float yStep = disXY->GetYaxis()->GetBinWidth(1);
01004   int   NbinY = disXY->GetYaxis()->GetNbins();
01005   for (int iy=0; iy<NbinY; ++iy){
01006     verX += dirX*yStep;
01007     verY += dirY*yStep;
01008     verZ += dirZ*yStep;
01009     float rXY = sqrt(verX*verX + verY*verY)*1000.; // [mm]
01010     float absZ = std::fabs(verZ)*1000.;                 // [mm]
01011     if (rXY < RadiusDet && absZ < Z_DistDet){
01012       disXY->Fill(verX,verY,log10(energy));
01013       disZY->Fill(verZ,verY,log10(energy));
01014       disC->cd(1); disXY->Draw("COLZ"); InteractionPoint->Draw("SAME"); r8m->Draw("SAME");
01015       logEaxis->DrawLatex((0.65*RadiusDet/1000.),(1.08*RadiusDet/1000.),"log_{10}E(#mu^{#pm})");
01016       disC->cd(2); disZY->Draw("COL"); InteractionPoint->Draw("SAME");
01017       gPad->Update();
01018     }
01019   }
01020 #endif
01021 }
01022 
01023 void CosmicMuonGenerator::setNumberOfEvents(unsigned int N){ if (NotInitialized) NumberOfEvents = N; }
01024 
01025 void CosmicMuonGenerator::setRanSeed(int N){ if (NotInitialized) RanSeed = N; }
01026 
01027 void CosmicMuonGenerator::setMinP(double P){ if (NotInitialized) MinP = P; }
01028 
01029 void CosmicMuonGenerator::setMinP_CMS(double P){ if (NotInitialized) MinP_CMS = P; }
01030 
01031 void CosmicMuonGenerator::setMaxP(double P){ if (NotInitialized) MaxP = P; }
01032 
01033 void CosmicMuonGenerator::setMinTheta(double Theta){ if (NotInitialized) MinTheta = Theta*Deg2Rad; }
01034 
01035 void CosmicMuonGenerator::setMaxTheta(double Theta){ if (NotInitialized) MaxTheta = Theta*Deg2Rad; } 
01036 
01037 void CosmicMuonGenerator::setMinPhi(double Phi){ if (NotInitialized) MinPhi = Phi*Deg2Rad; }
01038 
01039 void CosmicMuonGenerator::setMaxPhi(double Phi){ if (NotInitialized) MaxPhi = Phi*Deg2Rad; }
01040 
01041 void CosmicMuonGenerator::setMinT0(double T0){ if (NotInitialized) MinT0 = T0; }
01042 
01043 void CosmicMuonGenerator::setMaxT0(double T0){ if (NotInitialized) MaxT0 = T0; }
01044 
01045 void CosmicMuonGenerator::setElossScaleFactor(double ElossScaleFact){ if (NotInitialized) ElossScaleFactor = ElossScaleFact; }
01046 
01047 void CosmicMuonGenerator::setRadiusOfTarget(double R){ if (NotInitialized) RadiusOfTarget = R; }
01048 
01049 void CosmicMuonGenerator::setZDistOfTarget(double Z){ if (NotInitialized) ZDistOfTarget = Z; }
01050 
01051 void CosmicMuonGenerator::setZCentrOfTarget(double Z){ if (NotInitialized) ZCentrOfTarget = Z; }
01052 
01053 void CosmicMuonGenerator::setTrackerOnly(bool Tracker){ if (NotInitialized) TrackerOnly = Tracker; }
01054 
01055 void CosmicMuonGenerator::setMultiMuon(bool MultiMu){ if (NotInitialized) MultiMuon = MultiMu; }
01056 void CosmicMuonGenerator::setMultiMuonFileName(std::string MultiMuFile){ if (NotInitialized) MultiMuonFileName = MultiMuFile; }
01057 void CosmicMuonGenerator::setMultiMuonFileFirstEvent(int MultiMuFile1stEvt){ if (NotInitialized) MultiMuonFileFirstEvent = MultiMuFile1stEvt; }
01058 void CosmicMuonGenerator::setMultiMuonNmin(int MultiMuNmin){ if (NotInitialized) MultiMuonNmin = MultiMuNmin; }
01059 
01060 void CosmicMuonGenerator::setTIFOnly_constant(bool TIF){ if (NotInitialized) TIFOnly_constant = TIF; }
01061 
01062 void CosmicMuonGenerator::setTIFOnly_linear(bool TIF){ if (NotInitialized) TIFOnly_linear = TIF; }
01063 void CosmicMuonGenerator::setMTCCHalf(bool MTCC){ if (NotInitialized) MTCCHalf = MTCC; }
01064 
01065 void CosmicMuonGenerator::setPlugVx(double PlugVtx){ if (NotInitialized) PlugVx = PlugVtx; }
01066 void CosmicMuonGenerator::setPlugVz(double PlugVtz){ if (NotInitialized) PlugVz = PlugVtz; }
01067 void CosmicMuonGenerator::setRhoAir(double VarRhoAir){ if (NotInitialized) RhoAir = VarRhoAir; }
01068 void CosmicMuonGenerator::setRhoWall(double VarRhoWall){ if (NotInitialized) RhoWall = VarRhoWall; }
01069 void CosmicMuonGenerator::setRhoRock(double VarRhoRock){ if (NotInitialized) RhoRock = VarRhoRock; }
01070 void CosmicMuonGenerator::setRhoClay(double VarRhoClay){ if (NotInitialized) RhoClay = VarRhoClay; }
01071 void CosmicMuonGenerator::setRhoPlug(double VarRhoPlug){ if (NotInitialized) RhoPlug = VarRhoPlug; }
01072 void CosmicMuonGenerator::setClayWidth(double ClayLayerWidth){ if (NotInitialized) ClayWidth = ClayLayerWidth; }
01073 
01074 void CosmicMuonGenerator::setMinEnu(double MinEn){ if (NotInitialized) MinEnu = MinEn; }
01075 void CosmicMuonGenerator::setMaxEnu(double MaxEn){ if (NotInitialized) MaxEnu = MaxEn; }
01076 void CosmicMuonGenerator::setNuProdAlt(double NuPrdAlt){ if (NotInitialized) NuProdAlt = NuPrdAlt; }
01077 
01078 double CosmicMuonGenerator::getRate(){ return EventRate; }
01079 
01080 void CosmicMuonGenerator::setAcptAllMu(bool AllMu){ if (NotInitialized) AcptAllMu = AllMu; }