Skip to content

Commit

Permalink
Merge pull request #40415 from michael-pitt/CMSSW_10_6_X
Browse files Browse the repository at this point in the history
Add TaggedProton filter for analysis with tagged protons [10_6_X]
  • Loading branch information
cmsbuild committed Feb 7, 2023
2 parents 79ea630 + 0ba2b4d commit d3d76d3
Show file tree
Hide file tree
Showing 3 changed files with 53 additions and 2 deletions.
24 changes: 24 additions & 0 deletions GeneratorInterface/Core/interface/TaggedProtonHepMCFilter.h
@@ -0,0 +1,24 @@
#ifndef __TAGGEDPROTONHEPMCFILTER__
#define __TAGGEDPROTONHEPMCFILTER__

#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include "GeneratorInterface/Core/interface/BaseHepMCFilter.h"

class TaggedProtonHepMCFilter : public BaseHepMCFilter {
private:
const int proton_PDGID_ = 2212;
const int neutron_PDGID_ = 2112;
const double xiMin_;
const double xiMax_;
const double oneOverbeamEnergy_;
const int nProtons_;

public:
explicit TaggedProtonHepMCFilter(const edm::ParameterSet &);
~TaggedProtonHepMCFilter() override = default;

bool filter(const HepMC::GenEvent *evt) override;
};

#endif
7 changes: 5 additions & 2 deletions GeneratorInterface/Core/src/HepMCFilterDriver.cc
Expand Up @@ -3,6 +3,7 @@
#include "GeneratorInterface/Core/interface/PartonShowerBsHepMCFilter.h"
#include "GeneratorInterface/Core/interface/PartonShowerCsHepMCFilter.h"
#include "GeneratorInterface/Core/interface/EmbeddingHepMCFilter.h"
#include "GeneratorInterface/Core/interface/TaggedProtonHepMCFilter.h"
#include "GeneratorInterface/Core/interface/PythiaHepMCFilterGammaGamma.h"

#include "FWCore/MessageLogger/interface/MessageLogger.h"
Expand Down Expand Up @@ -31,10 +32,12 @@ HepMCFilterDriver::HepMCFilterDriver(const edm::ParameterSet& pset) :
else if (filterName=="PartonShowerCsHepMCFilter") {
filter_ = new PartonShowerCsHepMCFilter(filterParameters);
}

else if (filterName == "TaggedProtonHepMCFilter") {
filter_ = new TaggedProtonHepMCFilter(filterParameters);
}
else if (filterName=="EmbeddingHepMCFilter") {
filter_ = new EmbeddingHepMCFilter(filterParameters);
}
}
else if (filterName=="PythiaHepMCFilterGammaGamma") {
filter_ = new PythiaHepMCFilterGammaGamma(filterParameters);
}
Expand Down
24 changes: 24 additions & 0 deletions GeneratorInterface/Core/src/TaggedProtonHepMCFilter.cc
@@ -0,0 +1,24 @@
#include "GeneratorInterface/Core/interface/TaggedProtonHepMCFilter.h"

TaggedProtonHepMCFilter::TaggedProtonHepMCFilter(const edm::ParameterSet &iConfig)
: xiMin_(iConfig.getParameter<double>("xiMin")),
xiMax_(iConfig.getParameter<double>("xiMax")),
oneOverbeamEnergy_(2.0 / iConfig.getParameter<double>("comEnergy")),
nProtons_(iConfig.getParameter<int>("nProtons")) {}

bool TaggedProtonHepMCFilter::filter(const HepMC::GenEvent *evt) {
// Going through the particle list, and count good protons
int nGoodProtons = 0;
for (HepMC::GenEvent::particle_const_iterator particle = evt->particles_begin(); particle != evt->particles_end();
++particle) {
if ((*particle)->pdg_id() == proton_PDGID_ && 1 == (*particle)->status()) {
HepMC::FourVector p4 = (*particle)->momentum();
double xi = (1.0 - std::abs(p4.pz()) * oneOverbeamEnergy_);
if (xi > xiMin_ && xi < xiMax_)
nGoodProtons++;
if (nGoodProtons >= nProtons_)
return true;
}
}
return false;
}

0 comments on commit d3d76d3

Please sign in to comment.