00001
00002 #include <iostream>
00003
00004
00005
00006 #include "GeneratorInterface/ExternalDecays/interface/PhotosInterface.h"
00007
00008 #include "FWCore/ServiceRegistry/interface/Service.h"
00009
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
00022 extern "C"{
00023
00024 void phoini_( void );
00025 void photos_( int& );
00026
00027 double phoran_(int *idummy)
00028 {
00029 return decayRandomEngine->flat();
00030 }
00031
00032
00033
00034
00035
00036
00037
00038 extern struct {
00039
00040 bool qedrad[4000];
00041 } phoqed_;
00042
00043 }
00044
00045
00046 PhotosInterface::PhotosInterface()
00047 : fOnlyPDG(-1)
00048 {
00049 fSpecialSettings.push_back("QED-brem-off:all");
00050 fAvoidTauLeptonicDecays = false;
00051 fIsInitialized = false;
00052 }
00053
00054 PhotosInterface::PhotosInterface( const edm::ParameterSet& )
00055 : fOnlyPDG(-1)
00056 {
00057 fSpecialSettings.push_back("QED-brem-off:all");
00058 fIsInitialized = false;
00059 }
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076 void PhotosInterface::init()
00077 {
00078
00079 if ( fIsInitialized ) return;
00080
00081 phoini_();
00082
00083 fIsInitialized = true;
00084
00085 return;
00086 }
00087
00088 HepMC::GenEvent* PhotosInterface::apply( HepMC::GenEvent* evt )
00089 {
00090
00091
00092
00093
00094
00095
00096
00097 if ( !fIsInitialized ) return evt;
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111 fScaleFactors.clear();
00112
00113 for ( int ip=0; ip<evt->particles_size(); ip++ )
00114 {
00115 fScaleFactors.push_back( Scaling(HepMC::ThreeVector(1.,1.,1.),1) );
00116 phoqed_.qedrad[ip]=true;
00117 }
00118
00119
00120
00121
00122 bool tau_leptonic_decay = false;
00123 int iTauDescCounter = 0;
00124 int nTauDesc = 0;
00125
00126 for ( int iv=1; iv<=evt->vertices_size(); iv++ )
00127 {
00128
00129 HepMC::GenVertex* vtx = evt->barcode_to_vertex( -iv ) ;
00130 if ( vtx->particles_in_size() != 1 ) continue;
00131 if ( vtx->particles_out_size() <= 0 ) continue;
00132
00133
00134 bool legalVtx = false;
00135 for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00136 pitr != vtx->particles_end(HepMC::children); ++pitr)
00137 {
00138 if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
00139 {
00140
00141 legalVtx = true;
00142 break;
00143 }
00144 }
00145
00146 if ( !legalVtx ) continue;
00147
00148
00149
00150
00151
00152 HepMC::HEPEVT_Wrapper::zero_everything();
00153 fBarcodes.clear();
00154
00155
00156
00157 int index = 1;
00158 HepMC::HEPEVT_Wrapper::set_id( index, (*(vtx->particles_in_const_begin()))->pdg_id() );
00159 HepMC::FourVector vec4;
00160 vec4 = (*(vtx->particles_in_const_begin()))->momentum();
00161 HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
00162 HepMC::HEPEVT_Wrapper::set_mass( index, (*(vtx->particles_in_const_begin()))->generated_mass() );
00163 HepMC::HEPEVT_Wrapper::set_position( index, vtx->position().x(), vtx->position().y(),
00164 vtx->position().z(), vtx->position().t() );
00165 HepMC::HEPEVT_Wrapper::set_status( index, (*(vtx->particles_in_const_begin()))->status() );
00166 HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 );
00167 fBarcodes.push_back( (*(vtx->particles_in_const_begin()))->barcode() );
00168
00169
00170
00171
00172 if ( fAvoidTauLeptonicDecays && !tau_leptonic_decay && abs((*(vtx->particles_in_const_begin()))->pdg_id()) == 15 )
00173 {
00174 for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00175 pitr != vtx->particles_end(HepMC::children); ++pitr)
00176 {
00177 if ( abs((*pitr)->pdg_id()) == 11 || abs((*pitr)->pdg_id()) == 13 )
00178
00179 {
00180 tau_leptonic_decay = true;
00181 break;
00182 }
00183 }
00184 if ( vtx->particles_begin(HepMC::children) == vtx->particles_begin(HepMC::descendants) &&
00185 vtx->particles_end(HepMC::children) == vtx->particles_end(HepMC::descendants) )
00186 {
00187 nTauDesc = vtx->particles_out_size();
00188 }
00189 else
00190 {
00191 for ( HepMC::GenVertex::particle_iterator pitr1=vtx->particles_begin(HepMC::children);
00192 pitr1 != vtx->particles_end(HepMC::children); ++pitr1)
00193 {
00194 nTauDesc++;
00195 }
00196 }
00197
00198 phoqed_.qedrad[index-1]=true;
00199 iTauDescCounter = 0;
00200 }
00201
00202
00203
00204 int mbcode = (*(vtx->particles_in_const_begin()))->barcode();
00205
00206
00207
00208 int lastDau = 1;
00209 for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00210 pitr != vtx->particles_end(HepMC::children); ++pitr)
00211 {
00212 if ( fScaleFactors[mbcode-1].flag != 1. )
00213 {
00214
00215
00216 vec4 = (*pitr)->momentum();
00217 double mass2 = vec4.m2();
00218 double pxn = vec4.px() * fScaleFactors[mbcode-1].weights.x();
00219 double pyn = vec4.py() * fScaleFactors[mbcode-1].weights.y();
00220 double pzn = vec4.pz() * fScaleFactors[mbcode-1].weights.z();
00221 double en = sqrt( pxn*pxn + pyn*pyn + pzn*pzn + mass2 );
00222 (*pitr)->set_momentum( HepMC::FourVector(pxn,pyn,pzn,en) );
00223 int curbcode = (*pitr)->barcode();
00224 double scale = fScaleFactors[curbcode-1].weights.x();
00225 fScaleFactors[curbcode-1].weights.setX( scale*fScaleFactors[mbcode-1].weights.x() );
00226 scale = fScaleFactors[curbcode-1].weights.y();
00227 fScaleFactors[curbcode-1].weights.setY( scale*fScaleFactors[mbcode-1].weights.y() );
00228 scale = fScaleFactors[curbcode-1].weights.z();
00229 fScaleFactors[curbcode-1].weights.setZ( scale*fScaleFactors[mbcode-1].weights.z() );
00230 fScaleFactors[curbcode-1].flag = 0;
00231 }
00232
00233 if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
00234 {
00235 index++;
00236 vec4 = (*pitr)->momentum();
00237 HepMC::HEPEVT_Wrapper::set_id( index, (*pitr)->pdg_id() );
00238 HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
00239 HepMC::HEPEVT_Wrapper::set_mass( index, (*pitr)->generated_mass() );
00240 vec4 = (*pitr)->production_vertex()->position();
00241 HepMC::HEPEVT_Wrapper::set_position( index, vec4.x(), vec4.y(), vec4.z(), vec4.t() );
00242 HepMC::HEPEVT_Wrapper::set_status( index, (*pitr)->status() );
00243 HepMC::HEPEVT_Wrapper::set_parents( index, 1, 1 );
00244 fBarcodes.push_back( (*pitr)->barcode() );
00245 lastDau++;
00246 if ( fAvoidTauLeptonicDecays && tau_leptonic_decay )
00247 {
00248 phoqed_.qedrad[index-1]=false;
00249 iTauDescCounter++;
00250 }
00251 }
00252 }
00253
00254
00255
00256 int nentries = index;
00257
00258
00259 index = 1;
00260 HepMC::HEPEVT_Wrapper::set_children ( index, 2, lastDau );
00261
00262
00263
00264
00265 HepMC::HEPEVT_Wrapper::set_number_entries( nentries );
00266
00267
00268
00269
00270
00271
00272 photos_( index ) ;
00273
00274
00275
00276
00277
00278
00279
00280 attachParticles( evt, vtx, nentries );
00281
00282
00283
00284
00285
00286 if ( fAvoidTauLeptonicDecays && tau_leptonic_decay && iTauDescCounter == nTauDesc )
00287 {
00288 tau_leptonic_decay = false;
00289 }
00290
00291 }
00292
00293
00294 HepMC::HEPEVT_Wrapper::set_event_number( evt->event_number() );
00295
00296
00297
00298
00299
00300 return evt;
00301
00302 }
00303
00304 void PhotosInterface::attachParticles( HepMC::GenEvent* evt, HepMC::GenVertex* vtx, int nentries )
00305 {
00306
00307 if ( HepMC::HEPEVT_Wrapper::number_entries() > nentries )
00308 {
00309
00310
00311
00312
00313
00314
00315 int largestBarcode = -1;
00316 int Nbcodes = fBarcodes.size();
00317 for ( int ip=1; ip<Nbcodes; ip++ )
00318 {
00319 int bcode = fBarcodes[ip];
00320 HepMC::GenParticle* prt = evt->barcode_to_particle( bcode );
00321 if ( bcode > largestBarcode ) largestBarcode = bcode;
00322 double px = HepMC::HEPEVT_Wrapper::px(ip+1);
00323 double py = HepMC::HEPEVT_Wrapper::py(ip+1);
00324 double pz = HepMC::HEPEVT_Wrapper::pz(ip+1);
00325 double e = HepMC::HEPEVT_Wrapper::e(ip+1);
00326 HepMC::FourVector mom4 = prt->momentum();
00327
00328
00329
00330
00331 double scale = fScaleFactors[bcode-1].weights.x();
00332 fScaleFactors[bcode-1].weights.setX( scale*(px/mom4.px()) );
00333 scale = fScaleFactors[bcode-1].weights.y();
00334 fScaleFactors[bcode-1].weights.setY( scale*(py/mom4.py()) );
00335 scale = fScaleFactors[bcode-1].weights.z();
00336 fScaleFactors[bcode-1].weights.setZ( scale*(pz/mom4.pz()) );
00337 fScaleFactors[bcode-1].flag = 0;
00338
00339 prt->set_momentum( HepMC::FourVector(px,py,pz,e) );
00340
00341
00342
00343
00344
00345 }
00346
00347 int newlyGen = HepMC::HEPEVT_Wrapper::number_entries() - nentries;
00348
00349 if ( largestBarcode < evt->particles_size() )
00350 {
00351
00352
00353
00354
00355 for ( int ipp=evt->particles_size(); ipp>largestBarcode; ipp-- )
00356 {
00357 (evt->barcode_to_particle(ipp))->suggest_barcode( ipp+newlyGen );
00358 }
00359 }
00360
00361
00362
00363 for ( int ipnw=1; ipnw<=newlyGen; ipnw++ )
00364 {
00365 int nbcode = largestBarcode+ipnw;
00366 int pdg_id = HepMC::HEPEVT_Wrapper::id( nentries+ipnw );
00367 int status = HepMC::HEPEVT_Wrapper::status( nentries+ipnw );
00368 double px = HepMC::HEPEVT_Wrapper::px( nentries+ipnw );
00369 double py = HepMC::HEPEVT_Wrapper::py( nentries+ipnw );
00370 double pz = HepMC::HEPEVT_Wrapper::pz( nentries+ipnw );
00371 double e = HepMC::HEPEVT_Wrapper::e( nentries+ipnw );
00372 double m = HepMC::HEPEVT_Wrapper::m( nentries+ipnw );
00373
00374 HepMC::GenParticle* NewPart = new HepMC::GenParticle( HepMC::FourVector(px,py,pz,e),
00375 pdg_id, status);
00376 NewPart->set_generated_mass( m );
00377 NewPart->suggest_barcode( nbcode );
00378 vtx->add_particle_out( NewPart ) ;
00379
00380 fScaleFactors.push_back( Scaling(HepMC::ThreeVector(1.,1.,1.),1) );
00381 }
00382 }
00383
00384 return;
00385 }