CMS 3D CMS Logo

EcalGenEvtSelectorFrag.cc
Go to the documentation of this file.
1 #include <iostream>
4 #include "HepMC/GenVertex.h"
5 using namespace std;
6 
8  partonId_ = pset.getParameter<vector<int> >("partons");
9  partonStatus_ = pset.getParameter<vector<int> >("partonStatus");
10  partonPt_ = pset.getParameter<vector<double> >("partonPt");
11 
12  particleId_ = pset.getParameter<vector<int> >("particles");
13  particleStatus_ = pset.getParameter<vector<int> >("particleStatus");
14  particlePt_ = pset.getParameter<vector<double> >("particlePt");
15 
16  etaMax_ = pset.getParameter<double>("etaMax");
17 
18  int id = partonId_.size();
19  int st = partonStatus_.size();
20  int pt = partonPt_.size();
21 
22  if (partonId_.size() != partonStatus_.size() || partonId_.size() != partonPt_.size()) {
23  throw edm::Exception(edm::errors::LogicError) << id << st << pt << endl;
24  }
25 
26  id = particleId_.size();
27  st = particleStatus_.size();
28  pt = particlePt_.size();
29 
30  if (particleId_.size() != particleStatus_.size() || particleId_.size() != particlePt_.size()) {
31  throw edm::Exception(edm::errors::LogicError) << id << st << pt << endl;
32  }
33 }
34 
36  HepMC::GenEvent::particle_const_iterator begin = evt->particles_begin();
37  HepMC::GenEvent::particle_const_iterator end = evt->particles_end();
38 
39  bool foundParticle = false;
40  bool foundParton = false;
41 
42  HepMC::GenEvent::particle_const_iterator it = begin;
43  while ((!foundParton || !foundParticle) && it != end) {
44  bool isFrag = false;
45  bool isThisPhoton = false;
46 
47  foundParton = true;
48  /*for(unsigned i = 0; i < partonId_.size(); ++i){
49  if(selectParticle(*it, partonStatus_[i], partonId_[i], partonPt_[i], etaMax_)) foundParton = true;
50  }*/
51 
52  for (unsigned i = 0; i < particleId_.size(); ++i) {
54  isThisPhoton = true;
55  }
56 
57  // Now you have to requre the partcile is "prompt", meaning its mom is parton
58 
59  if (!((*it)->production_vertex())) {
60  isFrag = false;
61  } else {
62  const HepMC::GenVertex* productionVertex = (*it)->production_vertex();
63 
64  size_t numberOfMothers = productionVertex->particles_in_size();
65  if (numberOfMothers <= 0) {
66  isFrag = false;
67  // cout << "number of mothers = " << numberOfMothers << endl;
68  } else {
69  // cout << "number of mothers = " << numberOfMothers << endl;
70  HepMC::GenVertex::particles_in_const_iterator motherIt = productionVertex->particles_in_const_begin();
71  for (; motherIt != productionVertex->particles_in_const_end(); motherIt++) {
72  if (fabs((*motherIt)->pdg_id()) <= 22) {
73  isFrag = true;
74  }
75  }
76  }
77  }
78 
79  if ((isFrag == true) && (isThisPhoton == true)) {
80  //cout << "Found fragmentation photon!!" << endl ;
81  foundParticle = true;
82  }
83 
84  ++it;
85  }
86 
87  foundParton = true; // We don't care of the parton
88  return (foundParton && foundParticle);
89 }
EcalGenEvtSelectorFrag::particlePt_
std::vector< double > particlePt_
Definition: EcalGenEvtSelectorFrag.h:21
EcalGenEvtSelectorFrag::particleId_
std::vector< int > particleId_
Definition: EcalGenEvtSelectorFrag.h:19
mps_fire.i
i
Definition: mps_fire.py:355
edm::errors::LogicError
Definition: EDMException.h:37
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
EcalGenEvtSelectorFrag.h
EcalGenEvtSelectorFrag::particleStatus_
std::vector< int > particleStatus_
Definition: EcalGenEvtSelectorFrag.h:20
EcalGenEvtSelectorFrag::partonPt_
std::vector< double > partonPt_
Definition: EcalGenEvtSelectorFrag.h:17
end
#define end
Definition: vmac.h:39
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
EDMException.h
BaseHiGenEvtSelector::selectParticle
bool selectParticle(HepMC::GenParticle *par, int status, int pdg, double ptMin, double etaMax)
Definition: BaseHiGenEvtSelector.h:13
EcalGenEvtSelectorFrag::etaMax_
double etaMax_
Definition: EcalGenEvtSelectorFrag.h:23
edm::ParameterSet
Definition: ParameterSet.h:36
BaseHiGenEvtSelector
Definition: BaseHiGenEvtSelector.h:8
std
Definition: JetResolutionObject.h:76
Exception
Definition: hltDiff.cc:246
EcalGenEvtSelectorFrag::EcalGenEvtSelectorFrag
EcalGenEvtSelectorFrag(const edm::ParameterSet &pset)
Definition: EcalGenEvtSelectorFrag.cc:7
EcalGenEvtSelectorFrag::partonStatus_
std::vector< int > partonStatus_
Definition: EcalGenEvtSelectorFrag.h:16
EcalGenEvtSelectorFrag::partonId_
std::vector< int > partonId_
Definition: EcalGenEvtSelectorFrag.h:15
begin
#define begin
Definition: vmac.h:32
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
EcalGenEvtSelectorFrag::filter
bool filter(HepMC::GenEvent *) override
Definition: EcalGenEvtSelectorFrag.cc:35