CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
BeamDivergenceVtxGenerator Class Reference
Inheritance diagram for BeamDivergenceVtxGenerator:
edm::stream::EDProducer<>

Classes

struct  SmearingParameters
 

Public Member Functions

 BeamDivergenceVtxGenerator (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~BeamDivergenceVtxGenerator () override=default
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &)
 

Private Member Functions

void addSmearedGenParticle (const reco::GenParticle &gp, const SmearingParameters &sp, HepMC::GenEvent *genEvt)
 
void applySmearingHepMC (const SmearingParameters &sp, HepMC::GenEvent *genEvt)
 

Static Private Member Functions

template<typename T >
static HepMC::FourVector smearMomentum (const T &mom, const SmearingParameters &sp)
 

Private Attributes

edm::ESGetToken
< CTPPSBeamParameters,
CTPPSBeamParametersRcd
beamParametersToken_
 
std::vector< edm::EDGetTokenT
< reco::GenParticleCollection > > 
genParticleTokens_
 
bool simulateBeamDivergence_
 
bool simulateVertex_
 
edm::EDGetTokenT
< edm::HepMCProduct
sourceToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T...>
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T...>
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 29 of file BeamDivergenceVtxGenerator.cc.

Constructor & Destructor Documentation

BeamDivergenceVtxGenerator::BeamDivergenceVtxGenerator ( const edm::ParameterSet iConfig)
explicit

Definition at line 61 of file BeamDivergenceVtxGenerator.cc.

References genParticleTokens_, edm::ParameterSet::getParameter(), edm::Service< T >::isAvailable(), edm::InputTag::label(), and sourceToken_.

62  : beamParametersToken_(esConsumes<CTPPSBeamParameters, CTPPSBeamParametersRcd>()),
63  simulateVertex_(iConfig.getParameter<bool>("simulateVertex")),
64  simulateBeamDivergence_(iConfig.getParameter<bool>("simulateBeamDivergence")) {
65  const edm::InputTag tagSrcHepMC = iConfig.getParameter<edm::InputTag>("src");
66  if (!tagSrcHepMC.label().empty())
67  sourceToken_ = consumes<edm::HepMCProduct>(tagSrcHepMC);
68 
69  for (const auto &it : iConfig.getParameter<std::vector<edm::InputTag>>("srcGenParticle"))
70  genParticleTokens_.push_back(consumes<reco::GenParticleCollection>(it));
71 
73  if (!rng.isAvailable())
74  throw cms::Exception("Configuration")
75  << "The BeamDivergenceVtxGenerator requires the RandomNumberGeneratorService\n"
76  "which is not present in the configuration file. \n"
77  "You must add the service\n"
78  "in the configuration file or remove the modules that require it.";
79 
80  produces<edm::HepMCProduct>();
81 }
std::vector< edm::EDGetTokenT< reco::GenParticleCollection > > genParticleTokens_
edm::ESGetToken< CTPPSBeamParameters, CTPPSBeamParametersRcd > beamParametersToken_
bool isAvailable() const
Definition: Service.h:40
edm::EDGetTokenT< edm::HepMCProduct > sourceToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::string const & label() const
Definition: InputTag.h:36
BeamDivergenceVtxGenerator::~BeamDivergenceVtxGenerator ( )
overridedefault

Member Function Documentation

void BeamDivergenceVtxGenerator::addSmearedGenParticle ( const reco::GenParticle gp,
const SmearingParameters sp,
HepMC::GenEvent genEvt 
)
private

Definition at line 204 of file BeamDivergenceVtxGenerator.cc.

References GenParticle::GenParticle, reco::LeafCandidate::p4(), reco::LeafCandidate::pdgId(), smearMomentum(), reco::LeafCandidate::status(), BeamDivergenceVtxGenerator::SmearingParameters::vtx_t, BeamDivergenceVtxGenerator::SmearingParameters::vtx_x, BeamDivergenceVtxGenerator::SmearingParameters::vtx_y, BeamDivergenceVtxGenerator::SmearingParameters::vtx_z, reco::LeafCandidate::vx(), reco::LeafCandidate::vy(), and reco::LeafCandidate::vz().

Referenced by produce().

206  {
207  // add vertex of the particle
208  HepMC::GenVertex *vtx = new HepMC::GenVertex(HepMC::FourVector(
209  (gp.vx() + sp.vtx_x) * 1E1, // conversion: cm to mm
210  (gp.vy() + sp.vtx_y) * 1E1,
211  (gp.vz() + sp.vtx_z) * 1E1,
212  (/*gp.vt()*/ +sp.vtx_t) * 1E1)); // TODO: GenParticle doesn't seem to have time component of the vertex
213  genEvt->add_vertex(vtx);
214 
215  // add the particle itself
216  HepMC::GenParticle *particle = new HepMC::GenParticle(smearMomentum(gp.p4(), sp), gp.pdgId(), gp.status());
217  vtx->add_particle_out(particle);
218 }
double vz() const override
z coordinate of vertex position
static HepMC::FourVector smearMomentum(const T &mom, const SmearingParameters &sp)
double vy() const override
y coordinate of vertex position
int status() const final
status word
const LorentzVector & p4() const final
four-momentum Lorentz vector
int pdgId() const final
PDG identifier.
double vx() const override
x coordinate of vertex position
void BeamDivergenceVtxGenerator::applySmearingHepMC ( const SmearingParameters sp,
HepMC::GenEvent genEvt 
)
private

