CMS 3D CMS Logo

CTPPSFastTrackingProducer.h
Go to the documentation of this file.
1 #ifndef FastSimulation_CTPPSFastTrackingProducer_h
2 #define FastSimulation_CTPPSFastTrackingProducer_h
3 
4 // user include files
7 
10 
13 
15 
19 
22 
25 
28 //#include "FastSimulation/CTPPSSimHitProducer/plugins/FastCTPPSParameters.h"
29 
30 //CLHEP
31 #include "CLHEP/Units/GlobalSystemOfUnits.h"
32 #include "CLHEP/Units/GlobalPhysicalConstants.h"
33 #include <CLHEP/Vector/LorentzVector.h>
34 
35 //C++ library
36 #include <iostream>
37 #include <fstream>
38 #include <string>
39 #include <cmath>
40 #include <math.h>
41 #include <map>
42 #include <vector>
43 #include <utility>
44 #include <math.h>
45 
46 #include <TMatrixD.h>
47 
48 // hector includes
49 #include "H_Parameters.h"
50 #include "H_BeamLine.h"
51 #include "H_RecRPObject.h"
52 #include "H_BeamParticle.h"
53 
57 
58 
59 //
60 // class declaration
61 //
62 
64  public:
67  typedef CLHEP::HepLorentzVector LorentzVector;
68 
69  private:
70  virtual void beginStream(edm::StreamID) override;
71  virtual void produce(edm::Event&, const edm::EventSetup&) override;
72  virtual void endStream() override;
73  //this function will only be called once per event
74  virtual void beginEvent(edm::Event& event, const edm::EventSetup& eventSetup);
75  virtual void endEvent(edm::Event& event, const edm::EventSetup& eventSetup);
76 
77  // ----------member data ---------------------------
78 
79  typedef std::vector<CTPPSFastRecHit> CTPPSFastRecHitContainer;
82  void FastReco(int Direction,H_RecRPObject* station);
83  void Reconstruction();
84  void ReconstructArm(H_RecRPObject* pps_station, double x1,double y1,double x2,double y2, double& tx, double& ty,double& eloss);
85  void MatchCellId(int cellId, vector<int> vrecCellId, vector<double> vrecTof, bool& match, double& recTof);
86  bool SearchTrack(int ,int ,int Direction,double& xi,double& t,double& partP,double& pt,double& thx,double& thy,double& x0,double& y0, double& xt, double& yt, double& X1d, double& Y1d, double& X2d, double& Y2d);
87  void LorentzBoost(LorentzVector& p_out, const string& frame);
88  void Get_t_and_xi(const LorentzVector* p,double& t, double& xi);
89  void TrackerStationClear();
91  void ProjectToToF(const double x1, const double y1, const double x2, const double y2, double& xt, double& yt) {
92  xt = ((fz_timing-fz_tracker2)*(x2-x1)/(fz_tracker2-fz_tracker1)) + x2;
93  yt = ((fz_timing-fz_tracker2)*(y2-y1)/(fz_tracker2-fz_tracker1)) + y2;
94  };
95  // Hector objects
96 
97  std::map<unsigned int, H_BeamParticle*> m_beamPart;
98  std::unique_ptr<H_BeamLine> m_beamlineCTPPS1;
99  std::unique_ptr<H_BeamLine> m_beamlineCTPPS2;
100 
101  std::unique_ptr<H_RecRPObject> pps_stationF;
102  std::unique_ptr<H_RecRPObject> pps_stationB;
103 
106 
107  // Defaults
108  double lengthctpps ;
110  double fBeamEnergy;
115  std::unique_ptr<CTPPSTrkStation> TrkStation_F; // auxiliary object with the tracker geometry
116  std::unique_ptr<CTPPSTrkStation> TrkStation_B;
117  std::unique_ptr<CTPPSTrkDetector> det1F;
118  std::unique_ptr<CTPPSTrkDetector> det1B;
119  std::unique_ptr<CTPPSTrkDetector> det2F;
120  std::unique_ptr<CTPPSTrkDetector> det2B;
121  std::unique_ptr<CTPPSToFDetector> detToF_F;
122  std::unique_ptr<CTPPSToFDetector> detToF_B;
123 
124  std::vector<CTPPSFastTrack> theCTPPSFastTrack;
125 
127 
128  std::vector<int> recCellId_F, recCellId_B ;
129  std::vector<double> recTof_F, recTof_B ;
130 
133  std::vector<double> fToFCellWidth;
137 
138 
139 };
140 #endif
141 
virtual void beginStream(edm::StreamID) override
void MatchCellId(int cellId, vector< int > vrecCellId, vector< double > vrecTof, bool &match, double &recTof)
CLHEP::HepLorentzVector LorentzVector
std::unique_ptr< CTPPSTrkDetector > det2B
std::unique_ptr< H_BeamLine > m_beamlineCTPPS1
virtual void endEvent(edm::Event &event, const edm::EventSetup &eventSetup)
void FastReco(int Direction, H_RecRPObject *station)
void ReconstructArm(H_RecRPObject *pps_station, double x1, double y1, double x2, double y2, double &tx, double &ty, double &eloss)
std::unique_ptr< H_RecRPObject > pps_stationB
std::unique_ptr< CTPPSToFDetector > detToF_F
std::unique_ptr< CTPPSTrkDetector > det2F
std::unique_ptr< H_RecRPObject > pps_stationF
CTPPSFastTrackingProducer(const edm::ParameterSet &)
std::unique_ptr< CTPPSTrkDetector > det1B
virtual void beginEvent(edm::Event &event, const edm::EventSetup &eventSetup)
std::unique_ptr< CTPPSTrkStation > TrkStation_B
std::unique_ptr< CTPPSToFDetector > detToF_B
std::map< unsigned int, H_BeamParticle * > m_beamPart
std::vector< CTPPSFastRecHit > CTPPSFastRecHitContainer
std::unique_ptr< H_BeamLine > m_beamlineCTPPS2
void ProjectToToF(const double x1, const double y1, const double x2, const double y2, double &xt, double &yt)
void Get_t_and_xi(const LorentzVector *p, double &t, double &xi)
std::unique_ptr< CTPPSTrkStation > TrkStation_F
std::vector< CTPPSFastTrack > theCTPPSFastTrack
edm::EDGetTokenT< CTPPSFastRecHitContainer > _recHitToken
std::unique_ptr< CTPPSTrkDetector > det1F
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
void LorentzBoost(LorentzVector &p_out, const string &frame)
virtual void produce(edm::Event &, const edm::EventSetup &) override
void ReadRecHits(edm::Handle< CTPPSFastRecHitContainer > &)
Definition: event.py:1
bool SearchTrack(int, int, int Direction, double &xi, double &t, double &partP, double &pt, double &thx, double &thy, double &x0, double &y0, double &xt, double &yt, double &X1d, double &Y1d, double &X2d, double &Y2d)