CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/SimG4Core/Generators/src/Generator.cc

Go to the documentation of this file.
00001 
00002 #include "SimG4Core/Generators/interface/Generator.h"
00003 #include "SimG4Core/Generators/interface/HepMCParticle.h"
00004 
00005 // #include "FWCore/Utilities/interface/Exception.h"
00006 #include "SimG4Core/Notification/interface/SimG4Exception.h"
00007 
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 
00010 #include "G4Event.hh"
00011 #include "G4EventManager.hh"
00012 #include "G4HEPEvtParticle.hh"
00013 #include "G4ParticleDefinition.hh"
00014 #include "G4UnitsTable.hh"
00015 
00016 using namespace edm;
00017 using std::cout;
00018 using std::endl;
00019 using std::string;
00020 
00021 
00022 Generator::Generator(const ParameterSet & p) : 
00023   fPCuts(p.getParameter<bool>("ApplyPCuts")),
00024   fEtaCuts(p.getParameter<bool>("ApplyEtaCuts")), 
00025   fPhiCuts(p.getParameter<bool>("ApplyPhiCuts")),
00026   theMinPhiCut(p.getParameter<double>("MinPhiCut")),   // now operates in radians (CMS standard)
00027   theMaxPhiCut(p.getParameter<double>("MaxPhiCut")),
00028   theMinEtaCut(p.getParameter<double>("MinEtaCut")),
00029   theMaxEtaCut(p.getParameter<double>("MaxEtaCut")),
00030   theMinPCut(p.getParameter<double>("MinPCut")),    // now operates in GeV (CMS standard)
00031   theMaxPCut(p.getParameter<double>("MaxPCut")),   
00032   theRDecLenCut(p.getParameter<double>("RDecLenCut")*cm),
00033   theEtaCutForHector(p.getParameter<double>("EtaCutForHector")),
00034   verbose(p.getUntrackedParameter<int>("Verbosity",0)),
00035   evt_(0),
00036   vtx_(0),
00037   weight_(0),
00038   Z_lmin(0),
00039   Z_lmax(0),
00040   Z_hector(0)
00041 {
00042 
00043   if(fEtaCuts){
00044     Z_lmax = theRDecLenCut*( ( 1 - exp(-2*theMaxEtaCut) ) / ( 2*exp(-theMaxEtaCut) ) );
00045     Z_lmin = theRDecLenCut*( ( 1 - exp(-2*theMinEtaCut) ) / ( 2*exp(-theMinEtaCut) ) );
00046   }
00047 
00048   Z_hector = theRDecLenCut*( ( 1 - exp(-2*theEtaCutForHector) ) / ( 2*exp(-theEtaCutForHector) ) );
00049 
00050   if ( verbose > 0 ) LogDebug("SimG4CoreGenerator") << "Z_min = " << Z_lmin << " Z_max = " << Z_lmax << " Z_hector = " << Z_hector;
00051 
00052 }
00053 
00054 Generator::~Generator() 
00055 { 
00056 }
00057 
00058 void Generator::HepMC2G4(const HepMC::GenEvent * evt_orig, G4Event * g4evt)
00059 {
00060 
00061   if ( *(evt_orig->vertices_begin()) == 0 )
00062     {
00063       throw SimG4Exception( "SimG4CoreGenerator: Corrupted Event - GenEvent with no vertex" ) ;
00064     }
00065   
00066   
00067   HepMC::GenEvent* evt = new HepMC::GenEvent(*evt_orig) ;
00068   
00069   if ( evt->weights().size() > 0 )
00070     {
00071       weight_ = evt->weights()[0] ;
00072       for ( unsigned int iw=1; iw<evt->weights().size(); iw++ )
00073         {
00074           // terminate if the versot of weights contains a zero-weight
00075           if ( evt->weights()[iw] <= 0 ) break;
00076           weight_ *= evt->weights()[iw] ;
00077         }     
00078     }
00079   
00080   if (vtx_ != 0) delete vtx_;
00081   vtx_ = new math::XYZTLorentzVector((*(evt->vertices_begin()))->position().x(),
00082                                      (*(evt->vertices_begin()))->position().y(),
00083                                      (*(evt->vertices_begin()))->position().z(),
00084                                      (*(evt->vertices_begin()))->position().t());
00085   
00086   if(verbose >0){
00087     evt->print();
00088     LogDebug("SimG4CoreGenerator") << "Primary Vertex = (" << vtx_->x() << "," 
00089                                    << vtx_->y() << ","
00090                                    << vtx_->z() << ")";
00091   }
00092   
00093   //  double x0 = vtx_->x();
00094   //  double y0 = vtx_->y();
00095   
00096   unsigned int ng4vtx = 0;
00097 
00098   for(HepMC::GenEvent::vertex_const_iterator vitr= evt->vertices_begin();
00099       vitr != evt->vertices_end(); ++vitr ) { 
00100     // loop for vertex ...
00101     // real vertex?
00102     G4bool qvtx=false;
00103     
00104     for (HepMC::GenVertex::particle_iterator pitr= (*vitr)->particles_begin(HepMC::children);
00105          pitr != (*vitr)->particles_end(HepMC::children); ++pitr) {
00106       // Admit also status=1 && end_vertex for long vertex special decay treatment 
00107       if ((*pitr)->status()==1) {
00108         qvtx=true;
00109         if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenVertex barcode = " << (*vitr)->barcode() 
00110                                      << " selected for GenParticle barcode = " << (*pitr)->barcode() << std::endl;
00111         break;
00112       }  
00113       // The selection is made considering if the partcile with status = 2 have the end_vertex
00114       // with a radius (R) greater then the theRDecLenCut that means: the end_vertex is outside
00115       // the beampipe cilinder (no requirement on the Z of the vertex is applyed).
00116       else if ( (*pitr)->status()== 2 ) {
00117         if ( (*pitr)->end_vertex() != 0  ) { 
00118           //double xx = x0-(*pitr)->end_vertex()->position().x();
00119           //double yy = y0-(*pitr)->end_vertex()->position().y();
00120           double xx = (*pitr)->end_vertex()->position().x();
00121           double yy = (*pitr)->end_vertex()->position().y();
00122           double r_dd=std::sqrt(xx*xx+yy*yy);
00123           if (r_dd>theRDecLenCut){
00124             qvtx=true;
00125             if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenVertex barcode = " << (*vitr)->barcode() 
00126                                          << " selected for GenParticle barcode = " << (*pitr)->barcode() << " radius = " << r_dd << std::endl;
00127             break;
00128           }
00129         }
00130       }
00131     }
00132 
00133 
00134     if (!qvtx) {
00135       continue;
00136     }
00137     
00138     double x1 = (*vitr)->position().x();
00139     double y1 = (*vitr)->position().y();
00140     double z1 = (*vitr)->position().z();
00141     double t1 = (*vitr)->position().t();        
00142     G4PrimaryVertex* g4vtx= 
00143       new G4PrimaryVertex(x1*mm, y1*mm, z1*mm, t1*mm/c_light);
00144     
00145     for (HepMC::GenVertex::particle_iterator vpitr= (*vitr)->particles_begin(HepMC::children);
00146          vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr){
00147 
00148       // Special cases:
00149       // 1) import in Geant4 a full decay chain (e.g. also particles with status == 2) 
00150       //    from the generator in case the decay radius is larger than theRDecLenCut
00151       //    In this case no cuts will be applied to select the particles hat has to be 
00152       //    processed by geant
00153       // 2) import in Geant4 particles with status == 1 but a final end vertex. 
00154       //    The time of the vertex is used as the time of flight to be forced for the particle 
00155       
00156       double r_decay_length=-1;
00157       double decay_length=-1;
00158 
00159       if ( (*vpitr)->status() == 1 || (*vpitr)->status() == 2 ) {
00160         // this particle has decayed
00161         if ( (*vpitr)->end_vertex() != 0 ) { 
00162           // needed some particles have status 2 and no end_vertex 
00163           // Which are the particles with status 2 and not end_vertex, what are suppose to di such kind of particles
00164           // when propagated to geant?
00165           
00166           double x2 = (*vpitr)->end_vertex()->position().x();
00167           double y2 = (*vpitr)->end_vertex()->position().y();
00168           double z2 = (*vpitr)->end_vertex()->position().z();
00169           r_decay_length=std::sqrt(x2*x2+y2*y2);
00170           decay_length=std::sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2));
00171 
00172         }
00173       } 
00174       // end modification
00175 
00176       bool toBeAdded = false;
00177       math::XYZTLorentzVector p((*vpitr)->momentum().px(),
00178                                 (*vpitr)->momentum().py(),
00179                                 (*vpitr)->momentum().pz(),
00180                                 (*vpitr)->momentum().e());
00181 
00182       // protection against numerical problems for extremely low momenta
00183       const double minTan = 1.e-20;
00184       double zimpact = Z_hector+1.;
00185       if ( fabs(z1) < Z_hector && fabs(p.pt()/p.pz()) >= minTan ) { 
00186         // write tan(p.Theta()) as p.Pt()/p.Pz()
00187          zimpact = (theRDecLenCut-sqrt(x1*x1+y1*y1))*(p.pz()/p.pt())+z1;
00188       }
00189 
00190       if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "Processing GenParticle barcode = " << (*vpitr)->barcode() 
00191                                                         << " status = " << (*vpitr)->status() 
00192                                                         << " zimpact = " << zimpact;
00193 
00194       // Standard case: particles not decayed by the generator
00195       if( (*vpitr)->status() == 1 && fabs(zimpact) < Z_hector ) {
00196         if ( !particlePassesPrimaryCuts( p, zimpact ) ) {
00197           continue ;
00198         }
00199         toBeAdded = true;
00200         if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*vpitr)->barcode() << " passed case 1" << std::endl;
00201       }
00202       
00203       // Decay chain entering exiting the fiducial cylinder defined by theRDecLenCut
00204       else if((*vpitr)->status() == 2 && r_decay_length > theRDecLenCut  && 
00205          fabs(zimpact) < Z_hector ) {
00206         toBeAdded=true;
00207         if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*vpitr)->barcode() << " passed case 2" << std::endl;
00208       }
00209       
00210       // Particles trasnported along the beam pipe for forward detectors (HECTOR)
00211       // Always pass to Geant4 without cuts (to be checked)
00212       else if( (*vpitr)->status() == 1 && fabs(zimpact) >= Z_hector && fabs(z1) >= Z_hector) {
00213         toBeAdded = true;
00214         if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "GenParticle barcode = " << (*vpitr)->barcode() << " passed case 3" << std::endl;
00215       }
00216 
00217       if(toBeAdded){
00218         
00219         G4int pdgcode= (*vpitr)-> pdg_id();
00220         G4PrimaryParticle* g4prim= 
00221           new G4PrimaryParticle(pdgcode, p.Px()*GeV, p.Py()*GeV, p.Pz()*GeV);
00222         
00223         if ( g4prim->GetG4code() != 0 ){ 
00224           g4prim->SetMass( g4prim->GetG4code()->GetPDGMass() ) ;
00225           g4prim->SetCharge( g4prim->GetG4code()->GetPDGCharge() ) ;  
00226         }
00227         
00228         g4prim->SetWeight( 10000*(*vpitr)->barcode() ) ;
00229         setGenId( g4prim, (*vpitr)->barcode() ) ;
00230 
00231         if ( (*vpitr)->status() == 2 ) 
00232           particleAssignDaughters(g4prim,(HepMC::GenParticle *) *vpitr, decay_length);
00233         if ( verbose > 1 ) g4prim->Print();
00234         g4vtx->SetPrimary(g4prim);
00235         // impose also proper time for status=1 and available end_vertex
00236         if ( (*vpitr)->status()==1 && (*vpitr)->end_vertex()!=0) {
00237           double proper_time=decay_length/(p.Beta()*p.Gamma()*c_light);
00238           if ( verbose > 1 ) LogDebug("SimG4CoreGenerator") <<"Setting proper time for beta="<<p.Beta()<<" gamma="<<p.Gamma()<<" Proper time=" <<proper_time<<" ns" ;
00239           g4prim->SetProperTime(proper_time*ns);
00240         }
00241       }
00242     }
00243 
00244     if (verbose > 1 ) g4vtx->Print();
00245     g4evt->AddPrimaryVertex(g4vtx);
00246     ng4vtx++;
00247   }
00248   
00249   // Add a protection for completely empty events (produced by LHCTransport): add a dummy vertex with no particle attached to it
00250   if ( ng4vtx == 0 ) {
00251     double x1 = 0.;
00252     double y1 = 0.;
00253     double z1 = 0.;
00254     double t1 = 0.;
00255     G4PrimaryVertex* g4vtx= 
00256       new G4PrimaryVertex(x1*mm, y1*mm, z1*mm, t1*mm/c_light);
00257     if (verbose > 1 ) g4vtx->Print();
00258     g4evt->AddPrimaryVertex(g4vtx);
00259   }
00260   
00261   delete evt ;  
00262   
00263   return ;
00264 }
00265 
00266 void Generator::particleAssignDaughters( G4PrimaryParticle* g4p, HepMC::GenParticle* vp, double decaylength)
00267 {
00268  
00269   if ( !(vp->end_vertex())  ) return ;
00270    
00271   if ( verbose > 1 ) 
00272     LogDebug("SimG4CoreGenerator") << "Special case of long decay length \n" 
00273                                    << "Assign daughters with to mother with decaylength=" << decaylength << "mm";
00274   math::XYZTLorentzVector p(vp->momentum().px(), vp->momentum().py(), vp->momentum().pz(), vp->momentum().e());
00275   double proper_time=decaylength/(p.Beta()*p.Gamma()*c_light);
00276   if ( verbose > 2 ) {
00277     LogDebug("SimG4CoreGenerator") <<" px="<<vp->momentum().px()
00278                                    <<" py="<<vp->momentum().py()
00279                                    <<" pz="<<vp->momentum().pz()
00280                                    <<" e="<<vp->momentum().e()
00281                                    <<" beta="<<p.Beta()
00282                                    <<" gamma="<<p.Gamma()
00283                                    <<" Proper time=" <<proper_time<<" ns" ;
00284   }
00285   g4p->SetProperTime(proper_time*ns); // the particle will decay after the same length if it has not interacted before
00286   double x1 = vp->end_vertex()->position().x();
00287   double y1 = vp->end_vertex()->position().y();
00288   double z1 = vp->end_vertex()->position().z();
00289   for (HepMC::GenVertex::particle_iterator 
00290          vpdec= vp->end_vertex()->particles_begin(HepMC::children);
00291        vpdec != vp->end_vertex()->particles_end(HepMC::children); ++vpdec) {
00292     
00293     //transform decay products such that in the rest frame of mother
00294     math::XYZTLorentzVector pdec((*vpdec)->momentum().px(),
00295                                  (*vpdec)->momentum().py(),
00296                                  (*vpdec)->momentum().pz(),
00297                                  (*vpdec)->momentum().e());
00298     // children should only be taken into account once
00299     G4PrimaryParticle * g4daught= 
00300       new G4PrimaryParticle((*vpdec)->pdg_id(), pdec.x()*GeV, pdec.y()*GeV, pdec.z()*GeV);
00301     if ( g4daught->GetG4code() != 0 )
00302       { 
00303         g4daught->SetMass( g4daught->GetG4code()->GetPDGMass() ) ;
00304         g4daught->SetCharge( g4daught->GetG4code()->GetPDGCharge() ) ;  
00305       }
00306     g4daught->SetWeight( 10000*(*vpdec)->barcode() ) ;
00307     setGenId( g4daught, (*vpdec)->barcode() ) ;
00308     if ( verbose > 1 ) LogDebug("SimG4CoreGenerator") <<"Assigning a "<<(*vpdec)->pdg_id()
00309                                                       <<" as daughter of a " <<vp->pdg_id() ;
00310     if ( (*vpdec)->status() == 2 && (*vpdec)->end_vertex() != 0 ) 
00311       {
00312         double x2 = (*vpdec)->end_vertex()->position().x();
00313         double y2 = (*vpdec)->end_vertex()->position().y();
00314         double z2 = (*vpdec)->end_vertex()->position().z();
00315         double dd = std::sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)+(z1-z2)*(z1-z2));
00316         particleAssignDaughters(g4daught,*vpdec,dd);
00317       }
00318     (*vpdec)->set_status(1000+(*vpdec)->status()); 
00319     g4p->SetDaughter(g4daught);
00320   }
00321   return;
00322 }
00323 
00324 bool Generator::particlePassesPrimaryCuts( const math::XYZTLorentzVector& mom, const double zimp ) const 
00325 {
00326   
00327   double phi = mom.Phi() ;   
00328   double nrg  = mom.P() ;
00329   bool   flag = true;
00330 
00331   if ( (fEtaCuts) && (zimp < Z_lmin       || zimp > Z_lmax)      ) {
00332     flag = false;
00333   } else if ( (fPCuts)  && (nrg   < theMinPCut   || nrg  > theMaxPCut) ) {
00334     flag = false;
00335   } else if ( (fPhiCuts) && (phi  < theMinPhiCut || phi  > theMaxPhiCut) ) {
00336     flag = false;
00337   }
00338 
00339   if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "Generator p = " << nrg << " z_imp = " << zimp << " phi = " << mom.Phi() << " Flag = " << flag;
00340   
00341   return flag;
00342 }
00343 
00344 bool Generator::particlePassesPrimaryCuts(const G4PrimaryParticle * p) const 
00345 {
00346 
00347   G4ThreeVector mom = p->GetMomentum();
00348   double        phi = mom.phi() ;
00349   double        nrg  = sqrt(p->GetPx()*p->GetPx() + p->GetPy()*p->GetPy() + p->GetPz()*p->GetPz());
00350   nrg /= GeV ;  // need to convert, since Geant4 operates in MeV
00351   double        eta = -log(tan(mom.theta()/2));
00352   bool   flag = true;
00353 
00354   if ( (fEtaCuts) && (eta < theMinEtaCut || eta > theMaxEtaCut) ) {
00355     flag = false;
00356   } else if ( (fPCuts)  &&  (nrg  < theMinPCut  || nrg > theMaxPCut)   ) {
00357     flag = false;
00358   } else if ( (fPhiCuts) && (phi < theMinPhiCut || phi > theMaxPhiCut) ) {
00359     flag = false;
00360   }
00361   
00362   if ( verbose > 2 ) LogDebug("SimG4CoreGenerator") << "Generator p = " << nrg  << " eta = " << eta << " theta = " << mom.theta() << " phi = " << phi << " Flag = " << flag;
00363 
00364   return flag;
00365 
00366 }
00367 
00368 void Generator::nonBeamEvent2G4(const HepMC::GenEvent * evt, G4Event * g4evt) 
00369 {
00370   int i = 0; 
00371   for(HepMC::GenEvent::particle_const_iterator it = evt->particles_begin(); 
00372       it != evt->particles_end(); ++it ) 
00373     {
00374       i++;
00375       HepMC::GenParticle * g = (*it);   
00376       int g_status = g->status();
00377       // storing only particle with status == 1         
00378       if (g_status == 1) 
00379         {
00380           int g_id = g->pdg_id();           
00381           G4PrimaryParticle * g4p = 
00382             new G4PrimaryParticle(g_id,g->momentum().px()*GeV,g->momentum().py()*GeV,g->momentum().pz()*GeV);
00383           if (g4p->GetG4code() != 0) 
00384             { 
00385               g4p->SetMass(g4p->GetG4code()->GetPDGMass());
00386               g4p->SetCharge(g4p->GetG4code()->GetPDGCharge()) ;
00387             }
00388           g4p->SetWeight(i*10000);
00389           setGenId(g4p,i);
00390           if (particlePassesPrimaryCuts(g4p)) 
00391             {
00392               G4PrimaryVertex * v = new 
00393                 G4PrimaryVertex(g->production_vertex()->position().x()*mm,
00394                                 g->production_vertex()->position().y()*mm,
00395                                 g->production_vertex()->position().z()*mm,
00396                                 g->production_vertex()->position().t()*mm/c_light);
00397               v->SetPrimary(g4p);
00398               g4evt->AddPrimaryVertex(v);
00399               if(verbose >0) {
00400                 v->Print();
00401               }
00402             }
00403         }
00404     } // end loop on HepMC particles
00405 }