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 // $Id: BetaBoostEvtVtxGenerator.cc,v 1.3 2012/09/25 10:00:12 fabiocos Exp $
3 /*
4  ________________________________________________________________________
5 
6  BetaBoostEvtVtxGenerator
7 
8  Smear vertex according to the Beta function on the transverse plane
9  and a Gaussian on the z axis. It allows the beam to have a crossing
10  angle (slopes dxdz and dydz).
11 
12  Based on GaussEvtVtxGenerator
13  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
14 
15  FERMILAB
16  2006
17  ________________________________________________________________________
18 */
19 
20 //lingshan: add beta for z-axis boost
21 
22 //#include "IOMC/EventVertexGenerators/interface/BetafuncEvtVtxGenerator.h"
23 
27 
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 class RandGaussQ;
48 
49 
51 public:
53  virtual ~BetaBoostEvtVtxGenerator();
54 
56  //virtual CLHEP::Hep3Vector * newVertex();
57  virtual HepMC::FourVector* newVertex() ;
58  virtual void produce( edm::Event&, const edm::EventSetup& );
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  CLHEP::HepRandomEngine& getEngine();
86 
87 private:
91  BetaBoostEvtVtxGenerator& operator = (const BetaBoostEvtVtxGenerator & rhs );
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 
106  CLHEP::HepRandomEngine* fEngine;
108 
109  CLHEP::RandGaussQ* fRandom ;
110 
112 
113 };
114 
115 
117  fVertex(0), boost_(0), fTimeOffset(0), fEngine(0),
118  sourceLabel(p.getParameter<edm::InputTag>("src")),
119  verbosity_(p.getUntrackedParameter<bool>("verbosity",false))
120 {
121 
123 
124  if ( ! rng.isAvailable()) {
125 
126  throw cms::Exception("Configuration")
127  << "The BaseEvtVtxGenerator requires the RandomNumberGeneratorService\n"
128  "which is not present in the configuration file. You must add the service\n"
129  "in the configuration file or remove the modules that require it.";
130  }
131 
132  CLHEP::HepRandomEngine& engine = rng->getEngine();
133  fEngine = &engine;
134  fRandom = new CLHEP::RandGaussQ(getEngine());
135 
136  fX0 = p.getParameter<double>("X0")*cm;
137  fY0 = p.getParameter<double>("Y0")*cm;
138  fZ0 = p.getParameter<double>("Z0")*cm;
139  fSigmaZ = p.getParameter<double>("SigmaZ")*cm;
140  alpha_ = p.getParameter<double>("Alpha")*radian;
141  phi_ = p.getParameter<double>("Phi")*radian;
142  fbetastar = p.getParameter<double>("BetaStar")*cm;
143  femittance = p.getParameter<double>("Emittance")*cm; // this is not the normalized emittance
144  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light; // HepMC time units are mm
145  beta_=p.getParameter<double>("Beta");
146  if (fSigmaZ <= 0) {
147  throw cms::Exception("Configuration")
148  << "Error in BetaBoostEvtVtxGenerator: "
149  << "Illegal resolution in Z (SigmaZ is negative)";
150  }
151 
152  produces<bool>();
153 
154 }
155 
157 {
158  delete fVertex ;
159  if (boost_ != 0 ) delete boost_;
160  delete fRandom;
161 }
162 
163 CLHEP::HepRandomEngine& BetaBoostEvtVtxGenerator::getEngine(){
164  return *fEngine;
165 }
166 
167 
168 //Hep3Vector* BetaBoostEvtVtxGenerator::newVertex() {
170 
171 
172  double X,Y,Z;
173 
174  double tmp_sigz = fRandom->fire(0., fSigmaZ);
175  Z = tmp_sigz + fZ0;
176 
177  double tmp_sigx = BetaFunction(Z,fZ0);
178  // need sqrt(2) for beamspot width relative to single beam width
179  tmp_sigx /= sqrt(2.0);
180  X = fRandom->fire(0.,tmp_sigx) + fX0; // + Z*fdxdz ;
181 
182  double tmp_sigy = BetaFunction(Z,fZ0);
183  // need sqrt(2) for beamspot width relative to single beam width
184  tmp_sigy /= sqrt(2.0);
185  Y = fRandom->fire(0.,tmp_sigy) + fY0; // + Z*fdydz;
186 
187  double tmp_sigt = fRandom->fire(0., fSigmaZ);
188  double T = tmp_sigt + fTimeOffset;
189 
190  if ( fVertex == 0 ) fVertex = new HepMC::FourVector();
191  fVertex->set(X,Y,Z,T);
192 
193  return fVertex;
194 }
195 
196 double BetaBoostEvtVtxGenerator::BetaFunction(double z, double z0)
197 {
198  return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
199 
200 }
201 
202 
204 {
205  if (s>=0 ) {
206  fSigmaZ=s;
207  }
208  else {
209  throw cms::Exception("LogicError")
210  << "Error in BetaBoostEvtVtxGenerator::sigmaZ: "
211  << "Illegal resolution in Z (negative)";
212  }
213 }
214 
216 
217  //alpha_ = 0;
218  //phi_ = 142.e-6;
219  // if (boost_ != 0 ) return boost_;
220 
221  //boost_.ResizeTo(4,4);
222  //boost_ = new TMatrixD(4,4);
223  TMatrixD tmpboost(4,4);
224  TMatrixD tmpboostZ(4,4);
225  TMatrixD tmpboostXYZ(4,4);
226 
227  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
228 
229  // Lorentz boost to frame where the collision is head-on
230  // phi is the half crossing angle in the plane ZS
231  // alpha is the angle to the S axis from the X axis in the XY plane
232 
233  tmpboost(0,0) = 1./cos(phi_);
234  tmpboost(0,1) = - cos(alpha_)*sin(phi_);
235  tmpboost(0,2) = - tan(phi_)*sin(phi_);
236  tmpboost(0,3) = - sin(alpha_)*sin(phi_);
237  tmpboost(1,0) = - cos(alpha_)*tan(phi_);
238  tmpboost(1,1) = 1.;
239  tmpboost(1,2) = cos(alpha_)*tan(phi_);
240  tmpboost(1,3) = 0.;
241  tmpboost(2,0) = 0.;
242  tmpboost(2,1) = - cos(alpha_)*sin(phi_);
243  tmpboost(2,2) = cos(phi_);
244  tmpboost(2,3) = - sin(alpha_)*sin(phi_);
245  tmpboost(3,0) = - sin(alpha_)*tan(phi_);
246  tmpboost(3,1) = 0.;
247  tmpboost(3,2) = sin(alpha_)*tan(phi_);
248  tmpboost(3,3) = 1.;
249  //cout<<"beta "<<beta_;
250  double gama=1.0/sqrt(1-beta_*beta_);
251  tmpboostZ(0,0)=gama;
252  tmpboostZ(0,1)=0.;
253  tmpboostZ(0,2)=-1.0*beta_*gama;
254  tmpboostZ(0,3)=0.;
255  tmpboostZ(1,0)=0.;
256  tmpboostZ(1,1) = 1.;
257  tmpboostZ(1,2)=0.;
258  tmpboostZ(1,3)=0.;
259  tmpboostZ(2,0)=-1.0*beta_*gama;
260  tmpboostZ(2,1) = 0.;
261  tmpboostZ(2,2)=gama;
262  tmpboostZ(2,3) = 0.;
263  tmpboostZ(3,0)=0.;
264  tmpboostZ(3,1)=0.;
265  tmpboostZ(3,2)=0.;
266  tmpboostZ(3,3) = 1.;
267 
268  tmpboostXYZ=tmpboostZ*tmpboost;
269  tmpboostXYZ.Invert();
270 
271 
272 
273  boost_ = new TMatrixD(tmpboostXYZ);
274  if ( verbosity_ ) { boost_->Print(); }
275 
276  return boost_;
277 }
278 
280 {
281 
282 
283  Handle<HepMCProduct> HepMCEvt ;
284  evt.getByLabel( sourceLabel, HepMCEvt ) ;
285 
286  // generate new vertex & apply the shift
287  //
288  HepMCEvt->applyVtxGen( newVertex() ) ;
289 
290  //HepMCEvt->LorentzBoost( 0., 142.e-6 );
291  HepMCEvt->boostToLab( GetInvLorentzBoost(), "vertex" );
292  HepMCEvt->boostToLab( GetInvLorentzBoost(), "momentum" );
293  // OK, create a (pseudo)product and put in into edm::Event
294  //
295  auto_ptr<bool> NewProduct(new bool(true)) ;
296  evt.put( NewProduct ) ;
297  return ;
298 }
299 
300 
301 
302 
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:49
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:94
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:47
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
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:361
double BetaFunction(double z, double z0)
beta function
CLHEP::HepRandomEngine * fEngine
void Alpha(double m=0)
angle between crossing plane and horizontal plane
void Z0(double m=0)
set mean in Z in cm
long double T
virtual TMatrixD * GetInvLorentzBoost()
virtual void produce(edm::Event &, const edm::EventSetup &)