CMS 3D CMS Logo

MixBoostEvtVtxGenerator.cc
Go to the documentation of this file.
1 
2 /*
3 ________________________________________________________________________
4 
5  MixBoostEvtVtxGenerator
6 
7  Smear vertex according to the Beta function on the transverse plane
8  and a Gaussian on the z axis. It allows the beam to have a crossing
9  angle (slopes dxdz and dydz).
10 
11  Based on GaussEvtVtxGenerator
12  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
13 
14  FERMILAB
15  2006
16 ________________________________________________________________________
17 */
18 
19 //lingshan: add beta for z-axis boost
20 
21 //#include "IOMC/EventVertexGenerators/interface/BetafuncEvtVtxGenerator.h"
22 
26 
29 
34 
35 #include "CLHEP/Units/GlobalSystemOfUnits.h"
36 #include "CLHEP/Units/GlobalPhysicalConstants.h"
37 //#include "CLHEP/Vector/ThreeVector.h"
38 #include "HepMC/SimpleVector.h"
39 #include "TMatrixD.h"
40 
41 #include <iostream>
42 
43 using namespace edm;
44 using namespace std;
45 using namespace CLHEP;
46 
47 class RandGaussQ;
48 class FourVector;
49 
51 public:
56  MixBoostEvtVtxGenerator& operator=(const MixBoostEvtVtxGenerator& rhs) = delete;
57  ~MixBoostEvtVtxGenerator() override;
58 
60  void produce(edm::Event&, const edm::EventSetup&) override;
61  virtual TMatrixD* GetInvLorentzBoost();
62  virtual HepMC::FourVector* getVertex(edm::Event&);
63  virtual HepMC::FourVector* getRecVertex(edm::Event&);
64 
66  void sigmaZ(double s = 1.0);
67 
69  void X0(double m = 0) { fX0 = m; }
71  void Y0(double m = 0) { fY0 = m; }
73  void Z0(double m = 0) { fZ0 = m; }
74 
76  void Phi(double m = 0) { phi_ = m; }
78  void Alpha(double m = 0) { alpha_ = m; }
79  void Beta(double m = 0) { beta_ = m; }
80 
82  void betastar(double m = 0) { fbetastar = m; }
84  void emittance(double m = 0) { femittance = m; }
85 
87  double BetaFunction(double z, double z0);
88 
89 private:
90  double alpha_, phi_;
91  //TMatrixD boost_;
92  double beta_;
93  double fX0, fY0, fZ0;
94  double fSigmaZ;
95  //double fdxdz, fdydz;
96  double fbetastar, femittance;
97  double falpha;
98 
99  HepMC::FourVector* fVertex;
100  TMatrixD* boost_;
101  double fTimeOffset;
102 
106  const bool useRecVertex;
107  std::vector<double> vtxOffset;
108 };
109 
111  : fVertex(nullptr),
112  boost_(nullptr),
113  fTimeOffset(0),
114  vtxLabel(mayConsume<reco::VertexCollection>(pset.getParameter<edm::InputTag>("vtxLabel"))),
115  signalLabel(consumes<HepMCProduct>(pset.getParameter<edm::InputTag>("signalLabel"))),
116  mixLabel(consumes<CrossingFrame<HepMCProduct> >(pset.getParameter<edm::InputTag>("mixLabel"))),
117  useRecVertex(pset.exists("useRecVertex") ? pset.getParameter<bool>("useRecVertex") : false) {
118  beta_ = pset.getParameter<double>("Beta");
119  alpha_ = 0;
120  phi_ = 0;
121  if (pset.exists("Alpha")) {
122  alpha_ = pset.getParameter<double>("Alpha") * radian;
123  phi_ = pset.getParameter<double>("Phi") * radian;
124  }
125 
126  vtxOffset.resize(3);
127  if (pset.exists("vtxOffset"))
128  vtxOffset = pset.getParameter<std::vector<double> >("vtxOffset");
129 
130  produces<edm::HepMCProduct>();
131 }
132 
134  if (fVertex != nullptr)
135  delete fVertex;
136  if (boost_ != nullptr)
137  delete boost_;
138 }
139 
141  return sqrt(femittance * (fbetastar + (((z - z0) * (z - z0)) / fbetastar)));
142 }
143 
145  if (s >= 0) {
146  fSigmaZ = s;
147  } else {
148  throw cms::Exception("LogicError") << "Error in MixBoostEvtVtxGenerator::sigmaZ: "
149  << "Illegal resolution in Z (negative)";
150  }
151 }
152 
154  //alpha_ = 0;
155  //phi_ = 142.e-6;
156  // if (boost_ != 0 ) return boost_;
157 
158  //boost_.ResizeTo(4,4);
159  //boost_ = new TMatrixD(4,4);
160  TMatrixD tmpboost(4, 4);
161  TMatrixD tmpboostZ(4, 4);
162  TMatrixD tmpboostXYZ(4, 4);
163 
164  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
165 
166  // Lorentz boost to frame where the collision is head-on
167  // phi is the half crossing angle in the plane ZS
168  // alpha is the angle to the S axis from the X axis in the XY plane
169 
170  tmpboost(0, 0) = 1. / cos(phi_);
171  tmpboost(0, 1) = -cos(alpha_) * sin(phi_);
172  tmpboost(0, 2) = -tan(phi_) * sin(phi_);
173  tmpboost(0, 3) = -sin(alpha_) * sin(phi_);
174  tmpboost(1, 0) = -cos(alpha_) * tan(phi_);
175  tmpboost(1, 1) = 1.;
176  tmpboost(1, 2) = cos(alpha_) * tan(phi_);
177  tmpboost(1, 3) = 0.;
178  tmpboost(2, 0) = 0.;
179  tmpboost(2, 1) = -cos(alpha_) * sin(phi_);
180  tmpboost(2, 2) = cos(phi_);
181  tmpboost(2, 3) = -sin(alpha_) * sin(phi_);
182  tmpboost(3, 0) = -sin(alpha_) * tan(phi_);
183  tmpboost(3, 1) = 0.;
184  tmpboost(3, 2) = sin(alpha_) * tan(phi_);
185  tmpboost(3, 3) = 1.;
186  //cout<<"beta "<<beta_;
187  double gama = 1.0 / sqrt(1 - beta_ * beta_);
188  tmpboostZ(0, 0) = gama;
189  tmpboostZ(0, 1) = 0.;
190  tmpboostZ(0, 2) = -1.0 * beta_ * gama;
191  tmpboostZ(0, 3) = 0.;
192  tmpboostZ(1, 0) = 0.;
193  tmpboostZ(1, 1) = 1.;
194  tmpboostZ(1, 2) = 0.;
195  tmpboostZ(1, 3) = 0.;
196  tmpboostZ(2, 0) = -1.0 * beta_ * gama;
197  tmpboostZ(2, 1) = 0.;
198  tmpboostZ(2, 2) = gama;
199  tmpboostZ(2, 3) = 0.;
200  tmpboostZ(3, 0) = 0.;
201  tmpboostZ(3, 1) = 0.;
202  tmpboostZ(3, 2) = 0.;
203  tmpboostZ(3, 3) = 1.;
204 
205  tmpboostXYZ = tmpboostZ * tmpboost;
206  tmpboostXYZ.Invert();
207 
208  //cout<<"Boosting with beta : "<<beta_<<endl;
209 
210  boost_ = new TMatrixD(tmpboostXYZ);
211  boost_->Print();
212 
213  return boost_;
214 }
215 
216 HepMC::FourVector* MixBoostEvtVtxGenerator::getVertex(Event& evt) {
217  const HepMC::GenEvent* inev = nullptr;
218 
221 
222  const HepMCProduct& bkg = mix.getObject(1);
223  if (!(bkg.isVtxGenApplied())) {
224  throw cms::Exception("MatchVtx") << "Input background does not have smeared vertex!" << endl;
225  } else {
226  inev = bkg.GetEvent();
227  }
228 
229  HepMC::GenVertex* genvtx = inev->signal_process_vertex();
230  if (!genvtx) {
231  cout << "No Signal Process Vertex!" << endl;
232  HepMC::GenEvent::particle_const_iterator pt = inev->particles_begin();
233  HepMC::GenEvent::particle_const_iterator ptend = inev->particles_end();
234  while (!genvtx || (genvtx->particles_in_size() == 1 && pt != ptend)) {
235  if (!genvtx)
236  cout << "No Gen Vertex!" << endl;
237  if (pt == ptend)
238  cout << "End reached!" << endl;
239  genvtx = (*pt)->production_vertex();
240  ++pt;
241  }
242  }
243  double aX, aY, aZ, aT;
244 
245  aX = genvtx->position().x();
246  aY = genvtx->position().y();
247  aZ = genvtx->position().z();
248  aT = genvtx->position().t();
249 
250  if (!fVertex)
251  fVertex = new HepMC::FourVector();
252  fVertex->set(aX, aY, aZ, aT);
253 
254  return fVertex;
255 }
256 
259 
260  double aX, aY, aZ;
261 
262  aX = input->begin()->position().x() + vtxOffset[0];
263  aY = input->begin()->position().y() + vtxOffset[1];
264  aZ = input->begin()->position().z() + vtxOffset[2];
265 
266  if (!fVertex)
267  fVertex = new HepMC::FourVector();
268  fVertex->set(10.0 * aX, 10.0 * aY, 10.0 * aZ, 0.0); // HepMC positions in mm (RECO in cm)
269 
270  return fVertex;
271 }
272 
274  const edm::Handle<HepMCProduct>& HepUnsmearedMCEvt = evt.getHandle(signalLabel);
275 
276  // Copy the HepMC::GenEvent
277  HepMC::GenEvent* genevt = new HepMC::GenEvent(*HepUnsmearedMCEvt->GetEvent());
278  std::unique_ptr<edm::HepMCProduct> HepMCEvt(new edm::HepMCProduct(genevt));
279  // generate new vertex & apply the shift
280  //
281 
282  HepMCEvt->boostToLab(GetInvLorentzBoost(), "vertex");
283  HepMCEvt->boostToLab(GetInvLorentzBoost(), "momentum");
284 
285  HepMCEvt->applyVtxGen(useRecVertex ? getRecVertex(evt) : getVertex(evt));
286 
287  evt.put(std::move(HepMCEvt));
288  return;
289 }
290 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
void Z0(double m=0)
set mean in Z in cm
const edm::EDGetTokenT< CrossingFrame< HepMCProduct > > mixLabel
std::vector< double > vtxOffset
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void emittance(double m=0)
emittance (no the normalized)
T const * product() const
Definition: Handle.h:70
void X0(double m=0)
set mean in X in cm
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
float float float z
static std::string const input
Definition: EdmProvDump.cc:47
const edm::EDGetTokenT< reco::VertexCollection > vtxLabel
void sigmaZ(double s=1.0)
set resolution in Z in cm
T sqrt(T t)
Definition: SSEVec.h:19
void produce(edm::Event &, const edm::EventSetup &) override
return a new event vertex
const edm::EDGetTokenT< HepMCProduct > signalLabel
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual TMatrixD * GetInvLorentzBoost()
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
MixBoostEvtVtxGenerator(const edm::ParameterSet &p)
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
void Alpha(double m=0)
angle between crossing plane and horizontal plane
void Phi(double m=0)
set half crossing angle
void betastar(double m=0)
set beta_star
void Y0(double m=0)
set mean in Y in cm
fixed size matrix
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
HLT enums.
virtual HepMC::FourVector * getRecVertex(edm::Event &)
virtual HepMC::FourVector * getVertex(edm::Event &)
def move(src, dest)
Definition: eostools.py:511
double BetaFunction(double z, double z0)
beta function