9 #include "HepPDT/ParticleID.hh"
14 #define update(a, b) do { (a) = (a) | (b); } while(0)
20 tracer_(config,std::
move(collector)),
21 hepMCLabel_( config.getUntrackedParameter<edm::
InputTag>(
"hepMC") )
130 VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin();
131 ivertex != genVertexTrail.end();
136 HepMC::GenVertex * vertex =
const_cast<HepMC::GenVertex *
>(*ivertex);
140 HepMC::GenVertex::particle_iterator iparent = vertex->particles_begin(
HepMC::parents);
146 int pdgid =
std::abs((*iparent)->pdg_id());
148 HepPDT::ParticleID particleID(pdgid);
151 if (particleID.isValid())
187 VertexHistory::SimVertexTrail::const_iterator ivertex = simVertexTrail.begin();
188 ivertex != simVertexTrail.end();
199 for (its = (*ivertex)->sourceTracks_begin(); its != (*ivertex)->sourceTracks_end(); ++its)
201 for (itd = (*ivertex)->daughterTracks_begin(); itd != (*ivertex)->daughterTracks_end(); ++itd)
211 if ( its != (*ivertex)->sourceTracks_end() )
218 unsigned int processG4 = 0;
220 if((*ivertex)->nG4Vertices() > 0) {
221 processG4 = (*(*ivertex)->g4Vertices_begin()).processType();
256 iparticle != (*ivertex)->daughterTracks_end();
261 if ( (*iparticle)->numberOfTrackerLayers() )
268 HepPDT::ParticleID particleID(pdgid);
270 if (particleID.isValid())
306 typedef std::multimap<double, HepMC::ThreeVector> Clusters;
307 typedef std::pair<double, HepMC::ThreeVector> ClusterPair;
318 double const mm = 0.1;
322 VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin();
323 ivertex != genVertexTrail.end();
331 HepMC::ThreeVector
p = (*ivertex)->point3d();
332 double distance =
sqrt(
pow(p.x() * mm - genpv.
x, 2) +
pow(p.y() * mm - genpv.
y, 2) +
pow(p.z() * mm - genpv.
z, 2) );
335 if ( clusters.empty() )
337 clusters.insert( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
346 clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
350 bool cluster =
false;
358 double difference =
sqrt (
359 pow(p.x() * mm - icluster->second.x(), 2) +
360 pow(p.y() * mm - icluster->second.y(), 2) +
361 pow(p.z() * mm - icluster->second.z(), 2)
371 if (!cluster) clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
379 VertexHistory::SimVertexTrail::const_reverse_iterator ivertex = simVertexTrail.rbegin();
380 ivertex != simVertexTrail.rend();
390 if ( clusters.empty() )
392 clusters.insert( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
401 clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
405 bool cluster =
false;
413 double difference =
sqrt (
414 pow(p.x() - icluster->second.x(), 2) +
415 pow(p.y() - icluster->second.y(), 2) +
416 pow(p.z() - icluster->second.z(), 2)
426 if (!cluster) clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
429 if ( clusters.size() == 1 )
431 else if ( clusters.size() == 2 )
440 return !p->end_vertex() && p->status() == 1;
448 return part->charge()!=0;
468 for ( HepMC::GenEvent::vertex_const_iterator ivertex =
event->vertices_begin(); ivertex !=
event->vertices_end(); ++ivertex )
470 bool hasParentVertex =
false;
474 HepMC::GenVertex::particle_iterator iparent = (*ivertex)->particles_begin(
HepMC::parents);
478 if ( (*iparent)->production_vertex() )
480 hasParentVertex =
true;
485 if (hasParentVertex)
continue;
488 HepMC::FourVector pos = (*ivertex)->position();
490 double const mm = 0.1;
494 std::vector<GeneratedPrimaryVertex>::iterator ientry =
genpvs_.begin();
497 for (; ientry !=
genpvs_.end(); ++ientry)
509 ientry->genVertex.push_back((*ivertex)->barcode());
513 HepMC::GenVertex::particle_iterator idecendants = (*ivertex)->particles_begin(HepMC::descendants);
514 idecendants != (*ivertex)->particles_end(HepMC::descendants);
519 if (
find(ientry->finalstateParticles.begin(), ientry->finalstateParticles.end(), (*idecendants)->barcode()) == ientry->finalstateParticles.end() )
521 ientry->finalstateParticles.push_back((*idecendants)->barcode());
522 HepMC::FourVector
m = (*idecendants)->momentum();
524 ientry->ptot.setPx(ientry->ptot.px() + m.px());
525 ientry->ptot.setPy(ientry->ptot.py() + m.py());
526 ientry->ptot.setPz(ientry->ptot.pz() + m.pz());
527 ientry->ptot.setE(ientry->ptot.e() + m.e());
528 ientry->ptsq += m.perp() * m.perp();
530 if ( m.perp() > 0.8 &&
std::abs(m.pseudoRapidity()) < 2.5 &&
isCharged(*idecendants) ) ientry->nGenTrk++;
edm::Handle< edm::HepMCProduct > mcInformation_
T getUntrackedParameter(std::string const &, T const &) const
int event() const
get the contents of the subdetector field (should be protected?)
void newEvent(const edm::Event &, const edm::EventSetup &)
Pre-process event information (for accessing reconstruction information)
void vertexInformation()
Get geometrical information about the vertices.
std::vector< GeneratedPrimaryVertex > genpvs_
const edm::InputTag hepMCLabel_
SimVertexTrail const & simVertexTrail() const
Return all the simulated vertices in the history.
bool isCharged(const HepMC::GenParticle *)
Get track history and classify it in function of their .
void reset()
Reset the categories flags.
bool isNonnull() const
Checks for non-null.
const reco::VertexBaseRef & recoVertex() const
Return a reference to the reconstructed track.
bool isFinalstateParticle(const HepMC::GenParticle *)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
void genPrimaryVertices()
void getData(T &iHolder) const
math::XYZTLorentzVectorD LorentzVector
bool evaluate(TrackingVertexRef tvr)
Evaluate track history using a TrackingParticleRef.
Flags flags_
Flag containers.
int bunchCrossing() const
get the detector field from this detid
double longLivedDecayLength_
Abs< T >::type abs(const T &t)
std::vector< const HepMC::GenVertex * > GenVertexTrail
GenVertex trail type.
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
edm::ESHandle< ParticleDataTable > particleDataTable_
const G4toCMSLegacyProcTypeMap g4toCMSProcMap_
const unsigned int processId(unsigned int g4ProcessId) const
HepPDT::ParticleData ParticleData
const TrackingVertexRef & simVertex() const
Return the initial tracking vertex from the history.
VertexClassifier(edm::ParameterSet const &pset, edm::ConsumesCollector &&)
Constructor by ParameterSet.
void processesAtGenerator()
Get all the information related to decay process.
VertexClassifier const & evaluate(reco::VertexBaseRef const &)
Classify the RecoVertex in categories.
void depth(int d)
Set the depth of the history.
double vertexClusteringDistance_
virtual void newEvent(edm::Event const &, edm::EventSetup const &)
Pre-process event information (for accessing reconstraction information)
void processesAtSimulation()
Get information about conversion and other interactions.
GenVertexTrail const & genVertexTrail() const
Return all generated vertex in the history.
std::vector< TrackingVertexRef > SimVertexTrail
SimVertex trail type.
Auxiliary class holding simulated primary vertices.
Power< A, B >::type pow(const A &a, const B &b)
void simulationInformation()
Get all the information related to the simulation details.