LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
evgen::BaseRadioGen Class Referenceabstract

#include "BaseRadioGen.h"

Inheritance diagram for evgen::BaseRadioGen:
art::EDProducer art::detail::Producer art::detail::LegacyModule art::Modifier art::ModuleBase art::ProductRegistryHelper evgen::Decay0Gen evgen::SpectrumVolumeGen

Public Types

using ModuleType = EDProducer
 
template<typename UserConfig , typename KeysToIgnore = void>
using Table = Modifier::Table< UserConfig, KeysToIgnore >
 

Public Member Functions

 BaseRadioGen (fhicl::ParameterSet const &pset)
 
void produce (art::Event &evt)
 
void beginRun (art::Run &run)
 
void beginJob ()
 
void endJob ()
 
virtual void produce_radio (art::Event &evt)=0
 
virtual void beginRun_radio (art::Run &)
 
virtual void beginJob_radio ()
 
virtual void endJob_radio ()
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
void fillProductDescriptions ()
 
void registerProducts (ProductDescriptions &productsToRegister)
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
std::unique_ptr< Worker > makeWorker (WorkerParams const &wp)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Protected Member Functions

int GetNDecays ()
 
bool GetGoodPositionTime (TLorentzVector &position)
 
TLorentzVector dirCalc (double p, double m)
 
void FillHistos (simb::MCParticle &part)
 
CLHEP::HepRandomEngine & GetRandomEngine ()
 
int GetNEvents ()
 
void GetTs (double &T0, double &T1)
 
void GetXs (double &X0, double &X1)
 
void GetYs (double &Y0, double &Y1)
 
void GetZs (double &Z0, double &Z1)
 
ConsumesCollector & consumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Protected Attributes

std::vector< std::string > m_isotope
 isotope to simulate. Example: "Ar39" More...
 

Private Member Functions

void CalculateActiveVolumeFromAllNodes ()
 
void CalculateActiveVolumeFromXYZ ()
 
void SimplePDG (int pdg, int &simple, std::string &name)
 
void DeclareOutputHistos ()
 
void FillAllNodes (const TGeoNode *curnode)
 
void GetNodeXYZMinMax (const TGeoNode *from, const TGeoNode *to, double &x_min, double &x_max, double &y_min, double &y_max, double &z_min, double &z_max)
 
bool IsDaughterNode (const TGeoNode *daughter_node, const TGeoNode *mother_node)
 
bool findNode (const TGeoNode *curnode, std::string &tgtnname, const TGeoNode *&targetnode)
 Shamelessly stolen from here: https://cdcvs.fnal.gov/redmine/attachments/6719/calc_bbox.C. More...
 
bool findMotherNode (const TGeoNode *cur_node, std::string &daughter_name, const TGeoNode *&mother_node)
 Adapted from above. More...
 
bool NodeSupported (TGeoNode const *node) const
 

Private Attributes

std::set< const TGeoNode * > m_all_nodes
 
std::string m_material
 regex of materials in which to generate the decays. Example: "LAr" More...
 
std::string m_volume_rand
 The volume in which to generate the decays. More...
 
std::string m_volume_gen
 The volume in which to generate the decays. More...
 
std::regex m_regex_material
 
std::regex m_regex_volume
 
double m_Bq {signaling_NaN}
 Radioactivity in Becquerels (decay per sec) per cubic cm. More...
 
double m_rate
 Radioactivity in Becquerels (decay per sec) use either of this of Bq. More...
 
double m_T0 {signaling_NaN}
 Beginning of time window to simulate in ns. More...
 
double m_T1 {signaling_NaN}
 End of time window to simulate in ns. More...
 
double m_X0 {signaling_NaN}
 Bottom corner x position (cm) in world coordinates. More...
 
double m_Y0 {signaling_NaN}
 Bottom corner y position (cm) in world coordinates. More...
 
double m_Z0 {signaling_NaN}
 Bottom corner z position (cm) in world coordinates. More...
 
double m_X1 {signaling_NaN}
 Top corner x position (cm) in world coordinates. More...
 
double m_Y1 {signaling_NaN}
 Top corner y position (cm) in world coordinates. More...
 
double m_Z1 {signaling_NaN}
 Top corner z position (cm) in world coordinates. More...
 
double m_volume_cc {signaling_NaN}
 
art::ServiceHandle< geo::Geometry const > m_geo_service
 
TGeoManager * m_geo_manager
 
CLHEP::HepRandomEngine & m_engine
 
CLHEP::RandFlat m_random_flat
 
CLHEP::RandPoisson m_random_poisson
 
bool m_geo_volume_mode {false}
 
bool m_rate_mode
 
bool m_volume_rand_present
 
size_t m_max_tries_event
 
size_t m_max_tries_rate_calculation
 
size_t m_target_n_point_rate_calculation
 
std::map< const TGeoNode *, double > m_good_nodes {}
 
std::vector< const TGeoVolume * > m_good_volumes {}
 
std::vector< const TGeoMaterial * > m_good_materials {}
 
TF1 * m_distrib_xpos {nullptr}
 
TF1 * m_distrib_ypos {nullptr}
 
TF1 * m_distrib_zpos {nullptr}
 
int m_nevent {0}
 
std::map< int, TH2D * > m_pos_xy_TH2D {}
 
std::map< int, TH2D * > m_pos_xz_TH2D {}
 
std::map< int, TH1D * > m_dir_x_TH1D {}
 
std::map< int, TH1D * > m_dir_y_TH1D {}
 
std::map< int, TH1D * > m_dir_z_TH1D {}
 
std::map< int, TH1D * > m_mom_TH1D {}
 
std::map< int, TH1D * > m_ke_TH1D {}
 
std::map< int, TH1D * > m_time_TH1D {}
 
TH1D * m_pdg_TH1D {nullptr}
 

Static Private Attributes

static constexpr double signaling_NaN = std::numeric_limits<double>::signaling_NaN()
 

Detailed Description

Definition at line 71 of file BaseRadioGen.h.

Member Typedef Documentation

Definition at line 17 of file EDProducer.h.

template<typename UserConfig , typename KeysToIgnore = void>
using art::detail::Producer::Table = Modifier::Table<UserConfig, KeysToIgnore>
inherited

Definition at line 26 of file Producer.h.

Constructor & Destructor Documentation

evgen::BaseRadioGen::BaseRadioGen ( fhicl::ParameterSet const &  pset)
inlineexplicit

If we specified X, Y, Z we will throw with what the user provided This will throw an error if somebody specifies only one X0.

If we specified volume_rand we will throw with what the coordinates of the volume

Else we throw in the whole geometry but only in the nodes who have the specified volume/material

Definition at line 269 of file BaseRadioGen.h.

References art::detail::EngineCreator::createEngine(), m_engine, m_geo_manager, m_geo_service, m_material, m_max_tries_event, m_max_tries_rate_calculation, m_random_flat, m_random_poisson, m_rate_mode, m_regex_material, m_regex_volume, m_target_n_point_rate_calculation, m_volume_gen, m_volume_rand, m_volume_rand_present, and geo::GeometryCore::ROOTGeoManager().

