00001
00002
00004
00005
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);
00024 delRanGen = true;
00025 } else {
00026 RanGen = rng;
00027 delRanGen = false;
00028 }
00029 checkIn();
00030 if (NumberOfEvents > 0){
00031
00032 double RadiusTargetEff = RadiusOfTarget;
00033 double Z_DistTargetEff = ZDistOfTarget;
00034
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)
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
00056 if (MinTheta >= 90.*Deg2Rad)
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
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
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;
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
00107 bool notSelected = true;
00108 while (notSelected){
00109 bool badMomentumGenerated = true;
00110 while (badMomentumGenerated){
00111
00112 if (MinTheta > 90.*Deg2Rad)
00113 Cosmics->generateNuMu();
00114 else
00115 Cosmics->generate();
00116
00117 E = sqrt(Cosmics->momentum_times_charge()*Cosmics->momentum_times_charge() + MuonMass*MuonMass);
00118 Theta = TMath::ACos( Cosmics->cos_theta() ) ;
00119 Ngen+=1.;
00120 badMomentumGenerated = false;
00121 Phi = RanGen->flat()*(MaxPhi-MinPhi) + MinPhi;
00122 }
00123 Norm->events_n100cos(E, Theta);
00124 Ndiced += 1;
00125
00126
00127 double Nver = 0.;
00128 bool badVertexGenerated = true;
00129 while (badVertexGenerated){
00130 RxzV = sqrt(RanGen->flat())*SurfaceRadius;
00131 PhiV = RanGen->flat()*TwoPi;
00132
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.);
00140
00141
00142 int id = 13;
00143 if (Cosmics->momentum_times_charge() >0.) id = -13;
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);
00148 double Py = -verMom;
00149 double Pz = horMom*cos(Phi);
00150 double Vx = RxzV*sin(PhiV);
00151
00152 double Vy;
00153 if (MinTheta > 90.*Deg2Rad)
00154 Vy = -RadiusCMS;
00155 else
00156 Vy = SurfaceOfEarth + PlugWidth;
00157
00158 double Vz = RxzV*cos(PhiV);
00159 double T0 = (RanGen->flat()*(MaxT0-MinT0) + MinT0)*SpeedOfLight;
00160
00161 Id_at = id;
00162 Px_at = Px; Py_at = Py; Pz_at = Pz; E_at = E;
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
00167 if (goodOrientation()) {
00168 if (MinTheta > 90.*Deg2Rad)
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.;
00177 notSelected = false;
00178 }
00179 }
00180
00181 EventWeight = 1.;
00182
00183
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
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;
00197 Vx_sf[0] = Vx_at; Vy_sf[0] = Vy_at; Vz_sf[0] = Vz_at; T0_sf[0] = T0_at;
00198
00199
00200
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
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
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
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;
00249
00250 while (EvtRejected) {
00251
00252
00253
00254 Long64_t ientry = SimTree->GetEntry(SimTree_jentry);
00255 std::cout << "CosmicMuonGenerator::nextMultiEvent(): SimTree_jentry=" << SimTree_jentry
00256
00257 << " SimTreeEntries=" << SimTreeEntries << std::endl;
00258 if (ientry < 0) return false;
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;
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;
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;
00284 double MuMuDist;
00285 MuInMaxDist = false;
00286
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];
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.;
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;
00328 }
00329 }
00330
00331 if (NmuPmin < MultiMuonNmin && AcptAllMu==false) {
00332 NskippedMultiMuonEvents++;
00333 continue;
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
00344 Id_at = SimTree->shower_EventID;
00345
00346 double M_at = 0.;
00347
00348 Id_at = 2212;
00349 M_at = 938.272e-3;
00350
00351
00352 E_at = SimTree->shower_Energy;
00353 Theta_at = SimTree->shower_Theta;
00354 double phi_at = SimTree->shower_Phi - NorthCMSzDeltaPhi;
00355 if (phi_at < -Pi) phi_at +=TwoPi;
00356 else if (phi_at > Pi) phi_at -= TwoPi;
00357 double P_at = sqrt(E_at*E_at - M_at*M_at);
00358
00359
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
00365 double theta_mu_max = Theta_at;
00366 double theta_mu_min = Theta_at;
00367
00368 double phi_rel_min = 0.;
00369 double phi_rel_max = 0.;
00370
00371
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]);
00380 double phi_rel = phi_mu - phi_at;
00381 if (phi_rel < -Pi) phi_rel += TwoPi;
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;
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
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
00410
00411
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));
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
00426 double psRall = psR3;
00427
00428 double fR1=psR1/psRall, fR2=psR2/psRall, fR3=psR3/psRall;
00429 double pR1=0.25, pR2=0.7, pR3=0.05;
00430
00431
00432
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;
00439
00440 double pPh1=0.25, pPh2=0.7, pPh3=0.05;
00441
00442 Trials = 0;
00443 double trials = 0.;
00444 Vx_mu.resize(nmuons); Vy_mu.resize(nmuons); Vz_mu.resize(nmuons);
00445 int NmuHitTarget = 0;
00446 while (NmuHitTarget < MultiMuonNmin) {
00447
00448 NmuHitTarget = 0;
00449 double Nver = 0.;
00450
00451
00452
00453 double RxzV;
00454 double which_R_class = RanGen->flat();
00455 if (which_R_class < pR1) {
00456 RxzV = psR1min + psR1 * RanGen->flat();
00457 JdRxzV_dR_trans = fR1/pR1 * SurfaceRadius/psR1;
00458 }
00459 else if (which_R_class < pR1+pR2) {
00460 RxzV = psR2min + psR2 * RanGen->flat();
00461 JdRxzV_dR_trans = fR2/pR2 * SurfaceRadius/psR2;
00462 }
00463 else {
00464 RxzV = psR3min + psR3 * RanGen->flat();
00465 JdRxzV_dR_trans = fR3/pR3 * SurfaceRadius/psR3;
00466 }
00467
00468 JdR_trans_sqrt = 2.*RxzV/SurfaceRadius;
00469
00470
00471 double PhiV;
00472 double which_phi_class = RanGen->flat();
00473 if (which_phi_class < pPh1) {
00474 PhiV = phi_at + phi_rel_min + psPh1 * RanGen->flat();
00475 JdPhiV_dPhi_trans = fPh1/pPh1 * TwoPi/psPh1;
00476 }
00477 else if (which_phi_class < pPh1+pPh2) {
00478 PhiV = phi_at + phi_rel_min + psPh2 * RanGen->flat();
00479 JdPhiV_dPhi_trans = fPh2/pPh2 * TwoPi/psPh2;
00480 }
00481 else {
00482 PhiV = phi_at + phi_rel_min + psPh3 * RanGen->flat();
00483 JdPhiV_dPhi_trans = fPh3/pPh3 * TwoPi/psPh3;
00484 }
00485
00486
00487 if (PhiV < -Pi) PhiV+=TwoPi;
00488 else if (PhiV > Pi) PhiV-=TwoPi;
00489
00490
00491 Nver++;
00492 trials += JdR_trans_sqrt * JdRxzV_dR_trans * JdPhiV_dPhi_trans;
00493 Trials++;
00494 if (trials > max_Trials) break;
00495 Ngen += (Nver-1.);
00496
00497
00498 Vx_at = RxzV*sin(PhiV);
00499
00500 Vy_at = h_sf;
00501
00502
00503 Vz_at = RxzV*cos(PhiV);
00504
00505 int NmuHitTargetSphere = 0;
00506 for (int imu=0; imu<nmuons; ++imu) {
00507
00508 Vx_mu[imu] = Vx_at + (-SimTree->particle__x[imu]*sin(NorthCMSzDeltaPhi)
00509 +SimTree->particle__y[imu]*cos(NorthCMSzDeltaPhi) )*10;
00510 Vy_mu[imu] = h_sf;
00511 Vz_mu[imu] = Vz_at + ( SimTree->particle__x[imu]*cos(NorthCMSzDeltaPhi)
00512 +SimTree->particle__y[imu]*sin(NorthCMSzDeltaPhi) )*10;
00513
00514
00515
00516 double pt_sec = sqrt(Px_mu[imu]*Px_mu[imu]+Pz_mu[imu]*Pz_mu[imu]);
00517 double theta_sec = atan2(std::fabs(Py_mu[imu]),pt_sec);
00518 double r_sec = sqrt((Vx_mu[imu]-Vx_at)*(Vx_mu[imu]-Vx_at)
00519 +(Vz_mu[imu]-Vz_at)*(Vz_mu[imu]-Vz_at));
00520 double h_prod = r_sec * tan(theta_sec);
00521 if (h_prod + h_sf > Vy_at) Vy_at = h_prod + h_sf;
00522
00523
00524 if (SimTree->particle__ParticleID[imu] != 5 &&
00525 SimTree->particle__ParticleID[imu] != 6) continue;
00526
00527 if (P_mu[imu] < MinP_CMS && AcptAllMu==false) continue;
00528
00529
00530 double Vxz_mu = sqrt(Vx_mu[imu]*Vx_mu[imu] + Vz_mu[imu]*Vz_mu[imu]);
00531 theta_mu_max = atan((Vxz_mu+Target3dRadius)/(h_sf-Target3dRadius));
00532 theta_mu_min = atan((Vxz_mu-Target3dRadius)/(h_sf+Target3dRadius));
00533 if (Theta_mu[imu] > theta_mu_min && Theta_mu[imu] < theta_mu_max) {
00534
00535
00536 double dPhi = Pi; if (Vxz_mu > Target3dRadius) dPhi = asin(Target3dRadius/Vxz_mu);
00537 double PhiPmu = atan2(Px_mu[imu],Pz_mu[imu]);
00538 double PhiVmu = atan2(Vx_mu[imu],Vz_mu[imu]);
00539 double rotPhi = PhiVmu + Pi; if (rotPhi > Pi) rotPhi -= TwoPi;
00540 double disPhi = std::fabs(rotPhi - PhiPmu); if (disPhi > Pi) disPhi = TwoPi - disPhi;
00541 if (disPhi < dPhi) {
00542 NmuHitTargetSphere++;
00543 }
00544
00545 }
00546
00547 }
00548
00549
00550
00551
00552 if (NmuHitTargetSphere < MultiMuonNmin) continue;
00553
00554
00555 Id_sf.clear();
00556 Px_sf.clear();
00557 Py_sf.clear();
00558 Pz_sf.clear();
00559 E_sf.clear();
00560
00561 Vx_sf.clear();
00562 Vy_sf.clear();
00563 Vz_sf.clear();
00564 T0_sf.clear();
00565
00566
00567 Id_ug.clear();
00568 Px_ug.clear();
00569 Py_ug.clear();
00570 Pz_ug.clear();
00571 E_ug.clear();
00572
00573 Vx_ug.clear();
00574 Vy_ug.clear();
00575 Vz_ug.clear();
00576 T0_ug.clear();
00577
00578 int Id_sf_this =0;
00579 double Px_sf_this =0., Py_sf_this=0., Pz_sf_this=0.;
00580 double E_sf_this=0.;
00581
00582 double Vx_sf_this=0., Vy_sf_this=0., Vz_sf_this=0.;
00583 double T0_sf_this=0.;
00584
00585 T0_at = SimTree->shower_GH_t0 * SpeedOfLight;
00586
00587 for (int imu=0; imu<nmuons; ++imu) {
00588
00589 if (P_mu[imu] < MinP_CMS && AcptAllMu==false) continue;
00590
00591 if (SimTree->particle__ParticleID[imu] != 5 &&
00592 SimTree->particle__ParticleID[imu] != 6) continue;
00593
00594 Id_sf_this = SimTree->particle__ParticleID[imu];
00595 if (Id_sf_this == 5) Id_sf_this = -13;
00596 else if (Id_sf_this == 6) Id_sf_this = 13;
00597
00598 else if (Id_sf_this == 1) Id_sf_this = 22;
00599 else if (Id_sf_this == 2) Id_sf_this = -11;
00600 else if (Id_sf_this == 3) Id_sf_this = 11;
00601 else if (Id_sf_this == 7) Id_sf_this = 111;
00602 else if (Id_sf_this == 8) Id_sf_this = 211;
00603 else if (Id_sf_this == 9) Id_sf_this = -211;
00604 else if (Id_sf_this == 10) Id_sf_this = 130;
00605 else if (Id_sf_this == 11) Id_sf_this = 321;
00606 else if (Id_sf_this == 12) Id_sf_this = -321;
00607 else if (Id_sf_this == 13) Id_sf_this = 2112;
00608 else if (Id_sf_this == 14) Id_sf_this = 2212;
00609 else if (Id_sf_this == 15) Id_sf_this = -2212;
00610 else if (Id_sf_this == 16) Id_sf_this = 310;
00611 else if (Id_sf_this == 17) Id_sf_this = 221;
00612 else if (Id_sf_this == 18) Id_sf_this = 3122;
00613
00614 else {
00615 std::cout << "CosmicMuonGenerator.cc: Warning! Muon Id=" << Id_sf_this
00616 << " from file read in" <<std::endl;
00617 Id_sf_this = 99999;
00618 }
00619
00620 T0_sf_this = SimTree->particle__Time[imu] * SpeedOfLight;
00621
00622 Px_sf_this = Px_mu[imu];
00623 Py_sf_this = Py_mu[imu];
00624 Pz_sf_this = Pz_mu[imu];
00625 E_sf_this = sqrt(P_mu[imu]*P_mu[imu] + MuonMass*MuonMass);
00626 Vx_sf_this = Vx_mu[imu];
00627 Vy_sf_this = Vy_mu[imu];
00628 Vz_sf_this = Vz_mu[imu];
00629
00630
00631 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);
00632
00633 if (goodOrientation()) {
00634 OneMuoEvt.propagate(ElossScaleFactor, RadiusOfTarget, ZDistOfTarget, ZCentrOfTarget, TrackerOnly, MTCCHalf);
00635 }
00636
00637 if ( (OneMuoEvt.hitTarget()
00638 && sqrt(OneMuoEvt.e()*OneMuoEvt.e() - MuonMass*MuonMass) > MinP_CMS)
00639 || AcptAllMu==true ) {
00640
00641 Id_sf.push_back(Id_sf_this);
00642 Px_sf.push_back(Px_sf_this);
00643 Py_sf.push_back(Py_sf_this);
00644 Pz_sf.push_back(Pz_sf_this);
00645 E_sf.push_back(E_sf_this);
00646
00647 Vx_sf.push_back(Vx_sf_this);
00648 Vy_sf.push_back(Vy_sf_this);
00649 Vz_sf.push_back(Vz_sf_this);
00650 T0_sf.push_back(T0_sf_this);
00651
00652
00653 Id_ug.push_back(OneMuoEvt.id());
00654 Px_ug.push_back(OneMuoEvt.px());
00655 Py_ug.push_back(OneMuoEvt.py());
00656 Pz_ug.push_back(OneMuoEvt.pz());
00657 E_ug.push_back(OneMuoEvt.e());
00658
00659 Vx_ug.push_back(OneMuoEvt.vx());
00660 Vy_ug.push_back(OneMuoEvt.vy());
00661 Vz_ug.push_back(OneMuoEvt.vz());
00662 T0_ug.push_back(OneMuoEvt.t0());
00663
00664 NmuHitTarget++;
00665 }
00666 }
00667
00668
00669 }
00670
00671
00672 if (trials > max_Trials) {
00673 std::cout << "CosmicMuonGenerator.cc: Warning! trials reach max_trials=" << max_Trials
00674 << " without accepting event!" << std::endl;
00675 if (Debug) {
00676 std::cout << " N(mu)=" << Id_ug.size();
00677 if (Id_ug.size()>=1)
00678 std::cout << " E[0]=" << E_ug[0] << " theta="
00679 << acos(-Py_ug[0]/sqrt(Px_ug[0]*Px_ug[0]+Py_ug[0]*Py_ug[0]+Pz_ug[0]*Pz_ug[0]))
00680 << " R_xz=" << sqrt(Vx_sf[0]*Vx_sf[0]+Vy_sf[0]*Vy_sf[0]);
00681 std::cout << " Theta_at=" << Theta_at << std::endl;
00682 }
00683 std::cout << "Unweighted int num of Trials = " << Trials << std::endl;
00684 std::cout << "trying next event (" << SimTree_jentry << ") from file" << std::endl;
00685 NskippedMultiMuonEvents++;
00686 continue;
00687 }
00688 else {
00689 if (NmuHitTarget < MultiMuonNmin) {
00690 std::cout << "CosmicMuonGenerator.cc: Warning! less than " << MultiMuonNmin <<
00691 " muons hit target: N(mu=)" << NmuHitTarget << std::endl;
00692 std::cout << "trying next event (" << SimTree_jentry << ") from file" << std::endl;
00693 NskippedMultiMuonEvents++;
00694 continue;
00695 }
00696 else {
00697
00698
00699
00700 double T0_ug_min, T0_ug_max;
00701 T0_ug_min = T0_ug_max = T0_ug[0];
00702 double Tbox = (MaxT0 - MinT0) * SpeedOfLight;
00703 double minDeltaT0 = 2*Tbox;
00704 for (unsigned int imu=0; imu<Id_ug.size(); ++imu) {
00705 double T0_this = T0_ug[imu];
00706 if (T0_this < T0_ug_min) T0_ug_min = T0_this;
00707 if (T0_this > T0_ug_max) T0_ug_max = T0_this;
00708 if (Debug) std::cout << "imu=" << imu << " T0_this=" << T0_this
00709 << " P=" << sqrt(pow(Px_ug[imu],2) + pow(Py_ug[imu],2) + pow(Pz_ug[imu],2))
00710 << std::endl;
00711 for (unsigned int jmu=0; jmu<imu; ++jmu) {
00712 if (std::fabs(T0_ug[imu]-T0_ug[jmu]) < minDeltaT0) minDeltaT0 = std::fabs(T0_ug[imu]-T0_ug[jmu]);
00713 }
00714 }
00715
00716 if (int(Id_ug.size()) >= MultiMuonNmin && MultiMuonNmin>=2 && minDeltaT0 > Tbox)
00717 continue;
00718
00719 double T0_min = T0_ug_min +MinT0*SpeedOfLight;
00720 double T0_max = T0_ug_max +MaxT0*SpeedOfLight;
00721
00722
00723 int TboxTrials = 0;
00724 int NmuInTbox;
00725 double T0_offset, T0diff;
00726 for (int tboxtrial=0; tboxtrial<1000; ++tboxtrial) {
00727 T0_offset = RanGen->flat()*(T0_max -T0_min);
00728 TboxTrials++;
00729 T0diff = T0_offset - T0_max;
00730 NmuInTbox = 0;
00731 for (unsigned int imu=0; imu<Id_ug.size(); ++imu) {
00732 if (T0_ug[imu]+T0diff > MinT0*SpeedOfLight && T0_ug[imu]+T0diff < MaxT0*SpeedOfLight)
00733 NmuInTbox++;
00734 }
00735 if (NmuInTbox >= MultiMuonNmin) break;
00736
00737 }
00738 if (NmuInTbox < MultiMuonNmin) continue;
00739
00740
00741 if (Debug) std::cout << "initial T0_at=" << T0_at << " T0_min=" << T0_min << " T0_max=" << T0_max
00742 << " T0_offset=" << T0_offset;
00743 T0_at += T0diff;
00744 if (Debug) std::cout << " T0diff=" << T0diff << std::endl;
00745 for (unsigned int imu=0; imu<Id_ug.size(); ++imu) {
00746 if (Debug) std::cout << "before: T0_sf[" << imu << "]=" << T0_sf[imu] << " T0_ug=" << T0_ug[imu];
00747 T0_sf[imu] += T0diff;
00748 T0_ug[imu] += T0diff;
00749 if (Debug)
00750 std::cout << " after: T0_sf[" << imu << "]=" << T0_sf[imu] << " T0_ug=" << T0_ug[imu] << std::endl;
00751 }
00752 if (Debug) std::cout << "T0diff=" << T0diff << " T0_at=" << T0_at << std::endl;
00753
00754
00755
00756 Nsel += 1;
00757 EventWeight = JdR_trans_sqrt * JdRxzV_dR_trans * JdPhiV_dPhi_trans
00758 / (trials * TboxTrials);
00759 EvtRejected = false;
00760 if (Debug) std::cout << "CosmicMuonGenerator.cc: Theta_at=" << Theta_at << " phi_at=" << phi_at
00761 << " Px_at=" << Px_at << " Py_at=" << Py_at << " Pz_at=" << Pz_at
00762 << " T0_at=" << T0_at
00763 << " Vx_at=" << Vx_at << " Vy_at=" << Vy_at << " Vz_at=" << Vz_at
00764 << " EventWeight=" << EventWeight << " Nmuons=" << Id_sf.size() << std::endl;
00765 }
00766 }
00767
00768
00769 }
00770
00771
00772 return true;
00773
00774 }
00775
00776
00777
00778 void CosmicMuonGenerator::terminate(){
00779 if (NumberOfEvents > 0){
00780 std::cout << std::endl;
00781 std::cout << "*********************************************************" << std::endl;
00782 std::cout << "*********************************************************" << std::endl;
00783 std::cout << "*** ***" << std::endl;
00784 std::cout << "*** C O S M I C M U O N S T A T I S T I C S ***" << std::endl;
00785 std::cout << "*** ***" << std::endl;
00786 std::cout << "*********************************************************" << std::endl;
00787 std::cout << "*********************************************************" << std::endl;
00788 std::cout << std::endl;
00789 std::cout << " number of initial cosmic events: " << int(Ngen) << std::endl;
00790 std::cout << " number of actually diced events: " << int(Ndiced) << std::endl;
00791 std::cout << " number of generated and accepted events: " << int(Nsel) << std::endl;
00792 double selEff = Nsel/Ngen;
00793 std::cout << " event selection efficiency: " << selEff*100. << "%" << std::endl;
00794 int n100cos = Norm->events_n100cos(0., 0.);
00795 std::cout << " events with ~100 GeV and 1 - cos(theta) < 1/2pi: " << n100cos << std::endl;
00796 std::cout << std::endl;
00797 std::cout << " momentum range: " << MinP << " ... " << MaxP << " GeV" << std::endl;
00798 std::cout << " theta range: " << MinTheta*Rad2Deg << " ... " << MaxTheta*Rad2Deg << " deg" << std::endl;
00799 std::cout << " phi range: " << MinPhi*Rad2Deg << " ... " << MaxPhi*Rad2Deg << " deg" << std::endl;
00800 std::cout << " time range: " << MinT0 << " ... " << MaxT0 << " ns" << std::endl;
00801 std::cout << " energy loss: " << ElossScaleFactor*100. << "%" << std::endl;
00802 std::cout << std::endl;
00803 double area = 1.e-6*Pi*SurfaceRadius*SurfaceRadius;
00804 if (MinTheta > 90.*Deg2Rad)
00805 std::cout << " area of initial cosmics at CMS detector bottom surface: " << area << " m^2" << std::endl;
00806 else
00807 std::cout << " area of initial cosmics on Surface + PlugWidth: " << area << " m^2" << std::endl;
00808 std::cout << " depth of CMS detector (from Surface): " << SurfaceOfEarth/1000 << " m" << std::endl;
00809
00810
00811
00812
00813 if( (n100cos>0 && MaxTheta<84.26*Deg2Rad)
00814 || MinTheta>90.*Deg2Rad) {
00815
00816
00817
00818
00819
00820
00821
00822 if (MinTheta > 90.*Deg2Rad) {
00823 double Omega = (cos(MinTheta)-cos(MaxTheta)) * (MaxPhi-MinPhi);
00824
00825 EventRate = (Ndiced * 3.e-13) * Omega * 4.*RadiusOfTarget*ZDistOfTarget*1.e-2 * selEff;
00826 rateErr_stat = EventRate/sqrt( (double) Ndiced);
00827 rateErr_syst = EventRate/3.e-13 * 1.0e-13;
00828 }
00829 else {
00830 EventRate= (Ndiced * Norm->norm(n100cos)) * area * selEff;
00831 rateErr_stat = EventRate/sqrt( (double) n100cos);
00832 rateErr_syst = EventRate/2.63e-3 * 0.06e-3;
00833 }
00834
00835
00836 if(MaxTheta<0.572){
00837 double spacean = 2.*Pi*(1.-cos(MaxTheta));
00838 EventRate= (Ndiced * Norm->norm(n100cos)) * area * selEff * spacean;
00839 rateErr_stat = EventRate/sqrt( (double) n100cos);
00840 rateErr_syst = EventRate/2.63e-3 * 0.06e-3;
00841 }
00842
00843 }else{
00844 EventRate=Nsel;
00845 rateErr_stat =Nsel;
00846 rateErr_syst =Nsel;
00847 std::cout << std::endl;
00848 if (MinP > 100.)
00849 std::cout << " !!! MinP > 100 GeV. Cannot apply normalisation!" << std::endl;
00850 else if (MaxTheta > 84.26*Deg2Rad)
00851 std::cout << " !!! Note: generated cosmics exceed parameterisation. No flux calculated!" << std::endl;
00852
00853 else
00854 std::cout << " !!! Not enough statistics to apply normalisation (rate=1 +- 1) !!!" << std::endl;
00855 }
00856
00857 std::cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << std::endl;
00858 std::cout << " rate is " << EventRate << " +-" << rateErr_stat <<" (stat) " << "+-" <<
00859 rateErr_syst << " (syst) " <<" muons per second" << std::endl;
00860 if(EventRate!=0) std::cout << " number of events corresponds to " << Nsel/EventRate << " s" << std::endl;
00861 std::cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << std::endl;
00862 std::cout << std::endl;
00863 std::cout << "*********************************************************" << std::endl;
00864 std::cout << "*********************************************************" << std::endl;
00865 }
00866 }
00867
00868 void CosmicMuonGenerator::checkIn(){
00869 if (MinP < 0.){ NumberOfEvents = 0;
00870 std::cout << " CMG-ERR: min.energy is out of range (0 GeV ... inf]" << std::endl << std::endl; }
00871 if (MaxP < 0.){ NumberOfEvents = 0;
00872 std::cout << " CMG-ERR: max.energy is out of range (0 GeV ... inf]" << std::endl << std::endl; }
00873 if (MaxP <= MinP){ NumberOfEvents = 0;
00874 std::cout << " CMG-ERR: max.energy is not greater than min.energy" << std::endl << std::endl; }
00875 if (MinTheta < 0.){ NumberOfEvents = 0;
00876 std::cout << " CMG-ERR: min.theta is out of range [0 deg ... 90 deg)" << std::endl << std::endl; }
00877 if (MaxTheta < 0.){ NumberOfEvents = 0;
00878 std::cout << " CMG-ERR: max.theta is out of range [0 deg ... 90 deg)" << std::endl << std::endl; }
00879 if (MaxTheta <= MinTheta){ NumberOfEvents = 0;
00880 std::cout << " CMG-ERR: max.theta is not greater than min.theta" << std::endl << std::endl; }
00881 if (MinPhi < 0.){ NumberOfEvents = 0;
00882 std::cout << " CMG-ERR: min.phi is out of range [0 deg ... 360 deg]" << std::endl << std::endl; }
00883 if (MaxPhi < 0.){ NumberOfEvents = 0;
00884 std::cout << " CMG-ERR: max.phi is out of range [0 deg ... 360 deg]" << std::endl << std::endl; }
00885 if (MaxPhi <= MinPhi){ NumberOfEvents = 0;
00886 std::cout << " CMG-ERR: max.phi is not greater than min.phi" << std::endl << std::endl; }
00887 if (MaxT0 <= MinT0){ NumberOfEvents = 0;
00888 std::cout << " CMG-ERR: max.t0 is not greater than min.t0" << std::endl << std::endl; }
00889 if (ElossScaleFactor < 0.){ NumberOfEvents = 0;
00890 std::cout << " CMG-ERR: E-loss scale factor is out of range [0 ... inf)" << std::endl << std::endl; }
00891 if (MinEnu < 0.){ NumberOfEvents = 0;
00892 std::cout << " CMG-ERR: min.Enu is out of range [0 GeV ... inf]" << std::endl << std::endl; }
00893 if (MaxEnu < 0.){ NumberOfEvents = 0;
00894 std::cout << " CMG-ERR: max.Enu is out of range [0 GeV ... inf]" << std::endl << std::endl; }
00895 if (MaxEnu <= MinEnu){ NumberOfEvents = 0;
00896 std::cout << " CMG-ERR: max.Enu is not greater than min.Enu" << std::endl << std::endl; }
00897
00898 }
00899
00900 bool CosmicMuonGenerator::goodOrientation(){
00901
00902 bool goodAngles = false;
00903 bool phiaccepted = false;
00904 bool thetaaccepted = false;
00905 double RxzV = sqrt(OneMuoEvt.vx()*OneMuoEvt.vx() + OneMuoEvt.vz()*OneMuoEvt.vz());
00906
00907 double rVY;
00908 if (MinTheta > 90.*Deg2Rad)
00909 rVY = -sqrt(RxzV*RxzV + RadiusCMS*RadiusCMS);
00910 else
00911 rVY = sqrt(RxzV*RxzV + (SurfaceOfEarth+PlugWidth)*(SurfaceOfEarth+PlugWidth));
00912
00913 double Phi = OneMuoEvt.phi();
00914 double PhiV = atan2(OneMuoEvt.vx(),OneMuoEvt.vz()) + Pi; if (PhiV > TwoPi) PhiV -= TwoPi;
00915 double disPhi = std::fabs(PhiV - Phi); if (disPhi > Pi) disPhi = TwoPi - disPhi;
00916 double dPhi = Pi; if (RxzV > Target3dRadius) dPhi = asin(Target3dRadius/RxzV);
00917 if (disPhi < dPhi) phiaccepted = true;
00918 double Theta = OneMuoEvt.theta();
00919 double ThetaV = asin(RxzV/rVY);
00920 double dTheta = Pi; if (std::fabs(rVY) > Target3dRadius) dTheta = asin(Target3dRadius/std::fabs(rVY));
00921
00922
00923
00924 if (!phiaccepted && RxzV < Target3dRadius)
00925
00926 std::cout << "Rejected phi=" << Phi << " PhiV=" << PhiV
00927 << " dPhi=" << dPhi << " disPhi=" << disPhi
00928 << " RxzV=" << RxzV << " Target3dRadius=" << Target3dRadius
00929 << " Theta=" << Theta << std::endl;
00930
00931 if (std::fabs(Theta-ThetaV) < dTheta) thetaaccepted = true;
00932 if (phiaccepted && thetaaccepted) goodAngles = true;
00933 return goodAngles;
00934 }
00935
00936 void CosmicMuonGenerator::initEvDis(){
00937 #if ROOT_INTERACTIVE
00938 float rCMS = RadiusCMS/1000.;
00939 float zCMS = Z_DistCMS/1000.;
00940 if(TrackerOnly==true){
00941 rCMS = RadiusTracker/1000.;
00942 zCMS = Z_DistTracker/1000.;
00943 }
00944 TH2F* disXY = new TH2F("disXY","X-Y view",160,-rCMS,rCMS,160,-rCMS,rCMS);
00945 TH2F* disZY = new TH2F("disZY","Z-Y view",150,-zCMS,zCMS,160,-rCMS,rCMS);
00946 gStyle->SetPalette(1,0);
00947 gStyle->SetMarkerColor(1);
00948 gStyle->SetMarkerSize(1.5);
00949 TCanvas *disC = new TCanvas("disC","Cosmic Muon Event Display",0,0,800,410);
00950 disC->Divide(2,1);
00951 disC->cd(1);
00952 gPad->SetTicks(1,1);
00953 disXY->SetMinimum(log10(MinP));
00954 disXY->SetMaximum(log10(MaxP));
00955 disXY->GetXaxis()->SetLabelSize(0.05);
00956 disXY->GetXaxis()->SetTitleSize(0.05);
00957 disXY->GetXaxis()->SetTitleOffset(1.0);
00958 disXY->GetXaxis()->SetTitle("X [m]");
00959 disXY->GetYaxis()->SetLabelSize(0.05);
00960 disXY->GetYaxis()->SetTitleSize(0.05);
00961 disXY->GetYaxis()->SetTitleOffset(0.8);
00962 disXY->GetYaxis()->SetTitle("Y [m]");
00963 disC->cd(2);
00964 gPad->SetGrid(1,1);
00965 gPad->SetTicks(1,1);
00966 disZY->SetMinimum(log10(MinP));
00967 disZY->SetMaximum(log10(MaxP));
00968 disZY->GetXaxis()->SetLabelSize(0.05);
00969 disZY->GetXaxis()->SetTitleSize(0.05);
00970 disZY->GetXaxis()->SetTitleOffset(1.0);
00971 disZY->GetXaxis()->SetTitle("Z [m]");
00972 disZY->GetYaxis()->SetLabelSize(0.05);
00973 disZY->GetYaxis()->SetTitleSize(0.05);
00974 disZY->GetYaxis()->SetTitleOffset(0.8);
00975 disZY->GetYaxis()->SetTitle("Y [m]");
00976 #endif
00977 }
00978
00979 void CosmicMuonGenerator::displayEv(){
00980 #if ROOT_INTERACTIVE
00981 double RadiusDet=RadiusCMS;
00982 double Z_DistDet=Z_DistCMS;
00983 if(TrackerOnly==true){
00984 RadiusDet = RadiusTracker;
00985 Z_DistDet = Z_DistTracker;
00986 }
00987 disXY->Reset();
00988 disZY->Reset();
00989 TMarker* InteractionPoint = new TMarker(0.,0.,2);
00990 TArc* r8m = new TArc(0.,0.,(RadiusDet/1000.));
00991 TLatex* logEaxis = new TLatex(); logEaxis->SetTextSize(0.05);
00992 float energy = float(OneMuoEvt.e());
00993 float verX = float(OneMuoEvt.vx()/1000.);
00994 float verY = float(OneMuoEvt.vy()/1000.);
00995 float verZ = float(OneMuoEvt.vz()/1000.);
00996 float dirX = float(OneMuoEvt.px())/std::fabs(OneMuoEvt.py());
00997 float dirY = float(OneMuoEvt.py())/std::fabs(OneMuoEvt.py());
00998 float dirZ = float(OneMuoEvt.pz())/std::fabs(OneMuoEvt.py());
00999 float yStep = disXY->GetYaxis()->GetBinWidth(1);
01000 int NbinY = disXY->GetYaxis()->GetNbins();
01001 for (int iy=0; iy<NbinY; ++iy){
01002 verX += dirX*yStep;
01003 verY += dirY*yStep;
01004 verZ += dirZ*yStep;
01005 float rXY = sqrt(verX*verX + verY*verY)*1000.;
01006 float absZ = std::fabs(verZ)*1000.;
01007 if (rXY < RadiusDet && absZ < Z_DistDet){
01008 disXY->Fill(verX,verY,log10(energy));
01009 disZY->Fill(verZ,verY,log10(energy));
01010 disC->cd(1); disXY->Draw("COLZ"); InteractionPoint->Draw("SAME"); r8m->Draw("SAME");
01011 logEaxis->DrawLatex((0.65*RadiusDet/1000.),(1.08*RadiusDet/1000.),"log_{10}E(#mu^{#pm})");
01012 disC->cd(2); disZY->Draw("COL"); InteractionPoint->Draw("SAME");
01013 gPad->Update();
01014 }
01015 }
01016 #endif
01017 }
01018
01019 void CosmicMuonGenerator::setNumberOfEvents(unsigned int N){ if (NotInitialized) NumberOfEvents = N; }
01020
01021 void CosmicMuonGenerator::setRanSeed(int N){ if (NotInitialized) RanSeed = N; }
01022
01023 void CosmicMuonGenerator::setMinP(double P){ if (NotInitialized) MinP = P; }
01024
01025 void CosmicMuonGenerator::setMinP_CMS(double P){ if (NotInitialized) MinP_CMS = P; }
01026
01027 void CosmicMuonGenerator::setMaxP(double P){ if (NotInitialized) MaxP = P; }
01028
01029 void CosmicMuonGenerator::setMinTheta(double Theta){ if (NotInitialized) MinTheta = Theta*Deg2Rad; }
01030
01031 void CosmicMuonGenerator::setMaxTheta(double Theta){ if (NotInitialized) MaxTheta = Theta*Deg2Rad; }
01032
01033 void CosmicMuonGenerator::setMinPhi(double Phi){ if (NotInitialized) MinPhi = Phi*Deg2Rad; }
01034
01035 void CosmicMuonGenerator::setMaxPhi(double Phi){ if (NotInitialized) MaxPhi = Phi*Deg2Rad; }
01036
01037 void CosmicMuonGenerator::setMinT0(double T0){ if (NotInitialized) MinT0 = T0; }
01038
01039 void CosmicMuonGenerator::setMaxT0(double T0){ if (NotInitialized) MaxT0 = T0; }
01040
01041 void CosmicMuonGenerator::setElossScaleFactor(double ElossScaleFact){ if (NotInitialized) ElossScaleFactor = ElossScaleFact; }
01042
01043 void CosmicMuonGenerator::setRadiusOfTarget(double R){ if (NotInitialized) RadiusOfTarget = R; }
01044
01045 void CosmicMuonGenerator::setZDistOfTarget(double Z){ if (NotInitialized) ZDistOfTarget = Z; }
01046
01047 void CosmicMuonGenerator::setZCentrOfTarget(double Z){ if (NotInitialized) ZCentrOfTarget = Z; }
01048
01049 void CosmicMuonGenerator::setTrackerOnly(bool Tracker){ if (NotInitialized) TrackerOnly = Tracker; }
01050
01051 void CosmicMuonGenerator::setMultiMuon(bool MultiMu){ if (NotInitialized) MultiMuon = MultiMu; }
01052 void CosmicMuonGenerator::setMultiMuonFileName(std::string MultiMuFile){ if (NotInitialized) MultiMuonFileName = MultiMuFile; }
01053 void CosmicMuonGenerator::setMultiMuonFileFirstEvent(int MultiMuFile1stEvt){ if (NotInitialized) MultiMuonFileFirstEvent = MultiMuFile1stEvt; }
01054 void CosmicMuonGenerator::setMultiMuonNmin(int MultiMuNmin){ if (NotInitialized) MultiMuonNmin = MultiMuNmin; }
01055
01056 void CosmicMuonGenerator::setTIFOnly_constant(bool TIF){ if (NotInitialized) TIFOnly_constant = TIF; }
01057
01058 void CosmicMuonGenerator::setTIFOnly_linear(bool TIF){ if (NotInitialized) TIFOnly_linear = TIF; }
01059 void CosmicMuonGenerator::setMTCCHalf(bool MTCC){ if (NotInitialized) MTCCHalf = MTCC; }
01060
01061 void CosmicMuonGenerator::setPlugVx(double PlugVtx){ if (NotInitialized) PlugVx = PlugVtx; }
01062 void CosmicMuonGenerator::setPlugVz(double PlugVtz){ if (NotInitialized) PlugVz = PlugVtz; }
01063 void CosmicMuonGenerator::setRhoAir(double VarRhoAir){ if (NotInitialized) RhoAir = VarRhoAir; }
01064 void CosmicMuonGenerator::setRhoWall(double VarRhoWall){ if (NotInitialized) RhoWall = VarRhoWall; }
01065 void CosmicMuonGenerator::setRhoRock(double VarRhoRock){ if (NotInitialized) RhoRock = VarRhoRock; }
01066 void CosmicMuonGenerator::setRhoClay(double VarRhoClay){ if (NotInitialized) RhoClay = VarRhoClay; }
01067 void CosmicMuonGenerator::setRhoPlug(double VarRhoPlug){ if (NotInitialized) RhoPlug = VarRhoPlug; }
01068 void CosmicMuonGenerator::setClayWidth(double ClayLayerWidth){ if (NotInitialized) ClayWidth = ClayLayerWidth; }
01069
01070 void CosmicMuonGenerator::setMinEnu(double MinEn){ if (NotInitialized) MinEnu = MinEn; }
01071 void CosmicMuonGenerator::setMaxEnu(double MaxEn){ if (NotInitialized) MaxEnu = MaxEn; }
01072 void CosmicMuonGenerator::setNuProdAlt(double NuPrdAlt){ if (NotInitialized) NuProdAlt = NuPrdAlt; }
01073
01074 double CosmicMuonGenerator::getRate(){ return EventRate; }
01075
01076 void CosmicMuonGenerator::setAcptAllMu(bool AllMu){ if (NotInitialized) AcptAllMu = AllMu; }