CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/IOMC/EventVertexGenerators/src/BetafuncEvtVtxGenerator.cc

Go to the documentation of this file.
00001 
00002 // $Id: BetafuncEvtVtxGenerator.cc,v 1.14 2012/01/18 16:22:25 vlimant Exp $
00003 /*
00004 ________________________________________________________________________
00005 
00006  BetafuncEvtVtxGenerator
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 
00021 
00022 #include "IOMC/EventVertexGenerators/interface/BetafuncEvtVtxGenerator.h"
00023 #include "FWCore/Utilities/interface/Exception.h"
00024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00025 #include "FWCore/Framework/interface/EventSetup.h"
00026 #include "FWCore/Framework/interface/ESHandle.h"
00027 
00028 #include "CLHEP/Random/RandGaussQ.h"
00029 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00030 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00031 //#include "CLHEP/Vector/ThreeVector.h"
00032 #include "HepMC/SimpleVector.h"
00033 
00034 #include "CondFormats/DataRecord/interface/SimBeamSpotObjectsRcd.h"
00035 #include "CondFormats/BeamSpotObjects/interface/SimBeamSpotObjects.h"
00036 
00037 #include <iostream>
00038 
00039 
00040 
00041 BetafuncEvtVtxGenerator::BetafuncEvtVtxGenerator(const edm::ParameterSet & p )
00042 : BaseEvtVtxGenerator(p)
00043 { 
00044   
00045   fRandom = new CLHEP::RandGaussQ(getEngine());
00046 
00047   readDB_=p.getParameter<bool>("readDB");
00048   if (!readDB_){
00049     fX0 =        p.getParameter<double>("X0")*cm;
00050     fY0 =        p.getParameter<double>("Y0")*cm;
00051     fZ0 =        p.getParameter<double>("Z0")*cm;
00052     fSigmaZ =    p.getParameter<double>("SigmaZ")*cm;
00053     alpha_ =     p.getParameter<double>("Alpha")*radian;
00054     phi_ =       p.getParameter<double>("Phi")*radian;
00055     fbetastar =  p.getParameter<double>("BetaStar")*cm;
00056     femittance = p.getParameter<double>("Emittance")*cm; // this is not the normalized emittance
00057     fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light; // HepMC time units are mm
00058     
00059     if (fSigmaZ <= 0) {
00060       throw cms::Exception("Configuration")
00061         << "Error in BetafuncEvtVtxGenerator: "
00062         << "Illegal resolution in Z (SigmaZ is negative)";
00063     }
00064   }
00065 
00066   
00067 }
00068 
00069 BetafuncEvtVtxGenerator::~BetafuncEvtVtxGenerator() 
00070 {
00071     delete fRandom; 
00072 }
00073 
00074 void BetafuncEvtVtxGenerator::beginLuminosityBlock(edm::LuminosityBlock&, edm::EventSetup const& iEventSetup){
00075   update(iEventSetup);
00076 }
00077 void BetafuncEvtVtxGenerator::beginRun( edm::Run & , const edm::EventSetup& iEventSetup){
00078   update(iEventSetup);
00079 }
00080 
00081 void BetafuncEvtVtxGenerator::update(const edm::EventSetup& iEventSetup){
00082   if (readDB_ &&  parameterWatcher_.check(iEventSetup)){
00083     edm::ESHandle< SimBeamSpotObjects > beamhandle;
00084     iEventSetup.get<SimBeamSpotObjectsRcd>().get(beamhandle);
00085 
00086     fX0=beamhandle->fX0;
00087     fY0=beamhandle->fY0;
00088     fZ0=beamhandle->fZ0;
00089     //    falpha=beamhandle->fAlpha;
00090     alpha_=beamhandle->fAlpha;
00091     phi_=beamhandle->fPhi;
00092     fSigmaZ=beamhandle->fSigmaZ;
00093     fTimeOffset=beamhandle->fTimeOffset;
00094     fbetastar=beamhandle->fbetastar;
00095     femittance=beamhandle->femittance;
00096 
00097     //re-initialize the boost matrix
00098     delete boost_;
00099     boost_=0;
00100   }
00101 }
00102 
00103 //Hep3Vector* BetafuncEvtVtxGenerator::newVertex() {
00104 HepMC::FourVector* BetafuncEvtVtxGenerator::newVertex() {
00105 
00106         
00107         double X,Y,Z;
00108         
00109         double tmp_sigz = fRandom->fire(0., fSigmaZ);
00110         Z = tmp_sigz + fZ0;
00111 
00112         double tmp_sigx = BetaFunction(Z,fZ0); 
00113         // need sqrt(2) for beamspot width relative to single beam width
00114         tmp_sigx /= sqrt(2.0);
00115         X = fRandom->fire(0.,tmp_sigx) + fX0; // + Z*fdxdz ;
00116 
00117         double tmp_sigy = BetaFunction(Z,fZ0);
00118         // need sqrt(2) for beamspot width relative to single beam width
00119         tmp_sigy /= sqrt(2.0);
00120         Y = fRandom->fire(0.,tmp_sigy) + fY0; // + Z*fdydz;
00121 
00122         double tmp_sigt = fRandom->fire(0., fSigmaZ);
00123         double T = tmp_sigt + fTimeOffset; 
00124 
00125         if ( fVertex == 0 ) fVertex = new HepMC::FourVector();
00126         fVertex->set(X,Y,Z,T);
00127                 
00128         return fVertex;
00129 }
00130 
00131 double BetafuncEvtVtxGenerator::BetaFunction(double z, double z0)
00132 {
00133         return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
00134 
00135 }
00136 
00137 
00138 void BetafuncEvtVtxGenerator::sigmaZ(double s) 
00139 { 
00140         if (s>=0 ) {
00141                 fSigmaZ=s; 
00142         }
00143         else {
00144                 throw cms::Exception("LogicError")
00145                         << "Error in BetafuncEvtVtxGenerator::sigmaZ: "
00146                         << "Illegal resolution in Z (negative)";
00147         }
00148 }
00149 
00150 TMatrixD* BetafuncEvtVtxGenerator::GetInvLorentzBoost() {
00151 
00152         //alpha_ = 0;
00153         //phi_ = 142.e-6;
00154         
00155         if (boost_ != 0 ) return boost_;
00156         
00157         //boost_.ResizeTo(4,4);
00158         //boost_ = new TMatrixD(4,4);
00159         TMatrixD tmpboost(4,4);
00160         
00161         //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
00162         
00163         // Lorentz boost to frame where the collision is head-on
00164         // phi is the half crossing angle in the plane ZS
00165         // alpha is the angle to the S axis from the X axis in the XY plane
00166         
00167         tmpboost(0,0) = 1./cos(phi_);
00168         tmpboost(0,1) = - cos(alpha_)*sin(phi_);
00169         tmpboost(0,2) = - tan(phi_)*sin(phi_);
00170         tmpboost(0,3) = - sin(alpha_)*sin(phi_);
00171         tmpboost(1,0) = - cos(alpha_)*tan(phi_);
00172         tmpboost(1,1) = 1.;
00173         tmpboost(1,2) = cos(alpha_)*tan(phi_);
00174         tmpboost(1,3) = 0.;
00175         tmpboost(2,0) = 0.;
00176         tmpboost(2,1) = - cos(alpha_)*sin(phi_);
00177         tmpboost(2,2) = cos(phi_);
00178         tmpboost(2,3) = - sin(alpha_)*sin(phi_);
00179         tmpboost(3,0) = - sin(alpha_)*tan(phi_);
00180         tmpboost(3,1) = 0.;
00181         tmpboost(3,2) = sin(alpha_)*tan(phi_);
00182         tmpboost(3,3) = 1.;
00183 
00184         tmpboost.Invert();
00185         boost_ = new TMatrixD(tmpboost);
00186         //boost_->Print();
00187         
00188         return boost_;
00189 }