270  : EDProducer{pset}
271  , m_material{pset.get<std::string>("material", ".*")}
272  , m_volume_rand{""}
273  , m_volume_gen{pset.get<std::string>("volume_gen", ".*")}
274  , m_regex_material{(std::regex)m_material}
275  , m_regex_volume{(std::regex)m_volume_gen}
277  , m_engine(art::ServiceHandle<rndm::NuRandomService>()->registerAndSeedEngine(
278  createEngine(0, "HepJamesRandom", "BaseRadioGen"),
279  "HepJamesRandom",
280  "BaseRadioGen",
281  pset,
282  "SeedBaseRadioGen"))
285  , m_rate_mode{pset.has_key("rate")}
286  , m_volume_rand_present{pset.has_key("volume_rand")}
287  , m_max_tries_event{pset.get<size_t>("max_tries_event", 1'000'000)}
288  , m_max_tries_rate_calculation{pset.get<size_t>("max_tries_rate_calculation", 40'000'000)}
290  pset.get<size_t>("target_n_point_rate_calculation", 100'000)}
291  {
292  produces<std::vector<simb::MCTruth>>();
293  produces<sumdata::RunData, art::InRun>();
294 
295  m_rate_mode = pset.get_if_present<double>("rate", m_rate);
296  if (not m_rate_mode) m_Bq = pset.get<double>("BqPercc");
297 
298  bool timed_mode = pset.get_if_present<double>("T0", m_T0);
299  if (timed_mode) { m_T1 = pset.get<double>("T1"); }
300  else {
301  auto const clockData =
302  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
303  auto const detProp =
304  art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataForJob(clockData);
305  int nsample = detProp.NumberTimeSamples();
306  m_T0 = -nsample * sampling_rate(clockData);
307  m_T1 = -m_T0;
308  }
309 
310  if (m_material != ".*" || m_volume_gen != ".*") {
311  MF_LOG_INFO("BaseRadioGen") << "Calculating the volume of " << m_material
312  << " and the volume " << m_volume_gen << " in the geometry.\n";
313 
314  TObjArray* volumes = m_geo_manager->GetListOfVolumes();
315  TList* materials = m_geo_manager->GetListOfMaterials();
316 
317  for (int i = 0; i < volumes->GetEntries(); ++i) {
318  TGeoVolume* volume = (TGeoVolume*)volumes->At(i);
319  std::string volume_name = volume->GetName();
320  bool good = std::regex_match(volume_name, m_regex_volume);
321  if (good) {
322 
323  if (m_volume_gen != ".*") {
324  MF_LOG_INFO("BaseRadioGen")
325  << " Volume " << volume_name << " is an accepted volume of "
326  << volume->GetShape()->Capacity() << "cm^3.\n";
327  if (volume->GetShape()->TestShapeBits(TGeoShape::kGeoBox)) {
328  TGeoBBox* box = dynamic_cast<TGeoBBox*>(volume->GetShape());
329  if (box)
330  MF_LOG_INFO("BaseRadioGen")
331  << " It is a box of size " << 2. * box->GetDX() << " x " << 2. * box->GetDY()
332  << " x " << 2. * box->GetDZ() << " cm^3\n";
333  }
334  }
335 
336  m_good_volumes.push_back(volume);
337  }
338  }
339 
340  for (int i = 0; i < materials->GetEntries(); ++i) {
341  TGeoMaterial* material = (TGeoMaterial*)materials->At(i);
342  std::string material_name = material->GetName();
343  bool good = std::regex_match(material_name, m_regex_material);
344  if (good) m_good_materials.push_back(material);
345  }
346  }
347  MF_LOG_INFO("BaseRadioGen") << m_good_volumes.size() << " volumes correspond to the regex \""
348  << m_volume_gen << "\".\n";
349  MF_LOG_INFO("BaseRadioGen") << m_good_materials.size()
350  << " materials correspond to the regex \"" << m_material << "\".\n";
351 
352  if (pset.has_key("X0") or pset.has_key("X1") or pset.has_key("Y0") or pset.has_key("Y1") or
353  pset.has_key("Z0") or pset.has_key("Z1")) {
356 
357  m_geo_volume_mode = false;
358  m_volume_rand_present = false;
359 
360  m_X0 = pset.get<double>("X0");
361  m_Y0 = pset.get<double>("Y0");
362  m_Z0 = pset.get<double>("Z0");
363  m_X1 = pset.get<double>("X1");
364  m_Y1 = pset.get<double>("Y1");
365  m_Z1 = pset.get<double>("Z1");
366 
367  CalculateActiveVolumeFromXYZ();
368  }
369  else {
370  const TGeoNode* world = gGeoManager->GetNode(0);
371 
372  if (pset.get_if_present<std::string>("volume_rand", m_volume_rand)) {
373  m_geo_volume_mode = false;
374  m_volume_rand_present = true;
375 
377  const TGeoNode* node = nullptr;
378  findNode(world, m_volume_rand, node);
379  GetNodeXYZMinMax(node, world, m_X0, m_X1, m_Y0, m_Y1, m_Z0, m_Z1);
380 
381  CalculateActiveVolumeFromXYZ();
382  }
383  else {
385  m_geo_volume_mode = true;
386  m_volume_rand_present = false;
387  CalculateActiveVolumeFromAllNodes();
388  }
389  }
390 
391  auto rand = CLHEP::RandFlat(m_engine);
392  gRandom->SetSeed(rand.fireInt(std::numeric_limits<long>::max()));
393 
394  if (pset.has_key("distrib_x"))
395  m_distrib_xpos = new TF1("distrib_x", pset.get<std::string>("distrib_x").c_str(), m_X0, m_X1);
396  if (pset.has_key("distrib_y"))
397  m_distrib_ypos = new TF1("distrib_y", pset.get<std::string>("distrib_y").c_str(), m_Y0, m_Y1);
398  if (pset.has_key("distrib_z"))
399  m_distrib_zpos = new TF1("distrib_z", pset.get<std::string>("distrib_z").c_str(), m_Z0, m_Z1);
400  }
CLHEP::RandPoisson m_random_poisson
Definition: BaseRadioGen.h:162
base_engine_t & createEngine(seed_t seed)
size_t m_target_n_point_rate_calculation
Definition: BaseRadioGen.h:170
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.cc:6
std::string m_volume_gen
The volume in which to generate the decays.
Definition: BaseRadioGen.h:139
TGeoManager * ROOTGeoManager() const
Access to the ROOT geometry description manager.
std::string m_material
regex of materials in which to generate the decays. Example: "LAr"
Definition: BaseRadioGen.h:137
TGeoManager * m_geo_manager
Definition: BaseRadioGen.h:158
std::string m_volume_rand
The volume in which to generate the decays.
Definition: BaseRadioGen.h:138
std::regex m_regex_volume
Definition: BaseRadioGen.h:141
size_t m_max_tries_rate_calculation
Definition: BaseRadioGen.h:169
art::ServiceHandle< geo::Geometry const > m_geo_service
Definition: BaseRadioGen.h:157
CLHEP::RandFlat m_random_flat
Definition: BaseRadioGen.h:161
std::regex m_regex_material
Definition: BaseRadioGen.h:140
CLHEP::HepRandomEngine & m_engine
Definition: BaseRadioGen.h:160

Member Function Documentation

void evgen::BaseRadioGen::beginJob ( )
inlinevirtual

Reimplemented from art::EDProducer.

Definition at line 490 of file BaseRadioGen.h.

491  {
492  if (m_isotope.empty()) {
493  throw cet::exception("BaseRadioGen")
494  << "m_isotope is empty, you need to fill it yourself in the constructor of your module\n";
495  }
496 
498  beginJob_radio();
499  m_nevent = 0;
500  }
virtual void beginJob_radio()
Definition: BaseRadioGen.h:81
std::vector< std::string > m_isotope
isotope to simulate. Example: "Ar39"
Definition: BaseRadioGen.h:90
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
virtual void evgen::BaseRadioGen::beginJob_radio ( )
inlinevirtual

Reimplemented in evgen::Decay0Gen.

Definition at line 81 of file BaseRadioGen.h.

81 {}
void evgen::BaseRadioGen::beginRun ( art::Run run)
inlinevirtual

Reimplemented from art::EDProducer.

Definition at line 483 of file BaseRadioGen.h.

484  {
486  run.put(std::make_unique<sumdata::RunData>(m_geo_service->DetectorName()), art::fullRun());
487  beginRun_radio(run);
488  }
constexpr auto fullRun()
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Run.h:121
virtual void beginRun_radio(art::Run &)
Definition: BaseRadioGen.h:80
art::ServiceHandle< geo::Geometry const > m_geo_service
Definition: BaseRadioGen.h:157
std::string const & DetectorName() const
Returns a string with the name of the detector, as configured.
Definition: GeometryCore.h:203
Namespace collecting geometry-related classes utilities.
virtual void evgen::BaseRadioGen::beginRun_radio ( art::Run )
inlinevirtual

Definition at line 80 of file BaseRadioGen.h.

80 {}
void evgen::BaseRadioGen::CalculateActiveVolumeFromAllNodes ( )
inlineprivate

Definition at line 402 of file BaseRadioGen.h.

Referenced by GetZs().

403  {
404  FillAllNodes(gGeoManager->GetTopNode());
405  m_volume_cc = 0;
406 
407  for (auto const& node : m_all_nodes) {
408  if (not NodeSupported(node)) continue;
409 
410  double capa = node->GetVolume()->GetShape()->Capacity();
411  m_volume_cc += capa;
412  m_good_nodes[node] = m_volume_cc;
413  }
414 
415  MF_LOG_INFO("BaseRadioGen")
416  << m_good_nodes.size() << " nodes (i.e. instance of the volumes) satisfy both the regexes.\n";
417 
418  if (empty(m_good_nodes))
419  throw cet::exception("BaseRadioGen")
420  << "Didn't find an instance of material " << m_material << " and the volume "
421  << m_volume_gen << " in the geometry.\n";
422  }
std::set< const TGeoNode * > m_all_nodes
Definition: BaseRadioGen.h:120
std::map< const TGeoNode *, double > m_good_nodes
Definition: BaseRadioGen.h:172
bool NodeSupported(TGeoNode const *node) const
Definition: BaseRadioGen.h:192
std::string m_volume_gen
The volume in which to generate the decays.
Definition: BaseRadioGen.h:139
std::string m_material
regex of materials in which to generate the decays. Example: "LAr"
Definition: BaseRadioGen.h:137
#define MF_LOG_INFO(category)
void FillAllNodes(const TGeoNode *curnode)
Definition: BaseRadioGen.h:658
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
decltype(auto) constexpr empty(T &&obj)
ADL-aware version of std::empty.
Definition: StdUtils.h:109
void evgen::BaseRadioGen::CalculateActiveVolumeFromXYZ ( )
inlineprivate

Definition at line 424 of file BaseRadioGen.h.

Referenced by GetZs().

425  {
426  m_volume_cc = (m_X1 - m_X0) * (m_Y1 - m_Y0) * (m_Z1 - m_Z0);
427 
428  if (m_material != ".*" || m_volume_gen != ".*") {
429  MF_LOG_INFO("BaseRadioGen") << "Calculating the proportion of " << m_material
430  << " and the volume " << m_volume_gen
431  << " in the specified volume " << m_volume_rand << ".\n";
432  size_t nfound = 0;
433  size_t ntries = 0;
434 
435  double xyz[3];
436  TGeoNode* node = nullptr;
437 
439  ntries++;
440  node = nullptr;
441  xyz[0] = m_X0 + (m_X1 - m_X0) * m_random_flat.fire(0, 1.);
442  xyz[1] = m_Y0 + (m_Y1 - m_Y0) * m_random_flat.fire(0, 1.);
443  xyz[2] = m_Z0 + (m_Z1 - m_Z0) * m_random_flat.fire(0, 1.);
444  m_geo_manager->SetCurrentPoint(xyz);
445  node = m_geo_manager->FindNode();
446  if (!node) continue;
447  if (node->IsOverlapping()) continue;
448 
449  if (!NodeSupported(node)) continue;
450  nfound++;
451  }
452 
453  if (nfound == 0) {
454  throw cet::exception("BaseRadioGen")
455  << "Didn't find the material " << m_material << " and the volume " << m_volume_gen
456  << " in the specified volume " << m_volume_rand << ".\n"
457  << "Position of the box:\n"
458  << m_X0 << " " << m_X1 << "\n"
459  << m_Y0 << " " << m_Y1 << "\n"
460  << m_Z0 << " " << m_Z1 << "\n";
461  }
462 
463  double proportion = (double)nfound / ntries;
464  double proportion_error = proportion * sqrt(1. / nfound + 1. / ntries);
465  m_volume_cc *= proportion;
466 
467  MF_LOG_INFO("BaseRadioGen")
468  << "There is " << proportion * 100. << "% (+/- " << proportion_error * 100. << "%) of "
469  << m_material << " and " << m_volume_gen << " in the specified volume ("
470  << 100. * proportion_error / proportion << "% relat. uncert.).\n"
471  << "If the uncertainty is too big, crank up the parameters \"max_tries_rate_calculation\" "
472  "(default=40,000,000) and/or \"target_n_point_rate_calculation\" (default=100,000) in "
473  "your fhicl.\n\n\n";
474  }
475  }
double m_X1
Top corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:152
bool NodeSupported(TGeoNode const *node) const
Definition: BaseRadioGen.h:192
size_t m_target_n_point_rate_calculation
Definition: BaseRadioGen.h:170
std::string m_volume_gen
The volume in which to generate the decays.
Definition: BaseRadioGen.h:139
std::string m_material
regex of materials in which to generate the decays. Example: "LAr"
Definition: BaseRadioGen.h:137
TGeoManager * m_geo_manager
Definition: BaseRadioGen.h:158
double m_Y1
Top corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:153
std::string m_volume_rand
The volume in which to generate the decays.
Definition: BaseRadioGen.h:138
double m_Z1
Top corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:154
#define MF_LOG_INFO(category)
size_t m_max_tries_rate_calculation
Definition: BaseRadioGen.h:169
double m_Z0
Bottom corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:151
CLHEP::RandFlat m_random_flat
Definition: BaseRadioGen.h:161
double m_Y0
Bottom corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:150
double m_X0
Bottom corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:149
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
template<typename T , BranchType BT>
ProductToken< T > art::ModuleBase::consumes ( InputTag const &  tag)
protectedinherited