Definition at line 185 of file BeamDivergenceVtxGenerator.cc.

References simulateBeamDivergence_, simulateVertex_, smearMomentum(), BeamDivergenceVtxGenerator::SmearingParameters::vtx_t, BeamDivergenceVtxGenerator::SmearingParameters::vtx_x, BeamDivergenceVtxGenerator::SmearingParameters::vtx_y, and BeamDivergenceVtxGenerator::SmearingParameters::vtx_z.

Referenced by produce().

185  {
186  if (simulateVertex_) {
187  for (HepMC::GenEvent::vertex_iterator vit = genEvt->vertices_begin(); vit != genEvt->vertices_end(); ++vit) {
188  const auto &pos = (*vit)->position();
189  (*vit)->set_position(HepMC::FourVector(pos.x() + sp.vtx_x * 1E1, // conversion: cm to mm
190  pos.y() + sp.vtx_y * 1E1,
191  pos.z() + sp.vtx_z * 1E1,
192  pos.t() + sp.vtx_t * 1E1));
193  }
194  }
195 
197  for (HepMC::GenEvent::particle_iterator part = genEvt->particles_begin(); part != genEvt->particles_end(); ++part)
198  (*part)->set_momentum(smearMomentum((*part)->momentum(), sp));
199  }
200 }
static HepMC::FourVector smearMomentum(const T &mom, const SmearingParameters &sp)
part
Definition: HCALResponse.h:20
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
Definition: VertexFwd.h:19
void BeamDivergenceVtxGenerator::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 85 of file BeamDivergenceVtxGenerator.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), submitPVResolutionJobs::desc, HLT_FULL_cff::InputTag, and edm::ParameterDescriptionNode::setComment().

85  {
87 
88  desc.add<edm::InputTag>("src", edm::InputTag("generator", "unsmeared"))
89  ->setComment("input collection in HepMC format");
90 
91  desc.add<std::vector<edm::InputTag>>("srcGenParticle", std::vector<edm::InputTag>())
92  ->setComment("input collections in GenParticle format");
93 
94  desc.add<bool>("simulateBeamDivergence", true)->setComment("account for the beam angular divergence?");
95  desc.add<bool>("simulateVertex", true)->setComment("account for the vertex transverse smearing?");
96 
97  descriptions.add("beamDivergenceVtxGenerator", desc);
98 }
void setComment(std::string const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void BeamDivergenceVtxGenerator::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 102 of file BeamDivergenceVtxGenerator.cc.

References addSmearedGenParticle(), applySmearingHepMC(), BeamDivergenceVtxGenerator::SmearingParameters::bd_x_45, BeamDivergenceVtxGenerator::SmearingParameters::bd_x_56, BeamDivergenceVtxGenerator::SmearingParameters::bd_y_45, BeamDivergenceVtxGenerator::SmearingParameters::bd_y_56, beamParametersToken_, TtGenEvtProducer_cfi::genEvt, genParticleTokens_, edm::Event::getByToken(), edm::RandomNumberGenerator::getEngine(), edm::EventSetup::getHandle(), runTauDisplay::gp, edm::EDGetTokenT< T >::isUninitialized(), eostools::move(), convertSQLitetoXML_cfg::output, edm::Event::put(), simulateBeamDivergence_, simulateVertex_, sourceToken_, edm::Event::streamID(), unpackBuffers-CaloStage2::token, BeamDivergenceVtxGenerator::SmearingParameters::vtx_t, BeamDivergenceVtxGenerator::SmearingParameters::vtx_x, BeamDivergenceVtxGenerator::SmearingParameters::vtx_y, and BeamDivergenceVtxGenerator::SmearingParameters::vtx_z.

102  {
103  // get random engine
105  CLHEP::HepRandomEngine *rnd = &(rng->getEngine(iEvent.streamID()));
106 
107  // get conditions
109 
110  // get HepMC input (if given)
113  genEvt = new HepMC::GenEvent();
114  } else {
115  edm::Handle<edm::HepMCProduct> hepUnsmearedMCEvt;
116  iEvent.getByToken(sourceToken_, hepUnsmearedMCEvt);
117 
118  genEvt = new HepMC::GenEvent(*hepUnsmearedMCEvt->GetEvent());
119  }
120 
121  // prepare output
122  std::unique_ptr<edm::HepMCProduct> output(new edm::HepMCProduct(genEvt));
123 
124  // generate smearing parameters
125  SmearingParameters sp;
126 
127  if (simulateVertex_) {
128  // NB: the separtion between effective offsets in LHC sectors 45 and 56 cannot be applied, thus the values for 45 are used
129  sp.vtx_x = hBeamParameters->getVtxOffsetX45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevX();
130  sp.vtx_y = hBeamParameters->getVtxOffsetY45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevY();
131  sp.vtx_z = hBeamParameters->getVtxOffsetZ45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevZ();
132  sp.vtx_t = hBeamParameters->getVtxOffsetT45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevT();
133  }
134 
136  sp.bd_x_45 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceX45();
137  sp.bd_x_56 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceX56();
138  sp.bd_y_45 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceY45();
139  sp.bd_y_56 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceY56();
140  }
141 
142  // apply smearing
143  applySmearingHepMC(sp, genEvt);
144 
145  for (const auto &token : genParticleTokens_) {
147  iEvent.getByToken(token, hGPCollection);
148 
149  for (const auto &gp : *hGPCollection)
150  addSmearedGenParticle(gp, sp, genEvt);
151  }
152 
153  // save output
154  iEvent.put(std::move(output));
155 }
std::vector< edm::EDGetTokenT< reco::GenParticleCollection > > genParticleTokens_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:99
edm::ESGetToken< CTPPSBeamParameters, CTPPSBeamParametersRcd > beamParametersToken_
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
def move
Definition: eostools.py:511
edm::EDGetTokenT< edm::HepMCProduct > sourceToken_
void applySmearingHepMC(const SmearingParameters &sp, HepMC::GenEvent *genEvt)
void addSmearedGenParticle(const reco::GenParticle &gp, const SmearingParameters &sp, HepMC::GenEvent *genEvt)
StreamID streamID() const
Definition: Event.h:98
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
template<typename T >
HepMC::FourVector BeamDivergenceVtxGenerator::smearMomentum ( const T mom,
const SmearingParameters sp 
)
staticprivate

