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 class RandGaussQ;
47 
48 
50 public:
52  virtual ~BetaBoostEvtVtxGenerator();
53 
55  //virtual CLHEP::Hep3Vector * newVertex();
56  virtual HepMC::FourVector* newVertex() ;
57  virtual void produce( edm::Event&, const edm::EventSetup& ) override;
58  virtual TMatrixD* GetInvLorentzBoost();
59 
60 
62  void sigmaZ(double s=1.0);
63 
65  void X0(double m=0) { fX0=m; }
67  void Y0(double m=0) { fY0=m; }
69  void Z0(double m=0) { fZ0=m; }
70 
72  void Phi(double m=0) { phi_=m; }
74  void Alpha(double m=0) { alpha_=m; }
75  void Beta(double m=0) { beta_=m; }
76 
78  void betastar(double m=0) { fbetastar=m; }
80  void emittance(double m=0) { femittance=m; }
81 
83  double BetaFunction(double z, double z0);
84  CLHEP::HepRandomEngine& getEngine();
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 
105  CLHEP::HepRandomEngine* fEngine;
107 
108  CLHEP::RandGaussQ* fRandom ;
109 
111 
112 };
113 
114 
116  fVertex(0), boost_(0), fTimeOffset(0), fEngine(0),
117  sourceLabel(p.getParameter<edm::InputTag>("src")),
118  verbosity_(p.getUntrackedParameter<bool>("verbosity",false))
119 {
120 
122 
123  if ( ! rng.isAvailable()) {
124 
125  throw cms::Exception("Configuration")
126  << "The BaseEvtVtxGenerator requires the RandomNumberGeneratorService\n"
127  "which is not present in the configuration file. You must add the service\n"
128  "in the configuration file or remove the modules that require it.";
129  }
130 
131  CLHEP::HepRandomEngine& engine = rng->getEngine();
132  fEngine = &engine;
133  fRandom = new CLHEP::RandGaussQ(getEngine());
134 
135  fX0 = p.getParameter<double>("X0")*cm;
136  fY0 = p.getParameter<double>("Y0")*cm;
137  fZ0 = p.getParameter<double>("Z0")*cm;
138  fSigmaZ = p.getParameter<double>("SigmaZ")*cm;
139  alpha_ = p.getParameter<double>("Alpha")*radian;
140  phi_ = p.getParameter<double>("Phi")*radian;
141  fbetastar = p.getParameter<double>("BetaStar")*cm;
142  femittance = p.getParameter<double>("Emittance")*cm; // this is not the normalized emittance
143  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light; // HepMC time units are mm
144  beta_=p.getParameter<double>("Beta");
145  if (fSigmaZ <= 0) {
146  throw cms::Exception("Configuration")
147  << "Error in BetaBoostEvtVtxGenerator: "
148  << "Illegal resolution in Z (SigmaZ is negative)";
149  }
150 
151  produces<bool>();
152 
153 }
154 
156 {
157  delete fVertex ;
158  if (boost_ != 0 ) delete boost_;
159  delete fRandom;
160 }
161 
162 CLHEP::HepRandomEngine& BetaBoostEvtVtxGenerator::getEngine(){
163  return *fEngine;
164 }
165 
166 
167 //Hep3Vector* BetaBoostEvtVtxGenerator::newVertex() {
169 
170 
171  double X,Y,Z;
172 
173  double tmp_sigz = fRandom->fire(0., fSigmaZ);
174  Z = tmp_sigz + fZ0;
175 
176  double tmp_sigx = BetaFunction(Z,fZ0);
177  // need sqrt(2) for beamspot width relative to single beam width
178  tmp_sigx /= sqrt(2.0);
179  X = fRandom->fire(0.,tmp_sigx) + fX0; // + Z*fdxdz ;
180 
181  double tmp_sigy = BetaFunction(Z,fZ0);
182  // need sqrt(2) for beamspot width relative to single beam width
183  tmp_sigy /= sqrt(2.0);
184  Y = fRandom->fire(0.,tmp_sigy) + fY0; // + Z*fdydz;
185 
186  double tmp_sigt = fRandom->fire(0., fSigmaZ);
187  double T = tmp_sigt + fTimeOffset;
188 
189  if ( fVertex == 0 ) fVertex = new HepMC::FourVector();
190  fVertex->set(X,Y,Z,T);
191 
192  return fVertex;
193 }
194 
195 double BetaBoostEvtVtxGenerator::BetaFunction(double z, double z0)
196 {
197  return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
198 
199 }
200 
201 
203 {
204  if (s>=0 ) {
205  fSigmaZ=s;
206  }
207  else {
208  throw cms::Exception("LogicError")
209  << "Error in BetaBoostEvtVtxGenerator::sigmaZ: "
210  << "Illegal resolution in Z (negative)";
211  }
212 }
213 
215 
216  //alpha_ = 0;
217  //phi_ = 142.e-6;
218  // if (boost_ != 0 ) return boost_;
219 
220  //boost_.ResizeTo(4,4);
221  //boost_ = new TMatrixD(4,4);
222  TMatrixD tmpboost(4,4);
223  TMatrixD tmpboostZ(4,4);
224  TMatrixD tmpboostXYZ(4,4);
225 
226  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
227 
228  // Lorentz boost to frame where the collision is head-on
229  // phi is the half crossing angle in the plane ZS
230  // alpha is the angle to the S axis from the X axis in the XY plane
231 
232  tmpboost(0,0) = 1./cos(phi_);
233  tmpboost(0,1) = - cos(alpha_)*sin(phi_);
234  tmpboost(0,2) = - tan(phi_)*sin(phi_);
235  tmpboost(0,3) = - sin(alpha_)*sin(phi_);
236  tmpboost(1,0) = - cos(alpha_)*tan(phi_);
237  tmpboost(1,1) = 1.;
238  tmpboost(1,2) = cos(alpha_)*tan(phi_);
239  tmpboost(1,3) = 0.;
240  tmpboost(2,0) = 0.;
241  tmpboost(2,1) = - cos(alpha_)*sin(phi_);
242  tmpboost(2,2) = cos(phi_);
243  tmpboost(2,3) = - sin(alpha_)*sin(phi_);
244  tmpboost(3,0) = - sin(alpha_)*tan(phi_);
245  tmpboost(3,1) = 0.;
246  tmpboost(3,2) = sin(alpha_)*tan(phi_);
247  tmpboost(3,3) = 1.;
248  //cout<<"beta "<<beta_;
249  double gama=1.0/sqrt(1-beta_*beta_);
250  tmpboostZ(0,0)=gama;
251  tmpboostZ(0,1)=0.;
252  tmpboostZ(0,2)=-1.0*beta_*gama;
253  tmpboostZ(0,3)=0.;
254  tmpboostZ(1,0)=0.;
255  tmpboostZ(1,1) = 1.;
256  tmpboostZ(1,2)=0.;
257  tmpboostZ(1,3)=0.;
258  tmpboostZ(2,0)=-1.0*beta_*gama;
259  tmpboostZ(2,1) = 0.;
260  tmpboostZ(2,2)=gama;
261  tmpboostZ(2,3) = 0.;
262  tmpboostZ(3,0)=0.;
263  tmpboostZ(3,1)=0.;
264  tmpboostZ(3,2)=0.;
265  tmpboostZ(3,3) = 1.;
266 
267  tmpboostXYZ=tmpboostZ*tmpboost;
268  tmpboostXYZ.Invert();
269 
270 
271 
272  boost_ = new TMatrixD(tmpboostXYZ);
273  if ( verbosity_ ) { boost_->Print(); }
274 
275  return boost_;
276 }
277 
279 {
280 
281 
282  Handle<HepMCProduct> HepMCEvt ;
283  evt.getByLabel( sourceLabel, HepMCEvt ) ;
284 
285  // generate new vertex & apply the shift
286  //
287  HepMCEvt->applyVtxGen( newVertex() ) ;
288 
289  //HepMCEvt->LorentzBoost( 0., 142.e-6 );
290  HepMCEvt->boostToLab( GetInvLorentzBoost(), "vertex" );
291  HepMCEvt->boostToLab( GetInvLorentzBoost(), "momentum" );
292  // OK, create a (pseudo)product and put in into edm::Event
293  //
294  auto_ptr<bool> NewProduct(new bool(true)) ;
295  evt.put( NewProduct ) ;
296  return ;
297 }
298 
299 
300 
301 
const double Z[kNumberCalorimeter]
T getParameter(std::string const &) const
void Y0(double m=0)
set mean in Y in cm
virtual HepMC::FourVector * newVertex()
return a new event vertex
void emittance(double m=0)
emittance (no the normalized)
CLHEP::HepRandomEngine & getEngine()
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
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
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:116
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
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
virtual CLHEP::HepRandomEngine & getEngine() const =0
Use this to get the random number engine, this is the only function most users should call...
void betastar(double m=0)
set beta_star
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
double BetaFunction(double z, double z0)
beta function
CLHEP::HepRandomEngine * fEngine
void Alpha(double m=0)
angle between crossing plane and horizontal plane
return(e1-e2)*(e1-e2)+dp *dp
volatile std::atomic< bool > shutdown_flag false
void Z0(double m=0)
set mean in Z in cm
long double T
virtual TMatrixD * GetInvLorentzBoost()