CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/IOMC/ParticleGuns/src/MultiParticleInConeGunProducer.cc

Go to the documentation of this file.
00001 /*
00002  *  $Date: 2010/10/03 10:25:53 $
00003  *  $Revision: 1.6 $
00004  *  \author Jean-Roch Vlimant
00005  */
00006 
00007 #include <ostream>
00008 
00009 #include "IOMC/ParticleGuns/interface/MultiParticleInConeGunProducer.h"
00010 
00011 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00012 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00013 
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 
00017 // #include "FWCore/Utilities/interface/Exception.h"
00018 
00019 // #include "CLHEP/Random/RandFlat.h"
00020 
00021 using namespace edm;
00022 using namespace std;
00023 
00024 MultiParticleInConeGunProducer::MultiParticleInConeGunProducer(const ParameterSet& pset) : 
00025    BaseFlatGunProducer(pset)
00026 {
00027 
00028 
00029    ParameterSet defpset ;
00030    ParameterSet pgun_params = 
00031       pset.getParameter<ParameterSet>("PGunParameters") ;
00032   
00033    fMinPt = pgun_params.getParameter<double>("MinPt");
00034    fMaxPt = pgun_params.getParameter<double>("MaxPt");
00035 
00036    fInConeIds = pgun_params.getParameter< vector<int> >("InConeID");
00037    fMinDeltaR = pgun_params.getParameter<double>("MinDeltaR");
00038    fMaxDeltaR = pgun_params.getParameter<double>("MaxDeltaR");
00039    fMinMomRatio = pgun_params.getParameter<double>("MinMomRatio");
00040    fMaxMomRatio = pgun_params.getParameter<double>("MaxMomRatio");
00041 
00042    fInConeMinEta = pgun_params.getParameter<double>("InConeMinEta");
00043    fInConeMaxEta = pgun_params.getParameter<double>("InConeMaxEta");
00044    fInConeMinPhi = pgun_params.getParameter<double>("InConeMinPhi");
00045    fInConeMaxPhi = pgun_params.getParameter<double>("InConeMaxPhi");
00046    fInConeMaxTry = pgun_params.getParameter<unsigned int>("InConeMaxTry");
00047    
00048    produces<HepMCProduct>();
00049    produces<GenEventInfoProduct>();
00050 }
00051 
00052 MultiParticleInConeGunProducer::~MultiParticleInConeGunProducer()
00053 {
00054    // no need to cleanup GenEvent memory - done in HepMCProduct
00055 }
00056 
00057 void MultiParticleInConeGunProducer::produce(Event &e, const EventSetup& es) 
00058 {
00059 
00060    if ( fVerbosity > 0 )
00061    {
00062       cout << " MultiParticleInConeGunProducer : Begin New Event Generation" << endl ; 
00063    }
00064    // event loop (well, another step in it...)
00065           
00066    // no need to clean up GenEvent memory - done in HepMCProduct
00067    // 
00068    
00069    // here re-create fEvt (memory)
00070    //
00071    fEvt = new HepMC::GenEvent() ;
00072    
00073    // now actualy, cook up the event from PDGTable and gun parameters
00074    //
00075    // 1st, primary vertex
00076    //
00077    //HepMC::GenVertex* Vtx = new HepMC::GenVertex(CLHEP::HepLorentzVector(0.,0.,0.));
00078    HepMC::GenVertex* Vtx = new HepMC::GenVertex(HepMC::FourVector(0.,0.,0.));
00079 
00080    // loop over particles
00081    //
00082    int barcode = 1 ;
00083    for (unsigned int ip=0; ip<fPartIDs.size(); ++ip)
00084    {
00085 
00086        double pt     = fRandomGenerator->fire(fMinPt, fMaxPt) ;
00087        double eta    = fRandomGenerator->fire(fMinEta, fMaxEta) ;
00088        double phi    = fRandomGenerator->fire(fMinPhi, fMaxPhi) ;
00089        int PartID = fPartIDs[ip] ;
00090        const HepPDT::ParticleData* 
00091           PData = fPDGTable->particle(HepPDT::ParticleID(abs(PartID))) ;
00092        double mass   = PData->mass().value() ;
00093        double theta  = 2.*atan(exp(-eta)) ;
00094        double mom    = pt/sin(theta) ;
00095        double px     = pt*cos(phi) ;
00096        double py     = pt*sin(phi) ;
00097        double pz     = mom*cos(theta) ;
00098        double energy2= mom*mom + mass*mass ;
00099        double energy = sqrt(energy2) ; 
00100 
00101        HepMC::FourVector p(px,py,pz,energy) ;
00102        HepMC::GenParticle* Part = 
00103            new HepMC::GenParticle(p,PartID,1);
00104        Part->suggest_barcode( barcode ) ;
00105        barcode++ ;
00106        Vtx->add_particle_out(Part);
00107 
00108        if ( fAddAntiParticle ){}
00109 
00110        // now add the particles in the cone
00111        for (unsigned iPic=0; iPic!=fInConeIds.size();iPic++){
00112          unsigned int nTry=0;
00113          while(true){
00114            //shoot flat Deltar
00115            double dR = fRandomGenerator->fire(fMinDeltaR, fMaxDeltaR);
00116            //shoot flat eta/phi mixing
00117            double alpha = fRandomGenerator->fire(-3.14159265358979323846, 3.14159265358979323846);
00118            double dEta = dR*cos(alpha);
00119            double dPhi = dR*sin(alpha);
00120            
00121            /*
00122            //shoot Energy of associated particle         
00123            double energyIc = fRandomGenerator->fire(fMinEInCone, fMaxEInCone);
00124            if (mom2Ic>0){ momIC = sqrt(mom2Ic);}
00125            */
00126            //    get kinematics
00127            double etaIc = eta+dEta;
00128            double phiIc = phi+dPhi;
00129            //put it back in -Pi:Pi if necessary. multiple time might be necessary if dR > 3
00130            const unsigned int maxL=100;
00131            unsigned int iL=0;
00132            while(iL++<maxL){
00133              if (phiIc > 3.14159265358979323846) phiIc-=2*3.14159265358979323846;
00134              else if(phiIc <-3.14159265358979323846) phiIc+=2*3.14159265358979323846;
00135              
00136              if (abs(phiIc)<3.14159265358979323846) break;
00137            }
00138              
00139 
00140            //allow to skip it if you did not run out of possible drawing
00141            if (nTry++<=fInConeMaxTry){
00142              //draw another one if this one is not in acceptance
00143              if (etaIc<fInConeMinEta || etaIc > fInConeMaxEta) continue;
00144              if (phiIc<fInConeMinPhi || phiIc > fInConeMaxPhi) continue;
00145            }
00146            else{
00147              if ( fVerbosity > 0 )
00148                {
00149                  cout << " MultiParticleInConeGunProducer : could not produce a particle "
00150                       <<  fInConeIds[iPic]<<" in cone "
00151                       <<  fMinDeltaR<<" to "<<fMaxDeltaR<<" within the "<<fInConeMaxTry<<" allowed tries."<<endl;
00152                }
00153              /*      edm::LogError("MultiParticleInConeGunProducer")<< " MultiParticleInConeGunProducer : could not produce a particle "<<
00154                fInConeIds[iPic]<<" in cone "<<
00155                fMinDeltaR<<" to "<<fMaxDeltaR<<" within the "<<fInConeMaxTry<<" allowed tries.";*/
00156            }
00157            int PartIDIc=fInConeIds[iPic];
00158            const HepPDT::ParticleData* 
00159              PDataIc = fPDGTable->particle(HepPDT::ParticleID(abs(PartIDIc)));
00160            
00161            //shoot momentum ratio
00162            double momR = fRandomGenerator->fire(fMinMomRatio, fMaxMomRatio);
00163            double massIc= PDataIc->mass().value() ;
00164            double momIc = momR * mom;
00165            double energyIc = sqrt(momIc*momIc + massIc*massIc);
00166 
00167            double thetaIc = 2.*atan(exp(-etaIc)) ;
00168            double pxIc = momIc*sin(thetaIc)*cos(phiIc);
00169            double pyIc = momIc*sin(thetaIc)*sin(phiIc);
00170            double pzIc = momIc*cos(thetaIc);
00171 
00172            HepMC::FourVector pIc(pxIc,pyIc,pzIc,energyIc) ;
00173            HepMC::GenParticle* PartIc = new HepMC::GenParticle(pIc, PartIDIc, 1);
00174            PartIc->suggest_barcode( barcode ) ;
00175            barcode++ ;
00176            Vtx->add_particle_out(PartIc);
00177            break;
00178          }//try many times while not in acceptance
00179        }//loop over the particle Ids in the cone
00180    }
00181 
00182    fEvt->add_vertex(Vtx) ;
00183    fEvt->set_event_number(e.id().event()) ;
00184    fEvt->set_signal_process_id(20) ; 
00185         
00186    if ( fVerbosity > 0 )
00187    {
00188       fEvt->print() ;  
00189    }
00190 
00191    auto_ptr<HepMCProduct> BProduct(new HepMCProduct()) ;
00192    BProduct->addHepMCData( fEvt );
00193    e.put(BProduct);
00194 
00195    auto_ptr<GenEventInfoProduct> genEventInfo(new GenEventInfoProduct(fEvt));
00196    e.put(genEventInfo);
00197 
00198    if ( fVerbosity > 0 )
00199      {
00200        cout << " MultiParticleInConeGunProducer : Event Generation Done " << endl;
00201      }
00202 }
00203