CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/GeneratorInterface/HiGenCommon/plugins/BetaBoostEvtVtxGenerator.cc

Go to the documentation of this file.
00001 
00002 // $Id: BetaBoostEvtVtxGenerator.cc,v 1.2 2012/09/06 20:41:40 lixu Exp $
00003 /*
00004 ________________________________________________________________________
00005 
00006  BetaBoostEvtVtxGenerator
00007 
00008  Smear vertex according to the Beta function on the transverse plane
00009  and a Gaussian on the z axis. It allows the beam to have a crossing
00010  angle (slopes dxdz and dydz).
00011 
00012  Based on GaussEvtVtxGenerator
00013  implemented by Francisco Yumiceva (yumiceva@fnal.gov)
00014 
00015  FERMILAB
00016  2006
00017 ________________________________________________________________________
00018 */
00019 
00020 //lingshan: add beta for z-axis boost
00021 
00022 //#include "IOMC/EventVertexGenerators/interface/BetafuncEvtVtxGenerator.h"
00023 
00024 #include "FWCore/Framework/interface/Event.h"
00025 #include "FWCore/Utilities/interface/Exception.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00027 
00028 #include "FWCore/Framework/interface/EDProducer.h"
00029 #include "FWCore/Utilities/interface/InputTag.h"
00030 #include "FWCore/ServiceRegistry/interface/Service.h"
00031 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00032 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00033 
00034 #include "CLHEP/Random/RandGaussQ.h"
00035 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00036 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00037 //#include "CLHEP/Vector/ThreeVector.h"
00038 #include "HepMC/SimpleVector.h"
00039 #include "TMatrixD.h"
00040 
00041 #include <iostream>
00042 
00043 using namespace edm;
00044 using namespace std;
00045 using namespace CLHEP;
00046 
00047 class RandGaussQ;
00048 
00049 
00050 class BetaBoostEvtVtxGenerator : public edm::EDProducer{
00051 public:
00052   BetaBoostEvtVtxGenerator(const edm::ParameterSet & p);
00053   virtual ~BetaBoostEvtVtxGenerator();
00054 
00056   //virtual CLHEP::Hep3Vector * newVertex();
00057   virtual HepMC::FourVector* newVertex() ;
00058   virtual void produce( edm::Event&, const edm::EventSetup& );
00059   virtual TMatrixD* GetInvLorentzBoost();
00060 
00061 
00063   void sigmaZ(double s=1.0);
00064 
00066   void X0(double m=0) { fX0=m; }
00068   void Y0(double m=0) { fY0=m; }
00070   void Z0(double m=0) { fZ0=m; }
00071 
00073   void Phi(double m=0) { phi_=m; }
00075   void Alpha(double m=0) { alpha_=m; }
00076   void Beta(double m=0) { beta_=m; }
00077 
00079   void betastar(double m=0) { fbetastar=m; }
00081   void emittance(double m=0) { femittance=m; }
00082 
00084   double BetaFunction(double z, double z0);
00085   CLHEP::HepRandomEngine& getEngine();
00086 
00087 private:
00089   BetaBoostEvtVtxGenerator(const BetaBoostEvtVtxGenerator &p);
00091   BetaBoostEvtVtxGenerator&  operator = (const BetaBoostEvtVtxGenerator & rhs );
00092 
00093   double alpha_, phi_;
00094   //TMatrixD boost_;
00095   double beta_;
00096   double fX0, fY0, fZ0;
00097   double fSigmaZ;
00098   //double fdxdz, fdydz;
00099   double fbetastar, femittance;
00100   double falpha;
00101 
00102   HepMC::FourVector*       fVertex ;
00103   TMatrixD *boost_;
00104   double fTimeOffset;
00105   
00106   CLHEP::HepRandomEngine*  fEngine;
00107   edm::InputTag            sourceLabel;
00108 
00109   CLHEP::RandGaussQ*  fRandom ;
00110 
00111 };
00112 
00113 
00114 BetaBoostEvtVtxGenerator::BetaBoostEvtVtxGenerator(const edm::ParameterSet & p ):
00115   fVertex(0), boost_(0), fTimeOffset(0), fEngine(0),
00116   sourceLabel(p.getParameter<edm::InputTag>("src"))
00117 { 
00118   
00119   edm::Service<edm::RandomNumberGenerator> rng;
00120 
00121   if ( ! rng.isAvailable()) {
00122     
00123     throw cms::Exception("Configuration")
00124       << "The BaseEvtVtxGenerator requires the RandomNumberGeneratorService\n"
00125       "which is not present in the configuration file.  You must add the service\n"
00126       "in the configuration file or remove the modules that require it.";
00127   }
00128 
00129   CLHEP::HepRandomEngine& engine = rng->getEngine();
00130   fEngine = &engine;
00131   fRandom = new CLHEP::RandGaussQ(getEngine());
00132 
00133   fX0 =        p.getParameter<double>("X0")*cm;
00134   fY0 =        p.getParameter<double>("Y0")*cm;
00135   fZ0 =        p.getParameter<double>("Z0")*cm;
00136   fSigmaZ =    p.getParameter<double>("SigmaZ")*cm;
00137   alpha_ =     p.getParameter<double>("Alpha")*radian;
00138   phi_ =       p.getParameter<double>("Phi")*radian;
00139   fbetastar =  p.getParameter<double>("BetaStar")*cm;
00140   femittance = p.getParameter<double>("Emittance")*cm; // this is not the normalized emittance
00141   fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light; // HepMC time units are mm
00142   beta_=p.getParameter<double>("Beta"); 
00143   if (fSigmaZ <= 0) {
00144           throw cms::Exception("Configuration")
00145                   << "Error in BetaBoostEvtVtxGenerator: "
00146                   << "Illegal resolution in Z (SigmaZ is negative)";
00147   }
00148 
00149   produces<bool>(); 
00150   
00151 }
00152 
00153 BetaBoostEvtVtxGenerator::~BetaBoostEvtVtxGenerator() 
00154 {
00155   delete fVertex ;
00156   if (boost_ != 0 ) delete boost_;
00157   delete fRandom; 
00158 }
00159 
00160 CLHEP::HepRandomEngine& BetaBoostEvtVtxGenerator::getEngine(){
00161   return *fEngine;
00162 }
00163 
00164 
00165 //Hep3Vector* BetaBoostEvtVtxGenerator::newVertex() {
00166 HepMC::FourVector* BetaBoostEvtVtxGenerator::newVertex() {
00167 
00168         
00169         double X,Y,Z;
00170         
00171         double tmp_sigz = fRandom->fire(0., fSigmaZ);
00172         Z = tmp_sigz + fZ0;
00173 
00174         double tmp_sigx = BetaFunction(Z,fZ0); 
00175         // need sqrt(2) for beamspot width relative to single beam width
00176         tmp_sigx /= sqrt(2.0);
00177         X = fRandom->fire(0.,tmp_sigx) + fX0; // + Z*fdxdz ;
00178 
00179         double tmp_sigy = BetaFunction(Z,fZ0);
00180         // need sqrt(2) for beamspot width relative to single beam width
00181         tmp_sigy /= sqrt(2.0);
00182         Y = fRandom->fire(0.,tmp_sigy) + fY0; // + Z*fdydz;
00183 
00184         double tmp_sigt = fRandom->fire(0., fSigmaZ);
00185         double T = tmp_sigt + fTimeOffset; 
00186 
00187         if ( fVertex == 0 ) fVertex = new HepMC::FourVector();
00188         fVertex->set(X,Y,Z,T);
00189                 
00190         return fVertex;
00191 }
00192 
00193 double BetaBoostEvtVtxGenerator::BetaFunction(double z, double z0)
00194 {
00195         return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
00196 
00197 }
00198 
00199 
00200 void BetaBoostEvtVtxGenerator::sigmaZ(double s) 
00201 { 
00202         if (s>=0 ) {
00203                 fSigmaZ=s; 
00204         }
00205         else {
00206                 throw cms::Exception("LogicError")
00207                         << "Error in BetaBoostEvtVtxGenerator::sigmaZ: "
00208                         << "Illegal resolution in Z (negative)";
00209         }
00210 }
00211 
00212 TMatrixD* BetaBoostEvtVtxGenerator::GetInvLorentzBoost() {
00213 
00214         //alpha_ = 0;
00215         //phi_ = 142.e-6;
00216 //      if (boost_ != 0 ) return boost_;
00217         
00218         //boost_.ResizeTo(4,4);
00219         //boost_ = new TMatrixD(4,4);
00220         TMatrixD tmpboost(4,4);
00221         TMatrixD tmpboostZ(4,4);
00222         TMatrixD tmpboostXYZ(4,4);
00223 
00224         //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
00225         
00226         // Lorentz boost to frame where the collision is head-on
00227         // phi is the half crossing angle in the plane ZS
00228         // alpha is the angle to the S axis from the X axis in the XY plane
00229         
00230         tmpboost(0,0) = 1./cos(phi_);
00231         tmpboost(0,1) = - cos(alpha_)*sin(phi_);
00232         tmpboost(0,2) = - tan(phi_)*sin(phi_);
00233         tmpboost(0,3) = - sin(alpha_)*sin(phi_);
00234         tmpboost(1,0) = - cos(alpha_)*tan(phi_);
00235         tmpboost(1,1) = 1.;
00236         tmpboost(1,2) = cos(alpha_)*tan(phi_);
00237         tmpboost(1,3) = 0.;
00238         tmpboost(2,0) = 0.;
00239         tmpboost(2,1) = - cos(alpha_)*sin(phi_);
00240         tmpboost(2,2) = cos(phi_);
00241         tmpboost(2,3) = - sin(alpha_)*sin(phi_);
00242         tmpboost(3,0) = - sin(alpha_)*tan(phi_);
00243         tmpboost(3,1) = 0.;
00244         tmpboost(3,2) = sin(alpha_)*tan(phi_);
00245         tmpboost(3,3) = 1.;
00246        //cout<<"beta "<<beta_;
00247        double gama=1.0/sqrt(1-beta_*beta_);
00248        tmpboostZ(0,0)=gama;
00249        tmpboostZ(0,1)=0.;
00250        tmpboostZ(0,2)=-1.0*beta_*gama;
00251        tmpboostZ(0,3)=0.;
00252        tmpboostZ(1,0)=0.;
00253        tmpboostZ(1,1) = 1.;
00254        tmpboostZ(1,2)=0.;
00255        tmpboostZ(1,3)=0.;
00256        tmpboostZ(2,0)=-1.0*beta_*gama;
00257        tmpboostZ(2,1) = 0.;
00258        tmpboostZ(2,2)=gama;
00259        tmpboostZ(2,3) = 0.;
00260        tmpboostZ(3,0)=0.;
00261        tmpboostZ(3,1)=0.;
00262        tmpboostZ(3,2)=0.;
00263        tmpboostZ(3,3) = 1.;
00264 
00265        tmpboostXYZ=tmpboostZ*tmpboost;
00266        tmpboostXYZ.Invert();
00267 
00268 
00269 
00270        boost_ = new TMatrixD(tmpboostXYZ);
00271        boost_->Print();
00272         
00273         return boost_;
00274 }
00275 
00276 void BetaBoostEvtVtxGenerator::produce( Event& evt, const EventSetup& )
00277 {
00278   
00279   
00280   Handle<HepMCProduct> HepMCEvt ;    
00281   evt.getByLabel( sourceLabel, HepMCEvt ) ;
00282   
00283   // generate new vertex & apply the shift 
00284   //
00285   HepMCEvt->applyVtxGen( newVertex() ) ;
00286  
00287   //HepMCEvt->LorentzBoost( 0., 142.e-6 );
00288   HepMCEvt->boostToLab( GetInvLorentzBoost(), "vertex" );
00289  HepMCEvt->boostToLab( GetInvLorentzBoost(), "momentum" );    
00290   // OK, create a (pseudo)product and put in into edm::Event
00291   //
00292   auto_ptr<bool> NewProduct(new bool(true)) ;      
00293  evt.put( NewProduct ) ;       
00294   return ;
00295 }
00296 
00297 
00298 
00299 
00300 #include "FWCore/Framework/interface/MakerMacros.h"
00301 DEFINE_FWK_MODULE(BetaBoostEvtVtxGenerator);