CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/GeneratorInterface/ExternalDecays/src/PhotosInterface.cc

Go to the documentation of this file.
00001 
00002 #include <iostream>
00003 
00004 // #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Service.h"
00005 
00006 #include "GeneratorInterface/ExternalDecays/interface/PhotosInterface.h"
00007 
00008 #include "FWCore/ServiceRegistry/interface/Service.h"
00009 // #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00010 #include "GeneratorInterface/ExternalDecays/interface/DecayRandomEngine.h"
00011 
00012 #include "HepMC/GenEvent.h"
00013 #include "HepMC/IO_HEPEVT.h"
00014 #include "HepMC/HEPEVT_Wrapper.h"
00015 
00016 using namespace gen;
00017 using namespace edm;
00018 using namespace std;
00019 
00020 
00021 extern "C"{
00022 
00023    void phoini_( void );
00024    void photos_( int& );
00025 
00026    double phoran_(int *idummy)
00027    {
00028       return decayRandomEngine->flat();
00029    }
00030 /*
00031    double phoranc_(int *idummy)
00032    {
00033       return decayRandomEngine->flat();
00034    }
00035 */
00036 
00037    extern struct {
00038       // bool qedrad[NMXHEP];
00039       bool qedrad[4000]; // hardcoded for now...
00040    } phoqed_;
00041 
00042 }
00043 
00044 
00045 PhotosInterface::PhotosInterface()
00046    : fOnlyPDG(-1)
00047 {
00048    fSpecialSettings.push_back("QED-brem-off:all");
00049    fAvoidTauLeptonicDecays = false;
00050    fIsInitialized = false; 
00051 }
00052 
00053 PhotosInterface::PhotosInterface( const edm::ParameterSet& )
00054    : fOnlyPDG(-1)
00055 {
00056    fSpecialSettings.push_back("QED-brem-off:all");
00057    fIsInitialized = false;
00058 }
00059 
00060 /*  -->
00061 void PhotosInterface::configureOnlyFor( int ipdg )
00062 {
00063 
00064    fOnlyPDG = ipdg;
00065    std::ostringstream command;
00066    command << "QED-brem-off:" << fOnlyPDG ;
00067    fSpecialSettings.clear();
00068    fSpecialSettings.push_back( command.str() );
00069    
00070    return;
00071 
00072 }
00073 */
00074 
00075 void PhotosInterface::init()
00076 {
00077    
00078    if ( fIsInitialized ) return; // do init only once
00079    
00080    phoini_();
00081    
00082    fIsInitialized = true; 
00083 
00084    return;
00085 }
00086 
00087 HepMC::GenEvent* PhotosInterface::apply( HepMC::GenEvent* evt )
00088 {
00089    
00090    
00091    if ( !fIsInitialized ) return evt; // conv.read_next_event();
00092       
00093    // loop over HepMC::GenEvent, find vertices
00094       
00095    for ( int ip=0; ip<evt->particles_size(); ip++ )
00096    {
00097       phoqed_.qedrad[ip]=true;
00098    }
00099    
00100    
00101    // variables for special treatment of tau leptonic decays
00102    //
00103    bool tau_leptonic_decay = false;
00104    int iTauDescCounter = 0;
00105    int nTauDesc = 0;
00106    
00107    for ( int iv=1; iv<=evt->vertices_size(); iv++ )
00108    {
00109       
00110       bool legalVtx = false;
00111       
00112       HepMC::GenVertex* vtx = evt->barcode_to_vertex( -iv ) ;
00113       
00114       if ( vtx->particles_in_size() != 1 ) continue; // more complex than we need
00115       if ( vtx->particles_out_size() <= 1 ) continue; // no outcoming particles
00116       
00117       if ( (*(vtx->particles_in_const_begin()))->pdg_id() == 111 ) continue; // pi0 decay vtx - no point to try
00118       
00119       // --> if ( fOnlyPDG !=-1 && (*(vtx->particles_in_const_begin()))->pdg_id() == fOnlyPDG ) continue;
00120       
00121       for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00122             pitr != vtx->particles_end(HepMC::children); ++pitr) 
00123       {
00124 
00125          // quark or gluon out of this vertex - no good !
00126          if ( abs((*pitr)->pdg_id()) >=1 &&  abs((*pitr)->pdg_id()) <=8 ) break;
00127          if ( abs((*pitr)->pdg_id()) == 21 ) break;
00128 
00129          if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
00130          {
00131             // OK, legal already !
00132             legalVtx = true;
00133             break;
00134          }
00135       }
00136       
00137       if ( !legalVtx ) continue;
00138       
00139       // now do all the loops again
00140       //
00141       // first, flush out HEPEVT & tmp barcode storage
00142       //
00143       HepMC::HEPEVT_Wrapper::zero_everything();
00144       fBarcodes.clear();
00145       
00146       // add incoming particle
00147       //
00148       int index = 1;      
00149       HepMC::HEPEVT_Wrapper::set_id( index, (*(vtx->particles_in_const_begin()))->pdg_id() );
00150       HepMC::FourVector vec4;
00151       vec4 = (*(vtx->particles_in_const_begin()))->momentum();
00152       HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
00153       HepMC::HEPEVT_Wrapper::set_mass( index, (*(vtx->particles_in_const_begin()))->generated_mass() );
00154       HepMC::HEPEVT_Wrapper::set_position( index, vtx->position().x(), vtx->position().y(),
00155                                                   vtx->position().z(), vtx->position().t() );
00156       HepMC::HEPEVT_Wrapper::set_status( index, (*(vtx->particles_in_const_begin()))->status() );
00157       HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 );
00158       fBarcodes.push_back( (*(vtx->particles_in_const_begin()))->barcode() );
00159                   
00160       // special case: avoid tau leptonic decays
00161       //
00162 
00163       if ( fAvoidTauLeptonicDecays && !tau_leptonic_decay && abs((*(vtx->particles_in_const_begin()))->pdg_id()) == 15 )
00164       {      
00165          for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00166                pitr != vtx->particles_end(HepMC::children); ++pitr) 
00167          {
00168             if ( abs((*pitr)->pdg_id()) == 11 || abs((*pitr)->pdg_id()) == 13 ) // leptonic decay !!!
00169                                                                                 // do brem off tau but NOT off decay products
00170             {
00171                tau_leptonic_decay = true;
00172                break;
00173             }    
00174          }
00175          if ( vtx->particles_begin(HepMC::children) == vtx->particles_begin(HepMC::descendants) && 
00176               vtx->particles_end(HepMC::children) == vtx->particles_end(HepMC::descendants) ) // FIXME !!!!!
00177                                                                                               // Maybe better vtx nested loop(s) 
00178                                                                                               // instead of "descendants" ???
00179          {
00180             nTauDesc = vtx->particles_out_size();
00181          }
00182          else
00183          {
00184             for ( HepMC::GenVertex::particle_iterator pitr1=vtx->particles_begin(HepMC::children);
00185                   pitr1 != vtx->particles_end(HepMC::children); ++pitr1) 
00186             {
00187                nTauDesc++;
00188             }
00189          }
00190          // this is just the 1st tau in the branch, so it's allowed to emit
00191          phoqed_.qedrad[index-1]=true;
00192          iTauDescCounter = 0;
00193       }
00194            
00195       // add outcoming particles (decay products)
00196       //
00197       int lastDau = 1;
00198       for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00199             pitr != vtx->particles_end(HepMC::children); ++pitr) 
00200       {
00201 
00202          if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
00203          {
00204             index++;
00205             vec4 = (*pitr)->momentum();
00206             HepMC::HEPEVT_Wrapper::set_id( index, (*pitr)->pdg_id() );
00207             HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
00208             HepMC::HEPEVT_Wrapper::set_mass( index, (*pitr)->generated_mass() );
00209             vec4 = (*pitr)->production_vertex()->position();
00210             HepMC::HEPEVT_Wrapper::set_position( index, vec4.x(), vec4.y(), vec4.z(), vec4.t() );
00211             HepMC::HEPEVT_Wrapper::set_status( index, (*pitr)->status() );
00212             HepMC::HEPEVT_Wrapper::set_parents( index, 1, 1 );
00213             fBarcodes.push_back( (*pitr)->barcode() );
00214             lastDau++;
00215             if ( fAvoidTauLeptonicDecays && tau_leptonic_decay )
00216             {
00217                phoqed_.qedrad[index-1]=false;
00218                iTauDescCounter++;
00219             }
00220          }
00221       }
00222             
00223       // store, further to set NHEP in HEPEVT
00224       //
00225       int nentries = index;
00226       
00227       // reset master pointer to mother
00228       index = 1;
00229       HepMC::HEPEVT_Wrapper::set_children ( index, 2, lastDau ); // FIXME: need to check 
00230                                                                  // if last daughter>=2 !!!
00231       
00232       // finally, set number of entries (NHEP) in HEPEVT
00233       //
00234       HepMC::HEPEVT_Wrapper::set_number_entries( nentries );
00235 
00236       // cross-check printout HEPEVT
00237       // HepMC::HEPEVT_Wrapper::print_hepevt();
00238      
00239       // OK, 1-level vertex is formed - now, call PHOTOS
00240       //
00241       photos_( index ) ;
00242       
00243       // another cross-check printout HEPEVT - after photos
00244       // HepMC::HEPEVT_Wrapper::print_hepevt();
00245 
00246 
00247       // now check if something has been generated 
00248       // and make all adjustments to underlying vtx/parts
00249       //
00250       attachParticles( evt, vtx, nentries );
00251 
00252       // ugh, done with this vertex !
00253       
00254       // now some resets
00255       //
00256       if ( fAvoidTauLeptonicDecays && tau_leptonic_decay && iTauDescCounter == nTauDesc ) // treated tau leptonic decay and have come to the last descendent
00257       {
00258          tau_leptonic_decay = false;
00259       }
00260    
00261    }
00262    
00263    // restore event number in HEPEVT (safety measure, somehow needed by Hw6)
00264    HepMC::HEPEVT_Wrapper::set_event_number( evt->event_number() );
00265 
00266    return evt;
00267       
00268 }
00269 
00270 void PhotosInterface::attachParticles( HepMC::GenEvent* evt, HepMC::GenVertex* vtx, int nentries )
00271 {
00272 
00273    if ( HepMC::HEPEVT_Wrapper::number_entries() > nentries ) 
00274    {
00275          // yes, need all corrections and adjustments -
00276          // figure out how many photons and what particles in 
00277          // the decay branch have changes;
00278          // also, follow up each one and correct accordingly;
00279          // at the same time, add photon(s) to the GenVertex
00280          //      
00281          
00282          // vtx->print();
00283          
00284          int largestBarcode = -1;
00285          int Nbcodes = fBarcodes.size();
00286          
00287          for ( int ip=1; ip<Nbcodes; ip++ )
00288          {
00289 
00290             int bcode = fBarcodes[ip];
00291             HepMC::GenParticle* prt = evt->barcode_to_particle( bcode );
00292             if ( bcode > largestBarcode ) largestBarcode = bcode;
00293             double px = HepMC::HEPEVT_Wrapper::px(ip+1);
00294             double py = HepMC::HEPEVT_Wrapper::py(ip+1);
00295             double pz = HepMC::HEPEVT_Wrapper::pz(ip+1);
00296             double e  = HepMC::HEPEVT_Wrapper::e(ip+1);
00297             double m  = HepMC::HEPEVT_Wrapper::m(ip+1);   
00298                     
00299             if ( prt->end_vertex() )
00300             {
00301                
00302                HepMC::GenVertex* endVtx = prt->end_vertex();
00303                
00304                std::vector<int> secVtxStorage;
00305                secVtxStorage.clear();
00306                
00307                secVtxStorage.push_back( endVtx->barcode() );
00308             
00309                HepMC::FourVector mom4 = prt->momentum();
00310             
00311                // now rescale all descendants
00312                double bet1[3], bet2[3], gam1, gam2, pb;
00313                double mass = mom4.m();
00314                bet1[0] = -(mom4.px()/mass);
00315                bet1[1] = -(mom4.py()/mass);
00316                bet1[2] = -(mom4.pz()/mass);
00317                bet2[0] = px/m;
00318                bet2[1] = py/m;
00319                bet2[2] = pz/m;
00320                gam1 = mom4.e()/mass;
00321                gam2 = e/m;
00322             
00323                unsigned int vcounter = 0;
00324                     
00325                while ( vcounter < secVtxStorage.size() )
00326                {
00327                   
00328                   HepMC::GenVertex* theVtx = evt->barcode_to_vertex( secVtxStorage[vcounter] );
00329                                   
00330                   for ( HepMC::GenVertex::particle_iterator pitr=theVtx->particles_begin(HepMC::children);
00331                         pitr != theVtx->particles_end(HepMC::children); ++pitr) 
00332                   {
00333                
00334                      if ( (*pitr)->end_vertex() )
00335                      {
00336                         secVtxStorage.push_back( (*pitr)->end_vertex()->barcode() );
00337                      }
00338                      
00339                      if ( theVtx->particles_out_size() == 1 && (*pitr)->pdg_id() == prt->pdg_id() ) 
00340                      {
00341                         // carbon copy
00342                         (*pitr)->set_momentum( HepMC::FourVector(px,py,pz,e) );
00343                         continue;
00344                      }
00345                                
00346                      HepMC::FourVector dmom4 = (*pitr)->momentum();
00347                
00348                      // Boost vector to parent rest frame...
00349                      pb = bet1[0]*dmom4.px() + bet1[1]*dmom4.py() + bet1[2]*dmom4.pz();
00350                      double dpx = dmom4.px() + bet1[0] * (dmom4.e() + pb/(gam1+1.) );
00351                      double dpy = dmom4.py() + bet1[1] * (dmom4.e() + pb/(gam1+1.) );
00352                      double dpz = dmom4.pz() + bet1[2] * (dmom4.e() + pb/(gam1+1.) );
00353                      double de  = gam1*dmom4.e() + pb;
00354                      // ...and boost back to modified parent frame
00355                      pb = bet2[0]*dpx + bet2[1]*dpy + bet2[2]*dpz;
00356                      dpx += bet2[0] * ( de + pb/(gam2+1.) );
00357                      dpy += bet2[1] * ( de + pb/(gam2+1.) );
00358                      dpz += bet2[2] * ( de + pb/(gam2+1.) );
00359                      de *= gam2;
00360                      de += pb;
00361                
00362                      (*pitr)->set_momentum( HepMC::FourVector(dpx,dpy,dpz,de) );
00363                      
00364                   }                               
00365                   vcounter++;       
00366                }
00367                
00368                secVtxStorage.clear();
00369             }
00370             
00371             prt->set_momentum( HepMC::FourVector(px,py,pz,e) );
00372             
00373          } // ok, all affected particles update, but the photon(s) still not inserted
00374                  
00375          int newlyGen =  HepMC::HEPEVT_Wrapper::number_entries() - nentries;
00376          
00377          if ( largestBarcode < evt->particles_size() )
00378          {
00379             // need to adjust barcodes down from the affected vertex/particles
00380             // such that we "free up" barcodes for newly generated photons
00381             // in the middle of the event record
00382             //
00383             for ( int ipp=evt->particles_size(); ipp>largestBarcode; ipp-- )
00384             {
00385                (evt->barcode_to_particle(ipp))->suggest_barcode( ipp+newlyGen );
00386             }
00387          }
00388          
00389          // now attach new generated photons to the vertex
00390          //
00391          for ( int ipnw=1; ipnw<=newlyGen; ipnw++ )
00392          {
00393              int nbcode = largestBarcode+ipnw;
00394              int pdg_id = HepMC::HEPEVT_Wrapper::id( nentries+ipnw );
00395              int status = HepMC::HEPEVT_Wrapper::status( nentries+ipnw );
00396              double px  =  HepMC::HEPEVT_Wrapper::px( nentries+ipnw );
00397              double py  =  HepMC::HEPEVT_Wrapper::py( nentries+ipnw );
00398              double pz  =  HepMC::HEPEVT_Wrapper::pz( nentries+ipnw );
00399              double e   =  HepMC::HEPEVT_Wrapper::e(  nentries+ipnw );
00400              double m   =  HepMC::HEPEVT_Wrapper::m(  nentries+ipnw );
00401           
00402              HepMC::GenParticle* NewPart = new HepMC::GenParticle( HepMC::FourVector(px,py,pz,e),
00403                                                                     pdg_id, status);
00404              NewPart->set_generated_mass( m );
00405              NewPart->suggest_barcode( nbcode );
00406              vtx->add_particle_out( NewPart ) ;
00407          }  
00408    
00409       // vtx->print();
00410       // std::cout << " leaving attachParticles() " << std::endl;
00411    
00412    } // end of global if-statement 
00413 
00414    return;
00415 }