CMS 3D CMS Logo

BeamDivergenceVtxGenerator.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
16 
18 
20 
23 
24 #include <CLHEP/Random/RandGauss.h>
25 
26 //----------------------------------------------------------------------------------------------------
27 
29 public:
31  ~BeamDivergenceVtxGenerator() override = default;
32 
34 
35  void produce(edm::Event &, const edm::EventSetup &) override;
36 
37 private:
39  std::vector<edm::EDGetTokenT<reco::GenParticleCollection>> genParticleTokens_;
40 
41  bool simulateVertex_;
43 
44  struct SmearingParameters {
45  double vtx_x, vtx_y, vtx_z; // cm
46  double bd_x_45, bd_y_45, bd_x_56, bd_y_56; // rad
47  };
48 
49  template <typename T>
50  static HepMC::FourVector smearMomentum(const T &mom, const SmearingParameters &sp);
51 
53 
55 };
56 
57 //----------------------------------------------------------------------------------------------------
58 
60  : simulateVertex_(iConfig.getParameter<bool>("simulateVertex")),
61  simulateBeamDivergence_(iConfig.getParameter<bool>("simulateBeamDivergence")) {
62  const edm::InputTag tagSrcHepMC = iConfig.getParameter<edm::InputTag>("src");
63  if (!tagSrcHepMC.label().empty())
64  sourceToken_ = consumes<edm::HepMCProduct>(tagSrcHepMC);
65 
66  for (const auto &it : iConfig.getParameter<std::vector<edm::InputTag>>("srcGenParticle"))
67  genParticleTokens_.push_back(consumes<reco::GenParticleCollection>(it));
68 
70  if (!rng.isAvailable())
71  throw cms::Exception("Configuration")
72  << "The BeamDivergenceVtxGenerator requires the RandomNumberGeneratorService\n"
73  "which is not present in the configuration file. \n"
74  "You must add the service\n"
75  "in the configuration file or remove the modules that require it.";
76 
77  produces<edm::HepMCProduct>();
78 }
79 
80 //----------------------------------------------------------------------------------------------------
81 
84 
85  desc.add<edm::InputTag>("src", edm::InputTag("generator", "unsmeared"))
86  ->setComment("input collection in HepMC format");
87 
88  desc.add<std::vector<edm::InputTag>>("srcGenParticle", std::vector<edm::InputTag>())
89  ->setComment("input collections in GenParticle format");
90 
91  desc.add<bool>("simulateBeamDivergence", true)->setComment("account for the beam angular divergence?");
92  desc.add<bool>("simulateVertex", true)->setComment("account for the vertex transverse smearing?");
93 
94  descriptions.add("beamDivergenceVtxGenerator", desc);
95 }
96 
97 //----------------------------------------------------------------------------------------------------
98 
100  // get random engine
102  CLHEP::HepRandomEngine *rnd = &(rng->getEngine(iEvent.streamID()));
103 
104  // get conditions
105  edm::ESHandle<CTPPSBeamParameters> hBeamParameters;
106  iSetup.get<CTPPSBeamParametersRcd>().get(hBeamParameters);
107 
108  // get HepMC input (if given)
111  genEvt = new HepMC::GenEvent();
112  } else {
113  edm::Handle<edm::HepMCProduct> hepUnsmearedMCEvt;
114  iEvent.getByToken(sourceToken_, hepUnsmearedMCEvt);
115 
116  genEvt = new HepMC::GenEvent(*hepUnsmearedMCEvt->GetEvent());
117  }
118 
119  // prepare output
120  std::unique_ptr<edm::HepMCProduct> output(new edm::HepMCProduct(genEvt));
121 
122  // generate smearing parameters
124 
125  if (simulateVertex_) {
126  // NB: the separtion between effective offsets in LHC sectors 45 and 56 cannot be applied, thus the values for 45 are used
127  sp.vtx_x = hBeamParameters->getVtxOffsetX45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevX();
128  sp.vtx_y = hBeamParameters->getVtxOffsetY45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevY();
129  sp.vtx_z = hBeamParameters->getVtxOffsetZ45() + CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getVtxStddevZ();
130  }
131 
133  sp.bd_x_45 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceX45();
134  sp.bd_x_56 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceX56();
135  sp.bd_y_45 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceY45();
136  sp.bd_y_56 = CLHEP::RandGauss::shoot(rnd) * hBeamParameters->getBeamDivergenceY56();
137  }
138 
139  // apply smearing
141 
142  for (const auto &token : genParticleTokens_) {
144  iEvent.getByToken(token, hGPCollection);
145 
146  for (const auto &gp : *hGPCollection)
148  }
149 
150  // save output
151  iEvent.put(std::move(output));
152 }
153 
154 //----------------------------------------------------------------------------------------------------
155 
156 template <typename T>
157 HepMC::FourVector BeamDivergenceVtxGenerator::smearMomentum(const T &mom, const SmearingParameters &sp) {
158  // TODO: this is an oversimplified implemetation
159  // the TOTEM smearing module should be taken as reference
160 
161  double th_x = mom.x() / mom.z();
162  double th_y = mom.y() / mom.z();
163 
164  if (mom.z() > 0.) {
165  th_x += sp.bd_x_45;
166  th_y += sp.bd_y_45;
167  } else {
168  th_x += sp.bd_x_56;
169  th_y += sp.bd_y_56;
170  }
171 
172  // calculate consistent p_z component
173  const double sign = (mom.z() > 0.) ? 1. : -1.;
174  const double p = sqrt(mom.x() * mom.x() + mom.y() * mom.y() + mom.z() * mom.z());
175  const double p_z = sign * p / sqrt(1. + th_x * th_x + th_y * th_y);
176 
177  return HepMC::FourVector(p_z * th_x, p_z * th_y, p_z, mom.e());
178 }
179 
180 //----------------------------------------------------------------------------------------------------
181 
183  if (simulateVertex_) {
184  for (HepMC::GenEvent::vertex_iterator vit = genEvt->vertices_begin(); vit != genEvt->vertices_end(); ++vit) {
185  const auto &pos = (*vit)->position();
186  (*vit)->set_position(HepMC::FourVector(pos.x() + sp.vtx_x * 1E1, // conversion: cm to mm
187  pos.y() + sp.vtx_y * 1E1,
188  pos.z() + sp.vtx_z * 1E1,
189  pos.t()));
190  }
191  }
192 
194  for (HepMC::GenEvent::particle_iterator part = genEvt->particles_begin(); part != genEvt->particles_end(); ++part)
195  (*part)->set_momentum(smearMomentum((*part)->momentum(), sp));
196  }
197 }
198 
199 //----------------------------------------------------------------------------------------------------
200 
202  const SmearingParameters &sp,
204  // add vertex of the particle
205  HepMC::GenVertex *vtx = new HepMC::GenVertex(HepMC::FourVector((gp.vx() + sp.vtx_x) * 1E1, // conversion: cm to mm
206  (gp.vy() + sp.vtx_y) * 1E1,
207  (gp.vz() + sp.vtx_z) * 1E1,
208  0.));
209  genEvt->add_vertex(vtx);
210 
211  // add the particle itself
212  HepMC::GenParticle *particle = new HepMC::GenParticle(smearMomentum(gp.p4(), sp), gp.pdgId(), gp.status());
213  vtx->add_particle_out(particle);
214 }
215 
216 //----------------------------------------------------------------------------------------------------
217 
ConfigurationDescriptions.h
CTPPSBeamParameters::getVtxStddevY
double getVtxStddevY() const
Definition: CTPPSBeamParameters.cc:78
BeamDivergenceVtxGenerator
Definition: BeamDivergenceVtxGenerator.cc:27
edm::RandomNumberGenerator::getEngine
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
CTPPSBeamParameters.h
electrons_cff.bool
bool
Definition: electrons_cff.py:372
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
CTPPSBeamParametersRcd.h
ESHandle.h
BeamDivergenceVtxGenerator::SmearingParameters::bd_x_45
double bd_x_45
Definition: BeamDivergenceVtxGenerator.cc:47
reco::GenParticle
Definition: GenParticle.h:21
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:32
edm::EDGetTokenT< edm::HepMCProduct >
RandomNumberGenerator.h
BeamDivergenceVtxGenerator::SmearingParameters::bd_y_56
double bd_y_56
Definition: BeamDivergenceVtxGenerator.cc:47
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
pos
Definition: PixelAliasList.h:18
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
EDProducer.h
Validation_hcalonly_cfi.sign
sign
Definition: Validation_hcalonly_cfi.py:32
edm::Handle< edm::HepMCProduct >
edm::Service::isAvailable
bool isAvailable() const
Definition: Service.h:40
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
GenParticle.h
CTPPSBeamParameters::getBeamDivergenceY45
double getBeamDivergenceY45() const
Definition: CTPPSBeamParameters.cc:61
edm::InputTag::label
std::string const & label() const
Definition: InputTag.h:36
MakerMacros.h
CTPPSBeamParameters::getVtxOffsetZ45
double getVtxOffsetZ45() const
Definition: CTPPSBeamParameters.cc:72
part
part
Definition: HCALResponse.h:20
BeamDivergenceVtxGenerator::SmearingParameters
Definition: BeamDivergenceVtxGenerator.cc:45
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CTPPSBeamParameters::getBeamDivergenceY56
double getBeamDivergenceY56() const
Definition: CTPPSBeamParameters.cc:63
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
CTPPSBeamParameters::getBeamDivergenceX45
double getBeamDivergenceX45() const
Definition: CTPPSBeamParameters.cc:60
Service.h
edm::EDGetTokenT::isUninitialized
bool isUninitialized() const
Definition: EDGetToken.h:70
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
BeamDivergenceVtxGenerator::smearMomentum
static HepMC::FourVector smearMomentum(const T &mom, const SmearingParameters &sp)
Definition: BeamDivergenceVtxGenerator.cc:156
edm::ESHandle< CTPPSBeamParameters >
BeamDivergenceVtxGenerator::sourceToken_
edm::EDGetTokenT< edm::HepMCProduct > sourceToken_
Definition: BeamDivergenceVtxGenerator.cc:39
ParameterSetDescription.h
BeamDivergenceVtxGenerator::SmearingParameters::bd_y_45
double bd_y_45
Definition: BeamDivergenceVtxGenerator.cc:47
BeamDivergenceVtxGenerator::SmearingParameters::bd_x_56
double bd_x_56
Definition: BeamDivergenceVtxGenerator.cc:47
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
badGlobalMuonTaggersAOD_cff.vtx
vtx
Definition: badGlobalMuonTaggersAOD_cff.py:5
CTPPSBeamParameters::getBeamDivergenceX56
double getBeamDivergenceX56() const
Definition: CTPPSBeamParameters.cc:62
runTauDisplay.gp
gp
Definition: runTauDisplay.py:431
CTPPSBeamParameters::getVtxOffsetY45
double getVtxOffsetY45() const
Definition: CTPPSBeamParameters.cc:71
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
BeamDivergenceVtxGenerator::addSmearedGenParticle
void addSmearedGenParticle(const reco::GenParticle &gp, const SmearingParameters &sp, HepMC::GenEvent *genEvt)
Definition: BeamDivergenceVtxGenerator.cc:200
BeamDivergenceVtxGenerator::SmearingParameters::vtx_y
double vtx_y
Definition: BeamDivergenceVtxGenerator.cc:46
Event.h
edm::Service< edm::RandomNumberGenerator >
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::EventSetup
Definition: EventSetup.h:57
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
get
#define get
BeamDivergenceVtxGenerator::applySmearingHepMC
void applySmearingHepMC(const SmearingParameters &sp, HepMC::GenEvent *genEvt)
Definition: BeamDivergenceVtxGenerator.cc:181
CTPPSBeamParameters::getVtxOffsetX45
double getVtxOffsetX45() const
Definition: CTPPSBeamParameters.cc:70
CTPPSBeamParametersRcd
Definition: CTPPSBeamParametersRcd.h:14
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
BeamDivergenceVtxGenerator::BeamDivergenceVtxGenerator
BeamDivergenceVtxGenerator(const edm::ParameterSet &)
Definition: BeamDivergenceVtxGenerator.cc:58
BeamDivergenceVtxGenerator::~BeamDivergenceVtxGenerator
~BeamDivergenceVtxGenerator() override=default
GenParticle.GenParticle
GenParticle
Definition: GenParticle.py:18
eostools.move
def move(src, dest)
Definition: eostools.py:511
CTPPSBeamParameters::getVtxStddevZ
double getVtxStddevZ() const
Definition: CTPPSBeamParameters.cc:79
T
long double T
Definition: Basic3DVectorLD.h:48
EventSetup.h
BeamDivergenceVtxGenerator::SmearingParameters::vtx_x
double vtx_x
Definition: BeamDivergenceVtxGenerator.cc:46
BeamDivergenceVtxGenerator::SmearingParameters::vtx_z
double vtx_z
Definition: BeamDivergenceVtxGenerator.cc:46
cms::Exception
Definition: Exception.h:70
BeamDivergenceVtxGenerator::genParticleTokens_
std::vector< edm::EDGetTokenT< reco::GenParticleCollection > > genParticleTokens_
Definition: BeamDivergenceVtxGenerator.cc:40
edm::HepMCProduct
Definition: HepMCProduct.h:18
ParameterSet.h
HepMCProduct.h
CTPPSBeamParameters::getVtxStddevX
double getVtxStddevX() const
Definition: CTPPSBeamParameters.cc:77
edm::Event
Definition: Event.h:73
BeamDivergenceVtxGenerator::simulateVertex_
bool simulateVertex_
Definition: BeamDivergenceVtxGenerator.cc:42
edm::InputTag
Definition: InputTag.h:15
BeamDivergenceVtxGenerator::simulateBeamDivergence_
bool simulateBeamDivergence_
Definition: BeamDivergenceVtxGenerator.cc:43
TtGenEvtProducer_cfi.genEvt
genEvt
Definition: TtGenEvtProducer_cfi.py:7
BeamDivergenceVtxGenerator::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &)
Definition: BeamDivergenceVtxGenerator.cc:81
BeamDivergenceVtxGenerator::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: BeamDivergenceVtxGenerator.cc:98
reco::vertex_iterator
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
Definition: VertexFwd.h:19
unpackBuffers-CaloStage2.token
token
Definition: unpackBuffers-CaloStage2.py:316