CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 #include "FWCore/Framework/interface/EDAnalyzer.h"
9 #include "FWCore/Framework/interface/Frameworkfwd.h"
10 
11 #include "FWCore/Framework/interface/Event.h"
12 #include "FWCore/Framework/interface/MakerMacros.h"
13 
14 #include "FWCore/ParameterSet/interface/ParameterSet.h"
15 #include "FWCore/Utilities/interface/InputTag.h"
16 */
17 // HepMC headers
19 
20 #include "HepMC/GenEvent.h"
21 #include "HepMC/GenParticle.h"
22 #include "HepMC/GenVertex.h"
24 
26 #include <vector>
27 
28 // SimpleConfigurable replacement
30 
31 // Hector headers
32 #include "H_BeamLine.h"
33 #include "H_BeamParticle.h"
34 #include "H_RecRPObject.h"
35 #include <map>
36 #include <string>
37 
38 class TRandom3;
39 
40 class Hector {
41 public:
42  // Hector(const edm::ParameterSet & ps);
43  Hector(const edm::ParameterSet &ps,
45  bool verbosity,
46  bool FP420Transport,
47  bool ZDCTransport);
48  // Hector();
49  virtual ~Hector();
50 
52  void clearApertureFlags();
55  void clear();
57  void add(const HepMC::GenEvent *ev, const edm::EventSetup &es);
59  void filterFP420(TRandom3 *);
61  void filterZDC(TRandom3 *);
63  void filterD1(TRandom3 *);
64 
65  int getDirect(unsigned int part_n) const;
66 
68  void print() const;
71  // std::vector<unsigned int> part_list() const;
72 
73  // bool isCharged(const HepMC::GenParticle * p);
74 
76 
77  std::vector<LHCTransportLink> &getCorrespondenceMap() { return theCorrespondenceMap; }
78 
79  /*
80  private:
81  // edm::ParameterSet m_pBeamLine;
82  */
83 private:
85  // Defaults
86  double lengthfp420;
87  double lengthzdc;
88  double lengthd1;
89 
90  double etacut;
91  bool m_smearAng;
92  double m_sig_e;
93  bool m_smearE;
94  double m_sigmaSTX;
95  double m_sigmaSTY;
96 
97  float m_rpp420_f;
98  float m_rpp420_b;
99  float m_rppzdc;
100  float m_rppd1;
101 
103 
104  // Hector
105  H_BeamLine *m_beamlineFP4201;
106  H_BeamLine *m_beamlineFP4202;
107  H_BeamLine *m_beamlineZDC1;
108  H_BeamLine *m_beamlineZDC2;
109  H_BeamLine *m_beamlineD11;
110  H_BeamLine *m_beamlineD12;
111  //
112 
113  H_RecRPObject *m_rp420_f;
114  H_RecRPObject *m_rp420_b;
115 
116  std::map<unsigned int, H_BeamParticle *> m_beamPart;
117  std::map<unsigned int, int> m_direct;
118  std::map<unsigned int, bool> m_isStoppedfp420;
119  std::map<unsigned int, bool> m_isStoppedzdc;
120  std::map<unsigned int, bool> m_isStoppedd1;
121  std::map<unsigned int, double> m_xAtTrPoint;
122  std::map<unsigned int, double> m_yAtTrPoint;
123  std::map<unsigned int, double> m_TxAtTrPoint;
124  std::map<unsigned int, double> m_TyAtTrPoint;
125  std::map<unsigned int, double> m_eAtTrPoint;
126 
127  std::map<unsigned int, double> m_eta;
128  std::map<unsigned int, int> m_pdg;
129  std::map<unsigned int, double> m_pz;
130  std::map<unsigned int, bool> m_isCharged;
131 
134 
138 
139  std::vector<LHCTransportLink> theCorrespondenceMap;
140 };
141 #endif
std::map< unsigned int, int > m_direct
Definition: Hector.h:117
std::map< unsigned int, bool > m_isCharged
Definition: Hector.h:130
int getDirect(unsigned int part_n) const
Definition: Hector.cc:492
float m_rpp420_f
Definition: Hector.h:97
bool m_ZDCTransport
Definition: Hector.h:137
std::map< unsigned int, double > m_eta
Definition: Hector.h:127
float m_rpp420_b
Definition: Hector.h:98
H_RecRPObject * m_rp420_f
Definition: Hector.h:113
std::map< unsigned int, bool > m_isStoppedzdc
Definition: Hector.h:119
bool m_verbosity
Definition: Hector.h:135
bool ev
double m_sig_e
Definition: Hector.h:92
H_BeamLine * m_beamlineZDC1
Definition: Hector.h:107
edm::ESHandle< ParticleDataTable > pdt
Definition: Hector.h:102
H_BeamLine * m_beamlineD11
Definition: Hector.h:109
std::map< unsigned int, bool > m_isStoppedd1
Definition: Hector.h:120
std::map< unsigned int, bool > m_isStoppedfp420
Definition: Hector.h:118
Definition: Hector.h:40
void clear()
Definition: Hector.cc:155
double m_sigmaSTX
Definition: Hector.h:94
std::map< unsigned int, double > m_eAtTrPoint
Definition: Hector.h:125
float m_rppd1
Definition: Hector.h:100
H_BeamLine * m_beamlineZDC2
Definition: Hector.h:108
void filterD1(TRandom3 *)
Definition: Hector.cc:405
std::map< unsigned int, double > m_pz
Definition: Hector.h:129
double lengthfp420
Definition: Hector.h:86
string beam2filename
Definition: Hector.h:133
std::map< unsigned int, double > m_TyAtTrPoint
Definition: Hector.h:124
void filterFP420(TRandom3 *)
Definition: Hector.cc:230
std::vector< LHCTransportLink > theCorrespondenceMap
Definition: Hector.h:139
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:95
void print() const
Definition: Hector.cc:500
virtual ~Hector()
Definition: Hector.cc:136
void filterZDC(TRandom3 *)
Definition: Hector.cc:322
H_BeamLine * m_beamlineFP4201
Definition: Hector.h:105
bool m_smearE
Definition: Hector.h:93
double etacut
Definition: Hector.h:90
string beam1filename
Definition: Hector.h:132
const edm::ESGetToken< HepPDT::ParticleDataTable, PDTRecord > tok_pdt_
Definition: Hector.h:84
void add(const HepMC::GenEvent *ev, const edm::EventSetup &es)
Definition: Hector.cc:167
bool m_smearAng
Definition: Hector.h:91
std::vector< LHCTransportLink > & getCorrespondenceMap()
Definition: Hector.h:77
H_BeamLine * m_beamlineD12
Definition: Hector.h:110
H_BeamLine * m_beamlineFP4202
Definition: Hector.h:106
float m_rppzdc
Definition: Hector.h:99
HepMC::GenEvent * addPartToHepMC(HepMC::GenEvent *event)
Definition: Hector.cc:506
double lengthzdc
Definition: Hector.h:87
std::map< unsigned int, double > m_yAtTrPoint
Definition: Hector.h:122
std::map< unsigned int, int > m_pdg
Definition: Hector.h:128
std::map< unsigned int, double > m_TxAtTrPoint
Definition: Hector.h:123
std::map< unsigned int, double > m_xAtTrPoint
Definition: Hector.h:121
double lengthd1
Definition: Hector.h:88
std::map< unsigned int, H_BeamParticle * > m_beamPart
Definition: Hector.h:116
H_RecRPObject * m_rp420_b
Definition: Hector.h:114
bool m_FP420Transport
Definition: Hector.h:136
void clearApertureFlags()
Definition: Hector.cc:149