CMS 3D CMS Logo

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
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

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

std::vector< edm::EDGetTokenT< reco::GenParticleCollection > > genParticleTokens_
 
bool simulateBeamDivergence_
 
bool simulateVertex_
 
edm::EDGetTokenT< edm::HepMCProductsourceToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Definition at line 28 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_.

61  :
62  simulateVertex_ (iConfig.getParameter<bool>("simulateVertex")),
63  simulateBeamDivergence_(iConfig.getParameter<bool>("simulateBeamDivergence"))
64 {
65  const edm::InputTag tagSrcHepMC = iConfig.getParameter<edm::InputTag>("src");
66  if (tagSrcHepMC.label() != "")
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 }
T getParameter(std::string const &) const
std::vector< edm::EDGetTokenT< reco::GenParticleCollection > > genParticleTokens_
bool isAvailable() const
Definition: Service.h:40
edm::EDGetTokenT< edm::HepMCProduct > sourceToken_
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 218 of file BeamDivergenceVtxGenerator.cc.

References DEFINE_FWK_MODULE, GenParticle::GenParticle, reco::LeafCandidate::p4(), reco::LeafCandidate::pdgId(), smearMomentum(), reco::LeafCandidate::status(), extraflags_cff::vtx, 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().

219 {
220  // add vertex of the particle
221  HepMC::GenVertex *vtx = new HepMC::GenVertex(HepMC::FourVector(
222  (gp.vx() + sp.vtx_x)*1E1, // conversion: cm to mm
223  (gp.vy() + sp.vtx_y)*1E1,
224  (gp.vz() + sp.vtx_z)*1E1,
225  0.
226  ));
227  genEvt->add_vertex(vtx);
228 
229  // add the particle itself
230  HepMC::GenParticle* particle = new HepMC::GenParticle(smearMomentum(gp.p4(), sp), gp.pdgId(), gp.status());
231  vtx->add_particle_out(particle);
232 }
int pdgId() const final
PDG identifier.
double vy() const override
y coordinate of vertex position
static HepMC::FourVector smearMomentum(const T &mom, const SmearingParameters &sp)
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
double vz() const override
z coordinate of vertex position
int status() const final
status word
double vx() const override
x coordinate of vertex position
void BeamDivergenceVtxGenerator::applySmearingHepMC ( const SmearingParameters sp,
HepMC::GenEvent genEvt 
)
private

Definition at line 193 of file BeamDivergenceVtxGenerator.cc.

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

Referenced by produce().

194 {
195  if (simulateVertex_)
196  {
197  for (HepMC::GenEvent::vertex_iterator vit = genEvt->vertices_begin(); vit != genEvt->vertices_end(); ++vit)
198  {
199  const auto &pos = (*vit)->position();
200  (*vit)->set_position(HepMC::FourVector(
201  pos.x() + sp.vtx_x*1E1, // conversion: cm to mm
202  pos.y() + sp.vtx_y*1E1,
203  pos.z() + sp.vtx_z*1E1,
204  pos.t()
205  ));
206  }
207  }
208 
210  {
211  for (HepMC::GenEvent::particle_iterator part = genEvt->particles_begin(); part != genEvt->particles_end(); ++part)
212  (*part)->set_momentum(smearMomentum((*part)->momentum(), sp));
213  }
214 }
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(), and edm::ParameterDescriptionNode::setComment().

86 {
88 
89  desc.add<edm::InputTag>("src", edm::InputTag("generator", "unsmeared"))
90  ->setComment("input collection in HepMC format");
91 
92  desc.add<std::vector<edm::InputTag>>("srcGenParticle", std::vector<edm::InputTag>())
93  ->setComment("input collections in GenParticle format");
94 
95  desc.add<bool>("simulateBeamDivergence", true)->setComment("account for the beam angular divergence?");
96  desc.add<bool>("simulateVertex", true)->setComment("account for the vertex transverse smearing?");
97 
98  descriptions.add("beamDivergenceVtxGenerator", desc);
99 }
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 103 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, TtGenEvtProducer_cfi::genEvt, genParticleTokens_, edm::EventSetup::get(), CTPPSBeamParameters::getBeamDivergenceX45(), CTPPSBeamParameters::getBeamDivergenceX56(), CTPPSBeamParameters::getBeamDivergenceY45(), CTPPSBeamParameters::getBeamDivergenceY56(), edm::Event::getByToken(), edm::RandomNumberGenerator::getEngine(), edm::HepMCProduct::GetEvent(), CTPPSBeamParameters::getVtxOffsetX45(), CTPPSBeamParameters::getVtxOffsetY45(), CTPPSBeamParameters::getVtxOffsetZ45(), CTPPSBeamParameters::getVtxStddevX(), CTPPSBeamParameters::getVtxStddevY(), CTPPSBeamParameters::getVtxStddevZ(), runTauDisplay::gp, edm::EDGetTokenT< T >::isUninitialized(), eostools::move(), convertSQLitetoXML_cfg::output, edm::Event::put(), simulateBeamDivergence_, simulateVertex_, sourceToken_, edm::Event::streamID(), BeamDivergenceVtxGenerator::SmearingParameters::vtx_x, BeamDivergenceVtxGenerator::SmearingParameters::vtx_y, and BeamDivergenceVtxGenerator::SmearingParameters::vtx_z.

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

Definition at line 166 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, Validation_hcalonly_cfi::sign, and mathSSE::sqrt().

Referenced by addSmearedGenParticle(), and applySmearingHepMC().

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

Member Data Documentation

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

Definition at line 40 of file BeamDivergenceVtxGenerator.cc.

Referenced by BeamDivergenceVtxGenerator(), and produce().

bool BeamDivergenceVtxGenerator::simulateBeamDivergence_
private

Definition at line 43 of file BeamDivergenceVtxGenerator.cc.

Referenced by applySmearingHepMC(), and produce().

bool BeamDivergenceVtxGenerator::simulateVertex_
private

Definition at line 42 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().