38 #include <ext/algorithm>
87 std::vector<int> &genBarcodes,
88 bool &barcodesAreSorted)
const;
107 setStatus_(cfg.getParameter<int32_t>(
"setStatus")),
109 cfg.getParameter<
InputTag>(
"genParticles"))),
110 genBarcodesToken_(consumes<std::
vector<int>>(
111 cfg.getParameter<
InputTag>(
"genParticles")))
114 if (cfg.
exists(
"particleTypes")) {
116 if (!
pdts_.empty()) {
122 if (cfg.
existsAs<
string>(
"filter")) {
124 if (!filter.empty()) {
128 produces<GenParticleCollection>();
129 produces<vector<int>>();
139 std::vector<int> &genBarcodes,
140 bool &barcodesAreSorted)
const {
147 vtx = simvertices[stDau.
vertIndex()].position();
151 if (
filter_.get() !=
nullptr) {
157 genp.addMother(parentRef);
158 mergedGens.push_back(
genp);
160 unsigned int dauidx = mergedGens.size() - 1;
168 for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin(); isimtrk != simtracksSorted.end();
170 if (!isimtrk->noVertex()) {
172 const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
178 SimTrackContainer::const_iterator it =
180 if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
181 if (it->trackId() == stDau.
trackId()) {
184 stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
194 if (!
pdts_.empty()) {
197 for (vector<PdtEntry>::iterator itp =
pdts_.begin(), edp =
pdts_.end(); itp != edp; ++itp) {
211 std::unique_ptr<SimTrackContainer> simtracksTmp;
213 if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(),
LessById())) {
214 simtracksTmp = std::make_unique<SimTrackContainer>(*simtracks);
215 std::sort(simtracksTmp->begin(), simtracksTmp->end(),
LessById());
216 simtracksSorted = &*simtracksTmp;
226 bool barcodesAreSorted =
true;
229 if (gens->size() != genBarcodes->size())
230 throw cms::Exception(
"Corrupt data") <<
"Barcodes not of the same size as GenParticles!\n";
233 auto candsPtr = std::make_unique<GenParticleCollection>();
239 for (
size_t i = 0;
i < gens->size(); ++
i) {
241 cands.push_back(cand);
245 auto newGenBarcodes = std::make_unique<vector<int>>();
246 for (
unsigned int i = 0;
i < genBarcodes->size(); ++
i) {
247 newGenBarcodes->push_back((*genBarcodes)[
i]);
249 barcodesAreSorted = __gnu_cxx::is_sorted(newGenBarcodes->begin(), newGenBarcodes->end());
251 for (
size_t i = 0;
i < cands.size(); ++
i) {
256 for (
size_t d = 0;
d < nDaus; ++
d) {
263 for (
size_t m = 0;
m < nMoms; ++
m) {
268 for (SimTrackContainer::const_iterator isimtrk = simtracks->begin(); isimtrk != simtracks->end(); ++isimtrk) {
270 if (isimtrk->genpartIndex() != -1)
281 if (!isimtrk->noVertex()) {
283 const SimVertex &vtx = (*simvertices)[isimtrk->vertIndex()];
289 SimTrackContainer::const_iterator it =
291 if ((it != simtracksSorted->end()) && (it->trackId() == idx)) {
292 if (it->genpartIndex() != -1) {
293 std::vector<int>::const_iterator itIndex;
294 if (barcodesAreSorted) {
295 itIndex =
std::lower_bound(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
297 itIndex =
std::find(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
302 if ((itIndex != genBarcodes->end()) && (*itIndex == it->genpartIndex())) {
305 unsigned int momidx = itIndex - genBarcodes->begin();
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool operator()(const SimTrack &tk1, const SimTrack &tk2) const
bool existsAs(std::string const ¶meterName, bool trackiness=true) const
checks if a parameter exists as a given type
edm::EDGetTokenT< edm::SimVertexContainer > simverticesToken_
uint16_t *__restrict__ id
#define DEFINE_FWK_MODULE(type)
std::vector< PdtEntry > pdts_
edm::EDGetTokenT< std::vector< int > > genBarcodesToken_
size_t numberOfDaughters() const override
number of daughters
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_
bool exists(std::string const ¶meterName) const
checks if a parameter exists
~GenPlusSimParticleProducer() override
const daughters & daughterRefVector() const
references to daughtes
std::unique_ptr< StrFilter > filter_
size_t numberOfMothers() const override
number of mothers
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
float charge() const
charge
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
bool getData(T &iHolder) const
void resetDaughters(const edm::ProductID &id)
set daughters product ID
StringCutObjectSelector< reco::GenParticle > StrFilter
const mothers & motherRefVector() const
references to mothers
tuple key
prepare the HTCondor submission files and eventually submit them
Abs< T >::type abs(const T &t)
bool operator()(unsigned int id, const SimTrack &tk2) const
GenPlusSimParticleProducer(const edm::ParameterSet &)
Produces reco::GenParticle from SimTracks.
tuple genp
produce generated paricles in acceptance #
void addMother(const typename mothers::value_type &)
add a daughter via a reference
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
bool operator()(const SimTrack &tk1, unsigned int id) const
unsigned int trackId() const
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
edm::EDGetTokenT< reco::GenParticleCollection > gensToken_
Collection of GenParticles I need to make refs to. It must also have its associated vector<int> of ba...
std::vector< SimVertex > SimVertexContainer
T getParameter(std::string const &) const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
void produce(edm::Event &, const edm::EventSetup &) override
int type() const
particle type (HEP PDT convension)
__host__ __device__ constexpr RandomIt lower_bound(RandomIt first, RandomIt last, const T &value, Compare comp={})
const math::XYZTLorentzVectorD & momentum() const
ProductID id() const
Accessor for product ID.
void addGenParticle(const SimTrack &stMom, const SimTrack &stDau, unsigned int momidx, const edm::SimTrackContainer &simtks, const edm::SimVertexContainer &simvtxs, reco::GenParticleCollection &mergedGens, const reco::GenParticleRefProd ref, std::vector< int > &genBarcodes, bool &barcodesAreSorted) const
Try to link the GEANT particle to the generator particle it came from.
math::XYZPoint Point
point in the space
void resetMothers(const edm::ProductID &id)
set mother product ID
std::vector< SimTrack > SimTrackContainer
edm::EDGetTokenT< edm::SimTrackContainer > simtracksToken_