CMS 3D CMS Logo

CTPPSSimHitProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FastSimulation/CTPPSSimHitProducer
4 // Class: CTPPSSimHitProducer
5 //
13 //
14 // Original Author: Dilson De Jesus Damiao
15 // Created: Mon, 05 Sep 2016 18:49:10 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
25 
28 
31 
32 // SimHitContainer
34 
35 // STL headers
36 #include <vector>
37 #include <iostream>
38 
39 // HepMC headers
40 #include "HepMC/GenEvent.h"
41 #include "HepMC/GenVertex.h"
42 #include "HepMC/GenParticle.h"
45 
47 
50 
52 #include <CLHEP/Vector/LorentzVector.h>
53 
54 //
55 // class declaration
56 //
57 
59  public:
60  explicit CTPPSSimHitProducer(const edm::ParameterSet&);
62  typedef CLHEP::HepLorentzVector LorentzVector;
63 
64  private:
65  virtual void beginStream(edm::StreamID) override;
66  virtual void produce(edm::Event&, const edm::EventSetup&) override;
67  virtual void endStream() override;
68 
71  // ----------member data ---------------------------
73 
74 
75 
76 };
77 
78 //
79 // constants, enums and typedefs
80 //
81 
82 //
83 // static data member definitions
84 //
85 
86 //
87 // constructors and destructor
88 //
90 {
91 
92  produces<edm::PSimHitContainer>("CTPPSHits");
93  // consumes
94  mcEventToken = mayConsume<edm::HepMCProduct>(iConfig.getUntrackedParameter<edm::InputTag>("MCEvent",std::string("")));
95 
96  // Read the position of the trackers and timing detectors
97  fz_tracker1 = iConfig.getParameter<double>("Z_Tracker1");
98  fz_tracker2 = iConfig.getParameter<double>("Z_Tracker2");
99  fz_timing = iConfig.getParameter<double>("Z_Timing");
100 }
101 
103 {
104 
105 }
106 
107 
108 //
109 // member functions
110 //
111 
112 // ------------ method called to produce the data ------------
113  void
115 {
116  using namespace edm;
117 
118  std::vector<PSimHit> theCTPPSHits;
119  iEvent.getByToken( mcEventToken , EvtHandle );
120 
121  const HepMC::GenEvent* Evt = EvtHandle->GetEvent() ;
122  std::vector<math::XYZTLorentzVector> protonCTPPS;
123  protonCTPPS.clear();
124  for(HepMC::GenEvent::vertex_const_iterator ivtx = Evt->vertices_begin();ivtx!=Evt->vertices_end();ivtx++) {
125  if ((*ivtx)->id()!=0) continue;
126 
127  // Get the vertices at the entrance of CTPPS and get the protons coming out of them (propagated by Hector)
128  for(HepMC::GenVertex::particles_out_const_iterator i=(*ivtx)->particles_out_const_begin(); i != (*ivtx)->particles_out_const_end();++i) {
129  int pid = (*i)->pdg_id();
130  if(pid!=2212) continue;
131 
132  HepMC::GenVertex* pv = (*i)->production_vertex();
133  HepMC::FourVector vertex = pv->position();
134  const HepMC::FourVector p((*i)->momentum());
135  protonCTPPS.push_back(math::XYZTLorentzVector(p.x(),p.y(),p.z(),p.t()));
136 
137  LocalPoint initialPosition_tr1, initialPosition_tr2, initialPosition_tof;
138 
139  double t0_tr1 = 0., t0_tr2 = 0., t0_tof = 0.;
140  // t0_* has dimensions of mm
141  // Convert to ns for internal calculations.
142  const double c_light_s = 2.99792458e+11;// mm/s;
143  const double s_to_ns = 1.e9;
144  const double m_to_mm = 1.e3;
145  double x_tr1 = 0. , x_tr2 = 0. , x_tof = 0., y_tr1 = 0. , y_tr2 = 0. , y_tof = 0., z_tr1 = 0.;
146  double z_tr2 = fz_tracker2; //m
147  double z_tof = fz_timing; //m
148  int Direction =0;
149 
150  // Read the vertex made by SimTransport at the entrance of det1
151  if(std::abs(vertex.eta())>8. && (*i)->status() == 1){
152  if (vertex.z()>0) Direction = 1;
153  else if (vertex.z() < 0 ) Direction = -1;
154 
155  //Get the global coordinates at Tracker1, equal to those of the vertex at CTPPS
156  x_tr1 = vertex.x();
157  y_tr1 = vertex.y();
158  z_tr1 = vertex.z()*mm_to_m;//move z from mm to meters
159  Local3DPoint xyzzy_tr1(x_tr1, y_tr1, z_tr1);
160  initialPosition_tr1 = xyzzy_tr1;
161  t0_tr1 = vertex.t()/c_light_s*s_to_ns;
162  //Get the global coordinates at Tracker2 by propagating as a straight line from Tracker1
163  t0_tr2 = z_tr2*m_to_mm/c_light_s*s_to_ns;//discuss latter if needs to be corrected with vertex position
164  z_tr2 *= Direction;
165  x_tr2 = x_tr1 + (p.x()/p.z())*(z_tr2-z_tr1)*m_to_mm;
166  y_tr2 = y_tr1 + (p.y()/p.z())*(z_tr2-z_tr1)*m_to_mm;
167  Local3DPoint xyzzy_tr2(x_tr2, y_tr2, z_tr2);
168  initialPosition_tr2 = xyzzy_tr2;
169  //Propagate as a straight line from Tracker1
170  t0_tof = z_tof*m_to_mm/c_light_s*s_to_ns;//discuss latter if needs to be corrected with vertex position
171  z_tof *= Direction;
172  x_tof = x_tr1 + (p.x()/p.z())*(z_tof-z_tr1)*m_to_mm;
173  y_tof = y_tr1 + (p.y()/p.z())*(z_tof-z_tr1)*m_to_mm;
174  Local3DPoint xyzzy_tof(x_tof, y_tof, z_tof);
175  initialPosition_tof = xyzzy_tof;
176 
177  // DetId codification for PSimHit from CTPPSPixel- It will be replaced by CTPPSDetId
178  // 2014314496 -> Tracker1 zPositive
179  // 2014838784 -> Tracker2 zPositive
180  // 2046820352 -> Timing zPositive
181  // 2031091712 -> Tracker1 zNegative
182  // 2031616000 -> Tracker2 zNegative
183  // 2063597568 -> Timing zNegative
184 
185  if(Direction > 0.) {
186  PSimHit hit_tr1(xyzzy_tr1,xyzzy_tr1,0.,t0_tr1,0.,pid,2014314496,0,0.,0.,2);
187  PSimHit hit_tr2(xyzzy_tr2,xyzzy_tr2,0.,t0_tr2,0.,pid,2014838784,0,0.,0.,2);
188  PSimHit hit_tof(xyzzy_tof,xyzzy_tof,0.,t0_tof,0.,pid,2046820352,0,0.,0.,2);
189  theCTPPSHits.push_back(hit_tr1);
190  theCTPPSHits.push_back(hit_tr2);
191  theCTPPSHits.push_back(hit_tof);
192  }
193  if(Direction < 0.) {
194  PSimHit hit_tr1(xyzzy_tr1,xyzzy_tr1,0.,t0_tr1,0.,pid,2031091712,0,0.,0.,2);
195  PSimHit hit_tr2(xyzzy_tr2,xyzzy_tr2,0.,t0_tr2,0.,pid,2031616000,0,0.,0.,2);
196  PSimHit hit_tof(xyzzy_tof,xyzzy_tof,0.,t0_tof,0.,pid,2063597568,0,0.,0.,2);
197  theCTPPSHits.push_back(hit_tr1);
198  theCTPPSHits.push_back(hit_tr2);
199  theCTPPSHits.push_back(hit_tof);
200  }
201  }
202  }
203  }
204  std::unique_ptr<edm::PSimHitContainer> pctpps(new edm::PSimHitContainer);
205  int n = 0;
206  for ( std::vector<PSimHit>::const_iterator i = theCTPPSHits.begin();
207  i != theCTPPSHits.end(); i++ ) {
208  pctpps->push_back(*i);
209  n += 1;
210  }
211  iEvent.put(std::move(pctpps),"CTPPSHits");
212 }
213 
214 // ------------ method called once each stream before processing any runs, lumis or events ------------
215  void
217 {
218 }
219 
220 // ------------ method called once each stream after processing all runs, lumis and events ------------
221 void
223 }
224 
225 
226 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
virtual void produce(edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
virtual void beginStream(edm::StreamID) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::Handle< edm::HepMCProduct > EvtHandle
const double mm_to_m
const double c_light_s
virtual void endStream() override
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::HepMCProduct > mcEventToken
def pv(vc)
Definition: MetAnalyzer.py:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CLHEP::HepLorentzVector LorentzVector
const double s_to_ns
const double m_to_mm
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
HLT enums.
std::vector< PSimHit > PSimHitContainer
CTPPSSimHitProducer(const edm::ParameterSet &)
def move(src, dest)
Definition: eostools.py:510