CMS 3D CMS Logo

BetaBoostEvtVtxGenerator.cc
Go to the documentation of this file.
1 
2 /*
3  ________________________________________________________________________
4 
5  BetaBoostEvtVtxGenerator
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 
33 
34 #include "CLHEP/Random/RandGaussQ.h"
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 namespace CLHEP {
48  class HepRandomEngine;
49 }
50 
52 public:
58  ~BetaBoostEvtVtxGenerator() override = default;
59 
61  HepMC::FourVector newVertex(CLHEP::HepRandomEngine*) const;
62  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
63  TMatrixD GetInvLorentzBoost() const;
64 
66  double BetaFunction(double z, double z0) const;
67 
68 private:
69  const double alpha_; // angle between crossing plane and horizontal plane
70  const double phi_; // half crossing angle
71  const double beta_;
72  const double fX0; // mean in X in cm
73  const double fY0; // mean in Y in cm
74  const double fZ0; // mean in Z in cm
75  const double fSigmaZ; // resolution in Z in cm
76  //double fdxdz, fdydz;
77  const double fbetastar;
78  const double femittance; // emittance (no the normalized)a
79 
80  const double fTimeOffset;
81 
82  const TMatrixD boost_;
83 
85 
86  const bool verbosity_;
87 };
88 
90  : alpha_(p.getParameter<double>("Alpha") * radian),
91  phi_(p.getParameter<double>("Phi") * radian),
92  beta_(p.getParameter<double>("Beta")),
93  fX0(p.getParameter<double>("X0") * cm),
94  fY0(p.getParameter<double>("Y0") * cm),
95  fZ0(p.getParameter<double>("Z0") * cm),
96  fSigmaZ(p.getParameter<double>("SigmaZ") * cm),
97  fbetastar(p.getParameter<double>("BetaStar") * cm),
98  femittance(p.getParameter<double>("Emittance") * cm), // this is not the normalized emittance
99  fTimeOffset(p.getParameter<double>("TimeOffset") * ns * c_light), // HepMC time units are mm
100  boost_(GetInvLorentzBoost()),
101  sourceLabel(consumes<HepMCProduct>(p.getParameter<edm::InputTag>("src"))),
102  verbosity_(p.getUntrackedParameter<bool>("verbosity", false)) {
103  if (fSigmaZ <= 0) {
104  throw cms::Exception("Configuration") << "Error in BetaBoostEvtVtxGenerator: "
105  << "Illegal resolution in Z (SigmaZ is negative)";
106  }
107 
108  produces<edm::HepMCProduct>();
109 }
110 
111 HepMC::FourVector BetaBoostEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) const {
112  double X, Y, Z;
113 
114  double tmp_sigz = CLHEP::RandGaussQ::shoot(engine, 0.0, fSigmaZ);
115  Z = tmp_sigz + fZ0;
116 
117  double tmp_sigx = BetaFunction(Z, fZ0);
118  // need sqrt(2) for beamspot width relative to single beam width
119  tmp_sigx /= sqrt(2.0);
120  X = CLHEP::RandGaussQ::shoot(engine, 0.0, tmp_sigx) + fX0; // + Z*fdxdz ;
121 
122  double tmp_sigy = BetaFunction(Z, fZ0);
123  // need sqrt(2) for beamspot width relative to single beam width
124  tmp_sigy /= sqrt(2.0);
125  Y = CLHEP::RandGaussQ::shoot(engine, 0.0, tmp_sigy) + fY0; // + Z*fdydz;
126 
127  double tmp_sigt = CLHEP::RandGaussQ::shoot(engine, 0.0, fSigmaZ);
128  double T = tmp_sigt + fTimeOffset;
129 
130  return HepMC::FourVector(X, Y, Z, T);
131 }
132 
133 double BetaBoostEvtVtxGenerator::BetaFunction(double z, double z0) const {
134  return sqrt(femittance * (fbetastar + (((z - z0) * (z - z0)) / fbetastar)));
135 }
136 
138  //alpha_ = 0;
139  //phi_ = 142.e-6;
140  // if (boost_ != 0 ) return boost_;
141 
142  //boost_.ResizeTo(4,4);
143  //boost_ = new TMatrixD(4,4);
144  TMatrixD tmpboost(4, 4);
145  TMatrixD tmpboostZ(4, 4);
146  TMatrixD tmpboostXYZ(4, 4);
147 
148  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
149 
150  // Lorentz boost to frame where the collision is head-on
151  // phi is the half crossing angle in the plane ZS
152  // alpha is the angle to the S axis from the X axis in the XY plane
153 
154  tmpboost(0, 0) = 1. / cos(phi_);
155  tmpboost(0, 1) = -cos(alpha_) * sin(phi_);
156  tmpboost(0, 2) = -tan(phi_) * sin(phi_);
157  tmpboost(0, 3) = -sin(alpha_) * sin(phi_);
158  tmpboost(1, 0) = -cos(alpha_) * tan(phi_);
159  tmpboost(1, 1) = 1.;
160  tmpboost(1, 2) = cos(alpha_) * tan(phi_);
161  tmpboost(1, 3) = 0.;
162  tmpboost(2, 0) = 0.;
163  tmpboost(2, 1) = -cos(alpha_) * sin(phi_);
164  tmpboost(2, 2) = cos(phi_);
165  tmpboost(2, 3) = -sin(alpha_) * sin(phi_);
166  tmpboost(3, 0) = -sin(alpha_) * tan(phi_);
167  tmpboost(3, 1) = 0.;
168  tmpboost(3, 2) = sin(alpha_) * tan(phi_);
169  tmpboost(3, 3) = 1.;
170  //cout<<"beta "<<beta_;
171  double gama = 1.0 / sqrt(1 - beta_ * beta_);
172  tmpboostZ(0, 0) = gama;
173  tmpboostZ(0, 1) = 0.;
174  tmpboostZ(0, 2) = -1.0 * beta_ * gama;
175  tmpboostZ(0, 3) = 0.;
176  tmpboostZ(1, 0) = 0.;
177  tmpboostZ(1, 1) = 1.;
178  tmpboostZ(1, 2) = 0.;
179  tmpboostZ(1, 3) = 0.;
180  tmpboostZ(2, 0) = -1.0 * beta_ * gama;
181  tmpboostZ(2, 1) = 0.;
182  tmpboostZ(2, 2) = gama;
183  tmpboostZ(2, 3) = 0.;
184  tmpboostZ(3, 0) = 0.;
185  tmpboostZ(3, 1) = 0.;
186  tmpboostZ(3, 2) = 0.;
187  tmpboostZ(3, 3) = 1.;
188 
189  tmpboostXYZ = tmpboostZ * tmpboost;
190  tmpboostXYZ.Invert();
191 
192  if (verbosity_) {
193  tmpboostXYZ.Print();
194  }
195 
196  return tmpboostXYZ;
197 }
198 
201  if (!rng.isAvailable()) {
202  throw cms::Exception("Configuration")
203  << "Attempt to get a random engine when the RandomNumberGeneratorService is not configured.\n"
204  "You must configure the service if you want an engine.\n";
205  }
206  CLHEP::HepRandomEngine* engine = &rng->getEngine(evt.streamID());
207 
208  const auto& HepUnsmearedMCEvt = evt.get(sourceLabel);
209 
210  // Copy the HepMC::GenEvent
211  auto HepMCEvt = std::make_unique<edm::HepMCProduct>(new HepMC::GenEvent(*HepUnsmearedMCEvt.GetEvent()));
212 
213  // generate new vertex & apply the shift
214  //
215  auto vertex = newVertex(engine);
216  HepMCEvt->applyVtxGen(&vertex);
217 
218  //HepMCEvt->LorentzBoost( 0., 142.e-6 );
219  HepMCEvt->boostToLab(&boost_, "vertex");
220  HepMCEvt->boostToLab(&boost_, "momentum");
221  evt.put(std::move(HepMCEvt));
222 }
223 
Basic3DVector & operator=(const Basic3DVector &)=default
Assignment operator.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:344
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
HepMC::FourVector newVertex(CLHEP::HepRandomEngine *) const
return a new event vertex
double BetaFunction(double z, double z0) const
beta function
#define X(str)
Definition: MuonsGrabber.cc:38
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
float float float z
const edm::EDGetTokenT< HepMCProduct > sourceLabel
BetaBoostEvtVtxGenerator(const edm::ParameterSet &p)
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
StreamID streamID() const
Definition: Event.h:98
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
HLT enums.
bool isAvailable() const
Definition: Service.h:40
long double T
def move(src, dest)
Definition: eostools.py:511