CMS 3D CMS Logo

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