CMS 3D CMS Logo

Hector.h
Go to the documentation of this file.
1 #ifndef SimTransport_Hector_h
2 #define SimTransport_Hector_h
3 
4 // user include files
7 
8 // HepMC headers
10 
11 #include "HepMC/GenEvent.h"
12 #include "HepMC/GenParticle.h"
13 #include "HepMC/GenVertex.h"
15 
17 #include <vector>
18 
19 // SimpleConfigurable replacement
21 
22 // Hector headers
23 #include "H_BeamLine.h"
24 #include "H_BeamParticle.h"
25 #include "H_RecRPObject.h"
26 #include <map>
27 #include <string>
28 
29 class TRandom3;
30 
31 class Hector {
32 public:
33  // Hector(const edm::ParameterSet & ps);
34  Hector(const edm::ParameterSet &ps,
36  bool verbosity,
37  bool FP420Transport,
38  bool ZDCTransport);
39  // Hector();
40  virtual ~Hector();
41 
43  void clearApertureFlags();
46  void clear();
48  void add(const HepMC::GenEvent *ev, const edm::EventSetup &es);
50  void filterFP420(TRandom3 *);
52  void filterZDC(TRandom3 *);
54  void filterD1(TRandom3 *);
55 
56  int getDirect(unsigned int part_n) const;
57 
59  void print() const;
62  // std::vector<unsigned int> part_list() const;
63 
64  // bool isCharged(const HepMC::GenParticle * p);
65 
67 
68  std::vector<LHCTransportLink> &getCorrespondenceMap() { return theCorrespondenceMap; }
69 
70  /*
71  private:
72  // edm::ParameterSet m_pBeamLine;
73  */
74 private:
76  // Defaults
77  double lengthfp420;
78  double lengthzdc;
79  double lengthd1;
80 
81  double etacut;
82  bool m_smearAng;
83  double m_sig_e;
84  bool m_smearE;
85  double m_sigmaSTX;
86  double m_sigmaSTY;
87 
88  float m_rpp420_f;
89  float m_rpp420_b;
90  float m_rppzdc;
91  float m_rppd1;
92 
94 
95  // Hector
96  H_BeamLine *m_beamlineFP4201;
97  H_BeamLine *m_beamlineFP4202;
98  H_BeamLine *m_beamlineZDC1;
99  H_BeamLine *m_beamlineZDC2;
100  H_BeamLine *m_beamlineD11;
101  H_BeamLine *m_beamlineD12;
102  //
103 
104  H_RecRPObject *m_rp420_f;
105  H_RecRPObject *m_rp420_b;
106 
107  std::map<unsigned int, H_BeamParticle *> m_beamPart;
108  std::map<unsigned int, int> m_direct;
109  std::map<unsigned int, bool> m_isStoppedfp420;
110  std::map<unsigned int, bool> m_isStoppedzdc;
111  std::map<unsigned int, bool> m_isStoppedd1;
112  std::map<unsigned int, double> m_xAtTrPoint;
113  std::map<unsigned int, double> m_yAtTrPoint;
114  std::map<unsigned int, double> m_TxAtTrPoint;
115  std::map<unsigned int, double> m_TyAtTrPoint;
116  std::map<unsigned int, double> m_eAtTrPoint;
117 
118  std::map<unsigned int, double> m_eta;
119  std::map<unsigned int, int> m_pdg;
120  std::map<unsigned int, double> m_pz;
121  std::map<unsigned int, bool> m_isCharged;
122 
125 
129 
130  std::vector<LHCTransportLink> theCorrespondenceMap;
131 };
132 #endif
std::map< unsigned int, int > m_direct
Definition: Hector.h:108
std::map< unsigned int, bool > m_isCharged
Definition: Hector.h:121
ZDCTransport
HepMC source to be processed.
float m_rpp420_f
Definition: Hector.h:88
bool m_ZDCTransport
Definition: Hector.h:128
std::map< unsigned int, double > m_eta
Definition: Hector.h:118
float m_rpp420_b
Definition: Hector.h:89
H_RecRPObject * m_rp420_f
Definition: Hector.h:104
std::map< unsigned int, bool > m_isStoppedzdc
Definition: Hector.h:110
bool m_verbosity
Definition: Hector.h:126
double m_sig_e
Definition: Hector.h:83
H_BeamLine * m_beamlineZDC1
Definition: Hector.h:98
edm::ESHandle< ParticleDataTable > pdt
Definition: Hector.h:93
H_BeamLine * m_beamlineD11
Definition: Hector.h:100
std::map< unsigned int, bool > m_isStoppedd1
Definition: Hector.h:111
std::map< unsigned int, bool > m_isStoppedfp420
Definition: Hector.h:109
Definition: Hector.h:31
void clear()
Definition: Hector.cc:155
double m_sigmaSTX
Definition: Hector.h:85
FP420Transport
main flag to set transport for ZDC
std::map< unsigned int, double > m_eAtTrPoint
Definition: Hector.h:116
float m_rppd1
Definition: Hector.h:91
H_BeamLine * m_beamlineZDC2
Definition: Hector.h:99
void filterD1(TRandom3 *)
Definition: Hector.cc:405
std::map< unsigned int, double > m_pz
Definition: Hector.h:120
double lengthfp420
Definition: Hector.h:77
string beam2filename
Definition: Hector.h:124
std::map< unsigned int, double > m_TyAtTrPoint
Definition: Hector.h:115
void filterFP420(TRandom3 *)
Definition: Hector.cc:230
std::vector< LHCTransportLink > theCorrespondenceMap
Definition: Hector.h:130
Hector(const edm::ParameterSet &ps, const edm::ESGetToken< HepPDT::ParticleDataTable, PDTRecord > &, bool verbosity, bool FP420Transport, bool ZDCTransport)
Definition: Hector.cc:17
double m_sigmaSTY
Definition: Hector.h:86
void print() const
Definition: Hector.cc:500
const int verbosity
virtual ~Hector()
Definition: Hector.cc:136
int getDirect(unsigned int part_n) const
Definition: Hector.cc:492
void filterZDC(TRandom3 *)
Definition: Hector.cc:322
H_BeamLine * m_beamlineFP4201
Definition: Hector.h:96
bool m_smearE
Definition: Hector.h:84
double etacut
Definition: Hector.h:81
string beam1filename
Definition: Hector.h:123
const edm::ESGetToken< HepPDT::ParticleDataTable, PDTRecord > tok_pdt_
Definition: Hector.h:75
void add(const HepMC::GenEvent *ev, const edm::EventSetup &es)
Definition: Hector.cc:167
bool m_smearAng
Definition: Hector.h:82
std::vector< LHCTransportLink > & getCorrespondenceMap()
Definition: Hector.h:68
H_BeamLine * m_beamlineD12
Definition: Hector.h:101
H_BeamLine * m_beamlineFP4202
Definition: Hector.h:97
float m_rppzdc
Definition: Hector.h:90
HepMC::GenEvent * addPartToHepMC(HepMC::GenEvent *event)
Definition: Hector.cc:506
double lengthzdc
Definition: Hector.h:78
std::map< unsigned int, double > m_yAtTrPoint
Definition: Hector.h:113
std::map< unsigned int, int > m_pdg
Definition: Hector.h:119
std::map< unsigned int, double > m_TxAtTrPoint
Definition: Hector.h:114
std::map< unsigned int, double > m_xAtTrPoint
Definition: Hector.h:112
double lengthd1
Definition: Hector.h:79
std::map< unsigned int, H_BeamParticle * > m_beamPart
Definition: Hector.h:107
H_RecRPObject * m_rp420_b
Definition: Hector.h:105
Definition: event.py:1
bool m_FP420Transport
Definition: Hector.h:127
void clearApertureFlags()
Definition: Hector.cc:149