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 
34 
35 #include "CLHEP/Random/RandGaussQ.h"
36 #include "CLHEP/Units/GlobalSystemOfUnits.h"
37 #include "CLHEP/Units/GlobalPhysicalConstants.h"
38 //#include "CLHEP/Vector/ThreeVector.h"
39 #include "HepMC/SimpleVector.h"
40 #include "TMatrixD.h"
41 
42 #include <iostream>
43 
44 using namespace edm;
45 using namespace std;
46 using namespace CLHEP;
47 
48 namespace CLHEP {
49  class HepRandomEngine;
50 }
51 
53 public:
55  ~BetaBoostEvtVtxGenerator() override;
56 
58  //virtual CLHEP::Hep3Vector * newVertex();
59  virtual HepMC::FourVector* newVertex(CLHEP::HepRandomEngine*);
60  void produce(edm::Event&, const edm::EventSetup&) override;
61  virtual TMatrixD* GetInvLorentzBoost();
62 
64  void sigmaZ(double s = 1.0);
65 
67  void X0(double m = 0) { fX0 = m; }
69  void Y0(double m = 0) { fY0 = m; }
71  void Z0(double m = 0) { fZ0 = m; }
72 
74  void Phi(double m = 0) { phi_ = m; }
76  void Alpha(double m = 0) { alpha_ = m; }
77  void Beta(double m = 0) { beta_ = m; }
78 
80  void betastar(double m = 0) { fbetastar = m; }
82  void emittance(double m = 0) { femittance = m; }
83 
85  double BetaFunction(double z, double z0);
86 
87 private:
91  BetaBoostEvtVtxGenerator& operator=(const BetaBoostEvtVtxGenerator& rhs) = delete;
92 
93  double alpha_, phi_;
94  //TMatrixD boost_;
95  double beta_;
96  double fX0, fY0, fZ0;
97  double fSigmaZ;
98  //double fdxdz, fdydz;
99  double fbetastar, femittance;
100  double falpha;
101 
102  HepMC::FourVector* fVertex;
103  TMatrixD* boost_;
104  double fTimeOffset;
105 
107 
109 };
110 
112  : fVertex(nullptr),
113  boost_(nullptr),
114  fTimeOffset(0),
115  sourceLabel(consumes<HepMCProduct>(p.getParameter<edm::InputTag>("src"))),
116  verbosity_(p.getUntrackedParameter<bool>("verbosity", false)) {
117  fX0 = p.getParameter<double>("X0") * cm;
118  fY0 = p.getParameter<double>("Y0") * cm;
119  fZ0 = p.getParameter<double>("Z0") * cm;
120  fSigmaZ = p.getParameter<double>("SigmaZ") * cm;
121  alpha_ = p.getParameter<double>("Alpha") * radian;
122  phi_ = p.getParameter<double>("Phi") * radian;
123  fbetastar = p.getParameter<double>("BetaStar") * cm;
124  femittance = p.getParameter<double>("Emittance") * cm; // this is not the normalized emittance
125  fTimeOffset = p.getParameter<double>("TimeOffset") * ns * c_light; // HepMC time units are mm
126  beta_ = p.getParameter<double>("Beta");
127  if (fSigmaZ <= 0) {
128  throw cms::Exception("Configuration") << "Error in BetaBoostEvtVtxGenerator: "
129  << "Illegal resolution in Z (SigmaZ is negative)";
130  }
131 
132  produces<edm::HepMCProduct>();
133 }
134 
136  delete fVertex;
137  if (boost_ != nullptr)
138  delete boost_;
139 }
140 
141 //Hep3Vector* BetaBoostEvtVtxGenerator::newVertex() {
142 HepMC::FourVector* BetaBoostEvtVtxGenerator::newVertex(CLHEP::HepRandomEngine* engine) {
143  double X, Y, Z;
144 
145  double tmp_sigz = CLHEP::RandGaussQ::shoot(engine, 0.0, fSigmaZ);
146  Z = tmp_sigz + fZ0;
147 
148  double tmp_sigx = BetaFunction(Z, fZ0);
149  // need sqrt(2) for beamspot width relative to single beam width
150  tmp_sigx /= sqrt(2.0);
151  X = CLHEP::RandGaussQ::shoot(engine, 0.0, tmp_sigx) + fX0; // + Z*fdxdz ;
152 
153  double tmp_sigy = BetaFunction(Z, fZ0);
154  // need sqrt(2) for beamspot width relative to single beam width
155  tmp_sigy /= sqrt(2.0);
156  Y = CLHEP::RandGaussQ::shoot(engine, 0.0, tmp_sigy) + fY0; // + Z*fdydz;
157 
158  double tmp_sigt = CLHEP::RandGaussQ::shoot(engine, 0.0, fSigmaZ);
159  double T = tmp_sigt + fTimeOffset;
160 
161  if (fVertex == nullptr)
162  fVertex = new HepMC::FourVector();
163  fVertex->set(X, Y, Z, T);
164 
165  return fVertex;
166 }
167 
169  return sqrt(femittance * (fbetastar + (((z - z0) * (z - z0)) / fbetastar)));
170 }
171 
173  if (s >= 0) {
174  fSigmaZ = s;
175  } else {
176  throw cms::Exception("LogicError") << "Error in BetaBoostEvtVtxGenerator::sigmaZ: "
177  << "Illegal resolution in Z (negative)";
178  }
179 }
180 
182  //alpha_ = 0;
183  //phi_ = 142.e-6;
184  // if (boost_ != 0 ) return boost_;
185 
186  //boost_.ResizeTo(4,4);
187  //boost_ = new TMatrixD(4,4);
188  TMatrixD tmpboost(4, 4);
189  TMatrixD tmpboostZ(4, 4);
190  TMatrixD tmpboostXYZ(4, 4);
191 
192  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
193 
194  // Lorentz boost to frame where the collision is head-on
195  // phi is the half crossing angle in the plane ZS
196  // alpha is the angle to the S axis from the X axis in the XY plane
197 
198  tmpboost(0, 0) = 1. / cos(phi_);
199  tmpboost(0, 1) = -cos(alpha_) * sin(phi_);
200  tmpboost(0, 2) = -tan(phi_) * sin(phi_);
201  tmpboost(0, 3) = -sin(alpha_) * sin(phi_);
202  tmpboost(1, 0) = -cos(alpha_) * tan(phi_);
203  tmpboost(1, 1) = 1.;
204  tmpboost(1, 2) = cos(alpha_) * tan(phi_);
205  tmpboost(1, 3) = 0.;
206  tmpboost(2, 0) = 0.;
207  tmpboost(2, 1) = -cos(alpha_) * sin(phi_);
208  tmpboost(2, 2) = cos(phi_);
209  tmpboost(2, 3) = -sin(alpha_) * sin(phi_);
210  tmpboost(3, 0) = -sin(alpha_) * tan(phi_);
211  tmpboost(3, 1) = 0.;
212  tmpboost(3, 2) = sin(alpha_) * tan(phi_);
213  tmpboost(3, 3) = 1.;
214  //cout<<"beta "<<beta_;
215  double gama = 1.0 / sqrt(1 - beta_ * beta_);
216  tmpboostZ(0, 0) = gama;
217  tmpboostZ(0, 1) = 0.;
218  tmpboostZ(0, 2) = -1.0 * beta_ * gama;
219  tmpboostZ(0, 3) = 0.;
220  tmpboostZ(1, 0) = 0.;
221  tmpboostZ(1, 1) = 1.;
222  tmpboostZ(1, 2) = 0.;
223  tmpboostZ(1, 3) = 0.;
224  tmpboostZ(2, 0) = -1.0 * beta_ * gama;
225  tmpboostZ(2, 1) = 0.;
226  tmpboostZ(2, 2) = gama;
227  tmpboostZ(2, 3) = 0.;
228  tmpboostZ(3, 0) = 0.;
229  tmpboostZ(3, 1) = 0.;
230  tmpboostZ(3, 2) = 0.;
231  tmpboostZ(3, 3) = 1.;
232 
233  tmpboostXYZ = tmpboostZ * tmpboost;
234  tmpboostXYZ.Invert();
235 
236  boost_ = new TMatrixD(tmpboostXYZ);
237  if (verbosity_) {
238  boost_->Print();
239  }
240 
241  return boost_;
242 }
243 
246  if (!rng.isAvailable()) {
247  throw cms::Exception("Configuration")
248  << "Attempt to get a random engine when the RandomNumberGeneratorService is not configured.\n"
249  "You must configure the service if you want an engine.\n";
250  }
251  CLHEP::HepRandomEngine* engine = &rng->getEngine(evt.streamID());
252 
253  Handle<HepMCProduct> HepUnsmearedMCEvt;
254  evt.getByToken(sourceLabel, HepUnsmearedMCEvt);
255 
256  // Copy the HepMC::GenEvent
257  HepMC::GenEvent* genevt = new HepMC::GenEvent(*HepUnsmearedMCEvt->GetEvent());
258  std::unique_ptr<edm::HepMCProduct> HepMCEvt(new edm::HepMCProduct(genevt));
259 
260  // generate new vertex & apply the shift
261  //
262  HepMCEvt->applyVtxGen(newVertex(engine));
263 
264  //HepMCEvt->LorentzBoost( 0., 142.e-6 );
265  HepMCEvt->boostToLab(GetInvLorentzBoost(), "vertex");
266  HepMCEvt->boostToLab(GetInvLorentzBoost(), "momentum");
267  evt.put(std::move(HepMCEvt));
268 
269  return;
270 }
271 
edm::RandomNumberGenerator::getEngine
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
electrons_cff.bool
bool
Definition: electrons_cff.py:393
funct::false
false
Definition: Factorize.h:29
BetaBoostEvtVtxGenerator
Definition: BetaBoostEvtVtxGenerator.cc:52
BeamSpotPI::sigmaZ
Definition: BeamSpotPayloadInspectorHelper.h:35
EDProducer.h
BetaBoostEvtVtxGenerator::falpha
double falpha
Definition: BetaBoostEvtVtxGenerator.cc:100
detailsBasic3DVector::z
float float float z
Definition: extBasic3DVector.h:14
BetaBoostEvtVtxGenerator::Y0
void Y0(double m=0)
set mean in Y in cm
Definition: BetaBoostEvtVtxGenerator.cc:69
X
#define X(str)
Definition: MuonsGrabber.cc:38
edm::EDGetTokenT< HepMCProduct >
edm
HLT enums.
Definition: AlignableModifier.h:19
RandomNumberGenerator.h
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
BetaBoostEvtVtxGenerator::fVertex
HepMC::FourVector * fVertex
Definition: BetaBoostEvtVtxGenerator.cc:102
BetaBoostEvtVtxGenerator::beta_
double beta_
Definition: BetaBoostEvtVtxGenerator.cc:95
BetaBoostEvtVtxGenerator::alpha_
double alpha_
Definition: BetaBoostEvtVtxGenerator.cc:93
BetaBoostEvtVtxGenerator::~BetaBoostEvtVtxGenerator
~BetaBoostEvtVtxGenerator() override
Definition: BetaBoostEvtVtxGenerator.cc:135
edm::Handle
Definition: AssociativeIterator.h:50
edm::Service::isAvailable
bool isAvailable() const
Definition: Service.h:40
BetaBoostEvtVtxGenerator::fSigmaZ
double fSigmaZ
Definition: BetaBoostEvtVtxGenerator.cc:97
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
MakerMacros.h
alignCSCRings.s
s
Definition: alignCSCRings.py:92
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
BetaBoostEvtVtxGenerator::Z0
void Z0(double m=0)
set mean in Z in cm
Definition: BetaBoostEvtVtxGenerator.cc:71
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
BetaBoostEvtVtxGenerator::Phi
void Phi(double m=0)
set half crossing angle
Definition: BetaBoostEvtVtxGenerator.cc:74
Service.h
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
BetaBoostEvtVtxGenerator::fX0
double fX0
Definition: BetaBoostEvtVtxGenerator.cc:96
DDAxes::z
BetaBoostEvtVtxGenerator::sigmaZ
void sigmaZ(double s=1.0)
set resolution in Z in cm
Definition: BetaBoostEvtVtxGenerator.cc:172
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:531
HLTMuonOfflineAnalyzer_cfi.z0
z0
Definition: HLTMuonOfflineAnalyzer_cfi.py:98
EDGetToken.h
BetaBoostEvtVtxGenerator::boost_
TMatrixD * boost_
Definition: BetaBoostEvtVtxGenerator.cc:103
CLHEP
Definition: CocoaGlobals.h:27
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
BetaBoostEvtVtxGenerator::femittance
double femittance
Definition: BetaBoostEvtVtxGenerator.cc:99
funct::tan
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
edm::Service
Definition: Service.h:30
BetaBoostEvtVtxGenerator::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: BetaBoostEvtVtxGenerator.cc:244
edm::Event::put
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::EventSetup
Definition: EventSetup.h:57
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
BetaBoostEvtVtxGenerator::Beta
void Beta(double m=0)
Definition: BetaBoostEvtVtxGenerator.cc:77
InputTag.h
BetaBoostEvtVtxGenerator::fbetastar
double fbetastar
Definition: BetaBoostEvtVtxGenerator.cc:99
edm::Event::streamID
StreamID streamID() const
Definition: Event.h:98
BetaBoostEvtVtxGenerator::newVertex
virtual HepMC::FourVector * newVertex(CLHEP::HepRandomEngine *)
return a new event vertex
Definition: BetaBoostEvtVtxGenerator.cc:142
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
BetaBoostEvtVtxGenerator::sourceLabel
edm::EDGetTokenT< HepMCProduct > sourceLabel
Definition: BetaBoostEvtVtxGenerator.cc:106
T
long double T
Definition: Basic3DVectorLD.h:48
BetaBoostEvtVtxGenerator::X0
void X0(double m=0)
set mean in X in cm
Definition: BetaBoostEvtVtxGenerator.cc:67
BetaBoostEvtVtxGenerator::phi_
double phi_
Definition: BetaBoostEvtVtxGenerator.cc:93
Exception
Definition: hltDiff.cc:246
BetaBoostEvtVtxGenerator::verbosity_
bool verbosity_
Definition: BetaBoostEvtVtxGenerator.cc:108
BetaBoostEvtVtxGenerator::BetaFunction
double BetaFunction(double z, double z0)
beta function
Definition: BetaBoostEvtVtxGenerator.cc:168
edm::EDProducer
Definition: EDProducer.h:35
Exception.h
BeamSpotPI::Y
Definition: BeamSpotPayloadInspectorHelper.h:31
BetaBoostEvtVtxGenerator::emittance
void emittance(double m=0)
emittance (no the normalized)
Definition: BetaBoostEvtVtxGenerator.cc:82
BetaBoostEvtVtxGenerator::GetInvLorentzBoost
virtual TMatrixD * GetInvLorentzBoost()
Definition: BetaBoostEvtVtxGenerator.cc:181
ConsumesCollector.h
BetaBoostEvtVtxGenerator::fTimeOffset
double fTimeOffset
Definition: BetaBoostEvtVtxGenerator.cc:104
edm::HepMCProduct
Definition: HepMCProduct.h:18
ParameterSet.h
BeamSpotPI::Z
Definition: BeamSpotPayloadInspectorHelper.h:32
BetaBoostEvtVtxGenerator::betastar
void betastar(double m=0)
set beta_star
Definition: BetaBoostEvtVtxGenerator.cc:80
BetaBoostEvtVtxGenerator::BetaBoostEvtVtxGenerator
BetaBoostEvtVtxGenerator(const edm::ParameterSet &p)
Definition: BetaBoostEvtVtxGenerator.cc:111
HepMCProduct.h
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
BetaBoostEvtVtxGenerator::Alpha
void Alpha(double m=0)
angle between crossing plane and horizontal plane
Definition: BetaBoostEvtVtxGenerator.cc:76
BetaBoostEvtVtxGenerator::fY0
double fY0
Definition: BetaBoostEvtVtxGenerator.cc:96
BetaBoostEvtVtxGenerator::fZ0
double fZ0
Definition: BetaBoostEvtVtxGenerator.cc:96