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 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
00032
00033
00034
00035
00036
00037 extern struct {
00038
00039 bool qedrad[4000];
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 void PhotosInterface::configureOnlyFor( int ipdg )
00061 {
00062
00063 fOnlyPDG = ipdg;
00064
00065
00066 fSpecialSettings.clear();
00067
00068
00069 return;
00070
00071 }
00072
00073 void PhotosInterface::init()
00074 {
00075
00076 if ( fIsInitialized ) return;
00077
00078 phoini_();
00079
00080 fIsInitialized = true;
00081
00082 return;
00083 }
00084
00085 HepMC::GenEvent* PhotosInterface::apply( HepMC::GenEvent* evt )
00086 {
00087
00088 if ( !fIsInitialized ) return evt;
00089
00090
00091
00092
00093 for ( int ip=0; ip<4000; ip++ )
00094 {
00095 phoqed_.qedrad[ip]=true;
00096 }
00097
00098
00099
00100
00101
00102 for ( int iv=1; iv<=evt->vertices_size(); iv++ )
00103 {
00104
00105 bool legalVtx = false;
00106
00107 fSecVtxStore.clear();
00108
00109 HepMC::GenVertex* vtx = evt->barcode_to_vertex( -iv ) ;
00110
00111 if ( vtx->particles_in_size() != 1 ) continue;
00112 if ( vtx->particles_out_size() <= 1 ) continue;
00113
00114 if ( (*(vtx->particles_in_const_begin()))->pdg_id() == 111 ) continue;
00115
00116 if ( fOnlyPDG != 1 && (*(vtx->particles_in_const_begin()))->pdg_id() != fOnlyPDG )
00117 {
00118 continue;
00119 }
00120 else
00121 {
00122
00123
00124
00125
00126 bool same = false;
00127 for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00128 pitr != vtx->particles_end(HepMC::children); ++pitr)
00129 {
00130 if ( (*pitr)->pdg_id() == fOnlyPDG )
00131 {
00132 same = true;
00133 break;
00134 }
00135 }
00136 if ( same ) continue;
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149 if ( fOnlyPDG == 15 && fAvoidTauLeptonicDecays && isTauLeptonicDecay( vtx ) ) continue;
00150
00151 applyToBranch( evt, -iv );
00152 continue;
00153 }
00154
00155
00156
00157 for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00158 pitr != vtx->particles_end(HepMC::children); ++pitr)
00159 {
00160
00161
00162 if ( abs((*pitr)->pdg_id()) >=1 && abs((*pitr)->pdg_id()) <=8 ) break;
00163 if ( abs((*pitr)->pdg_id()) == 21 ) break;
00164
00165 if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
00166 {
00167
00168 legalVtx = true;
00169 break;
00170 }
00171 }
00172
00173 if ( !legalVtx ) continue;
00174
00175 applyToVertex( evt, -iv );
00176
00177 }
00178
00179
00180 HepMC::HEPEVT_Wrapper::set_event_number( evt->event_number() );
00181
00182 return evt;
00183
00184 }
00185
00186 void PhotosInterface::applyToVertex( HepMC::GenEvent* evt, int vtxbcode )
00187 {
00188
00189 HepMC::GenVertex* vtx = evt->barcode_to_vertex( vtxbcode );
00190
00191 if ( fAvoidTauLeptonicDecays && isTauLeptonicDecay( vtx ) ) return;
00192
00193
00194
00195
00196
00197
00198
00199 HepMC::HEPEVT_Wrapper::zero_everything();
00200 fBarcodes.clear();
00201
00202
00203
00204 int index = 1;
00205 HepMC::HEPEVT_Wrapper::set_id( index, (*(vtx->particles_in_const_begin()))->pdg_id() );
00206 HepMC::FourVector vec4;
00207 vec4 = (*(vtx->particles_in_const_begin()))->momentum();
00208 HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
00209 HepMC::HEPEVT_Wrapper::set_mass( index, (*(vtx->particles_in_const_begin()))->generated_mass() );
00210 HepMC::HEPEVT_Wrapper::set_position( index, vtx->position().x(), vtx->position().y(),
00211 vtx->position().z(), vtx->position().t() );
00212 HepMC::HEPEVT_Wrapper::set_status( index, (*(vtx->particles_in_const_begin()))->status() );
00213 HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 );
00214 fBarcodes.push_back( (*(vtx->particles_in_const_begin()))->barcode() );
00215
00216
00217
00218 int lastDau = 1;
00219 for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00220 pitr != vtx->particles_end(HepMC::children); ++pitr)
00221 {
00222
00223 if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
00224 {
00225 index++;
00226 vec4 = (*pitr)->momentum();
00227 HepMC::HEPEVT_Wrapper::set_id( index, (*pitr)->pdg_id() );
00228 HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
00229 HepMC::HEPEVT_Wrapper::set_mass( index, (*pitr)->generated_mass() );
00230 vec4 = (*pitr)->production_vertex()->position();
00231 HepMC::HEPEVT_Wrapper::set_position( index, vec4.x(), vec4.y(), vec4.z(), vec4.t() );
00232 HepMC::HEPEVT_Wrapper::set_status( index, (*pitr)->status() );
00233 HepMC::HEPEVT_Wrapper::set_parents( index, 1, 1 );
00234 fBarcodes.push_back( (*pitr)->barcode() );
00235 lastDau++;
00236 }
00237 if ( (*pitr)->end_vertex() )
00238 {
00239 fSecVtxStore.push_back( (*pitr)->end_vertex()->barcode() );
00240 }
00241 }
00242
00243
00244
00245 int nentries = index;
00246
00247
00248 index = 1;
00249 HepMC::HEPEVT_Wrapper::set_children ( index, 2, lastDau );
00250
00251
00252
00253
00254 HepMC::HEPEVT_Wrapper::set_number_entries( nentries );
00255
00256
00257
00258
00259
00260
00261 photos_( index ) ;
00262
00263
00264
00265
00266
00267
00268
00269
00270 attachParticles( evt, vtx, nentries );
00271
00272
00273
00274 return;
00275
00276 }
00277
00278 void PhotosInterface::applyToBranch( HepMC::GenEvent* evt, int vtxbcode )
00279 {
00280
00281
00282 fSecVtxStore.clear();
00283
00284
00285
00286 applyToVertex( evt, vtxbcode );
00287
00288
00289
00290
00291
00292 unsigned int vcounter = 0;
00293
00294 while ( vcounter < fSecVtxStore.size() )
00295 {
00296 applyToVertex( evt, fSecVtxStore[vcounter] );
00297 vcounter++;
00298 }
00299
00300 fSecVtxStore.clear();
00301
00302 return;
00303
00304 }
00305
00306 void PhotosInterface::attachParticles( HepMC::GenEvent* evt, HepMC::GenVertex* vtx, int nentries )
00307 {
00308
00309 if ( HepMC::HEPEVT_Wrapper::number_entries() > nentries )
00310 {
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320 int largestBarcode = -1;
00321 int Nbcodes = fBarcodes.size();
00322
00323 for ( int ip=1; ip<Nbcodes; ip++ )
00324 {
00325
00326 int bcode = fBarcodes[ip];
00327 HepMC::GenParticle* prt = evt->barcode_to_particle( bcode );
00328 if ( bcode > largestBarcode ) largestBarcode = bcode;
00329 double px = HepMC::HEPEVT_Wrapper::px(ip+1);
00330 double py = HepMC::HEPEVT_Wrapper::py(ip+1);
00331 double pz = HepMC::HEPEVT_Wrapper::pz(ip+1);
00332 double e = HepMC::HEPEVT_Wrapper::e(ip+1);
00333 double m = HepMC::HEPEVT_Wrapper::m(ip+1);
00334
00335 if ( prt->end_vertex() )
00336 {
00337
00338 HepMC::GenVertex* endVtx = prt->end_vertex();
00339
00340 std::vector<int> secVtxStorage;
00341 secVtxStorage.clear();
00342
00343 secVtxStorage.push_back( endVtx->barcode() );
00344
00345 HepMC::FourVector mom4 = prt->momentum();
00346
00347
00348 double bet1[3], bet2[3], gam1, gam2, pb;
00349 double mass = mom4.m();
00350 bet1[0] = -(mom4.px()/mass);
00351 bet1[1] = -(mom4.py()/mass);
00352 bet1[2] = -(mom4.pz()/mass);
00353 bet2[0] = px/m;
00354 bet2[1] = py/m;
00355 bet2[2] = pz/m;
00356 gam1 = mom4.e()/mass;
00357 gam2 = e/m;
00358
00359 unsigned int vcounter = 0;
00360
00361 while ( vcounter < secVtxStorage.size() )
00362 {
00363
00364 HepMC::GenVertex* theVtx = evt->barcode_to_vertex( secVtxStorage[vcounter] );
00365
00366 for ( HepMC::GenVertex::particle_iterator pitr=theVtx->particles_begin(HepMC::children);
00367 pitr != theVtx->particles_end(HepMC::children); ++pitr)
00368 {
00369
00370 if ( (*pitr)->end_vertex() )
00371 {
00372 secVtxStorage.push_back( (*pitr)->end_vertex()->barcode() );
00373 }
00374
00375 if ( theVtx->particles_out_size() == 1 && (*pitr)->pdg_id() == prt->pdg_id() )
00376 {
00377
00378 (*pitr)->set_momentum( HepMC::FourVector(px,py,pz,e) );
00379 continue;
00380 }
00381
00382 HepMC::FourVector dmom4 = (*pitr)->momentum();
00383
00384
00385 pb = bet1[0]*dmom4.px() + bet1[1]*dmom4.py() + bet1[2]*dmom4.pz();
00386 double dpx = dmom4.px() + bet1[0] * (dmom4.e() + pb/(gam1+1.) );
00387 double dpy = dmom4.py() + bet1[1] * (dmom4.e() + pb/(gam1+1.) );
00388 double dpz = dmom4.pz() + bet1[2] * (dmom4.e() + pb/(gam1+1.) );
00389 double de = gam1*dmom4.e() + pb;
00390
00391 pb = bet2[0]*dpx + bet2[1]*dpy + bet2[2]*dpz;
00392 dpx += bet2[0] * ( de + pb/(gam2+1.) );
00393 dpy += bet2[1] * ( de + pb/(gam2+1.) );
00394 dpz += bet2[2] * ( de + pb/(gam2+1.) );
00395 de *= gam2;
00396 de += pb;
00397
00398 (*pitr)->set_momentum( HepMC::FourVector(dpx,dpy,dpz,de) );
00399
00400 }
00401 vcounter++;
00402 }
00403
00404 secVtxStorage.clear();
00405 }
00406
00407 prt->set_momentum( HepMC::FourVector(px,py,pz,e) );
00408
00409 }
00410
00411 int newlyGen = HepMC::HEPEVT_Wrapper::number_entries() - nentries;
00412
00413 if ( largestBarcode < evt->particles_size() )
00414 {
00415
00416
00417
00418
00419 for ( int ipp=evt->particles_size(); ipp>largestBarcode; ipp-- )
00420 {
00421 (evt->barcode_to_particle(ipp))->suggest_barcode( ipp+newlyGen );
00422 }
00423 }
00424
00425
00426
00427 for ( int ipnw=1; ipnw<=newlyGen; ipnw++ )
00428 {
00429 int nbcode = largestBarcode+ipnw;
00430 int pdg_id = HepMC::HEPEVT_Wrapper::id( nentries+ipnw );
00431 int status = HepMC::HEPEVT_Wrapper::status( nentries+ipnw );
00432 double px = HepMC::HEPEVT_Wrapper::px( nentries+ipnw );
00433 double py = HepMC::HEPEVT_Wrapper::py( nentries+ipnw );
00434 double pz = HepMC::HEPEVT_Wrapper::pz( nentries+ipnw );
00435 double e = HepMC::HEPEVT_Wrapper::e( nentries+ipnw );
00436 double m = HepMC::HEPEVT_Wrapper::m( nentries+ipnw );
00437
00438 HepMC::GenParticle* NewPart = new HepMC::GenParticle( HepMC::FourVector(px,py,pz,e),
00439 pdg_id, status);
00440 NewPart->set_generated_mass( m );
00441 NewPart->suggest_barcode( nbcode );
00442 vtx->add_particle_out( NewPart ) ;
00443 }
00444
00445
00446
00447
00448 }
00449
00450 return;
00451 }
00452
00453 bool PhotosInterface::isTauLeptonicDecay( HepMC::GenVertex* vtx )
00454 {
00455
00456 if ( abs((*(vtx->particles_in_const_begin()))->pdg_id()) != 15 ) return false;
00457
00458 for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
00459 pitr != vtx->particles_end(HepMC::children); ++pitr)
00460 {
00461 if ( abs((*pitr)->pdg_id()) == 11 || abs((*pitr)->pdg_id()) == 13 )
00462 {
00463 return true;
00464 }
00465 }
00466
00467 return false;
00468
00469 }
00470
00471
00472
00473
00474
00475
00476
00477
00478
00479
00480
00481
00482
00483
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500
00501
00502
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515
00516
00517
00518
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529
00530
00531
00532
00533
00534
00535
00536
00537
00538
00539
00540
00541
00542
00543
00544
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554
00555
00556
00557
00558
00559
00560
00561
00562
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576
00577
00578
00579
00580
00581
00582
00583
00584
00585
00586
00587
00588
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598
00599
00600
00601
00602
00603