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/GenVertex.h"
13 #include "HepMC/GenParticle.h"
15 
17 #include <vector>
18 
19 // SimpleConfigurable replacement
21 
22 //Hector headers
23 #include "H_BeamLine.h"
24 #include "H_RecRPObject.h"
25 #include "H_BeamParticle.h"
26 #include <string>
27 #include <map>
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();
45  void clear();
47  void add(const HepMC::GenEvent *ev, const edm::EventSetup &es);
49  void filterFP420(TRandom3 *);
51  void filterZDC(TRandom3 *);
53  void filterD1(TRandom3 *);
54 
55  int getDirect(unsigned int part_n) const;
56 
58  void print() const;
60  // std::vector<unsigned int> part_list() const;
61 
62  // bool isCharged(const HepMC::GenParticle * p);
63 
65 
66  std::vector<LHCTransportLink> &getCorrespondenceMap() { return theCorrespondenceMap; }
67 
68 private:
70 
71  // Defaults
72  double lengthfp420;
73  double lengthzdc;
74  double lengthd1;
75 
76  double etacut;
77  bool m_smearAng;
78  double m_sig_e;
79  bool m_smearE;
80  double m_sigmaSTX;
81  double m_sigmaSTY;
82 
83  float m_rpp420_f;
84  float m_rpp420_b;
85  float m_rppzdc;
86  float m_rppd1;
87 
89 
90  // Hector
91  H_BeamLine *m_beamlineFP4201;
92  H_BeamLine *m_beamlineFP4202;
93  H_BeamLine *m_beamlineZDC1;
94  H_BeamLine *m_beamlineZDC2;
95  H_BeamLine *m_beamlineD11;
96  H_BeamLine *m_beamlineD12;
97  //
98 
99  H_RecRPObject *m_rp420_f;
100  H_RecRPObject *m_rp420_b;
101 
102  std::map<unsigned int, H_BeamParticle *> m_beamPart;
103  std::map<unsigned int, int> m_direct;
104  std::map<unsigned int, bool> m_isStoppedfp420;
105  std::map<unsigned int, bool> m_isStoppedzdc;
106  std::map<unsigned int, bool> m_isStoppedd1;
107  std::map<unsigned int, double> m_xAtTrPoint;
108  std::map<unsigned int, double> m_yAtTrPoint;
109  std::map<unsigned int, double> m_TxAtTrPoint;
110  std::map<unsigned int, double> m_TyAtTrPoint;
111  std::map<unsigned int, double> m_eAtTrPoint;
112 
113  std::map<unsigned int, double> m_eta;
114  std::map<unsigned int, int> m_pdg;
115  std::map<unsigned int, double> m_pz;
116  std::map<unsigned int, bool> m_isCharged;
117 
120 
124 
125  std::vector<LHCTransportLink> theCorrespondenceMap;
126 };
127 #endif
std::map< unsigned int, int > m_direct
Definition: Hector.h:103
std::map< unsigned int, bool > m_isCharged
Definition: Hector.h:116
ZDCTransport
HepMC source to be processed.
float m_rpp420_f
Definition: Hector.h:83
bool m_ZDCTransport
Definition: Hector.h:123
std::map< unsigned int, double > m_eta
Definition: Hector.h:113
float m_rpp420_b
Definition: Hector.h:84
H_RecRPObject * m_rp420_f
Definition: Hector.h:99
std::map< unsigned int, bool > m_isStoppedzdc
Definition: Hector.h:105
bool m_verbosity
Definition: Hector.h:121
double m_sig_e
Definition: Hector.h:78
H_BeamLine * m_beamlineZDC1
Definition: Hector.h:93
edm::ESHandle< ParticleDataTable > pdt
Definition: Hector.h:88
H_BeamLine * m_beamlineD11
Definition: Hector.h:95
std::map< unsigned int, bool > m_isStoppedd1
Definition: Hector.h:106
std::map< unsigned int, bool > m_isStoppedfp420
Definition: Hector.h:104
Definition: Hector.h:31
void clear()
Definition: Hector.cc:137
double m_sigmaSTX
Definition: Hector.h:80
FP420Transport
main flag to set transport for ZDC
std::map< unsigned int, double > m_eAtTrPoint
Definition: Hector.h:111
std::string beam2filename
Definition: Hector.h:119
float m_rppd1
Definition: Hector.h:86
H_BeamLine * m_beamlineZDC2
Definition: Hector.h:94
void filterD1(TRandom3 *)
Definition: Hector.cc:378
std::map< unsigned int, double > m_pz
Definition: Hector.h:115
double lengthfp420
Definition: Hector.h:72
std::map< unsigned int, double > m_TyAtTrPoint
Definition: Hector.h:110
void filterFP420(TRandom3 *)
Definition: Hector.cc:211
std::vector< LHCTransportLink > theCorrespondenceMap
Definition: Hector.h:125
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:81
void print() const
Definition: Hector.cc:471
const int verbosity
virtual ~Hector()
Definition: Hector.cc:118
int getDirect(unsigned int part_n) const
Definition: Hector.cc:463
void filterZDC(TRandom3 *)
Definition: Hector.cc:302
H_BeamLine * m_beamlineFP4201
Definition: Hector.h:91
bool m_smearE
Definition: Hector.h:79
double etacut
Definition: Hector.h:76
const edm::ESGetToken< HepPDT::ParticleDataTable, PDTRecord > tok_pdt_
Definition: Hector.h:69
std::string beam1filename
Definition: Hector.h:118
void add(const HepMC::GenEvent *ev, const edm::EventSetup &es)
Definition: Hector.cc:149
bool m_smearAng
Definition: Hector.h:77
std::vector< LHCTransportLink > & getCorrespondenceMap()
Definition: Hector.h:66
H_BeamLine * m_beamlineD12
Definition: Hector.h:96
H_BeamLine * m_beamlineFP4202
Definition: Hector.h:92
float m_rppzdc
Definition: Hector.h:85
HepMC::GenEvent * addPartToHepMC(HepMC::GenEvent *event)
Definition: Hector.cc:477
double lengthzdc
Definition: Hector.h:73
std::map< unsigned int, double > m_yAtTrPoint
Definition: Hector.h:108
std::map< unsigned int, int > m_pdg
Definition: Hector.h:114
std::map< unsigned int, double > m_TxAtTrPoint
Definition: Hector.h:109
std::map< unsigned int, double > m_xAtTrPoint
Definition: Hector.h:107
double lengthd1
Definition: Hector.h:74
std::map< unsigned int, H_BeamParticle * > m_beamPart
Definition: Hector.h:102
H_RecRPObject * m_rp420_b
Definition: Hector.h:100
Definition: event.py:1
bool m_FP420Transport
Definition: Hector.h:122
void clearApertureFlags()
Definition: Hector.cc:131