1 |
|
---|
2 | /** \class Hector
|
---|
3 | *
|
---|
4 | * Propagates candidates using Hector library.
|
---|
5 | *
|
---|
6 | * $Date: 2014-04-14 14:53:45 +0000 (Mon, 14 Apr 2014) $
|
---|
7 | * $Revision: 1362 $
|
---|
8 | *
|
---|
9 | *
|
---|
10 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
11 | *
|
---|
12 | */
|
---|
13 |
|
---|
14 | #include "modules/Hector.h"
|
---|
15 |
|
---|
16 | #include "classes/DelphesClasses.h"
|
---|
17 | #include "classes/DelphesFactory.h"
|
---|
18 | #include "classes/DelphesFormula.h"
|
---|
19 |
|
---|
20 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
21 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
22 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
23 |
|
---|
24 | #include "TMath.h"
|
---|
25 | #include "TString.h"
|
---|
26 | #include "TFormula.h"
|
---|
27 | #include "TRandom3.h"
|
---|
28 | #include "TObjArray.h"
|
---|
29 | #include "TDatabasePDG.h"
|
---|
30 | #include "TLorentzVector.h"
|
---|
31 |
|
---|
32 | #include <algorithm>
|
---|
33 | #include <stdexcept>
|
---|
34 | #include <iostream>
|
---|
35 | #include <sstream>
|
---|
36 |
|
---|
37 | #include "Hector/H_BeamLine.h"
|
---|
38 | #include "Hector/H_RecRPObject.h"
|
---|
39 | #include "Hector/H_BeamParticle.h"
|
---|
40 |
|
---|
41 | using namespace std;
|
---|
42 |
|
---|
43 | //------------------------------------------------------------------------------
|
---|
44 |
|
---|
45 | Hector::Hector() :
|
---|
46 | fBeamLine(0), fItInputArray(0)
|
---|
47 | {
|
---|
48 | }
|
---|
49 |
|
---|
50 | //------------------------------------------------------------------------------
|
---|
51 |
|
---|
52 | Hector::~Hector()
|
---|
53 | {
|
---|
54 | }
|
---|
55 |
|
---|
56 | //------------------------------------------------------------------------------
|
---|
57 |
|
---|
58 | void Hector::Init()
|
---|
59 | {
|
---|
60 | // read Hector parameters
|
---|
61 |
|
---|
62 | fDirection = GetInt("Direction", 1);
|
---|
63 | fBeamLineLength = GetDouble("BeamLineLength", 430.0);
|
---|
64 | fDistance = GetDouble("Distance", 420.0);
|
---|
65 | fOffsetX = GetDouble("OffsetX", 0.0);
|
---|
66 | fOffsetS = GetDouble("OffsetS", 120.0);
|
---|
67 | fSigmaE = GetDouble("SigmaE", 0.0);
|
---|
68 | fSigmaX = GetDouble("SigmaX", 0.0);
|
---|
69 | fSigmaY = GetDouble("SigmaY", 0.0);
|
---|
70 | fSigmaT = GetDouble("SigmaT", 0.0);
|
---|
71 | fEtaMin = GetDouble("EtaMin", 5.0);
|
---|
72 |
|
---|
73 | fBeamLine = new H_BeamLine(fDirection, fBeamLineLength + 0.1);
|
---|
74 | fBeamLine->fill(GetString("BeamLineFile", "examples/LHCB1IR5_5TeV.tfs"), fDirection, GetString("IPName", "IP5"));
|
---|
75 | fBeamLine->offsetElements(fOffsetS, fOffsetX);
|
---|
76 | fBeamLine->calcMatrix();
|
---|
77 |
|
---|
78 | // import input array
|
---|
79 |
|
---|
80 | fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
|
---|
81 | fItInputArray = fInputArray->MakeIterator();
|
---|
82 |
|
---|
83 | // create output array
|
---|
84 |
|
---|
85 | fOutputArray = ExportArray(GetString("OutputArray", "hits"));
|
---|
86 | }
|
---|
87 |
|
---|
88 | //------------------------------------------------------------------------------
|
---|
89 |
|
---|
90 | void Hector::Finish()
|
---|
91 | {
|
---|
92 | if(fItInputArray) delete fItInputArray;
|
---|
93 | if(fBeamLine) delete fBeamLine;
|
---|
94 | }
|
---|
95 |
|
---|
96 | //------------------------------------------------------------------------------
|
---|
97 |
|
---|
98 | void Hector::Process()
|
---|
99 | {
|
---|
100 | Candidate *candidate, *mother;
|
---|
101 | Double_t pz;
|
---|
102 | Double_t x, y, z, tx, ty, theta;
|
---|
103 | Double_t distance, time;
|
---|
104 |
|
---|
105 | const Double_t c_light = 2.99792458E8;
|
---|
106 |
|
---|
107 | fItInputArray->Reset();
|
---|
108 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
109 | {
|
---|
110 | const TLorentzVector &candidatePosition = candidate->Position;
|
---|
111 | const TLorentzVector &candidateMomentum = candidate->Momentum;
|
---|
112 | pz = candidateMomentum.Pz();
|
---|
113 |
|
---|
114 | if(candidateMomentum.Eta() <= fEtaMin || fDirection*pz <= 0.0) continue;
|
---|
115 |
|
---|
116 | x = 1.0E3 * candidatePosition.X();
|
---|
117 | y = 1.0E3 * candidatePosition.Y();
|
---|
118 | z = 1.0E-3 * candidatePosition.Z();
|
---|
119 |
|
---|
120 | // tx = 1.0E6 * TMath::ATan(candidateMomentum.Px()/pz);
|
---|
121 | // ty = 1.0E6 * TMath::ATan(candidateMomentum.Py()/pz);
|
---|
122 |
|
---|
123 | tx = 0.0;
|
---|
124 | ty = 0.0;
|
---|
125 |
|
---|
126 | theta = TMath::Hypot(TMath::ATan(candidateMomentum.Px()/pz), TMath::ATan(candidateMomentum.Py()/pz));
|
---|
127 | distance = (fDistance - 1.0E-3 * candidatePosition.Z())/TMath::Cos(theta);
|
---|
128 | time = gRandom->Gaus((distance + 1.0E-3 * candidatePosition.T())/c_light, fSigmaT);
|
---|
129 |
|
---|
130 | H_BeamParticle particle(candidate->Mass, candidate->Charge);
|
---|
131 | particle.set4Momentum(candidateMomentum);
|
---|
132 | particle.setPosition(x, y, tx, ty, z);
|
---|
133 |
|
---|
134 | particle.smearAng(fSigmaX, fSigmaY, gRandom);
|
---|
135 | particle.smearE(fSigmaE, gRandom);
|
---|
136 |
|
---|
137 | particle.computePath(fBeamLine);
|
---|
138 |
|
---|
139 | if(particle.stopped(fBeamLine)) continue;
|
---|
140 |
|
---|
141 | particle.propagate(fDistance);
|
---|
142 |
|
---|
143 | mother = candidate;
|
---|
144 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
145 | candidate->Position.SetXYZT(particle.getX(), particle.getY(), particle.getS(), time);
|
---|
146 | candidate->Momentum.SetPxPyPzE(particle.getTX(), particle.getTY(), 0.0, particle.getE());
|
---|
147 | candidate->AddCandidate(mother);
|
---|
148 |
|
---|
149 | fOutputArray->Add(candidate);
|
---|
150 | }
|
---|
151 | }
|
---|
152 |
|
---|
153 | //------------------------------------------------------------------------------
|
---|