Fork me on GitHub

source: git/modules/Hector.cc@ 5b822e5

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5b822e5 was 5b822e5, checked in by pavel <pavel@…>, 10 years ago

add Hector module

  • Property mode set to 100644
File size: 3.7 KB
Line 
1
2/** \class Hector
3 *
4 * Propagates candidates using Hector library.
5 *
6 * $Date$
7 * $Revision$
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
41using namespace std;
42
43//------------------------------------------------------------------------------
44
45Hector::Hector() :
46 fBeamLine(0), fItInputArray(0)
47{
48}
49
50//------------------------------------------------------------------------------
51
52Hector::~Hector()
53{
54}
55
56//------------------------------------------------------------------------------
57
58void 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 fSigmaE = GetDouble("SigmaE", 0.0);
66 fSigmaX = GetDouble("SigmaX", 0.0);
67 fSigmaY = GetDouble("SigmaY", 0.0);
68 fEtaMin = GetDouble("EtaMin", 5.0);
69
70 fBeamLine = new H_BeamLine(fDirection, fBeamLineLength + 0.1);
71 fBeamLine->fill(GetString("BeamLineFile", "examples/LHCB1IR5_5TeV.tfs"), fDirection, "IP5" );
72 fBeamLine->offsetElements(120, -0.097*fDirection);
73 fBeamLine->calcMatrix();
74
75 // import input array
76
77 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
78 fItInputArray = fInputArray->MakeIterator();
79
80 // create output array
81
82 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
83}
84
85//------------------------------------------------------------------------------
86
87void Hector::Finish()
88{
89 if(fItInputArray) delete fItInputArray;
90 if(fBeamLine) delete fBeamLine;
91}
92
93//------------------------------------------------------------------------------
94
95void Hector::Process()
96{
97 Candidate *candidate, *mother;
98 Double_t pz;
99 Double_t x, y, z, tx, ty;
100
101 fItInputArray->Reset();
102 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
103 {
104 const TLorentzVector &candidatePosition = candidate->Position;
105 const TLorentzVector &candidateMomentum = candidate->Momentum;
106 pz = candidateMomentum.Pz();
107
108 if(candidateMomentum.Eta() <= fEtaMin || fDirection*pz <= 0.0) continue;
109
110 x = 1.0E3 * candidatePosition.X();
111 y = 1.0E3 * candidatePosition.Y();
112 z = 1.0E-2 * candidatePosition.Z();
113
114// tx = 1.0E6 * TMath::ATan(candidateMomentum.Px()/pz);
115// ty = 1.0E6 * TMath::ATan(candidateMomentum.Py()/pz);
116
117 tx = 0.0;
118 ty = 0.0;
119
120 H_BeamParticle particle(candidate->Mass, candidate->Charge);
121 particle.set4Momentum(candidateMomentum);
122 particle.setPosition(x, y, tx, ty, z);
123
124 particle.smearAng(fSigmaX, fSigmaY, gRandom);
125 particle.smearE(fSigmaE, gRandom);
126
127 particle.computePath(fBeamLine);
128
129 if(particle.stopped(fBeamLine)) continue;
130
131 particle.propagate(fDistance);
132
133 mother = candidate;
134 candidate = static_cast<Candidate*>(candidate->Clone());
135 candidate->Position.SetXYZT(particle.getX(), particle.getY(), fDistance, 0.0);
136 candidate->Momentum.SetPxPyPzE(particle.getTX(), particle.getTY(), 0.0, particle.getE());
137 candidate->AddCandidate(mother);
138
139 fOutputArray->Add(candidate);
140 }
141}
142
143//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.