CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
BeamDivergenceVtxGenerator Class Reference
Inheritance diagram for BeamDivergenceVtxGenerator:
edm::stream::EDProducer<>

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 hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &)
 

Private Attributes

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 26 of file BeamDivergenceVtxGenerator.cc.

Constructor & Destructor Documentation

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

Definition at line 44 of file BeamDivergenceVtxGenerator.cc.

References edm::Service< T >::isAvailable().

45  : sourceToken_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("src"))),
46  simulateVertex_(iConfig.getParameter<bool>("simulateVertex")),
47  simulateBeamDivergence_(iConfig.getParameter<bool>("simulateBeamDivergence")) {
49  if (!rng.isAvailable())
50  throw cms::Exception("Configuration")
51  << "The BeamDivergenceVtxGenerator requires the RandomNumberGeneratorService\n"
52  "which is not present in the configuration file. \n"
53  "You must add the service\n"
54  "in the configuration file or remove the modules that require it.";
55 
56  produces<edm::HepMCProduct>();
57 }
T getParameter(std::string const &) const
bool isAvailable() const
Definition: Service.h:40
edm::EDGetTokenT< edm::HepMCProduct > sourceToken_
BeamDivergenceVtxGenerator::~BeamDivergenceVtxGenerator ( )
overridedefault

Member Function Documentation

void BeamDivergenceVtxGenerator::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 132 of file BeamDivergenceVtxGenerator.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), DEFINE_FWK_MODULE, and HLT_2018_cff::InputTag.

132  {
134  desc.add<edm::InputTag>("src", edm::InputTag("generator", "unsmeared"))
135  ->setComment("input collection where to retrieve outgoing particles kinematics to be smeared");
136  desc.add<bool>("simulateBeamDivergence", true)->setComment("account for the beam angular divergence?");
137  desc.add<bool>("simulateVertex", true)->setComment("account for the vertex transverse smearing?");
138 
139  descriptions.add("beamDivergenceVtxGenerator", desc);
140 }
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 61 of file BeamDivergenceVtxGenerator.cc.

References 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(), eostools::move(), edm::Event::put(), edm::shift, Validation_hcalonly_cfi::sign, simulateBeamDivergence_, simulateVertex_, sourceToken_, mathSSE::sqrt(), and edm::Event::streamID().

61  {
62  // get random engine
64  CLHEP::HepRandomEngine* rnd = &(rng->getEngine(iEvent.streamID()));
65 
66  // get conditions
67  edm::ESHandle<CTPPSBeamParameters> hBeamParameters;
68  iSetup.get<CTPPSBeamParametersRcd>().get(hBeamParameters);
69 
70  // get input
71  edm::Handle<edm::HepMCProduct> hepUnsmearedMCEvt;
72  iEvent.getByToken(sourceToken_, hepUnsmearedMCEvt);
73 
74  // prepare output
75  HepMC::GenEvent* genevt = new HepMC::GenEvent(*hepUnsmearedMCEvt->GetEvent());
76  std::unique_ptr<edm::HepMCProduct> pEvent(new edm::HepMCProduct(genevt));
77 
78  // apply vertex smearing
79  if (simulateVertex_) {
80  // NB: the separtion between effective offsets in LHC sectors 45 and 56 cannot be applied, thus the values for 45 are used
81  const double vtx_x =
82  hBeamParameters->getVtxOffsetX45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevX();
83  const double vtx_y =
84  hBeamParameters->getVtxOffsetY45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevY();
85  const double vtx_z =
86  hBeamParameters->getVtxOffsetZ45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevZ();
87 
88  HepMC::FourVector shift(vtx_x * 1E1, vtx_y * 1E1, vtx_z * 1E1, 0.); // conversions: cm to mm
89  pEvent->applyVtxGen(&shift);
90  }
91 
92  // apply beam divergence
94  const double bd_x_45 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceX45();
95  const double bd_x_56 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceX56();
96 
97  const double bd_y_45 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceY45();
98  const double bd_y_56 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceY56();
99 
100  for (HepMC::GenEvent::particle_iterator part = genevt->particles_begin(); part != genevt->particles_end(); ++part) {
101  const HepMC::FourVector mom = (*part)->momentum();
102 
103  // TODO: this is an oversimplified implemetation
104  // the TOTEM smearing module should be taken as reference
105 
106  double th_x = mom.x() / mom.z();
107  double th_y = mom.y() / mom.z();
108 
109  if (mom.z() > 0.) {
110  th_x += bd_x_45;
111  th_y += bd_y_45;
112  } else {
113  th_x += bd_x_56;
114  th_y += bd_y_56;
115  }
116 
117  // calculate consistent p_z component
118  const double sign = (mom.z() > 0.) ? 1. : -1.;
119  const double p_z = sign * mom.rho() / sqrt(1. + th_x * th_x + th_y * th_y);
120 
121  // set smeared momentum
122  (*part)->set_momentum(HepMC::FourVector(p_z * th_x, p_z * th_y, p_z, mom.e()));
123  }
124  }
125 
126  // save output
127  iEvent.put(std::move(pEvent));
128 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
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
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< edm::HepMCProduct > sourceToken_
double getVtxStddevX() const
double getBeamDivergenceY45() const
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
part
Definition: HCALResponse.h:20
double getVtxOffsetZ45() const
double getVtxOffsetY45() const
T get() const
Definition: EventSetup.h:73
StreamID streamID() const
Definition: Event.h:96
static unsigned int const shift
double getVtxStddevY() const
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

bool BeamDivergenceVtxGenerator::simulateBeamDivergence_
private

Definition at line 39 of file BeamDivergenceVtxGenerator.cc.

Referenced by produce().

bool BeamDivergenceVtxGenerator::simulateVertex_
private

Definition at line 38 of file BeamDivergenceVtxGenerator.cc.

Referenced by produce().

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

Definition at line 36 of file BeamDivergenceVtxGenerator.cc.

Referenced by produce().