Definition at line 160 of file BeamDivergenceVtxGenerator.cc.

References BeamDivergenceVtxGenerator::SmearingParameters::bd_x_45, BeamDivergenceVtxGenerator::SmearingParameters::bd_x_56, BeamDivergenceVtxGenerator::SmearingParameters::bd_y_45, BeamDivergenceVtxGenerator::SmearingParameters::bd_y_56, AlCaHLTBitMon_ParallelJobs::p, jetcorrextractor::sign(), and mathSSE::sqrt().

Referenced by addSmearedGenParticle(), and applySmearingHepMC().

160  {
161  // TODO: this is an oversimplified implemetation
162  // the TOTEM smearing module should be taken as reference
163 
164  double th_x = mom.x() / mom.z();
165  double th_y = mom.y() / mom.z();
166 
167  if (mom.z() > 0.) {
168  th_x += sp.bd_x_45;
169  th_y += sp.bd_y_45;
170  } else {
171  th_x += sp.bd_x_56;
172  th_y += sp.bd_y_56;
173  }
174 
175  // calculate consistent p_z component
176  const double sign = (mom.z() > 0.) ? 1. : -1.;
177  const double p = sqrt(mom.x() * mom.x() + mom.y() * mom.y() + mom.z() * mom.z());
178  const double p_z = sign * p / sqrt(1. + th_x * th_x + th_y * th_y);
179 
180  return HepMC::FourVector(p_z * th_x, p_z * th_y, p_z, mom.e());
181 }
double sign(double x)
T sqrt(T t)
Definition: SSEVec.h:19

Member Data Documentation

edm::ESGetToken<CTPPSBeamParameters, CTPPSBeamParametersRcd> BeamDivergenceVtxGenerator::beamParametersToken_
private

Definition at line 40 of file BeamDivergenceVtxGenerator.cc.

Referenced by produce().

std::vector<edm::EDGetTokenT<reco::GenParticleCollection> > BeamDivergenceVtxGenerator::genParticleTokens_
private

Definition at line 41 of file BeamDivergenceVtxGenerator.cc.

Referenced by BeamDivergenceVtxGenerator(), and produce().

bool BeamDivergenceVtxGenerator::simulateBeamDivergence_
private

Definition at line 44 of file BeamDivergenceVtxGenerator.cc.

Referenced by applySmearingHepMC(), and produce().

bool BeamDivergenceVtxGenerator::simulateVertex_
private

Definition at line 43 of file BeamDivergenceVtxGenerator.cc.

Referenced by applySmearingHepMC(), and produce().

edm::EDGetTokenT<edm::HepMCProduct> BeamDivergenceVtxGenerator::sourceToken_
private

Definition at line 39 of file BeamDivergenceVtxGenerator.cc.

Referenced by BeamDivergenceVtxGenerator(), and produce().