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 
53 //
54 // class declaration
55 //
56 
58 public:
59  explicit CTPPSSimHitProducer(const edm::ParameterSet&);
60  ~CTPPSSimHitProducer() override;
61 
62 private:
63  void beginStream(edm::StreamID) override;
64  void produce(edm::Event&, const edm::EventSetup&) override;
65  void endStream() override;
66 
69  // ----------member data ---------------------------
71 };
72 
73 //
74 // constants, enums and typedefs
75 //
76 
77 //
78 // static data member definitions
79 //
80 
81 //
82 // constructors and destructor
83 //
85  produces<edm::PSimHitContainer>("CTPPSHits");
86  // consumes
87  mcEventToken =
88  mayConsume<edm::HepMCProduct>(iConfig.getUntrackedParameter<edm::InputTag>("MCEvent", std::string("")));
89 
90  // Read the position of the trackers and timing detectors
91  fz_tracker1 = iConfig.getParameter<double>("Z_Tracker1");
92  fz_tracker2 = iConfig.getParameter<double>("Z_Tracker2");
93  fz_timing = iConfig.getParameter<double>("Z_Timing");
94 }
95 
97 
98 //
99 // member functions
100 //
101 
102 // ------------ method called to produce the data ------------
104  using namespace edm;
105 
106  std::vector<PSimHit> theCTPPSHits;
107  iEvent.getByToken(mcEventToken, EvtHandle);
108 
109  const HepMC::GenEvent* Evt = EvtHandle->GetEvent();
110  std::vector<math::XYZTLorentzVector> protonCTPPS;
111  protonCTPPS.clear();
112  for (HepMC::GenEvent::vertex_const_iterator ivtx = Evt->vertices_begin(); ivtx != Evt->vertices_end(); ivtx++) {
113  if ((*ivtx)->id() != 0)
114  continue;
115  double prim_vtxZ = (*ivtx)->position().z() * mm_to_m; //in meters
116  // Get the vertices at the entrance of CTPPS and get the protons coming out of them (propagated by Hector)
117  for (HepMC::GenVertex::particles_out_const_iterator i = (*ivtx)->particles_out_const_begin();
118  i != (*ivtx)->particles_out_const_end();
119  ++i) {
120  int pid = (*i)->pdg_id();
121  if (pid != 2212)
122  continue;
123 
124  HepMC::GenVertex* pv = (*i)->production_vertex();
125  const HepMC::FourVector& vertex = pv->position();
126  const HepMC::FourVector p((*i)->momentum());
127  protonCTPPS.push_back(math::XYZTLorentzVector(p.x(), p.y(), p.z(), p.t()));
128 
129  LocalPoint initialPosition_tr1, initialPosition_tr2, initialPosition_tof;
130 
131  double t0_tr1 = 0., t0_tr2 = 0., t0_tof = 0.;
132  // t0_* has dimensions of mm
133  // Convert to ns for internal calculations.
134  const double c_light_s = 2.99792458e+11; // mm/s;
135  const double s_to_ns = 1.e9;
136  const double m_to_mm = 1.e3;
137  double x_tr1 = 0., x_tr2 = 0., x_tof = 0., y_tr1 = 0., y_tr2 = 0., y_tof = 0., z_tr1 = 0.;
138  double z_tr2 = fz_tracker2; //m
139  double z_tof = fz_timing; //m
140  int Direction = 0;
141 
142  // Read the vertex made by SimTransport at the entrance of det1
143  if (std::abs(vertex.eta()) > 8. && (*i)->status() == 1) {
144  if (vertex.z() > 0)
145  Direction = 1;
146  else if (vertex.z() < 0)
147  Direction = -1;
148 
149  //Get the global coordinates at Tracker1, equal to those of the vertex at CTPPS
150  x_tr1 = vertex.x();
151  y_tr1 = vertex.y();
152  z_tr1 = vertex.z() * mm_to_m; //move z from mm to meters
153  Local3DPoint xyzzy_tr1(x_tr1, y_tr1, z_tr1);
154  initialPosition_tr1 = xyzzy_tr1;
155  t0_tr1 = vertex.t() / c_light_s * s_to_ns;
156  //Get the global coordinates at Tracker2 by propagating as a straight line from Tracker1
157  t0_tr2 = z_tr2 * m_to_mm / c_light_s * s_to_ns; //discuss latter if needs to be corrected with vertex position
158  t0_tr2 = t0_tr1 + (z_tr2 - z_tr1) * m_to_mm / c_light_s * s_to_ns; //corrected with vertex position
159  z_tr2 *= Direction;
160  x_tr2 = x_tr1 + (p.x() / p.z()) * (z_tr2 - z_tr1) * m_to_mm;
161  y_tr2 = y_tr1 + (p.y() / p.z()) * (z_tr2 - z_tr1) * m_to_mm;
162  Local3DPoint xyzzy_tr2(x_tr2, y_tr2, z_tr2);
163  initialPosition_tr2 = xyzzy_tr2;
164  //Propagate as a straight line from Tracker1
165  t0_tof = z_tof * m_to_mm / c_light_s * s_to_ns; //discuss latter if needs to be corrected with vertex position
166  t0_tof = (z_tof - prim_vtxZ) * m_to_mm / c_light_s * s_to_ns; //corrected with vertex position
167  z_tof *= Direction;
168  x_tof = x_tr1 + (p.x() / p.z()) * (z_tof - z_tr1) * m_to_mm;
169  y_tof = y_tr1 + (p.y() / p.z()) * (z_tof - z_tr1) * m_to_mm;
170  Local3DPoint xyzzy_tof(x_tof, y_tof, z_tof);
171  initialPosition_tof = xyzzy_tof;
172 
173  // DetId codification for PSimHit from CTPPSPixel- It will be replaced by CTPPSDetId
174  // 2014314496 -> Tracker1 zPositive
175  // 2014838784 -> Tracker2 zPositive
176  // 2046820352 -> Timing zPositive
177  // 2031091712 -> Tracker1 zNegative
178  // 2031616000 -> Tracker2 zNegative
179  // 2063597568 -> Timing zNegative
180 
181  if (Direction > 0.) {
182  PSimHit hit_tr1(xyzzy_tr1, xyzzy_tr1, 0., t0_tr1, 0., pid, 2014314496, 0, 0., 0., 2);
183  PSimHit hit_tr2(xyzzy_tr2, xyzzy_tr2, 0., t0_tr2, 0., pid, 2014838784, 0, 0., 0., 2);
184  PSimHit hit_tof(xyzzy_tof, xyzzy_tof, 0., t0_tof, 0., pid, 2046820352, 0, 0., 0., 2);
185  theCTPPSHits.push_back(hit_tr1);
186  theCTPPSHits.push_back(hit_tr2);
187  theCTPPSHits.push_back(hit_tof);
188  }
189  if (Direction < 0.) {
190  PSimHit hit_tr1(xyzzy_tr1, xyzzy_tr1, 0., t0_tr1, 0., pid, 2031091712, 0, 0., 0., 2);
191  PSimHit hit_tr2(xyzzy_tr2, xyzzy_tr2, 0., t0_tr2, 0., pid, 2031616000, 0, 0., 0., 2);
192  PSimHit hit_tof(xyzzy_tof, xyzzy_tof, 0., t0_tof, 0., pid, 2063597568, 0, 0., 0., 2);
193  theCTPPSHits.push_back(hit_tr1);
194  theCTPPSHits.push_back(hit_tr2);
195  theCTPPSHits.push_back(hit_tof);
196  }
197  }
198  }
199  }
200  std::unique_ptr<edm::PSimHitContainer> pctpps(new edm::PSimHitContainer);
201  int n = 0;
202  for (std::vector<PSimHit>::const_iterator i = theCTPPSHits.begin(); i != theCTPPSHits.end(); i++) {
203  pctpps->push_back(*i);
204  n += 1;
205  }
206  iEvent.put(std::move(pctpps), "CTPPSHits");
207 }
208 
209 // ------------ method called once each stream before processing any runs, lumis or events ------------
211 
212 // ------------ method called once each stream after processing all runs, lumis and events ------------
214 
215 //define this as a plug-in
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void produce(edm::Event &, const edm::EventSetup &) override
static const double c_light_s
edm::Handle< edm::HepMCProduct > EvtHandle
void beginStream(edm::StreamID) override
static const double s_to_ns
T getUntrackedParameter(std::string const &, T const &) const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
static const double m_to_mm
HLT enums.
static const double mm_to_m
edm::EDGetTokenT< edm::HepMCProduct > mcEventToken
std::vector< PSimHit > PSimHitContainer
CTPPSSimHitProducer(const edm::ParameterSet &)
def move(src, dest)
Definition: eostools.py:511