Definition at line 61 of file ModuleBase.h.

References art::ModuleBase::collector_, and art::ConsumesCollector::consumes().

62  {
63  return collector_.consumes<T, BT>(tag);
64  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
ProductToken< T > consumes(InputTag const &)
ConsumesCollector & art::ModuleBase::consumesCollector ( )
protectedinherited

Definition at line 57 of file ModuleBase.cc.

References art::ModuleBase::collector_.

58  {
59  return collector_;
60  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
template<typename T , BranchType BT>
void art::ModuleBase::consumesMany ( )
protectedinherited

Definition at line 75 of file ModuleBase.h.

References art::ModuleBase::collector_, and art::ConsumesCollector::consumesMany().

76  {
77  collector_.consumesMany<T, BT>();
78  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
template<typename Element , BranchType = InEvent>
ViewToken<Element> art::ModuleBase::consumesView ( InputTag const &  )
protectedinherited
template<typename T , BranchType BT>
ViewToken<T> art::ModuleBase::consumesView ( InputTag const &  tag)
inherited

Definition at line 68 of file ModuleBase.h.

References art::ModuleBase::collector_, and art::ConsumesCollector::consumesView().

69  {
70  return collector_.consumesView<T, BT>(tag);
71  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
ViewToken< Element > consumesView(InputTag const &)
void evgen::BaseRadioGen::DeclareOutputHistos ( )
inlineprivate

Definition at line 724 of file BaseRadioGen.h.

Referenced by GetZs().

725  {
727  auto pdgs = {1000020040, 11, -11, 22, 2112, 2212, 9999};
728  m_pdg_TH1D = tfs->make<TH1D>("PDG", ";PDG;n particles/event", pdgs.size(), 0, pdgs.size());
729 
730  for (auto pdg : pdgs) {
731 
732  std::string part = "";
733  int simple_pdg;
734  SimplePDG(pdg, simple_pdg, part);
735  art::TFileDirectory dir = tfs->mkdir(part.c_str());
736 
737  m_pos_xy_TH2D[simple_pdg] =
738  dir.make<TH2D>("posXY", ";X [cm];Y [cm]", 100, m_X0, m_X1, 100, m_Y0, m_Y1);
739  m_pos_xz_TH2D[simple_pdg] =
740  dir.make<TH2D>("posXZ", ";X [cm];Z [cm]", 100, m_X0, m_X1, 100, m_Z0, m_Z1);
741  m_dir_x_TH1D[simple_pdg] = dir.make<TH1D>("dirX", ";X momentum projection", 100, -1, 1);
742  m_dir_y_TH1D[simple_pdg] = dir.make<TH1D>("dirY", ";Y momentum projection", 100, -1, 1);
743  m_dir_z_TH1D[simple_pdg] = dir.make<TH1D>("dirZ", ";Z momentum projection", 100, -1, 1);
744  m_mom_TH1D[simple_pdg] =
745  dir.make<TH1D>("Momentum", ";Momentum [MeV];n particles/event", 5000, 0, 500);
746  m_ke_TH1D[simple_pdg] = dir.make<TH1D>("KE", ";KE [MeV];n particles/event", 5000, 0, 500);
747  m_time_TH1D[simple_pdg] =
748  dir.make<TH1D>("Time", ";Time[ns];n particles/event", 100, m_T0, m_T1);
749 
750  m_pdg_TH1D->GetXaxis()->SetBinLabel(simple_pdg + 1, part.c_str());
751  }
752  }
void SimplePDG(int pdg, int &simple, std::string &name)
Definition: BaseRadioGen.h:754
double m_X1
Top corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:152
std::map< int, TH1D * > m_dir_y_TH1D
Definition: BaseRadioGen.h:184
std::map< int, TH1D * > m_dir_x_TH1D
Definition: BaseRadioGen.h:183
std::map< int, TH1D * > m_mom_TH1D
Definition: BaseRadioGen.h:186
TString part[npart]
Definition: Style.C:32
double m_Y1
Top corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:153
double m_Z1
Top corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:154
std::map< int, TH1D * > m_ke_TH1D
Definition: BaseRadioGen.h:187
double m_T0
Beginning of time window to simulate in ns.
Definition: BaseRadioGen.h:147
std::map< int, TH1D * > m_dir_z_TH1D
Definition: BaseRadioGen.h:185
double m_Z0
Bottom corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:151
TDirectory * dir
Definition: macro.C:5
std::map< int, TH2D * > m_pos_xz_TH2D
Definition: BaseRadioGen.h:182
std::map< int, TH1D * > m_time_TH1D
Definition: BaseRadioGen.h:188
double m_Y0
Bottom corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:150
double m_X0
Bottom corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:149
std::map< int, TH2D * > m_pos_xy_TH2D
Definition: BaseRadioGen.h:181
double m_T1
End of time window to simulate in ns.
Definition: BaseRadioGen.h:148
TLorentzVector evgen::BaseRadioGen::dirCalc ( double  p,
double  m 
)
inlineprotected

Definition at line 810 of file BaseRadioGen.h.

Referenced by endJob_radio(), and evgen::SpectrumVolumeGen::produce_radio().

811  {
812  // isotropic production angle for the decay product
813  double costheta = (2.0 * m_random_flat.fire() - 1.0);
814 
815  if (costheta < -1.0) costheta = -1.0;
816  if (costheta > 1.0) costheta = 1.0;
817 
818  double const sintheta = sqrt(1.0 - costheta * costheta);
819  double const phi = 2.0 * M_PI * m_random_flat.fire();
820 
821  return TLorentzVector{p * sintheta * std::cos(phi),
822  p * sintheta * std::sin(phi),
823  p * costheta,
824  std::sqrt(p * p + m * m)};
825  }
CLHEP::RandFlat m_random_flat
Definition: BaseRadioGen.h:161
void art::detail::Producer::doBeginJob ( SharedResources const &  resources)
inherited

Definition at line 22 of file Producer.cc.

References art::detail::Producer::beginJobWithFrame(), and art::detail::Producer::setupQueues().

23  {
24  setupQueues(resources);
25  ProcessingFrame const frame{ScheduleID{}};
26  beginJobWithFrame(frame);
27  }
virtual void setupQueues(SharedResources const &)=0
virtual void beginJobWithFrame(ProcessingFrame const &)=0
bool art::detail::Producer::doBeginRun ( RunPrincipal rp,
ModuleContext const &  mc 
)
inherited

Definition at line 65 of file Producer.cc.

References art::detail::Producer::beginRunWithFrame(), art::RangeSet::forRun(), art::RunPrincipal::makeRun(), r, art::RunPrincipal::runID(), and art::ModuleContext::scheduleID().

66  {
67  auto r = rp.makeRun(mc, RangeSet::forRun(rp.runID()));
68  ProcessingFrame const frame{mc.scheduleID()};
69  beginRunWithFrame(r, frame);
70  r.commitProducts();
71  return true;
72  }
TRandom r
Definition: spectrum.C:23
virtual void beginRunWithFrame(Run &, ProcessingFrame const &)=0
static RangeSet forRun(RunID)
Definition: RangeSet.cc:51
bool art::detail::Producer::doBeginSubRun ( SubRunPrincipal srp,
ModuleContext const &  mc 
)
inherited

Definition at line 85 of file Producer.cc.

References art::detail::Producer::beginSubRunWithFrame(), art::RangeSet::forSubRun(), art::SubRunPrincipal::makeSubRun(), art::ModuleContext::scheduleID(), and art::SubRunPrincipal::subRunID().

86  {
87  auto sr = srp.makeSubRun(mc, RangeSet::forSubRun(srp.subRunID()));
88  ProcessingFrame const frame{mc.scheduleID()};
89  beginSubRunWithFrame(sr, frame);
90  sr.commitProducts();
91  return true;
92  }
virtual void beginSubRunWithFrame(SubRun &, ProcessingFrame const &)=0
static RangeSet forSubRun(SubRunID)
Definition: RangeSet.cc:57
void art::detail::Producer::doEndJob ( )
inherited

Definition at line 30 of file Producer.cc.

References art::detail::Producer::endJobWithFrame().

31  {
32  ProcessingFrame const frame{ScheduleID{}};
33  endJobWithFrame(frame);
34  }
virtual void endJobWithFrame(ProcessingFrame const &)=0
bool art::detail::Producer::doEndRun ( RunPrincipal rp,
ModuleContext const &  mc 
)
inherited

Definition at line 75 of file Producer.cc.

References art::detail::Producer::endRunWithFrame(), art::RunPrincipal::makeRun(), r, art::ModuleContext::scheduleID(), and art::Principal::seenRanges().

76  {
77  auto r = rp.makeRun(mc, rp.seenRanges());
78  ProcessingFrame const frame{mc.scheduleID()};
79  endRunWithFrame(r, frame);
80  r.commitProducts();
81  return true;
82  }
TRandom r
Definition: spectrum.C:23
virtual void endRunWithFrame(Run &, ProcessingFrame const &)=0
bool art::detail::Producer::doEndSubRun ( SubRunPrincipal srp,
ModuleContext const &  mc 
)
inherited

Definition at line 95 of file Producer.cc.

References art::detail::Producer::endSubRunWithFrame(), art::SubRunPrincipal::makeSubRun(), art::ModuleContext::scheduleID(), and art::Principal::seenRanges().

96  {
97  auto sr = srp.makeSubRun(mc, srp.seenRanges());
98  ProcessingFrame const frame{mc.scheduleID()};
99  endSubRunWithFrame(sr, frame);
100  sr.commitProducts();
101  return true;
102  }
virtual void endSubRunWithFrame(SubRun &, ProcessingFrame const &)=0
bool art::detail::Producer::doEvent ( EventPrincipal ep,
ModuleContext const &  mc,
std::atomic< std::size_t > &  counts_run,
std::atomic< std::size_t > &  counts_passed,
std::atomic< std::size_t > &  counts_failed 
)
inherited

Definition at line 105 of file Producer.cc.

References art::detail::Producer::checkPutProducts_, e, art::EventPrincipal::makeEvent(), art::detail::Producer::produceWithFrame(), and art::ModuleContext::scheduleID().

110  {
111  auto e = ep.makeEvent(mc);
112  ++counts_run;
113  ProcessingFrame const frame{mc.scheduleID()};
114  produceWithFrame(e, frame);
115  e.commitProducts(checkPutProducts_, &expectedProducts<InEvent>());
116  ++counts_passed;
117  return true;
118  }
bool const checkPutProducts_
Definition: Producer.h:70
Float_t e
Definition: plot.C:35
virtual void produceWithFrame(Event &, ProcessingFrame const &)=0
void art::detail::Producer::doRespondToCloseInputFile ( FileBlock const &  fb)
inherited

Definition at line 44 of file Producer.cc.

References art::detail::Producer::respondToCloseInputFileWithFrame().

45  {
46  ProcessingFrame const frame{ScheduleID{}};
48  }
virtual void respondToCloseInputFileWithFrame(FileBlock const &, ProcessingFrame const &)=0
TFile fb("Li6.root")
void art::detail::Producer::doRespondToCloseOutputFiles ( FileBlock const &  fb)
inherited

Definition at line 58 of file Producer.cc.

References art::detail::Producer::respondToCloseOutputFilesWithFrame().

59  {
60  ProcessingFrame const frame{ScheduleID{}};
62  }
virtual void respondToCloseOutputFilesWithFrame(FileBlock const &, ProcessingFrame const &)=0
TFile fb("Li6.root")
void art::detail::Producer::doRespondToOpenInputFile ( FileBlock const &  fb)
inherited

Definition at line 37 of file Producer.cc.

References art::detail::Producer::respondToOpenInputFileWithFrame().

38  {
39  ProcessingFrame const frame{ScheduleID{}};
41  }
virtual void respondToOpenInputFileWithFrame(FileBlock const &, ProcessingFrame const &)=0
TFile fb("Li6.root")
void art::detail::Producer::doRespondToOpenOutputFiles ( FileBlock const &  fb)
inherited

Definition at line 51 of file Producer.cc.

References art::detail::Producer::respondToOpenOutputFilesWithFrame().

52  {
53  ProcessingFrame const frame{ScheduleID{}};
55  }
virtual void respondToOpenOutputFilesWithFrame(FileBlock const &, ProcessingFrame const &)=0
TFile fb("Li6.root")
void evgen::BaseRadioGen::endJob ( )
inlinevirtual

Reimplemented from art::EDProducer.

Definition at line 502 of file BaseRadioGen.h.

503  {
504  if (m_nevent) {
505  m_pdg_TH1D->Scale(1. / m_nevent);
506 
507  for (auto histo : m_pos_xy_TH2D) {
508  m_pos_xy_TH2D[histo.first]->Scale(1. / m_nevent);
509  m_pos_xz_TH2D[histo.first]->Scale(1. / m_nevent);
510  m_dir_x_TH1D[histo.first]->Scale(1. / m_nevent);
511  m_dir_y_TH1D[histo.first]->Scale(1. / m_nevent);
512  m_dir_z_TH1D[histo.first]->Scale(1. / m_nevent);
513  m_mom_TH1D[histo.first]->Scale(1. / m_nevent);
514  m_ke_TH1D[histo.first]->Scale(1. / m_nevent);
515  m_time_TH1D[histo.first]->Scale(1. / m_nevent);
516  }
517  }
518  endJob_radio();
519  }
std::map< int, TH1D * > m_dir_y_TH1D
Definition: BaseRadioGen.h:184
std::map< int, TH1D * > m_dir_x_TH1D
Definition: BaseRadioGen.h:183
virtual void endJob_radio()
Definition: BaseRadioGen.h:82
std::map< int, TH1D * > m_mom_TH1D
Definition: BaseRadioGen.h:186
std::map< int, TH1D * > m_ke_TH1D
Definition: BaseRadioGen.h:187
std::map< int, TH1D * > m_dir_z_TH1D
Definition: BaseRadioGen.h:185
std::map< int, TH2D * > m_pos_xz_TH2D
Definition: BaseRadioGen.h:182
std::map< int, TH1D * > m_time_TH1D
Definition: BaseRadioGen.h:188
std::map< int, TH2D * > m_pos_xy_TH2D
Definition: BaseRadioGen.h:181
virtual void evgen::BaseRadioGen::endJob_radio ( )
inlinevirtual

Reimplemented in evgen::Decay0Gen.

Definition at line 82 of file BaseRadioGen.h.

References dirCalc(), FillHistos(), GetGoodPositionTime(), GetNDecays(), and part.

82 {}
void evgen::BaseRadioGen::FillAllNodes ( const TGeoNode *  curnode)
inlineprivate

Definition at line 658 of file BaseRadioGen.h.

659  {
660  if (!curnode) return;
661  TObjArray* daunodes = curnode->GetNodes();
662  if (!daunodes) return;
663 
664  TIter next(daunodes);
665 
666  const TGeoNode* anode = 0;
667  while ((anode = (const TGeoNode*)next())) {
668  m_all_nodes.insert(anode);
669  FillAllNodes(anode);
670  }
671  }
std::set< const TGeoNode * > m_all_nodes
Definition: BaseRadioGen.h:120
void FillAllNodes(const TGeoNode *curnode)
Definition: BaseRadioGen.h:658
void evgen::BaseRadioGen::FillHistos ( simb::MCParticle part)
inlineprotected

Definition at line 788 of file BaseRadioGen.h.

Referenced by endJob_radio(), evgen::SpectrumVolumeGen::produce_radio(), and evgen::Decay0Gen::produce_radio().

789  {
790  int pdg = part.PdgCode();
791  int simple_pdg = 0;
792  std::string dummy = "";
793  SimplePDG(pdg, simple_pdg, dummy);
794 
795  TLorentzVector position = part.Position();
796  TLorentzVector mom = part.Momentum();
797  double mass = mom.M();
798  m_pos_xy_TH2D[simple_pdg]->Fill(position.X(), position.Y());
799  m_pos_xz_TH2D[simple_pdg]->Fill(position.X(), position.Z());
800  m_dir_x_TH1D[simple_pdg]->Fill(mom.Px() / mom.P());
801  m_dir_y_TH1D[simple_pdg]->Fill(mom.Py() / mom.P());
802  m_dir_z_TH1D[simple_pdg]->Fill(mom.Pz() / mom.P());
803  double ke = (sqrt(mom.P() * mom.P() + mass * mass) - mass) * 1000.;
804  m_ke_TH1D[simple_pdg]->Fill(ke);
805  m_mom_TH1D[simple_pdg]->Fill(mom.P() * 1000.);
806  m_time_TH1D[simple_pdg]->Fill(position.T());
807  m_pdg_TH1D->Fill(simple_pdg);
808  }
void SimplePDG(int pdg, int &simple, std::string &name)
Definition: BaseRadioGen.h:754
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:220
int PdgCode() const
Definition: MCParticle.h:213
std::map< int, TH1D * > m_dir_y_TH1D
Definition: BaseRadioGen.h:184
std::map< int, TH1D * > m_dir_x_TH1D
Definition: BaseRadioGen.h:183
std::map< int, TH1D * > m_mom_TH1D
Definition: BaseRadioGen.h:186
std::map< int, TH1D * > m_ke_TH1D
Definition: BaseRadioGen.h:187
std::map< int, TH1D * > m_dir_z_TH1D
Definition: BaseRadioGen.h:185
const TLorentzVector & Momentum(const int i=0) const
Definition: MCParticle.h:221
std::map< int, TH2D * > m_pos_xz_TH2D
Definition: BaseRadioGen.h:182
std::map< int, TH1D * > m_time_TH1D
Definition: BaseRadioGen.h:188
std::map< int, TH2D * > m_pos_xy_TH2D
Definition: BaseRadioGen.h:181
void art::Modifier::fillProductDescriptions ( )
inherited

Definition at line 10 of file Modifier.cc.

References art::ProductRegistryHelper::fillDescriptions(), and art::ModuleBase::moduleDescription().

11  {
13  }
void fillDescriptions(ModuleDescription const &md)
ModuleDescription const & moduleDescription() const
Definition: ModuleBase.cc:13
bool evgen::BaseRadioGen::findMotherNode ( const TGeoNode *  cur_node,
std::string &  daughter_name,
const TGeoNode *&  mother_node 
)
inlineprivate

Adapted from above.

Definition at line 697 of file BaseRadioGen.h.

Referenced by GetNodeXYZMinMax().

701  {
702  TObjArray* daunodes = cur_node->GetNodes();
703 
704  if (!daunodes) return false;
705 
706  TIter next(daunodes);
707 
708  const TGeoNode* anode = 0;
709  bool found = 0;
710  while ((anode = (const TGeoNode*)next())) {
711  std::string nname = anode->GetName();
712  std::string vname = anode->GetVolume()->GetName();
713 
714  if (nname == daughter_name || vname == daughter_name) {
715  mother_node = cur_node;
716  return true;
717  }
718  found = findMotherNode(anode, daughter_name, mother_node);
719  }
720 
721  return found;
722  }
bool findMotherNode(const TGeoNode *cur_node, std::string &daughter_name, const TGeoNode *&mother_node)
Adapted from above.
Definition: BaseRadioGen.h:697
bool evgen::BaseRadioGen::findNode ( const TGeoNode *  curnode,
std::string &  tgtnname,
const TGeoNode *&  targetnode 
)
inlineprivate

Shamelessly stolen from here: https://cdcvs.fnal.gov/redmine/attachments/6719/calc_bbox.C.

Definition at line 673 of file BaseRadioGen.h.

677  {
678  std::string nname = curnode->GetName();
679  std::string vname = curnode->GetVolume()->GetName();
680  if (nname == tgtnname || vname == tgtnname) {
681  targetnode = curnode;
682  return true;
683  }
684 
685  TObjArray* daunodes = curnode->GetNodes();
686  if (!daunodes) return false;
687  TIter next(daunodes);
688  const TGeoNode* anode = 0;
689  while ((anode = (const TGeoNode*)next())) {
690  bool found = findNode(anode, tgtnname, targetnode);
691  if (found) return true;
692  }
693 
694  return false;
695  }
bool findNode(const TGeoNode *curnode, std::string &tgtnname, const TGeoNode *&targetnode)
Shamelessly stolen from here: https://cdcvs.fnal.gov/redmine/attachments/6719/calc_bbox.C.
Definition: BaseRadioGen.h:673
std::array< std::vector< ProductInfo >, NumBranchTypes > const & art::ModuleBase::getConsumables ( ) const
inherited

Definition at line 43 of file ModuleBase.cc.

References art::ModuleBase::collector_, and art::ConsumesCollector::getConsumables().

44  {
45  return collector_.getConsumables();
46  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables() const
bool evgen::BaseRadioGen::GetGoodPositionTime ( TLorentzVector &  position)
inlineprotected

Deal with the time first

We use all the XYZ that have been set earlier, until we find the correct volume and material

We use the list of nodes

Definition at line 521 of file BaseRadioGen.h.

Referenced by endJob_radio(), evgen::SpectrumVolumeGen::produce_radio(), and evgen::Decay0Gen::produce_radio().

522  {
523 
525  double time = m_T0 + m_random_flat.fire() * (m_T1 - m_T0);
526 
527  if (not m_geo_volume_mode) {
529 
530  size_t n_tries = 0;
531 
532  double xpos = std::numeric_limits<double>::signaling_NaN();
533  double ypos = std::numeric_limits<double>::signaling_NaN();
534  double zpos = std::numeric_limits<double>::signaling_NaN();
535 
536  while (n_tries++ < m_max_tries_event) {
537 
538  if (m_distrib_xpos)
539  xpos = m_distrib_xpos->GetRandom(m_X0, m_X1);
540  else
541  xpos = m_X0 + m_random_flat.fire() * (m_X1 - m_X0);
542 
543  if (m_distrib_ypos)
544  ypos = m_distrib_ypos->GetRandom(m_Y0, m_Y1);
545  else
546  ypos = m_Y0 + m_random_flat.fire() * (m_Y1 - m_Y0);
547 
548  if (m_distrib_zpos)
549  zpos = m_distrib_zpos->GetRandom(m_Z0, m_Z1);
550  else
551  zpos = m_Z0 + m_random_flat.fire() * (m_Z1 - m_Z0);
552 
553  auto node = gGeoManager->FindNode(xpos, ypos, zpos);
554 
555  if (NodeSupported(node)) break;
556  }
557 
558  if (std::isnan(xpos) or std::isnan(ypos) or std::isnan(zpos)) {
559  MF_LOG_ERROR("BaseRadioGen") << "Error in generation of random position!";
560  }
561 
562  position.SetXYZT(xpos, ypos, zpos, time);
563  return true;
564  }
565  else {
567 
568  if (m_good_nodes.empty() or m_volume_cc == 0)
569  MF_LOG_ERROR("BaseRadioGen") << "There is no node to throw events in!";
570 
571  double which_vol = m_random_flat.fire() * m_volume_cc;
572 
573  const TGeoNode* node = nullptr;
574  size_t i = 0;
575 
576  for (auto const& nd : m_good_nodes) {
577  if (which_vol < nd.second) {
578  node = nd.first;
579  break;
580  }
581  i++;
582  }
583 
584  const TGeoNode* world = gGeoManager->GetNode(0);
585  double x_min, x_max;
586  double y_min, y_max;
587  double z_min, z_max;
588  GetNodeXYZMinMax(node, world, x_min, x_max, y_min, y_max, z_min, z_max);
589 
590  size_t n_tries = 0;
591 
592  while (n_tries++ < m_max_tries_event) {
593  double xpos = std::numeric_limits<double>::signaling_NaN();
594  double ypos = std::numeric_limits<double>::signaling_NaN();
595  double zpos = std::numeric_limits<double>::signaling_NaN();
596 
597  if (m_distrib_xpos)
598  xpos = m_distrib_xpos->GetRandom(x_min, x_max);
599  else
600  xpos = x_min + m_random_flat.fire() * (x_max - x_min);
601 
602  if (m_distrib_ypos)
603  ypos = m_distrib_ypos->GetRandom(y_min, y_max);
604  else
605  ypos = y_min + m_random_flat.fire() * (y_max - y_min);
606 
607  if (m_distrib_zpos)
608  zpos = m_distrib_zpos->GetRandom(z_min, z_max);
609  else
610  zpos = z_min + m_random_flat.fire() * (z_max - z_min);
611 
612  if (std::isnan(xpos) or std::isnan(ypos) or std::isnan(zpos)) {
613  MF_LOG_ERROR("BaseRadioGen") << "Error in generation of random position!";
614  }
615  position.SetXYZT(xpos, ypos, zpos, time);
616 
617  const TGeoNode* node_generated = gGeoManager->FindNode(xpos, ypos, zpos);
618  if (node_generated == node) return true;
619 
620  if (node->GetNdaughters() and IsDaughterNode(node_generated, node)) return true;
621  }
622  }
623 
624  return false;
625  }
double m_X1
Top corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:152
std::map< const TGeoNode *, double > m_good_nodes
Definition: BaseRadioGen.h:172
bool NodeSupported(TGeoNode const *node) const
Definition: BaseRadioGen.h:192
#define MF_LOG_ERROR(category)
STL namespace.
double x_min
Definition: berger.C:15
double m_Y1
Top corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:153
double m_Z1
Top corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:154
double m_T0
Beginning of time window to simulate in ns.
Definition: BaseRadioGen.h:147
double x_max
Definition: berger.C:16
double m_Z0
Bottom corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:151
bool IsDaughterNode(const TGeoNode *daughter_node, const TGeoNode *mother_node)
Definition: BaseRadioGen.h:639
CLHEP::RandFlat m_random_flat
Definition: BaseRadioGen.h:161
double m_Y0
Bottom corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:150
double m_X0
Bottom corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:149
double m_T1
End of time window to simulate in ns.
Definition: BaseRadioGen.h:148
void GetNodeXYZMinMax(const TGeoNode *from, const TGeoNode *to, double &x_min, double &x_max, double &y_min, double &y_max, double &z_min, double &z_max)
Definition: BaseRadioGen.h:199
int evgen::BaseRadioGen::GetNDecays ( )
inlineprotected

Definition at line 629 of file BaseRadioGen.h.

Referenced by endJob_radio(), evgen::SpectrumVolumeGen::produce_radio(), and evgen::Decay0Gen::produce_radio().

630  {
631 
632  if (m_rate_mode) { return m_random_poisson.fire(m_rate); }
633  else {
634  double rate = abs(m_Bq * (m_T1 - m_T0) * m_volume_cc / 1.0E9);
635  return m_random_poisson.fire(rate);
636  }
637  }
CLHEP::RandPoisson m_random_poisson
Definition: BaseRadioGen.h:162
constexpr auto abs(T v)
Returns the absolute value of the argument.
double m_rate
Radioactivity in Becquerels (decay per sec) use either of this of Bq.
Definition: BaseRadioGen.h:145
double m_T0
Beginning of time window to simulate in ns.
Definition: BaseRadioGen.h:147
double m_Bq
Radioactivity in Becquerels (decay per sec) per cubic cm.
Definition: BaseRadioGen.h:144
double m_T1
End of time window to simulate in ns.
Definition: BaseRadioGen.h:148
int evgen::BaseRadioGen::GetNEvents ( )
inlineprotected

Definition at line 93 of file BaseRadioGen.h.

References m_nevent.

Referenced by evgen::Decay0Gen::endJob_radio().

93 { return m_nevent; }
void evgen::BaseRadioGen::GetNodeXYZMinMax ( const TGeoNode *  from,
const TGeoNode *  to,
double &  x_min,
double &  x_max,
double &  y_min,
double &  y_max,
double &  z_min,
double &  z_max 
)
inlineprivate

Definition at line 199 of file BaseRadioGen.h.

References findMotherNode(), and geo::origin().

207  {
208  std::vector<const TGeoNode*> mother_nodes;
209  const TGeoNode* current_node = from;
210  std::string daughter_name = from->GetName();
211  int nmax = 20;
212  int iter = 0;
213  while (current_node != to and iter++ < nmax) {
214  const TGeoNode* mother_node = nullptr;
215  daughter_name = current_node->GetName();
216  bool found_mum = findMotherNode(to, daughter_name, mother_node);
217  if (not found_mum) {
218  throw cet::exception("BaseRadioGen")
219  << "Didn't find the mum of the following node: " << daughter_name;
220  }
221  mother_nodes.push_back(mother_node);
222  current_node = mother_node;
223  }
224 
225  TGeoVolume* vol = from->GetVolume();
226  TGeoShape* shape = vol->GetShape();
227  TGeoBBox* bbox = (TGeoBBox*)shape;
228 
229  double dx = bbox->GetDX();
230  double dy = bbox->GetDY();
231  double dz = bbox->GetDZ();
232 
233  double halfs[3] = {dx, dy, dz};
234  double posmin[3] = {1.0e30, 1.0e30, 1.0e30};
235  double posmax[3] = {-1.0e30, -1.0e30, -1.0e30};
236 
237  const double* origin = bbox->GetOrigin();
238  for (int ix = -1; ix <= 1; ix += 2) {
239  for (int iy = -1; iy <= 1; iy += 2) {
240  for (int iz = -1; iz <= 1; iz += 2) {
241  double local[3];
242  local[0] = origin[0] + (double)ix * halfs[0];
243  local[1] = origin[1] + (double)iy * halfs[1];
244  local[2] = origin[2] + (double)iz * halfs[2];
245  double master[3];
246  from->LocalToMaster(local, master);
247  for (auto const& mum : mother_nodes) {
248  local[0] = master[0];
249  local[1] = master[1];
250  local[2] = master[2];
251  mum->LocalToMaster(local, master);
252  }
253  for (int j = 0; j < 3; ++j) {
254  posmin[j] = TMath::Min(posmin[j], master[j]);
255  posmax[j] = TMath::Max(posmax[j], master[j]);
256  }
257  }
258  }
259  }
260 
261  x_min = posmin[0];
262  y_min = posmin[1];
263  z_min = posmin[2];
264  x_max = posmax[0];
265  y_max = posmax[1];
266  z_max = posmax[2];
267  }
double x_min
Definition: berger.C:15
bool findMotherNode(const TGeoNode *cur_node, std::string &daughter_name, const TGeoNode *&mother_node)
Adapted from above.
Definition: BaseRadioGen.h:697
double x_max
Definition: berger.C:16
constexpr Point origin()
Returns a origin position with a point of the specified type.
Definition: geo_vectors.h:229
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
CLHEP::HepRandomEngine& evgen::BaseRadioGen::GetRandomEngine ( )
inlineprotected

Definition at line 92 of file BaseRadioGen.h.

References m_engine.

Referenced by evgen::Decay0Gen::Decay0Gen().

92 { return m_engine; }
CLHEP::HepRandomEngine & m_engine
Definition: BaseRadioGen.h:160
void evgen::BaseRadioGen::GetTs ( double &  T0,
double &  T1 
)
inlineprotected

Definition at line 94 of file BaseRadioGen.h.

References m_T0, and m_T1.

Referenced by evgen::Decay0Gen::beginJob_radio().

95  {
96  T0 = m_T0;
97  T1 = m_T1;
98  }
double m_T0
Beginning of time window to simulate in ns.
Definition: BaseRadioGen.h:147
double m_T1
End of time window to simulate in ns.
Definition: BaseRadioGen.h:148
void evgen::BaseRadioGen::GetXs ( double &  X0,
double &  X1 
)
inlineprotected

Definition at line 99 of file BaseRadioGen.h.

References m_X0, and m_X1.

100  {
101  X0 = m_X0;
102  X1 = m_X1;
103  }
double m_X1
Top corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:152
Double_t X1
Definition: plot.C:263
double m_X0
Bottom corner x position (cm) in world coordinates.
Definition: BaseRadioGen.h:149
void evgen::BaseRadioGen::GetYs ( double &  Y0,
double &  Y1 
)
inlineprotected

Definition at line 104 of file BaseRadioGen.h.

References m_Y0, and m_Y1.

105  {
106  Y0 = m_Y0;
107  Y1 = m_Y1;
108  }
double m_Y1
Top corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:153
Double_t Y1
Definition: plot.C:263
double m_Y0
Bottom corner y position (cm) in world coordinates.
Definition: BaseRadioGen.h:150
void evgen::BaseRadioGen::GetZs ( double &  Z0,
double &  Z1 
)
inlineprotected

Definition at line 109 of file BaseRadioGen.h.

References CalculateActiveVolumeFromAllNodes(), CalculateActiveVolumeFromXYZ(), DeclareOutputHistos(), m_Z0, m_Z1, and SimplePDG().

110  {
111  Z0 = m_Z0;
112  Z1 = m_Z1;
113  }
double m_Z1
Top corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:154
Double_t Z1
Definition: plot.C:263
double m_Z0
Bottom corner z position (cm) in world coordinates.
Definition: BaseRadioGen.h:151
bool evgen::BaseRadioGen::IsDaughterNode ( const TGeoNode *  daughter_node,
const TGeoNode *  mother_node 
)
inlineprivate

Definition at line 639 of file BaseRadioGen.h.

641  {
642  if (mother_node == daughter_node) return true;
643 
644  TObjArray* daunodes = mother_node->GetNodes();
645  if (!daunodes) return false;
646 
647  TIter next(daunodes);
648  const TGeoNode* anode = 0;
649 
650  while ((anode = (const TGeoNode*)next())) {
651  bool found = IsDaughterNode(daughter_node, anode);
652  if (found) return true;
653  }
654 
655  return false;
656  }
bool IsDaughterNode(const TGeoNode *daughter_node, const TGeoNode *mother_node)
Definition: BaseRadioGen.h:639
std::unique_ptr< Worker > art::ModuleBase::makeWorker ( WorkerParams const &  wp)
inherited

Definition at line 37 of file ModuleBase.cc.

References art::ModuleBase::doMakeWorker(), and art::NumBranchTypes.

38  {
39  return doMakeWorker(wp);
40  }
virtual std::unique_ptr< Worker > doMakeWorker(WorkerParams const &wp)=0
template<typename T , BranchType BT>
ProductToken< T > art::ModuleBase::mayConsume ( InputTag const &  tag)
protectedinherited

Definition at line 82 of file ModuleBase.h.

References art::ModuleBase::collector_, and art::ConsumesCollector::mayConsume().

83  {
84  return collector_.mayConsume<T, BT>(tag);
85  }
ProductToken< T > mayConsume(InputTag const &)
ConsumesCollector collector_
Definition: ModuleBase.h:56
template<typename T , BranchType BT>
void art::ModuleBase::mayConsumeMany ( )
protectedinherited

Definition at line 96 of file ModuleBase.h.

References art::ModuleBase::collector_, and art::ConsumesCollector::mayConsumeMany().

97  {
98  collector_.mayConsumeMany<T, BT>();
99  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
template<typename Element , BranchType = InEvent>
ViewToken<Element> art::ModuleBase::mayConsumeView ( InputTag const &  )
protectedinherited
template<typename T , BranchType BT>
ViewToken<T> art::ModuleBase::mayConsumeView ( InputTag const &  tag)
inherited

Definition at line 89 of file ModuleBase.h.

References art::ModuleBase::collector_, and art::ConsumesCollector::mayConsumeView().

90  {
91  return collector_.mayConsumeView<T, BT>(tag);
92  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
ViewToken< Element > mayConsumeView(InputTag const &)
ModuleDescription const & art::ModuleBase::moduleDescription ( ) const
inherited

Definition at line 13 of file ModuleBase.cc.

References art::errors::LogicError.

Referenced by art::OutputModule::doRespondToOpenInputFile(), art::OutputModule::doWriteEvent(), art::Modifier::fillProductDescriptions(), art::OutputModule::makePlugins_(), art::OutputWorker::OutputWorker(), reco::shower::LArPandoraModularShowerCreation::produce(), art::Modifier::registerProducts(), and art::OutputModule::registerProducts().

14  {
15  if (md_.has_value()) {
16  return *md_;
17  }
18 
20  "There was an error while calling moduleDescription().\n"}
21  << "The moduleDescription() base-class member function cannot be called\n"
22  "during module construction. To determine which module is "
23  "responsible\n"
24  "for calling it, find the '<module type>:<module "
25  "label>@Construction'\n"
26  "tag in the message prefix above. Please contact artists@fnal.gov\n"
27  "for guidance.\n";
28  }
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
std::optional< ModuleDescription > md_
Definition: ModuleBase.h:55
bool evgen::BaseRadioGen::NodeSupported ( TGeoNode const *  node) const
inlineprivate

Definition at line 192 of file BaseRadioGen.h.

References m_good_materials, and m_good_volumes.

193  {
194  assert(node != nullptr);
195  return cet::search_all(m_good_volumes, node->GetVolume()) and
196  cet::search_all(m_good_materials, node->GetMedium()->GetMaterial());
197  }
std::vector< const TGeoVolume * > m_good_volumes
Definition: BaseRadioGen.h:173
std::vector< const TGeoMaterial * > m_good_materials
Definition: BaseRadioGen.h:174
void evgen::BaseRadioGen::produce ( art::Event evt)
inlinevirtual

Implements art::EDProducer.

Definition at line 477 of file BaseRadioGen.h.

478  {
479  m_nevent++;
480  produce_radio(evt);
481  }
virtual void produce_radio(art::Event &evt)=0
virtual void evgen::BaseRadioGen::produce_radio ( art::Event evt)
pure virtual
void art::Modifier::registerProducts ( ProductDescriptions productsToRegister)
inherited

Definition at line 16 of file Modifier.cc.

References art::ModuleBase::moduleDescription(), and art::ProductRegistryHelper::registerProducts().

17  {
18  ProductRegistryHelper::registerProducts(productsToRegister,
20  }
void registerProducts(ProductDescriptions &productsToRegister, ModuleDescription const &md)
ModuleDescription const & moduleDescription() const
Definition: ModuleBase.cc:13
void art::ModuleBase::setModuleDescription ( ModuleDescription const &  md)
inherited

Definition at line 31 of file ModuleBase.cc.

References art::ModuleBase::md_.

32  {
33  md_ = md;
34  }
std::optional< ModuleDescription > md_
Definition: ModuleBase.h:55
void evgen::BaseRadioGen::SimplePDG ( int  pdg,
int &  simple,
std::string &  name 
)
inlineprivate

Definition at line 754 of file BaseRadioGen.h.

Referenced by GetZs().

755  {
756  switch (pdg) {
757  case 1000020040:
758  name = "alpha";
759  simple = 0;
760  return;
761  case 22:
762  name = "gamma";
763  simple = 1;
764  return;
765  case -11:
766  name = "positron";
767  simple = 2;
768  return;
769  case 11:
770  name = "electron";
771  simple = 3;
772  return;
773  case 2112:
774  name = "neutron";
775  simple = 4;
776  return;
777  case 2212:
778  name = "proton";
779  simple = 5;
780  return;
781  default:
782  name = "other";
783  simple = 6;
784  return;
785  }
786  }
void art::ModuleBase::sortConsumables ( std::string const &  current_process_name)
inherited

Definition at line 49 of file ModuleBase.cc.

References art::ModuleBase::collector_, and art::ConsumesCollector::sortConsumables().

50  {
51  // Now that we know we have seen all the consumes declarations,
52  // sort the results for fast lookup later.
53  collector_.sortConsumables(current_process_name);
54  }
ConsumesCollector collector_
Definition: ModuleBase.h:56
void sortConsumables(std::string const &current_process_name)

Member Data Documentation

std::set<const TGeoNode*> evgen::BaseRadioGen::m_all_nodes
private

Definition at line 120 of file BaseRadioGen.h.

double evgen::BaseRadioGen::m_Bq {signaling_NaN}
private

Radioactivity in Becquerels (decay per sec) per cubic cm.

Definition at line 144 of file BaseRadioGen.h.

std::map<int, TH1D*> evgen::BaseRadioGen::m_dir_x_TH1D {}
private

Definition at line 183 of file BaseRadioGen.h.

std::map<int, TH1D*> evgen::BaseRadioGen::m_dir_y_TH1D {}
private

Definition at line 184 of file BaseRadioGen.h.

std::map<int, TH1D*> evgen::BaseRadioGen::m_dir_z_TH1D {}
private

Definition at line 185 of file BaseRadioGen.h.

TF1* evgen::BaseRadioGen::m_distrib_xpos {nullptr}
private

Definition at line 176 of file BaseRadioGen.h.

TF1* evgen::BaseRadioGen::m_distrib_ypos {nullptr}
private

Definition at line 177 of file BaseRadioGen.h.

TF1* evgen::BaseRadioGen::m_distrib_zpos {nullptr}
private

Definition at line 178 of file BaseRadioGen.h.

CLHEP::HepRandomEngine& evgen::BaseRadioGen::m_engine
private

Definition at line 160 of file BaseRadioGen.h.

Referenced by BaseRadioGen(), and GetRandomEngine().

TGeoManager* evgen::BaseRadioGen::m_geo_manager
private

Definition at line 158 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

art::ServiceHandle<geo::Geometry const> evgen::BaseRadioGen::m_geo_service
private

Definition at line 157 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

bool evgen::BaseRadioGen::m_geo_volume_mode {false}
private

Definition at line 164 of file BaseRadioGen.h.

std::vector<const TGeoMaterial*> evgen::BaseRadioGen::m_good_materials {}
private

Definition at line 174 of file BaseRadioGen.h.

Referenced by NodeSupported().

std::map<const TGeoNode*, double> evgen::BaseRadioGen::m_good_nodes {}
private

Definition at line 172 of file BaseRadioGen.h.

std::vector<const TGeoVolume*> evgen::BaseRadioGen::m_good_volumes {}
private

Definition at line 173 of file BaseRadioGen.h.

Referenced by NodeSupported().

std::vector<std::string> evgen::BaseRadioGen::m_isotope
protected

isotope to simulate. Example: "Ar39"

Definition at line 90 of file BaseRadioGen.h.

Referenced by evgen::Decay0Gen::Decay0Gen(), and evgen::SpectrumVolumeGen::SpectrumVolumeGen().

std::map<int, TH1D*> evgen::BaseRadioGen::m_ke_TH1D {}
private

Definition at line 187 of file BaseRadioGen.h.

std::string evgen::BaseRadioGen::m_material
private

regex of materials in which to generate the decays. Example: "LAr"

Definition at line 137 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

size_t evgen::BaseRadioGen::m_max_tries_event
private

Definition at line 168 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

size_t evgen::BaseRadioGen::m_max_tries_rate_calculation
private

Definition at line 169 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

std::map<int, TH1D*> evgen::BaseRadioGen::m_mom_TH1D {}
private

Definition at line 186 of file BaseRadioGen.h.

int evgen::BaseRadioGen::m_nevent {0}
private

Definition at line 180 of file BaseRadioGen.h.

Referenced by GetNEvents().

TH1D* evgen::BaseRadioGen::m_pdg_TH1D {nullptr}
private

Definition at line 189 of file BaseRadioGen.h.

std::map<int, TH2D*> evgen::BaseRadioGen::m_pos_xy_TH2D {}
private

Definition at line 181 of file BaseRadioGen.h.

std::map<int, TH2D*> evgen::BaseRadioGen::m_pos_xz_TH2D {}
private

Definition at line 182 of file BaseRadioGen.h.

CLHEP::RandFlat evgen::BaseRadioGen::m_random_flat
private

Definition at line 161 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

CLHEP::RandPoisson evgen::BaseRadioGen::m_random_poisson
private

Definition at line 162 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

double evgen::BaseRadioGen::m_rate
private
Initial value:

Radioactivity in Becquerels (decay per sec) use either of this of Bq.

Definition at line 145 of file BaseRadioGen.h.

bool evgen::BaseRadioGen::m_rate_mode
private

Definition at line 165 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

std::regex evgen::BaseRadioGen::m_regex_material
private

Definition at line 140 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

std::regex evgen::BaseRadioGen::m_regex_volume
private

Definition at line 141 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

double evgen::BaseRadioGen::m_T0 {signaling_NaN}
private

Beginning of time window to simulate in ns.

Definition at line 147 of file BaseRadioGen.h.

Referenced by GetTs().

double evgen::BaseRadioGen::m_T1 {signaling_NaN}
private

End of time window to simulate in ns.

Definition at line 148 of file BaseRadioGen.h.

Referenced by GetTs().

size_t evgen::BaseRadioGen::m_target_n_point_rate_calculation
private

Definition at line 170 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

std::map<int, TH1D*> evgen::BaseRadioGen::m_time_TH1D {}
private

Definition at line 188 of file BaseRadioGen.h.

double evgen::BaseRadioGen::m_volume_cc {signaling_NaN}
private

Definition at line 155 of file BaseRadioGen.h.

std::string evgen::BaseRadioGen::m_volume_gen
private

The volume in which to generate the decays.

Definition at line 139 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

std::string evgen::BaseRadioGen::m_volume_rand
private

The volume in which to generate the decays.

Definition at line 138 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

bool evgen::BaseRadioGen::m_volume_rand_present
private

Definition at line 166 of file BaseRadioGen.h.

Referenced by BaseRadioGen().

double evgen::BaseRadioGen::m_X0 {signaling_NaN}
private

Bottom corner x position (cm) in world coordinates.

Definition at line 149 of file BaseRadioGen.h.

Referenced by GetXs().

double evgen::BaseRadioGen::m_X1 {signaling_NaN}
private

Top corner x position (cm) in world coordinates.

Definition at line 152 of file BaseRadioGen.h.

Referenced by GetXs().

double evgen::BaseRadioGen::m_Y0 {signaling_NaN}
private

Bottom corner y position (cm) in world coordinates.

Definition at line 150 of file BaseRadioGen.h.

Referenced by GetYs().

double evgen::BaseRadioGen::m_Y1 {signaling_NaN}
private

Top corner y position (cm) in world coordinates.

Definition at line 153 of file BaseRadioGen.h.

Referenced by GetYs().

double evgen::BaseRadioGen::m_Z0 {signaling_NaN}
private

Bottom corner z position (cm) in world coordinates.

Definition at line 151 of file BaseRadioGen.h.

Referenced by GetZs().

double evgen::BaseRadioGen::m_Z1 {signaling_NaN}
private

Top corner z position (cm) in world coordinates.

Definition at line 154 of file BaseRadioGen.h.

Referenced by GetZs().

constexpr double evgen::BaseRadioGen::signaling_NaN = std::numeric_limits<double>::signaling_NaN()
staticprivate

Definition at line 143 of file BaseRadioGen.h.


The documentation for this class was generated from the following file: