CMS 3D CMS Logo

CTPPSHector.h
Go to the documentation of this file.
1 #ifndef SimTransport_CTPPSHector_h
2 #define SimTransport_CTPPSHector_h
3 
4 // user include files
7 // HepMC headers
9 
10 #include "HepMC/GenEvent.h"
11 #include "HepMC/GenVertex.h"
12 #include "HepMC/GenParticle.h"
14 
16 #include <vector>
17 
18 // SimpleConfigurable replacement
20 
21 //CTPPSHector headers
22 #include "H_BeamLine.h"
23 #include "H_BeamParticle.h"
24 
25 #include <string>
26 #include <map>
27 #include <cmath>
28 #include <math.h>
29 
31 
35 #include <CLHEP/Vector/LorentzVector.h>
36 
37 class TRandom3;
38 
39 class CTPPSHector {
40 
41  public:
42  // CTPPSHector(const edm::ParameterSet & ps);
44  // CTPPSHector();
45  virtual ~CTPPSHector();
46  typedef CLHEP::HepLorentzVector LorentzVector;
48  void clearApertureFlags();
50  void clear();
51 
53  void add( const HepMC::GenEvent * ev , const edm::EventSetup & es, CLHEP::HepRandomEngine * engine);
54 
56  void filterCTPPS(TRandom3*);
57 
58  // New function to calculate the LorentzBoost
59  void LorentzBoost(LorentzVector& p_out, const string& frame);
60 
62 
63  double get_BeamEnergy() {return fBeamEnergy;};
64 
65  double get_BeamMomentum() {return fBeamMomentum;};
66 
67  void ApplyBeamCorrection(LorentzVector&, CLHEP::HepRandomEngine * engine);
68 
69  int getDirect( unsigned int part_n ) const;
70 
72  void print() const;
73 
74  HepMC::GenEvent * addPartToHepMC( HepMC::GenEvent * event );
75 
76  std::vector<LHCTransportLink> & getCorrespondenceMap() { return theCorrespondenceMap; }
77 
78  private:
79  // Defaults
80  double lengthctpps ;
81 
82  double etacut;
83  bool m_smearAng;
84  double m_sig_e;
85  bool m_smearE;
86  double m_sigmaSTX;
87  double m_sigmaSTY;
88  float m_f_ctpps_f;
89  float m_b_ctpps_b;
90 
91  //HECTOR CTPPS Parameters
94  double fBeamMomentum;
95  double fBeamEnergy;
96  double fVtxMeanX;
97  double fVtxMeanY;
98  double fVtxMeanZ;
99  double fMomentumMin;
100 
102 
103  // CTPPSHector
104  std::unique_ptr<H_BeamLine> m_beamlineCTPPS1;
105  std::unique_ptr<H_BeamLine> m_beamlineCTPPS2;
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_isStoppedctpps;
110  std::map<unsigned int, double> m_xAtTrPoint;
111  std::map<unsigned int, double> m_yAtTrPoint;
112  std::map<unsigned int, double> m_TxAtTrPoint;
113  std::map<unsigned int, double> m_TyAtTrPoint;
114  std::map<unsigned int, double> m_eAtTrPoint;
115 
116  std::map<unsigned int, double> m_eta;
117  std::map<unsigned int, int> m_pdg;
118  std::map<unsigned int, double> m_pz;
119  std::map<unsigned int, bool> m_isCharged;
120 
123 
126 
127  std::vector<LHCTransportLink> theCorrespondenceMap;
128  int NEvent;
129 };
130 #endif
HepMC::GenEvent * addPartToHepMC(HepMC::GenEvent *event)
Definition: CTPPSHector.cc:345
double get_BeamMomentum()
Definition: CTPPSHector.h:65
CTPPSTransport
HepMC source to be processed.
void clearApertureFlags()
Definition: CTPPSHector.cc:107
double fBeamMomentum
Definition: CTPPSHector.h:94
bool m_verbosity
Definition: CTPPSHector.h:124
std::map< unsigned int, double > m_eAtTrPoint
Definition: CTPPSHector.h:114
double get_BeamEnergy()
Definition: CTPPSHector.h:63
double m_sigmaSTY
Definition: CTPPSHector.h:87
bool m_smearAng
Definition: CTPPSHector.h:83
std::map< unsigned int, double > m_yAtTrPoint
Definition: CTPPSHector.h:111
std::map< unsigned int, int > m_pdg
Definition: CTPPSHector.h:117
void print() const
Definition: CTPPSHector.cc:275
bool ev
string beam1filename
Definition: CTPPSHector.h:121
std::map< unsigned int, double > m_pz
Definition: CTPPSHector.h:118
std::map< unsigned int, int > m_direct
Definition: CTPPSHector.h:108
std::map< unsigned int, H_BeamParticle * > m_beamPart
Definition: CTPPSHector.h:107
bool m_CTPPSTransport
Definition: CTPPSHector.h:125
double etacut
Definition: CTPPSHector.h:82
edm::ESHandle< ParticleDataTable > pdt
Definition: CTPPSHector.h:101
double fVtxMeanX
Definition: CTPPSHector.h:96
std::map< unsigned int, bool > m_isStoppedctpps
Definition: CTPPSHector.h:109
T sqrt(T t)
Definition: SSEVec.h:18
void add(const HepMC::GenEvent *ev, const edm::EventSetup &es, CLHEP::HepRandomEngine *engine)
Definition: CTPPSHector.cc:123
double fVtxMeanY
Definition: CTPPSHector.h:97
double fVtxMeanZ
Definition: CTPPSHector.h:98
double lengthctpps
Definition: CTPPSHector.h:80
float m_f_ctpps_f
Definition: CTPPSHector.h:88
virtual ~CTPPSHector()
Definition: CTPPSHector.cc:99
double fMomentumMin
Definition: CTPPSHector.h:99
bool m_smearE
Definition: CTPPSHector.h:85
const double ProtonMassSQ
std::map< unsigned int, double > m_eta
Definition: CTPPSHector.h:116
double fCrossingAngle
Definition: CTPPSHector.h:93
std::vector< LHCTransportLink > theCorrespondenceMap
Definition: CTPPSHector.h:127
void clear()
Definition: CTPPSHector.cc:111
std::map< unsigned int, bool > m_isCharged
Definition: CTPPSHector.h:119
int getDirect(unsigned int part_n) const
Definition: CTPPSHector.cc:267
float m_b_ctpps_b
Definition: CTPPSHector.h:89
double m_sigmaSTX
Definition: CTPPSHector.h:86
double fBeamEnergy
Definition: CTPPSHector.h:95
std::map< unsigned int, double > m_TxAtTrPoint
Definition: CTPPSHector.h:112
std::map< unsigned int, double > m_xAtTrPoint
Definition: CTPPSHector.h:110
std::unique_ptr< H_BeamLine > m_beamlineCTPPS1
Definition: CTPPSHector.h:104
std::map< unsigned int, double > m_TyAtTrPoint
Definition: CTPPSHector.h:113
CTPPSHector(const edm::ParameterSet &ps, bool verbosity, bool CTPPSTransport)
Definition: CTPPSHector.cc:20
std::unique_ptr< H_BeamLine > m_beamlineCTPPS2
Definition: CTPPSHector.h:105
void filterCTPPS(TRandom3 *)
Definition: CTPPSHector.cc:197
bool fCrossAngleCorr
Definition: CTPPSHector.h:92
void set_BeamEnergy(double e)
Definition: CTPPSHector.h:61
string beam2filename
Definition: CTPPSHector.h:122
double m_sig_e
Definition: CTPPSHector.h:84
CLHEP::HepLorentzVector LorentzVector
Definition: CTPPSHector.h:46
void LorentzBoost(LorentzVector &p_out, const string &frame)
Definition: CTPPSHector.cc:305
std::vector< LHCTransportLink > & getCorrespondenceMap()
Definition: CTPPSHector.h:76
void ApplyBeamCorrection(LorentzVector &, CLHEP::HepRandomEngine *engine)
Definition: CTPPSHector.cc:281
Definition: event.py:1