LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
PDFastSimPAR_module.cc
Go to the documentation of this file.
1 // Class: PDFastSimPAR
3 // Plugin Type: producer
4 // File: PDFastSimPAR_module.cc
5 // Description:
6 // - acts on sim::SimEnergyDeposit from LArG4Main,
7 // - simulate (fast, photon visibility service) the OpDet response to optical
8 // photons Input: 'sim::SimEnergyDeposit' Output: 'sim::OpDetBacktrackerRecord'
9 // Fast simulation of propagating the photons created from SimEnergyDeposits.
10 
11 // This module does a fast simulation of propagating the photons created from
12 // SimEnergyDeposits, This simulation is done using the Semi-Analytical model, which
13 // stores the visibilities of each optical channel with respect to each optical
14 // voxel in the TPC volume, to avoid propagating single photons using Geant4. At
15 // the end of this module a collection of the propagated photons either as
16 //'sim::OpDetBacktrackerRecord' are placed into the art event.
17 
18 // The steps this module takes are:
19 // - to take number of photon and the vertex information from
20 // 'sim::SimEnergyDeposits',
21 // - use the Semi-Analytical model (visibilities) to determine the amount of visible
22 // photons at each optical channel,
23 // - visible photons: the number of photons times the visibility at the middle
24 // of the Geant4 step for a given optical channel.
25 // - other photon information is got from 'sim::SimEnergyDeposits'
26 // - add 'sim::OpDetBacktrackerRecord' to event
27 // Aug. 19 by Mu Wei
28 // Restructured Nov. 21 by P. Green
30 
31 // LArSoft libraries
46 
48 
49 // Art libraries
58 #include "cetlib_except/exception.h"
59 #include "fhiclcpp/ParameterSet.h"
60 #include "fhiclcpp/types/Atom.h"
61 #include "fhiclcpp/types/Comment.h"
63 #include "fhiclcpp/types/Name.h"
66 
67 // Random numbers
68 #include "CLHEP/Random/RandPoissonQ.h"
69 
70 #include <cmath>
71 #include <ctime>
72 #include <map>
73 #include <memory>
74 #include <vector>
75 
76 namespace phot {
77  class PDFastSimPAR : public art::EDProducer {
78  public:
79  // Define the fhicl configuration
80  struct Config {
81  using Name = fhicl::Name;
85 
87  Comment("SimEnergyDeposit label.")};
89  Comment("Simulate slow scintillation light, default true"),
90  true};
92  Comment("Simulate slow scintillation light")};
94  Comment("Simulate reflected visible light")};
96  Name("IncludeAnodeReflections"),
97  Comment("Simulate anode reflections, default false"),
98  false};
100  Comment("Simulate light propagation time")};
102  Name("GeoPropTimeOnly"),
103  Comment("Simulate light propagation time geometric approximation, default false"),
104  false};
106  Name("UseLitePhotons"),
107  Comment("Store SimPhotonsLite/OpDetBTRs instead of SimPhotons")};
109  Comment("Photons cannot cross the cathode")};
111  Name("OnlyActiveVolume"),
112  Comment("PAR fast sim usually only for active volume, default true"),
113  true};
115  Comment("Set to true if light is only supported in C:1")};
116  DP ScintTimeTool{Name("ScintTimeTool"),
117  Comment("Tool describing scintillation time structure")};
119  Name("UseXeAbsorption"),
120  Comment("Use Xe absorption length instead of Ar, default false"),
121  false};
122  ODP VUVTiming{Name("VUVTiming"), Comment("Configuration for UV timing parameterization")};
123  ODP VISTiming{Name("VISTiming"),
124  Comment("Configuration for visible timing parameterization")};
125  DP VUVHits{Name("VUVHits"), Comment("Configuration for UV visibility parameterization")};
126  ODP VISHits{Name("VISHits"),
127  Comment("Configuration for visibile visibility parameterization")};
128  };
130 
131  explicit PDFastSimPAR(Parameters const& config);
132  void produce(art::Event&) override;
133 
134  private:
135  void detectedNumPhotons(std::vector<int>& DetectedNumPhotons,
136  const std::vector<double>& OpDetVisibilities,
137  const int NumPhotons) const;
138 
139  void AddOpDetBTR(std::vector<sim::OpDetBacktrackerRecord>& opbtr,
140  std::vector<int>& ChannelMap,
141  const sim::OpDetBacktrackerRecord& btr) const;
142 
143  bool isOpDetInSameTPC(geo::Point_t const& ScintPoint, geo::Point_t const& OpDetPoint) const;
144  std::vector<geo::Point_t> opDetCenters() const;
145 
146  // semi-analytical model
147  std::unique_ptr<SemiAnalyticalModel> fVisibilityModel;
148 
149  // propagation time model
150  std::unique_ptr<PropagationTimeModel> fPropTimeModel;
151 
152  // random numbers
153  CLHEP::HepRandomEngine& fPhotonEngine;
154  std::unique_ptr<CLHEP::RandPoissonQ> fRandPoissPhot;
155  CLHEP::HepRandomEngine& fScintTimeEngine;
156  std::unique_ptr<ScintTime> fScintTime; // Tool to retrieve timinig of scintillation
157 
158  // geometry properties
161  const size_t fNOpChannels;
162  const std::vector<geo::BoxBoundedGeo> fActiveVolumes;
163  const int fNTPC;
164  const std::vector<geo::Point_t> fOpDetCenter;
165 
166  // Module behavior
168  const bool fDoFastComponent;
169  const bool fDoSlowComponent;
170  const bool fDoReflectedLight;
172  const bool fIncludePropTime;
173  const bool fGeoPropTimeOnly;
174  const bool fUseLitePhotons;
175  const bool fOpaqueCathode;
176  const bool fOnlyActiveVolume;
177  const bool fOnlyOneCryostat;
178  const bool fUseXeAbsorption;
179  };
180 
181  //......................................................................
183  : art::EDProducer{config}
184  , fPhotonEngine(art::ServiceHandle<rndm::NuRandomService>()->registerAndSeedEngine(
185  createEngine(0, "HepJamesRandom", "photon"),
186  "HepJamesRandom",
187  "photon",
188  config.get_PSet(),
189  "SeedPhoton"))
190  , fRandPoissPhot(std::make_unique<CLHEP::RandPoissonQ>(fPhotonEngine))
192  createEngine(0, "HepJamesRandom", "scinttime"),
193  "HepJamesRandom",
194  "scinttime",
195  config.get_PSet(),
196  "SeedScintTime"))
197  , fScintTime{art::make_tool<phot::ScintTime>(config().ScintTimeTool.get<fhicl::ParameterSet>())}
198  , fGeom(*(lar::providerFrom<geo::Geometry>()))
199  , fISTPC{fGeom}
202  , fNTPC(fGeom.NTPC())
204  , fSimTag(config().SimulationLabel())
205  , fDoFastComponent(config().DoFastComponent())
206  , fDoSlowComponent(config().DoSlowComponent())
207  , fDoReflectedLight(config().DoReflectedLight())
208  , fIncludeAnodeReflections(config().IncludeAnodeReflections())
209  , fIncludePropTime(config().IncludePropTime())
210  , fGeoPropTimeOnly(config().GeoPropTimeOnly())
211  , fUseLitePhotons(config().UseLitePhotons())
212  , fOpaqueCathode(config().OpaqueCathode())
213  , fOnlyActiveVolume(config().OnlyActiveVolume())
214  , fOnlyOneCryostat(config().OnlyOneCryostat())
215  , fUseXeAbsorption(config().UseXeAbsorption())
216  {
217  mf::LogInfo("PDFastSimPAR") << "Initializing PDFastSimPAR." << std::endl;
218 
219  // Parameterized Simulation
220  fhicl::ParameterSet VUVHitsParams = config().VUVHits.get<fhicl::ParameterSet>();
221  fhicl::ParameterSet VUVTimingParams;
222  fhicl::ParameterSet VISHitsParams;
223  fhicl::ParameterSet VISTimingParams;
224 
225  // Validate configuration options
226  if (fIncludePropTime &&
227  !config().VUVTiming.get_if_present<fhicl::ParameterSet>(VUVTimingParams)) {
229  << "Propagation time simulation requested, but VUVTiming not specified."
230  << "\n";
231  }
233  !config().VISHits.get_if_present<fhicl::ParameterSet>(VISHitsParams)) {
235  << "Reflected light or anode reflections simulation requested, but VisHits not specified."
236  << "\n";
237  }
239  !config().VISTiming.get_if_present<fhicl::ParameterSet>(VISTimingParams)) {
241  << "Reflected light propagation time simulation requested, but VISTiming not specified."
242  << "\n";
243  }
244  if (fGeom.Ncryostats() > 1U) {
245  if (fOnlyOneCryostat) {
246  mf::LogWarning("PDFastSimPAR")
247  << std::string(80, '=') << "\nA detector with " << fGeom.Ncryostats()
248  << " cryostats is configured"
249  << " , and semi-analytic model is requested for scintillation photon propagation."
250  << " THIS CONFIGURATION IS NOT SUPPORTED and it is open to bugs"
251  << " (e.g. scintillation may be detected only in cryostat #0)."
252  << "\nThis would be normally a fatal error, but it has been forcibly overridden."
253  << "\n"
254  << std::string(80, '=');
255  }
256  else {
258  << "Photon propagation via semi-analytic model is not supported yet"
259  << " on detectors with more than one cryostat.";
260  }
261  }
262 
263  // Initialise the Scintillation Time
264  fScintTime->initRand(fScintTimeEngine);
265 
266  // photo-detector visibility model (semi-analytical model)
267  fVisibilityModel = std::make_unique<SemiAnalyticalModel>(
268  VUVHitsParams, VISHitsParams, fDoReflectedLight, fIncludeAnodeReflections, fUseXeAbsorption);
269 
270  // propagation time model
271  if (fIncludePropTime)
272  fPropTimeModel = std::make_unique<PropagationTimeModel>(
273  VUVTimingParams, VISTimingParams, fScintTimeEngine, fDoReflectedLight, fGeoPropTimeOnly);
274 
275  {
276  auto log = mf::LogTrace("PDFastSimPAR") << "PDFastSimPAR: active volume boundaries from "
277  << fActiveVolumes.size() << " volumes:";
278  for (auto const& [iCryo, box] : util::enumerate(fActiveVolumes)) {
279  log << "\n - C:" << iCryo << ": " << box.Min() << " -- " << box.Max() << " cm";
280  }
281  }
282 
283  if (fUseLitePhotons) {
284  mf::LogInfo("PDFastSimPAR") << "Using Lite Photons";
285  produces<std::vector<sim::SimPhotonsLite>>();
286  produces<std::vector<sim::OpDetBacktrackerRecord>>();
287  if (fDoReflectedLight) {
288  mf::LogInfo("PDFastSimPAR") << "Storing Reflected Photons";
289  produces<std::vector<sim::SimPhotonsLite>>("Reflected");
290  produces<std::vector<sim::OpDetBacktrackerRecord>>("Reflected");
291  }
292  }
293  else {
294  mf::LogInfo("PDFastSimPAR") << "Using Sim Photons";
295  produces<std::vector<sim::SimPhotons>>();
296  if (fDoReflectedLight) {
297  mf::LogInfo("PDFastSimPAR") << "Storing Reflected Photons";
298  produces<std::vector<sim::SimPhotons>>("Reflected");
299  }
300  }
301  mf::LogInfo("PDFastSimPAR") << "PDFastSimPAR Initialization finish.\n"
302  << "Simulate using semi-analytic model for number of hits."
303  << std::endl;
304  }
305 
306  //......................................................................
308  {
309  mf::LogTrace("PDFastSimPAR") << "PDFastSimPAR Module Producer"
310  << "EventID: " << event.event();
311 
312  std::vector<int> PDChannelToSOCMapDirect(fNOpChannels, -1); // Where each OpChan is.
313  std::vector<int> PDChannelToSOCMapReflect(fNOpChannels, -1); // Where each OpChan is.
314 
315  // SimPhotonsLite
316  auto phlit = std::make_unique<std::vector<sim::SimPhotonsLite>>();
317  auto opbtr = std::make_unique<std::vector<sim::OpDetBacktrackerRecord>>();
318  auto phlit_ref = std::make_unique<std::vector<sim::SimPhotonsLite>>();
319  auto opbtr_ref = std::make_unique<std::vector<sim::OpDetBacktrackerRecord>>();
320  auto& dir_phlitcol(*phlit);
321  auto& ref_phlitcol(*phlit_ref);
322  // SimPhotons
323  auto phot = std::make_unique<std::vector<sim::SimPhotons>>();
324  auto phot_ref = std::make_unique<std::vector<sim::SimPhotons>>();
325  auto& dir_photcol(*phot);
326  auto& ref_photcol(*phot_ref);
327  if (fUseLitePhotons) {
328  dir_phlitcol.resize(fNOpChannels);
329  ref_phlitcol.resize(fNOpChannels);
330  for (unsigned int i = 0; i < fNOpChannels; i++) {
331  dir_phlitcol[i].OpChannel = i;
332  ref_phlitcol[i].OpChannel = i;
333  }
334  }
335  else { // SimPhotons
336  dir_photcol.resize(fNOpChannels);
337  ref_photcol.resize(fNOpChannels);
338  for (unsigned int i = 0; i < fNOpChannels; i++) {
339  dir_photcol[i].fOpChannel = i;
340  ref_photcol[i].fOpChannel = i;
341  }
342  }
343 
345  if (!event.getByLabel(fSimTag, edepHandle)) {
346  mf::LogError("PDFastSimPAR") << "PDFastSimPAR Module Cannot getByLabel: " << fSimTag;
347  return;
348  }
349 
350  auto const& edeps = edepHandle;
351 
352  int num_points = 0;
353  int num_fastph = 0;
354  int num_slowph = 0;
355  int num_fastdp = 0;
356  int num_slowdp = 0;
357 
358  for (auto const& edepi : *edeps) {
359  num_points++;
360 
361  int nphot_fast = edepi.NumFPhotons();
362  int nphot_slow = edepi.NumSPhotons();
363 
364  num_fastph += nphot_fast;
365  num_slowph += nphot_slow;
366 
367  if (!((nphot_fast > 0 && fDoFastComponent) || (nphot_slow > 0 && fDoSlowComponent))) continue;
368 
369  int trackID = edepi.TrackID();
370  int nphot = edepi.NumPhotons();
371  double edeposit = edepi.Energy() / nphot;
372  double pos[3] = {edepi.MidPointX(), edepi.MidPointY(), edepi.MidPointZ()};
373  geo::Point_t const ScintPoint = {pos[0], pos[1], pos[2]};
374 
375  if (fOnlyActiveVolume && !fISTPC.isScintInActiveVolume(ScintPoint)) continue;
376 
377  // direct light
378  std::vector<int> DetectedNumFast(fNOpChannels);
379  std::vector<int> DetectedNumSlow(fNOpChannels);
380 
381  std::vector<double> OpDetVisibilities;
382  fVisibilityModel->detectedDirectVisibilities(OpDetVisibilities, ScintPoint);
383  detectedNumPhotons(DetectedNumFast, OpDetVisibilities, nphot_fast);
384  detectedNumPhotons(DetectedNumSlow, OpDetVisibilities, nphot_slow);
385 
387  std::vector<int> AnodeDetectedNumFast(fNOpChannels);
388  std::vector<int> AnodeDetectedNumSlow(fNOpChannels);
389 
390  std::vector<double> OpDetVisibilitiesAnode;
391  fVisibilityModel->detectedReflectedVisibilities(OpDetVisibilitiesAnode, ScintPoint, true);
392  detectedNumPhotons(AnodeDetectedNumFast, OpDetVisibilitiesAnode, nphot_fast);
393  detectedNumPhotons(AnodeDetectedNumSlow, OpDetVisibilitiesAnode, nphot_slow);
394 
395  // add to existing count
396  for (size_t i = 0; i < AnodeDetectedNumFast.size(); ++i) {
397  DetectedNumFast[i] += AnodeDetectedNumFast[i];
398  }
399  for (size_t i = 0; i < AnodeDetectedNumSlow.size(); ++i) {
400  DetectedNumSlow[i] += AnodeDetectedNumSlow[i];
401  }
402  }
403 
404  // reflected light, if enabled
405  std::vector<int> ReflDetectedNumFast(fNOpChannels);
406  std::vector<int> ReflDetectedNumSlow(fNOpChannels);
407  if (fDoReflectedLight) {
408  std::vector<double> OpDetVisibilitiesRefl;
409  fVisibilityModel->detectedReflectedVisibilities(OpDetVisibilitiesRefl, ScintPoint, false);
410  detectedNumPhotons(ReflDetectedNumFast, OpDetVisibilitiesRefl, nphot_fast);
411  detectedNumPhotons(ReflDetectedNumSlow, OpDetVisibilitiesRefl, nphot_slow);
412  }
413 
414  // loop through direct photons then reflected photons cases
415  size_t DoReflected = (fDoReflectedLight) ? 1 : 0;
416  for (size_t Reflected = 0; Reflected <= DoReflected; ++Reflected) {
417  for (size_t channel = 0; channel < fNOpChannels; channel++) {
418 
419  if (fOpaqueCathode && !isOpDetInSameTPC(ScintPoint, fOpDetCenter[channel])) continue;
420 
421  int ndetected_fast = DetectedNumFast[channel];
422  int ndetected_slow = DetectedNumSlow[channel];
423  if (Reflected) {
424  ndetected_fast = ReflDetectedNumFast[channel];
425  ndetected_slow = ReflDetectedNumSlow[channel];
426  }
427  if (!((ndetected_fast > 0 && fDoFastComponent) ||
428  (ndetected_slow > 0 && fDoSlowComponent)))
429  continue;
430 
431  // calculate propagation time, does not matter whether fast or slow photon
432  std::vector<double> transport_time;
433  if (fIncludePropTime) {
434  transport_time.resize(ndetected_fast + ndetected_slow);
435  fPropTimeModel->propagationTime(transport_time, ScintPoint, channel, Reflected);
436  }
437 
438  // SimPhotonsLite case
439  if (fUseLitePhotons) {
440  sim::OpDetBacktrackerRecord tmpbtr(channel);
441  if (ndetected_fast > 0 && fDoFastComponent) {
442  int n = ndetected_fast;
443  num_fastdp += n;
444  for (int i = 0; i < n; ++i) {
445  // calculates the time at which the photon was produced
446  double dtime = edepi.StartT() + fScintTime->fastScintTime();
447  if (fIncludePropTime) dtime += transport_time[i];
448  int time = static_cast<int>(std::round(dtime));
449  if (Reflected)
450  ++ref_phlitcol[channel].DetectedPhotons[time];
451  else
452  ++dir_phlitcol[channel].DetectedPhotons[time];
453  tmpbtr.AddScintillationPhotons(trackID, time, 1, pos, edeposit);
454  }
455  }
456  if (ndetected_slow > 0 && fDoSlowComponent) {
457  int n = ndetected_slow;
458  num_slowdp += n;
459  for (int i = 0; i < n; ++i) {
460  double dtime = edepi.StartT() + fScintTime->slowScintTime();
461  if (fIncludePropTime) dtime += transport_time[ndetected_fast + i];
462  int time = static_cast<int>(std::round(dtime));
463  if (Reflected)
464  ++ref_phlitcol[channel].DetectedPhotons[time];
465  else
466  ++dir_phlitcol[channel].DetectedPhotons[time];
467  tmpbtr.AddScintillationPhotons(trackID, time, 1, pos, edeposit);
468  }
469  }
470  if (Reflected)
471  AddOpDetBTR(*opbtr_ref, PDChannelToSOCMapReflect, tmpbtr);
472  else
473  AddOpDetBTR(*opbtr, PDChannelToSOCMapDirect, tmpbtr);
474  }
475  // SimPhotons case
476  else {
477  sim::OnePhoton photon;
478  photon.SetInSD = false;
479  photon.InitialPosition = edepi.End();
480  photon.MotherTrackID = edepi.TrackID();
481  if (Reflected)
482  photon.Energy = 2.9 * CLHEP::eV; // 430 nm
483  else
484  photon.Energy = 9.7 * CLHEP::eV; // 128 nm
485  // TODO: un-hardcode and add another energy for Xe scintillation
486  if (ndetected_fast > 0 && fDoFastComponent) {
487  int n = ndetected_fast;
488  num_fastdp += n;
489  for (int i = 0; i < n; ++i) {
490  // calculates the time at which the photon was produced
491  double dtime = edepi.StartT() + fScintTime->fastScintTime();
492  if (fIncludePropTime) dtime += transport_time[i];
493  int time = static_cast<int>(std::round(dtime));
494  photon.Time = time;
495  if (Reflected)
496  ref_photcol[channel].insert(ref_photcol[channel].end(), 1, photon);
497  else
498  dir_photcol[channel].insert(dir_photcol[channel].end(), 1, photon);
499  }
500  }
501  if (ndetected_slow > 0 && fDoSlowComponent) {
502  int n = ndetected_slow;
503  num_slowdp += n;
504  for (int i = 0; i < n; ++i) {
505  double dtime = edepi.StartT() + fScintTime->slowScintTime();
506  if (fIncludePropTime) dtime += transport_time[ndetected_fast + i];
507  int time = static_cast<int>(std::round(dtime));
508  photon.Time = time;
509  if (Reflected)
510  ref_photcol[channel].insert(ref_photcol[channel].end(), 1, photon);
511  else
512  dir_photcol[channel].insert(dir_photcol[channel].end(), 1, photon);
513  }
514  }
515  }
516  }
517  }
518  }
519 
520  mf::LogTrace("PDFastSimPAR") << "Total points: " << num_points
521  << ", total fast photons: " << num_fastph
522  << ", total slow photons: " << num_slowph
523  << "\ndetected fast photons: " << num_fastdp
524  << ", detected slow photons: " << num_slowdp;
525 
526  if (fUseLitePhotons) {
527  event.put(move(phlit));
528  event.put(move(opbtr));
529  if (fDoReflectedLight) {
530  event.put(move(phlit_ref), "Reflected");
531  event.put(move(opbtr_ref), "Reflected");
532  }
533  }
534  else {
535  event.put(move(phot));
536  if (fDoReflectedLight) { event.put(move(phot_ref), "Reflected"); }
537  }
538 
539  return;
540  }
541 
542  //......................................................................
543  void PDFastSimPAR::AddOpDetBTR(std::vector<sim::OpDetBacktrackerRecord>& opbtr,
544  std::vector<int>& ChannelMap,
545  const sim::OpDetBacktrackerRecord& btr) const
546  {
547  int iChan = btr.OpDetNum();
548  if (ChannelMap[iChan] < 0) {
549  ChannelMap[iChan] = opbtr.size();
550  opbtr.emplace_back(std::move(btr));
551  }
552  else {
553  size_t idtest = ChannelMap[iChan];
554  auto const& timePDclockSDPsMap = btr.timePDclockSDPsMap();
555  for (auto const& timePDclockSDP : timePDclockSDPsMap) {
556  for (auto const& sdp : timePDclockSDP.second) {
557  double xyz[3] = {sdp.x, sdp.y, sdp.z};
558  opbtr.at(idtest).AddScintillationPhotons(
559  sdp.trackID, timePDclockSDP.first, sdp.numPhotons, xyz, sdp.energy);
560  }
561  }
562  }
563  }
564 
565  //......................................................................
566  // calculates number of photons detected given visibility and emitted number of photons
567  void PDFastSimPAR::detectedNumPhotons(std::vector<int>& DetectedNumPhotons,
568  const std::vector<double>& OpDetVisibilities,
569  const int NumPhotons) const
570  {
571  for (size_t i = 0; i < OpDetVisibilities.size(); ++i) {
572  DetectedNumPhotons[i] = fRandPoissPhot->fire(OpDetVisibilities[i] * NumPhotons);
573  }
574  }
575 
576  //......................................................................
577  // checks whether photo-detector is able to see the emitted light scintillation
579  geo::Point_t const& OpDetPoint) const
580  {
581  // check optical channel is in same TPC as scintillation light, if not doesn't see light
582  // temporary method working for SBND, uBooNE, DUNE 1x2x6; to be replaced to work in full DUNE geometry
583  // check x coordinate has same sign or is close to zero
584  if (((ScintPoint.X() < 0.) != (OpDetPoint.X() < 0.)) && std::abs(OpDetPoint.X()) > 10. &&
585  fNTPC == 2) { // TODO: replace with geometry service method
586  return false;
587  }
588  return true;
589  }
590 
591  std::vector<geo::Point_t> PDFastSimPAR::opDetCenters() const
592  {
593  std::vector<geo::Point_t> opDetCenter;
594  for (size_t const i : util::counter(fNOpChannels)) {
595  geo::OpDetGeo const& opDet = fGeom.OpDetGeoFromOpDet(i);
596  opDetCenter.push_back(opDet.GetCenter());
597  }
598  return opDetCenter;
599  }
600  // ---------------------------------------------------------------------------
601 
602 } // namespace phot
603 
const std::vector< geo::BoxBoundedGeo > fActiveVolumes
base_engine_t & createEngine(seed_t seed)
bool isOpDetInSameTPC(geo::Point_t const &ScintPoint, geo::Point_t const &OpDetPoint) const
unsigned int NTPC(CryostatID const &cryoid=cryostat_zero) const
Returns the total number of TPCs in the specified cryostat.
Definition: GeometryCore.h:686
void AddOpDetBTR(std::vector< sim::OpDetBacktrackerRecord > &opbtr, std::vector< int > &ChannelMap, const sim::OpDetBacktrackerRecord &btr) const
OpDetGeo const & OpDetGeoFromOpDet(unsigned int OpDet) const
Returns the geo::OpDetGeo object for the given detector number.
static std::vector< geo::BoxBoundedGeo > extractActiveLArVolume(geo::GeometryCore const &geom)
Definition: ISTPC.cxx:45
Utilities related to art service access.
Definition of util::enumerate().
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
std::unique_ptr< PropagationTimeModel > fPropTimeModel
fhicl::Atom< bool > DoSlowComponent
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.cc:6
PDFastSimPAR(Parameters const &config)
All information of a photon entering the sensitive optical detector volume.
Definition: SimPhotons.h:60
fhicl::Atom< bool > DoReflectedLight
constexpr auto abs(T v)
Returns the absolute value of the argument.
std::unique_ptr< ScintTime > fScintTime
MaybeLogger_< ELseverityLevel::ELsev_error, false > LogError
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Definition: GeometryCore.h:430
Energy deposited on a readout Optical Detector by simulated tracks.
geo::Point_t InitialPosition
Scintillation position in world coordinates [cm].
Definition: SimPhotons.h:63
auto enumerate(Iterables &&...iterables)
Range-for loop helper tracking the number of iteration.
Definition: enumerate.h:65
fhicl::Atom< bool > UseLitePhotons
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
std::unique_ptr< CLHEP::RandPoissonQ > fRandPoissPhot
CLHEP::HepRandomEngine & fScintTimeEngine
Simulation objects for optical detectors.
int OpDetNum() const
Returns the readout Optical Detector this object describes.
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:295
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
const std::vector< geo::Point_t > fOpDetCenter
Definitions of geometry vector data types.
std::unique_ptr< SemiAnalyticalModel > fVisibilityModel
fhicl::Atom< bool > OpaqueCathode
T get(std::string const &key) const
Definition: ParameterSet.h:314
fhicl::Atom< bool > IncludeAnodeReflections
void produce(art::Event &) override
Test of util::counter and support utilities.
Description of geometry of one entire detector.
Definition: GeometryCore.h:119
const art::InputTag fSimTag
An art service to assist in the distribution of guaranteed unique seeds to all engines within an art ...
unsigned int NOpDets() const
Number of OpDets in the whole detector.
Provides a base class aware of world box coordinates.
void AddScintillationPhotons(TrackID_t trackID, timePDclock_t timePDclock, double numberPhotons, double const *xyz, double energy)
Add scintillation photons and energy to this OpticalDetector.
Encapsulate the geometry of an optical detector.
geo::Point_t const & GetCenter() const
Definition: OpDetGeo.h:72
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:180
fhicl::Atom< bool > UseXeAbsorption
General LArSoft Utilities.
fhicl::Atom< bool > GeoPropTimeOnly
CLHEP::HepRandomEngine & fPhotonEngine
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
std::vector< geo::Point_t > opDetCenters() const
bool SetInSD
Whether the photon reaches the sensitive detector.
Definition: SimPhotons.h:84
contains information for a single step in the detector simulation
void detectedNumPhotons(std::vector< int > &DetectedNumPhotons, const std::vector< double > &OpDetVisibilities, const int NumPhotons) const
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
fhicl::Atom< art::InputTag > SimulationLabel
Definition: MVAAlg.h:12
Char_t n[5]
MaybeLogger_< ELseverityLevel::ELsev_success, true > LogTrace
fhicl::Atom< bool > IncludePropTime
float Energy
Scintillation photon energy [GeV].
Definition: SimPhotons.h:78
fhicl::Atom< bool > OnlyActiveVolume
bool isScintInActiveVolume(geo::Point_t const &ScintPoint)
Definition: ISTPC.cxx:36
const bool fIncludeAnodeReflections
timePDclockSDPs_t const & timePDclockSDPsMap() const
Returns all the deposited energy information as stored.
fhicl::Atom< bool > OnlyOneCryostat
art framework interface to geometry description
geo::GeometryCore const & fGeom
fhicl::Atom< bool > DoFastComponent
Event finding and building.
int MotherTrackID
ID of the GEANT4 track causing the scintillation.
Definition: SimPhotons.h:81