00001
00002 #include "HepMC/GenEvent.h"
00003 #include "HepMC/GenVertex.h"
00004 #include "HepMC/GenParticle.h"
00005
00006
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008
00009
00010 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00011 #include "DataFormats/Candidate/interface/Candidate.h"
00012
00013
00014 #include "FastSimulation/Event/interface/FBaseSimEvent.h"
00015 #include "FastSimulation/Event/interface/FSimTrack.h"
00016 #include "FastSimulation/Event/interface/FSimVertex.h"
00017 #include "FastSimulation/Event/interface/KineParticleFilter.h"
00018 #include "FastSimulation/BaseParticlePropagator/interface/BaseParticlePropagator.h"
00019 #include "FastSimulation/Event/interface/BetaFuncPrimaryVertexGenerator.h"
00020 #include "FastSimulation/Event/interface/GaussianPrimaryVertexGenerator.h"
00021 #include "FastSimulation/Event/interface/FlatPrimaryVertexGenerator.h"
00022 #include "FastSimulation/Event/interface/NoPrimaryVertexGenerator.h"
00023
00024
00025 using namespace HepPDT;
00026
00027
00028 #include <iostream>
00029 #include <iomanip>
00030 #include <map>
00031 #include <string>
00032
00033 FBaseSimEvent::FBaseSimEvent(const edm::ParameterSet& kine)
00034 :
00035 nSimTracks(0),
00036 nSimVertices(0),
00037 nGenParticles(0),
00038 nChargedParticleTracks(0),
00039 initialSize(5000),
00040 random(0)
00041 {
00042
00043 theVertexGenerator = new NoPrimaryVertexGenerator();
00044 theBeamSpot = math::XYZPoint(0.0,0.0,0.0);
00045
00046
00047 theGenParticles = new std::vector<HepMC::GenParticle*>();
00048 theSimTracks = new std::vector<FSimTrack>;
00049 theSimVertices = new std::vector<FSimVertex>;
00050 theChargedTracks = new std::vector<unsigned>();
00051
00052
00053
00054 theSimTracks->resize(initialSize);
00055 theSimVertices->resize(initialSize);
00056 theGenParticles->resize(initialSize);
00057 theChargedTracks->resize(initialSize);
00058 theTrackSize = initialSize;
00059 theVertexSize = initialSize;
00060 theGenSize = initialSize;
00061 theChargedSize = initialSize;
00062
00063
00064
00065 myFilter = new KineParticleFilter(kine);
00066
00067 }
00068
00069 FBaseSimEvent::FBaseSimEvent(const edm::ParameterSet& vtx,
00070 const edm::ParameterSet& kine,
00071 const RandomEngine* engine)
00072 :
00073 nSimTracks(0),
00074 nSimVertices(0),
00075 nGenParticles(0),
00076 nChargedParticleTracks(0),
00077 initialSize(5000),
00078 theVertexGenerator(0),
00079 random(engine)
00080 {
00081
00082
00083 std::string vtxType = vtx.getParameter<std::string>("type");
00084 if ( vtxType == "Gaussian" )
00085 theVertexGenerator = new GaussianPrimaryVertexGenerator(vtx,random);
00086 else if ( vtxType == "Flat" )
00087 theVertexGenerator = new FlatPrimaryVertexGenerator(vtx,random);
00088 else if ( vtxType == "BetaFunc" )
00089 theVertexGenerator = new BetaFuncPrimaryVertexGenerator(vtx,random);
00090 else
00091 theVertexGenerator = new NoPrimaryVertexGenerator();
00092
00093 theBeamSpot = math::XYZPoint(0.0,0.0,0.0);
00094
00095
00096
00097
00098 lateVertexPosition = 10.;
00099
00100
00101 theGenParticles = new std::vector<HepMC::GenParticle*>();
00102 theSimTracks = new std::vector<FSimTrack>;
00103 theSimVertices = new std::vector<FSimVertex>;
00104 theChargedTracks = new std::vector<unsigned>();
00105
00106
00107
00108 theSimTracks->resize(initialSize);
00109 theSimVertices->resize(initialSize);
00110 theGenParticles->resize(initialSize);
00111 theChargedTracks->resize(initialSize);
00112 theTrackSize = initialSize;
00113 theVertexSize = initialSize;
00114 theGenSize = initialSize;
00115 theChargedSize = initialSize;
00116
00117
00118
00119 myFilter = new KineParticleFilter(kine);
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130 }
00131
00132 FBaseSimEvent::~FBaseSimEvent(){
00133
00134
00135 theGenParticles->clear();
00136 theSimTracks->clear();
00137 theSimVertices->clear();
00138 theChargedTracks->clear();
00139
00140
00141 delete theGenParticles;
00142 delete theSimTracks;
00143 delete theSimVertices;
00144 delete theChargedTracks;
00145 delete myFilter;
00146
00147
00148
00149
00150 }
00151
00152 void
00153 FBaseSimEvent::initializePdt(const HepPDT::ParticleDataTable* aPdt) {
00154
00155 pdt = aPdt;
00156
00157 }
00158
00159
00160
00161
00162
00163
00164
00165
00166 void
00167 FBaseSimEvent::fill(const HepMC::GenEvent& myGenEvent) {
00168
00169
00170 clear();
00171
00172
00173 addParticles(myGenEvent);
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183 }
00184
00185 void
00186 FBaseSimEvent::fill(const reco::GenParticleCollection& myGenParticles) {
00187
00188
00189 clear();
00190
00191
00192 addParticles(myGenParticles);
00193
00194 }
00195
00196 void
00197 FBaseSimEvent::fill(const std::vector<SimTrack>& simTracks,
00198 const std::vector<SimVertex>& simVertices) {
00199
00200
00201
00202
00203 clear();
00204
00205 unsigned nVtx = simVertices.size();
00206 unsigned nTks = simTracks.size();
00207
00208
00209 if ( nVtx == 0 ) return;
00210
00211
00212 std::vector<int> myVertices(nVtx,-1);
00213 std::vector<int> myTracks(nTks,-1);
00214
00215
00216
00217
00218 std::map<unsigned, unsigned> geantToIndex;
00219 for( unsigned it=0; it<simTracks.size(); ++it ) {
00220 geantToIndex[ simTracks[it].trackId() ] = it;
00221 }
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238 XYZTLorentzVector primaryVertex(simVertices[0].position().x(),
00239 simVertices[0].position().y(),
00240 simVertices[0].position().z(),
00241 simVertices[0].position().t());
00242
00243 myFilter->setMainVertex(primaryVertex);
00244
00245 addSimVertex(myFilter->vertex());
00246 myVertices[0] = 0;
00247
00248 for( unsigned trackId=0; trackId<nTks; ++trackId ) {
00249
00250
00251 const SimTrack& track = simTracks[trackId];
00252
00253
00254
00255 int vertexId = track.vertIndex();
00256 const SimVertex& vertex = simVertices[vertexId];
00257
00258
00259
00260 int motherId = -1;
00261 if( !vertex.noParent() ) {
00262
00263 unsigned motherGeantId = vertex.parentIndex();
00264 std::map<unsigned, unsigned >::iterator association
00265 = geantToIndex.find( motherGeantId );
00266 if(association != geantToIndex.end() )
00267 motherId = association->second;
00268 }
00269 int originId = motherId == - 1 ? -1 : myTracks[motherId];
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281 XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
00282 vertex.position().pz(),vertex.position().e());
00283 if ( myVertices[vertexId] == -1 )
00284
00285
00286
00287
00288 myVertices[vertexId] = addSimVertex(position,originId);
00289
00290
00291 int motherType = motherId == -1 ? 0 : simTracks[motherId].type();
00292
00293 if ( abs(motherType) != 11 || motherType != track.type() ) {
00294
00295
00296
00297
00298 XYZTLorentzVector momentum(track.momentum().px(),track.momentum().py(),
00299 track.momentum().pz(),track.momentum().e());
00300 RawParticle part(momentum,position);
00301
00302 part.setID(track.type());
00303
00304
00305
00306
00307 myTracks[trackId] = addSimTrack(&part,myVertices[vertexId],track.genpartIndex());
00308 if ( track.trackerSurfacePosition().perp2() > 22500. || fabs(track.trackerSurfacePosition().z()) > 400. )
00309 (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition()/10.);
00310 else
00311 (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
00312 (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
00313 } else {
00314 myTracks[trackId] = myTracks[motherId];
00315 if ( track.trackerSurfacePosition().perp2() > 22500. || fabs(track.trackerSurfacePosition().z()) > 400. )
00316 (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition()/10.);
00317 else
00318 (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
00319 (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
00320 }
00321
00322 }
00323
00324
00325 for( unsigned vertexId=0; vertexId<nVtx; ++vertexId ) {
00326
00327
00328 if ( myVertices[vertexId] != -1 ) continue;
00329
00330
00331 const SimVertex& vertex = simVertices[vertexId];
00332
00333
00334 int motherId = -1;
00335 if( !vertex.noParent() ) {
00336
00337
00338 unsigned motherGeantId = vertex.parentIndex();
00339 std::map<unsigned, unsigned >::iterator association
00340 = geantToIndex.find( motherGeantId );
00341 if(association != geantToIndex.end() )
00342 motherId = association->second;
00343 }
00344 int originId = motherId == - 1 ? -1 : myTracks[motherId];
00345
00346
00347
00348
00349
00350
00351 XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
00352 vertex.position().pz(),vertex.position().e());
00353 myVertices[vertexId] = addSimVertex(position,originId);
00354 }
00355
00356
00357 BaseParticlePropagator myPart;
00358 XYZTLorentzVector mom;
00359 XYZTLorentzVector pos;
00360
00361
00362 for( int fsimi=0; fsimi < (int)nTracks() ; ++fsimi) {
00363
00364
00365 FSimTrack& myTrack = track(fsimi);
00366 double trackerSurfaceTime = myTrack.vertex().position().t()
00367 + myTrack.momentum().e()/myTrack.momentum().pz()
00368 * ( myTrack.trackerSurfacePosition().z()
00369 - myTrack.vertex().position().z() );
00370 pos = XYZTLorentzVector(myTrack.trackerSurfacePosition().x(),
00371 myTrack.trackerSurfacePosition().y(),
00372 myTrack.trackerSurfacePosition().z(),
00373 trackerSurfaceTime);
00374 mom = XYZTLorentzVector(myTrack.trackerSurfaceMomentum().x(),
00375 myTrack.trackerSurfaceMomentum().y(),
00376 myTrack.trackerSurfaceMomentum().z(),
00377 myTrack.trackerSurfaceMomentum().t());
00378
00379 if ( mom.T() > 0. ) {
00380
00381 myPart = BaseParticlePropagator(RawParticle(mom,pos),0.,0.,4.);
00382 myPart.setCharge(myTrack.charge());
00383
00384
00385 myPart.propagateToPreshowerLayer1(false);
00386 if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
00387 myTrack.setLayer1(myPart,myPart.getSuccess());
00388
00389
00390 myPart.propagateToPreshowerLayer2(false);
00391 if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
00392 myTrack.setLayer2(myPart,myPart.getSuccess());
00393
00394
00395 myPart.propagateToEcalEntrance(false);
00396 if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00397 myTrack.setEcal(myPart,myPart.getSuccess());
00398
00399
00400 myPart.propagateToHcalEntrance(false);
00401 if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00402 myTrack.setHcal(myPart,myPart.getSuccess());
00403
00404
00405 myPart.propagateToVFcalEntrance(false);
00406 if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00407 myTrack.setVFcal(myPart,myPart.getSuccess());
00408 }
00409
00410 }
00411
00412 }
00413
00414
00415 void
00416 FBaseSimEvent::addParticles(const HepMC::GenEvent& myGenEvent) {
00417
00419 int genEventSize = myGenEvent.particles_size();
00420 std::vector<int> myGenVertices(genEventSize, static_cast<int>(0));
00421
00422
00423 if ( myGenEvent.particles_empty() ) return;
00424
00425
00426 int offset = nGenParts();
00427
00428
00429 HepMC::GenVertex* primaryVertex = *(myGenEvent.vertices_begin());
00430
00431
00432 unsigned primaryMother = primaryVertex->particles_in_size();
00433 if ( primaryMother ) {
00434 unsigned partId = (*(primaryVertex->particles_in_const_begin()))->pdg_id();
00435 if ( abs(partId) == 2212 ) primaryMother = 0;
00436 }
00437
00438
00439 XYZTLorentzVector primaryVertexPosition(primaryVertex->position().x()/10.,
00440 primaryVertex->position().y()/10.,
00441 primaryVertex->position().z()/10.,
00442 primaryVertex->position().t()/10.);
00443
00444 primaryVertexPosition *= (1-primaryMother);
00445
00446
00447
00448
00449 XYZTLorentzVector smearedVertex;
00450 if ( primaryVertexPosition.Vect().Mag2() < 1E-16 ) {
00451 theVertexGenerator->generate();
00452 smearedVertex = XYZTLorentzVector(
00453 theVertexGenerator->X()-theVertexGenerator->beamSpot().X()+theBeamSpot.X(),
00454 theVertexGenerator->Y()-theVertexGenerator->beamSpot().Y()+theBeamSpot.Y(),
00455 theVertexGenerator->Z()-theVertexGenerator->beamSpot().Z()+theBeamSpot.Z(),
00456 0.);
00457 }
00458
00459
00460 myFilter->setMainVertex(primaryVertexPosition+smearedVertex);
00461
00462
00463 int mainVertex = addSimVertex(myFilter->vertex());
00464
00465 HepMC::GenEvent::particle_const_iterator piter;
00466 HepMC::GenEvent::particle_const_iterator pbegin = myGenEvent.particles_begin();
00467 HepMC::GenEvent::particle_const_iterator pend = myGenEvent.particles_end();
00468
00469 int initialBarcode = 0;
00470 if ( pbegin != pend ) initialBarcode = (*pbegin)->barcode();
00471
00472 for ( piter = pbegin; piter != pend; ++piter ) {
00473
00474
00475 HepMC::GenParticle* p = *piter;
00476
00477 if ( !offset ) {
00478 (*theGenParticles)[nGenParticles++] = p;
00479 if ( nGenParticles/theGenSize*theGenSize == nGenParticles ) {
00480 theGenSize *= 2;
00481 theGenParticles->resize(theGenSize);
00482 }
00483
00484 }
00485
00486
00487
00488
00489 XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
00490 HepMC::GenVertex* productionVertex = p->production_vertex();
00491 if ( productionVertex ) {
00492 unsigned productionMother = productionVertex->particles_in_size();
00493 if ( productionMother ) {
00494 int productionMotherId = (*(productionVertex->particles_in_const_begin()))->pdg_id();
00495 if ( abs(productionMotherId) < 1000000 )
00496 productionVertexPosition =
00497 XYZTLorentzVector(productionVertex->position().x()/10.,
00498 productionVertex->position().y()/10.,
00499 productionVertex->position().z()/10.,
00500 productionVertex->position().t()/10.) + smearedVertex;
00501 }
00502 }
00503
00504 if ( productionVertexPosition.Perp2() > lateVertexPosition ) continue;
00505
00506
00507
00508
00509 bool testStable = p->status()%1000==1;
00510
00511
00512 if ( p->status() == 2 && abs(p->pdg_id()) < 1000000) {
00513 HepMC::GenVertex* endVertex = p->end_vertex();
00514 if ( endVertex ) {
00515 XYZTLorentzVector decayPosition =
00516 XYZTLorentzVector(endVertex->position().x()/10.,
00517 endVertex->position().y()/10.,
00518 endVertex->position().z()/10.,
00519 endVertex->position().t()/10.) + smearedVertex;
00520
00521 if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
00522 }
00523 }
00524
00525
00526 bool testDaugh = false;
00527 if ( !testStable &&
00528 p->status() == 2 &&
00529 p->end_vertex() &&
00530 p->end_vertex()->particles_out_size() ) {
00531 HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
00532 p->end_vertex()->particles_out_const_begin();
00533 HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
00534 p->end_vertex()->particles_out_const_end();
00535 for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
00536 HepMC::GenParticle* daugh = *firstDaughterIt;
00537 if ( daugh->status()%1000==1 ) {
00538 testDaugh=true;
00539 break;
00540 }
00541 }
00542 }
00543
00544
00545 double dist = 0.;
00546 if ( !testStable && !testDaugh && p->production_vertex() ) {
00547 XYZTLorentzVector
00548 productionVertexPosition(p->production_vertex()->position().x()/10.,
00549 p->production_vertex()->position().y()/10.,
00550 p->production_vertex()->position().z()/10.,
00551 p->production_vertex()->position().t()/10.);
00552 dist = (primaryVertexPosition-productionVertexPosition).Vect().Mag2();
00553 }
00554 bool testDecay = ( dist > 1e-8 ) ? true : false;
00555
00556
00557 if ( testStable || testDaugh || testDecay ) {
00558
00559
00560
00561
00562
00563
00564 int motherBarcode = p->production_vertex() &&
00565 p->production_vertex()->particles_in_const_begin() !=
00566 p->production_vertex()->particles_in_const_end() ?
00567 (*(p->production_vertex()->particles_in_const_begin()))->barcode() : 0;
00568
00569 int originVertex =
00570 motherBarcode && myGenVertices[motherBarcode-initialBarcode] ?
00571 myGenVertices[motherBarcode-initialBarcode] : mainVertex;
00572
00573 XYZTLorentzVector momentum(p->momentum().px(),
00574 p->momentum().py(),
00575 p->momentum().pz(),
00576 p->momentum().e());
00577 RawParticle part(momentum, vertex(originVertex).position());
00578 part.setID(p->pdg_id());
00579
00580
00581
00582 int theTrack = testStable && p->end_vertex() ?
00583
00584 addSimTrack(&part,originVertex, nGenParts()-offset,p->end_vertex()) :
00585
00586 addSimTrack(&part,originVertex, nGenParts()-offset);
00587
00588 if (
00589
00590 !p->end_vertex() ||
00591
00592
00593 ( testStable && p->end_vertex() )
00594
00595 ) continue;
00596
00597
00598 XYZTLorentzVector decayVertex =
00599 XYZTLorentzVector(p->end_vertex()->position().x()/10.,
00600 p->end_vertex()->position().y()/10.,
00601 p->end_vertex()->position().z()/10.,
00602 p->end_vertex()->position().t()/10.) +
00603 smearedVertex;
00604
00605 int theVertex = addSimVertex(decayVertex,theTrack);
00606
00607 if ( theVertex != -1 ) myGenVertices[p->barcode()-initialBarcode] = theVertex;
00608
00609
00610 }
00611 }
00612
00613 }
00614
00615 void
00616 FBaseSimEvent::addParticles(const reco::GenParticleCollection& myGenParticles) {
00617
00618
00619 unsigned int nParticles = myGenParticles.size();
00620 nGenParticles = nParticles;
00621
00622 if ( !nParticles ) return;
00623
00625 std::map<const reco::Candidate*,int> myGenVertices;
00626
00627
00628 int offset = nTracks();
00629
00630
00631 nGenParticles = 0;
00632 unsigned int ip = 0;
00633 if ( nParticles > 1 &&
00634 myGenParticles[0].pdgId() == 2212 &&
00635 myGenParticles[1].pdgId() == 2212 ) {
00636 ip = 2;
00637 nGenParticles = 2;
00638 }
00639
00640
00641 XYZTLorentzVector primaryVertex (myGenParticles[ip].vx(),
00642 myGenParticles[ip].vy(),
00643 myGenParticles[ip].vz(),
00644 0.);
00645
00646
00647 XYZTLorentzVector smearedVertex;
00648 if ( primaryVertex.mag() < 1E-8 ) {
00649 theVertexGenerator->generate();
00650 smearedVertex = XYZTLorentzVector(
00651 theVertexGenerator->X()-theVertexGenerator->beamSpot().X()+theBeamSpot.X(),
00652 theVertexGenerator->Y()-theVertexGenerator->beamSpot().Y()+theBeamSpot.Y(),
00653 theVertexGenerator->Z()-theVertexGenerator->beamSpot().Z()+theBeamSpot.Z(),
00654 0.);
00655 }
00656
00657
00658 myFilter->setMainVertex(primaryVertex+smearedVertex);
00659
00660
00661 int mainVertex = addSimVertex(myFilter->vertex());
00662
00663
00664 for ( ; ip<nParticles; ++ip ) {
00665
00666
00667
00668 nGenParticles++;
00669 const reco::GenParticle& p = myGenParticles[ip];
00670
00671
00672
00673
00674 XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
00675 const reco::Candidate* productionMother = p.numberOfMothers() ? p.mother(0) : 0;
00676 if ( productionMother ) {
00677 int productionMotherId = productionMother->pdgId();
00678 if ( abs(productionMotherId) < 1000000 )
00679 productionVertexPosition = XYZTLorentzVector(p.vx(), p.vy(), p.vz(), 0.) + smearedVertex;
00680 }
00681 if ( productionVertexPosition.Perp2() > lateVertexPosition ) continue;
00682
00683
00684
00685 bool testStable = p.status()%1000==1;
00686
00687
00688 if ( p.status() == 2 && abs(p.pdgId()) < 1000000 ) {
00689 unsigned int nDaughters = p.numberOfDaughters();
00690 if ( nDaughters ) {
00691 const reco::Candidate* daughter = p.daughter(0);
00692 XYZTLorentzVector decayPosition =
00693 XYZTLorentzVector(daughter->vx(), daughter->vy(), daughter->vz(), 0.) + smearedVertex;
00694
00695 if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
00696 }
00697 }
00698
00699
00700 bool testDaugh = false;
00701 unsigned int nDaughters = p.numberOfDaughters();
00702 if ( !testStable &&
00703
00704 nDaughters ) {
00705 for ( unsigned iDaughter=0; iDaughter<nDaughters; ++iDaughter ) {
00706 const reco::Candidate* daughter = p.daughter(iDaughter);
00707 if ( daughter->status()%1000==1 ) {
00708 testDaugh=true;
00709 break;
00710 }
00711 }
00712 }
00713
00714
00715 double dist = 0.;
00716 if ( !testStable && !testDaugh ) {
00717 XYZTLorentzVector productionVertex(p.vx(),p.vy(),p.vz(),0.);
00718 dist = (primaryVertex-productionVertex).Vect().Mag2();
00719 }
00720 bool testDecay = ( dist > 1e-8 ) ? true : false;
00721
00722
00723 if ( testStable || testDaugh || testDecay ) {
00724
00725 const reco::Candidate* mother = p.numberOfMothers() ? p.mother(0) : 0;
00726
00727 int originVertex =
00728 mother &&
00729 myGenVertices.find(mother) != myGenVertices.end() ?
00730 myGenVertices[mother] : mainVertex;
00731
00732 XYZTLorentzVector momentum(p.px(),p.py(),p.pz(),p.energy());
00733 RawParticle part(momentum, vertex(originVertex).position());
00734 part.setID(p.pdgId());
00735
00736
00737 int theTrack = addSimTrack(&part,originVertex, nGenParts()-offset);
00738
00739
00740 if ( !nDaughters ) continue;
00741 const reco::Candidate* daughter = p.daughter(0);
00742
00743
00744 XYZTLorentzVector decayVertex =
00745 XYZTLorentzVector(daughter->vx(), daughter->vy(),
00746 daughter->vz(), 0.) + smearedVertex;
00747 int theVertex = addSimVertex(decayVertex,theTrack);
00748
00749 if ( theVertex != -1 ) myGenVertices[&p] = theVertex;
00750
00751
00752 }
00753 }
00754
00755
00756
00757
00758 }
00759
00760 int
00761 FBaseSimEvent::addSimTrack(const RawParticle* p, int iv, int ig,
00762 const HepMC::GenVertex* ev) {
00763
00764
00765
00766 if ( !myFilter->accept(p) && ig >= -1 ) return -1;
00767
00768
00769 int trackId = nSimTracks++;
00770 if ( nSimTracks/theTrackSize*theTrackSize == nSimTracks ) {
00771 theTrackSize *= 2;
00772 theSimTracks->resize(theTrackSize);
00773 }
00774
00775
00776 vertex(iv).addDaughter(trackId);
00777 if ( !vertex(iv).noParent() ) {
00778 track(vertex(iv).parent().id()).addDaughter(trackId);
00779
00780 if ( ig == -1 ) {
00781 int motherId = track(vertex(iv).parent().id()).genpartIndex();
00782 if ( motherId < -1 ) ig = motherId;
00783 }
00784 }
00785
00786
00787 (*theSimTracks)[trackId] = ev ?
00788
00789 FSimTrack(p,iv,ig,trackId,this,
00790 ev->position().t()/10.
00791 * p->PDGmass()
00792 / std::sqrt(p->momentum().Vect().Mag2())) :
00793
00794 FSimTrack(p,iv,ig,trackId,this);
00795
00796 return trackId;
00797
00798 }
00799
00800 int
00801 FBaseSimEvent::addSimVertex(const XYZTLorentzVector& v,int im) {
00802
00803
00804 if ( !myFilter->accept(v) ) return -1;
00805
00806
00807 int vertexId = nSimVertices++;
00808 if ( nSimVertices/theVertexSize*theVertexSize == nSimVertices ) {
00809 theVertexSize *= 2;
00810 theSimVertices->resize(theVertexSize);
00811 }
00812
00813
00814 if ( im !=-1 ) track(im).setEndVertex(vertexId);
00815
00816
00817 (*theSimVertices)[vertexId] = FSimVertex(v,im,vertexId,this);
00818
00819 return vertexId;
00820
00821 }
00822
00823 void
00824 FBaseSimEvent::printMCTruth(const HepMC::GenEvent& myGenEvent) {
00825
00826 std::cout << "Id Gen Name eta phi pT E Vtx1 "
00827 << " x y z "
00828 << "Moth Vtx2 eta phi R Z Da1 Da2 Ecal?" << std::endl;
00829
00830 for ( HepMC::GenEvent::particle_const_iterator
00831 piter = myGenEvent.particles_begin();
00832 piter != myGenEvent.particles_end();
00833 ++piter ) {
00834
00835 HepMC::GenParticle* p = *piter;
00836
00837 int partId = p->pdg_id();
00838 std::string name;
00839
00840 if ( pdt->particle(ParticleID(partId)) !=0 ) {
00841 name = (pdt->particle(ParticleID(partId)))->name();
00842 } else {
00843 name = "none";
00844 }
00845
00846 XYZTLorentzVector momentum1(p->momentum().px(),
00847 p->momentum().py(),
00848 p->momentum().pz(),
00849 p->momentum().e());
00850
00851 int vertexId1 = 0;
00852
00853 if ( !p->production_vertex() ) continue;
00854
00855 XYZVector vertex1 (p->production_vertex()->position().x()/10.,
00856 p->production_vertex()->position().y()/10.,
00857 p->production_vertex()->position().z()/10.);
00858 vertexId1 = p->production_vertex()->barcode();
00859
00860 std::cout.setf(std::ios::fixed, std::ios::floatfield);
00861 std::cout.setf(std::ios::right, std::ios::adjustfield);
00862
00863 std::cout << std::setw(4) << p->barcode() << " "
00864 << name;
00865
00866 for(unsigned int k=0;k<11-name.length() && k<12; k++) std::cout << " ";
00867
00868 double eta = momentum1.eta();
00869 if ( eta > +10. ) eta = +10.;
00870 if ( eta < -10. ) eta = -10.;
00871 std::cout << std::setw(6) << std::setprecision(2) << eta << " "
00872 << std::setw(6) << std::setprecision(2) << momentum1.phi() << " "
00873 << std::setw(7) << std::setprecision(2) << momentum1.pt() << " "
00874 << std::setw(7) << std::setprecision(2) << momentum1.e() << " "
00875 << std::setw(4) << vertexId1 << " "
00876 << std::setw(6) << std::setprecision(1) << vertex1.x() << " "
00877 << std::setw(6) << std::setprecision(1) << vertex1.y() << " "
00878 << std::setw(6) << std::setprecision(1) << vertex1.z() << " ";
00879
00880 const HepMC::GenParticle* mother =
00881 *(p->production_vertex()->particles_in_const_begin());
00882
00883 if ( mother )
00884 std::cout << std::setw(4) << mother->barcode() << " ";
00885 else
00886 std::cout << " " ;
00887
00888 if ( p->end_vertex() ) {
00889 XYZTLorentzVector vertex2(p->end_vertex()->position().x()/10.,
00890 p->end_vertex()->position().y()/10.,
00891 p->end_vertex()->position().z()/10.,
00892 p->end_vertex()->position().t()/10.);
00893 int vertexId2 = p->end_vertex()->barcode();
00894
00895 std::vector<const HepMC::GenParticle*> children;
00896 HepMC::GenVertex::particles_out_const_iterator firstDaughterIt =
00897 p->end_vertex()->particles_out_const_begin();
00898 HepMC::GenVertex::particles_out_const_iterator lastDaughterIt =
00899 p->end_vertex()->particles_out_const_end();
00900 for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
00901 children.push_back(*firstDaughterIt);
00902 }
00903
00904 std::cout << std::setw(4) << vertexId2 << " "
00905 << std::setw(6) << std::setprecision(2) << vertex2.eta() << " "
00906 << std::setw(6) << std::setprecision(2) << vertex2.phi() << " "
00907 << std::setw(5) << std::setprecision(1) << vertex2.pt() << " "
00908 << std::setw(6) << std::setprecision(1) << vertex2.z() << " ";
00909 for ( unsigned id=0; id<children.size(); ++id )
00910 std::cout << std::setw(4) << children[id]->barcode() << " ";
00911 }
00912 std::cout << std::endl;
00913
00914 }
00915
00916 }
00917
00918 void
00919 FBaseSimEvent::print() const {
00920
00921 std::cout << " Id Gen Name eta phi pT E Vtx1 "
00922 << " x y z "
00923 << "Moth Vtx2 eta phi R Z Daughters Ecal?" << std::endl;
00924
00925 for( int i=0; i<(int)nTracks(); i++ )
00926 std::cout << track(i) << std::endl;
00927 }
00928
00929 void
00930 FBaseSimEvent::clear() {
00931
00932 nSimTracks = 0;
00933 nSimVertices = 0;
00934 nGenParticles = 0;
00935 nChargedParticleTracks = 0;
00936
00937 }
00938
00939 void
00940 FBaseSimEvent::addChargedTrack(int id) {
00941 (*theChargedTracks)[nChargedParticleTracks++] = id;
00942 if ( nChargedParticleTracks/theChargedSize*theChargedSize
00943 == nChargedParticleTracks ) {
00944 theChargedSize *= 2;
00945 theChargedTracks->resize(theChargedSize);
00946 }
00947 }
00948
00949 int
00950 FBaseSimEvent::chargedTrack(int id) const {
00951 if (id>=0 && id<(int)nChargedParticleTracks)
00952 return (*theChargedTracks)[id];
00953 else
00954 return -1;
00955 }
00956
00957
00958
00959
00960
00961
00962
00963
00964
00965
00966
00967
00968
00969 const HepMC::GenParticle*
00970 FBaseSimEvent::embdGenpart(int i) const {
00971 return (*theGenParticles)[i];
00972 }
00973
00974
00975
00976
00977
00978
00979
00980
00981
00982
00983
00984
00985