- Timestamp:
- Apr 8, 2014, 6:43:52 PM (11 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- b67be31
- Parents:
- 5b822e5
- Location:
- modules
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
modules/Hector.cc
r5b822e5 r8f7db23 63 63 fBeamLineLength = GetDouble("BeamLineLength", 430.0); 64 64 fDistance = GetDouble("Distance", 420.0); 65 fOffsetX = GetDouble("OffsetX", 0.0); 66 fOffsetS = GetDouble("OffsetS", 120.0); 65 67 fSigmaE = GetDouble("SigmaE", 0.0); 66 68 fSigmaX = GetDouble("SigmaX", 0.0); 67 69 fSigmaY = GetDouble("SigmaY", 0.0); 70 fSigmaT = GetDouble("SigmaT", 0.0); 68 71 fEtaMin = GetDouble("EtaMin", 5.0); 69 72 70 73 fBeamLine = new H_BeamLine(fDirection, fBeamLineLength + 0.1); 71 fBeamLine->fill(GetString("BeamLineFile", "examples/LHCB1IR5_5TeV.tfs"), fDirection, "IP5");74 fBeamLine->fill(GetString("BeamLineFile", "examples/LHCB1IR5_5TeV.tfs"), fDirection, GetString("IPName", "IP5")); 72 75 fBeamLine->offsetElements(120, -0.097*fDirection); 76 fBeamLine->offsetElements(fOffsetS, fOffsetX); 73 77 fBeamLine->calcMatrix(); 74 78 … … 80 84 // create output array 81 85 82 fOutputArray = ExportArray(GetString("OutputArray", " stableParticles"));86 fOutputArray = ExportArray(GetString("OutputArray", "hits")); 83 87 } 84 88 … … 97 101 Candidate *candidate, *mother; 98 102 Double_t pz; 99 Double_t x, y, z, tx, ty; 103 Double_t x, y, z, tx, ty, theta; 104 Double_t distance, time; 105 106 const Double_t c_light = 2.99792458E8; 100 107 101 108 fItInputArray->Reset(); … … 118 125 ty = 0.0; 119 126 127 theta = TMath::Hypot(TMath::ATan(candidateMomentum.Px()/pz), TMath::ATan(candidateMomentum.Py()/pz)); 128 distance = (fDistance - 1.0E-3 * candidatePosition.Z())/TMath::Cos(theta); 129 time = gRandom->Gaus((distance + 1.0E-3 * candidatePosition.T())/c_light, fSigmaT); 130 120 131 H_BeamParticle particle(candidate->Mass, candidate->Charge); 121 132 particle.set4Momentum(candidateMomentum); … … 133 144 mother = candidate; 134 145 candidate = static_cast<Candidate*>(candidate->Clone()); 135 candidate->Position.SetXYZT(particle.getX(), particle.getY(), fDistance, 0.0);146 candidate->Position.SetXYZT(particle.getX(), particle.getY(), particle.getS(), time); 136 147 candidate->Momentum.SetPxPyPzE(particle.getTX(), particle.getTY(), 0.0, particle.getE()); 137 148 candidate->AddCandidate(mother); -
modules/Hector.h
r5b822e5 r8f7db23 34 34 35 35 Int_t fDirection; 36 36 37 37 Double_t fBeamLineLength, fDistance; 38 Double_t fSigmaE, fSigmaX, fSigmaY; 38 Double_t fOffsetX, fOffsetS; 39 Double_t fSigmaE, fSigmaX, fSigmaY, fSigmaT; 39 40 Double_t fEtaMin; 40 41 41 42 H_BeamLine *fBeamLine; 42 43 -
modules/TreeWriter.cc
r5b822e5 r8f7db23 67 67 fClassMap[Rho::Class()] = &TreeWriter::ProcessRho; 68 68 fClassMap[Weight::Class()] = &TreeWriter::ProcessWeight; 69 fClassMap[HectorHit::Class()] = &TreeWriter::ProcessHectorHit; 69 70 70 71 TBranchMap::iterator itBranchMap; … … 95 96 continue; 96 97 } 97 98 98 99 itClassMap = fClassMap.find(branchClass); 99 100 if(itClassMap == fClassMap.end()) … … 160 161 GenParticle *entry = 0; 161 162 Double_t pt, signPz, cosTheta, eta, rapidity; 162 163 const Double_t c_light = 2.99792458E8; 164 163 164 const Double_t c_light = 2.99792458E8; 165 165 166 // loop over all particles 166 167 iterator.Reset(); … … 220 221 Candidate *candidate = 0; 221 222 Vertex *entry = 0; 222 223 223 224 const Double_t c_light = 2.99792458E8; 224 225 … … 248 249 Double_t pt, signz, cosTheta, eta, rapidity; 249 250 const Double_t c_light = 2.99792458E8; 250 251 251 252 // loop over all tracks 252 253 iterator.Reset(); … … 310 311 Double_t pt, signPz, cosTheta, eta, rapidity; 311 312 const Double_t c_light = 2.99792458E8; 312 313 313 314 // loop over all towers 314 315 iterator.Reset(); … … 317 318 const TLorentzVector &momentum = candidate->Momentum; 318 319 const TLorentzVector &position = candidate->Position; 319 320 320 321 pt = momentum.Pt(); 321 322 cosTheta = TMath::Abs(momentum.CosTheta()); … … 339 340 entry->Edges[2] = candidate->Edges[2]; 340 341 entry->Edges[3] = candidate->Edges[3]; 341 342 342 343 entry->T = position.T()*1.0E-3/c_light; 343 344 344 345 FillParticles(candidate, &entry->Particles); 345 346 } … … 355 356 Double_t pt, signPz, cosTheta, eta, rapidity; 356 357 const Double_t c_light = 2.99792458E8; 357 358 358 359 array->Sort(); 359 360 … … 365 366 const TLorentzVector &momentum = candidate->Momentum; 366 367 const TLorentzVector &position = candidate->Position; 367 368 368 369 369 370 pt = momentum.Pt(); … … 379 380 entry->PT = pt; 380 381 entry->E = momentum.E(); 381 382 382 383 entry->T = position.T()*1.0E-3/c_light; 383 384 384 385 entry->EhadOverEem = candidate->Eem > 0.0 ? candidate->Ehad/candidate->Eem : 999.9; 385 386 … … 397 398 Double_t pt, signPz, cosTheta, eta, rapidity; 398 399 const Double_t c_light = 2.99792458E8; 399 400 400 401 array->Sort(); 401 402 … … 406 407 const TLorentzVector &momentum = candidate->Momentum; 407 408 const TLorentzVector &position = candidate->Position; 408 409 409 410 pt = momentum.Pt(); 410 411 cosTheta = TMath::Abs(momentum.CosTheta()); … … 418 419 entry->Phi = momentum.Phi(); 419 420 entry->PT = pt; 420 421 421 422 entry->T = position.T()*1.0E-3/c_light; 422 423 423 424 entry->Charge = candidate->Charge; 424 425 … … 437 438 Muon *entry = 0; 438 439 Double_t pt, signPz, cosTheta, eta, rapidity; 439 440 const Double_t c_light = 2.99792458E8; 441 440 441 const Double_t c_light = 2.99792458E8; 442 442 443 array->Sort(); 443 444 … … 448 449 const TLorentzVector &momentum = candidate->Momentum; 449 450 const TLorentzVector &position = candidate->Position; 450 451 451 452 452 453 pt = momentum.Pt(); … … 466 467 467 468 entry->T = position.T()*1.0E-3/c_light; 468 469 469 470 // cout<<entry->PT<<","<<entry->T<<endl; 470 471 … … 485 486 Double_t ecalEnergy, hcalEnergy; 486 487 const Double_t c_light = 2.99792458E8; 487 488 488 489 array->Sort(); 489 490 … … 493 494 { 494 495 TIter itConstituents(candidate->GetCandidates()); 495 496 const TLorentzVector &momentum = candidate->Momentum; 497 const TLorentzVector &position = candidate->Position; 498 496 497 const TLorentzVector &momentum = candidate->Momentum; 498 const TLorentzVector &position = candidate->Position; 499 499 500 pt = momentum.Pt(); 500 501 cosTheta = TMath::Abs(momentum.CosTheta()); … … 510 511 511 512 entry->T = position.T()*1.0E-3/c_light; 512 513 513 514 entry->Mass = momentum.M(); 514 515 … … 533 534 534 535 entry->EhadOverEem = ecalEnergy > 0.0 ? hcalEnergy/ecalEnergy : 999.9; 535 536 536 537 //--- Pile-Up Jet ID variables ---- 537 538 … … 632 633 //------------------------------------------------------------------------------ 633 634 635 void TreeWriter::ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array) 636 { 637 TIter iterator(array); 638 Candidate *candidate = 0; 639 HectorHit *entry = 0; 640 641 // loop over all roman pot hits 642 iterator.Reset(); 643 while((candidate = static_cast<Candidate*>(iterator.Next()))) 644 { 645 const TLorentzVector &position = candidate->Position; 646 const TLorentzVector &momentum = candidate->Momentum; 647 648 entry = static_cast<HectorHit*>(branch->NewEntry()); 649 650 entry->E = momentum.E(); 651 652 entry->Tx = momentum.Px(); 653 entry->Ty = momentum.Py(); 654 655 entry->T = position.T(); 656 657 entry->X = position.X(); 658 entry->Y = position.Y(); 659 entry->S = position.Z(); 660 } 661 } 662 663 //------------------------------------------------------------------------------ 664 634 665 void TreeWriter::Process() 635 666 { -
modules/TreeWriter.h
r5b822e5 r8f7db23 53 53 void ProcessRho(ExRootTreeBranch *branch, TObjArray *array); 54 54 void ProcessWeight(ExRootTreeBranch *branch, TObjArray *array); 55 void ProcessHectorHit(ExRootTreeBranch *branch, TObjArray *array); 55 56 56 57 #ifndef __CINT__
Note:
See TracChangeset
for help on using the changeset viewer.