All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PhotonVisibilityService.cc
Go to the documentation of this file.
1 // ////////////////////////////////////////////////////////////////////////
2 //
3 // \file PhotonVisibilityService_service.cc
4 //
5 ////////////////////////////////////////////////////////////////////////
6 //
7 // Ben Jones, MIT 2012
8 //
9 // This service reports the visibility of a particular point in
10 // the detector to each OpDet. This is used by the fast
11 // optical simulation and by track-light association algorithms.
12 //
13 // Visibility is defined as the fraction of isotropically produced
14 // photons from a detector voxel which are expected to reach the
15 // OpDet in question.
16 //
17 // This information is lookup up from a previously generated
18 // optical library file, whose path is specified to this service.
19 //
20 // Note that it is important that the voxelization schemes match
21 // between the library and the service instance for sensible results.
22 //
23 //
24 
25 // LArSoft includes
27 
33 
34 // framework libraries
35 #include "art_root_io/TFileService.h"
36 #include "art/Framework/Services/Registry/ServiceHandle.h"
37 #include "art/Utilities/make_tool.h"
38 #include "canvas/Utilities/Exception.h"
39 #include "messagefacility/MessageLogger/MessageLogger.h"
40 #include "fhiclcpp/ParameterSet.h"
41 #include "cetlib/search_path.h"
42 
43 // ROOT libraries
44 #include "TF1.h"
45 
46 // C/C++ standard libraries
47 
48 namespace phot {
49 
51  {
52  delete fparslogNorm;
53  delete fparslogNorm_far;
54  delete fparsMPV;
55  delete fparsMPV_far;
56  delete fparsWidth;
57  delete fparsCte;
58  delete fparsCte_far;
59  delete fparsSlope;
60  delete fparslogNorm_refl;
61  delete fparsMPV_refl;
62  delete fparsWidth_refl;
63  delete fparsCte_refl;
64  delete fparsSlope_refl;
65  delete fTheLibrary;
66  }
67 
68  //--------------------------------------------------------------------
69  PhotonVisibilityService::PhotonVisibilityService(fhicl::ParameterSet const& pset)
70  :
71 
72  fCurrentVoxel(0)
73  , fCurrentValue(0.)
74  , fXmin(0.)
75  , fXmax(0.)
76  , fYmin(0.)
77  , fYmax(0.)
78  , fZmin(0.)
79  , fZmax(0.)
80  , fNx(0)
81  , fNy(0)
82  , fNz(0)
83  , fUseCryoBoundary(false)
84  , fLibraryBuildJob(false)
85  , fDoNotLoadLibrary(false)
86  , fParameterization(false)
87  , fHybrid(false)
88  , fStoreReflected(false)
89  , fStoreReflT0(false)
90  , fIncludePropTime(false)
91  , fUseNhitsModel(false)
92  , fParPropTime(false)
93  , fParPropTime_npar(0)
94  , fParPropTime_formula()
95  , fParPropTime_MaxRange()
96  , fInterpolate(false)
97  , fReflectOverZeroX(false)
98  , fparslogNorm(nullptr)
99  , fparslogNorm_far(nullptr)
100  , fparsMPV(nullptr)
101  , fparsMPV_far(nullptr)
102  , fparsWidth(nullptr)
103  , fparsCte(nullptr)
104  , fparsCte_far(nullptr)
105  , fparsSlope(nullptr)
106  , fD_break(0.0)
107  , fD_max(0.0)
108  , fTF1_sampling_factor(0.0)
109  , fparslogNorm_refl(nullptr)
110  , fparsMPV_refl(nullptr)
111  , fparsWidth_refl(nullptr)
112  , fparsCte_refl(nullptr)
113  , fparsSlope_refl(nullptr)
114  , fT0_max(0.0)
115  , fT0_break_point(0.0)
116  , fLibraryFile()
117  , fTheLibrary(nullptr)
118  , fVoxelDef()
119  {
120  this->reconfigure(pset);
121 
122  if (pset.has_key("ReflectOverZeroX")) { // legacy parameter warning
123  if (pset.has_key("Mapping")) {
124  throw art::Exception(art::errors::Configuration)
125  << "`PhotonVisbilityService` configuration specifies both `Mapping` and "
126  "`ReflectOverZeroX`."
127  " Please remove the latter (and use `PhotonMappingXMirrorTransformations` tool).";
128  }
129  else {
130  mf::LogWarning("PhotonVisbilityService")
131  << "Please update the configuration of `PhotonVisbilityService` service"
132  " replacing `ReflectOverZeroX` with tool configuration:"
133  "\n Mapping: { tool_type: \"PhotonMappingXMirrorTransformations\" }";
134  }
135  } // if
136  fhicl::ParameterSet mapDefaultSet;
137  mapDefaultSet.put("tool_type",
138  fReflectOverZeroX ? "PhotonMappingXMirrorTransformations" :
139  "PhotonMappingIdentityTransformations");
140  fMapping = art::make_tool<phot::IPhotonMappingTransformations>(
141  pset.get<fhicl::ParameterSet>("Mapping", mapDefaultSet));
142 
143  mf::LogInfo("PhotonVisibilityService") << "PhotonVisbilityService initializing" << std::endl;
144  }
145 
146  //--------------------------------------------------------------------
147  void
149  {
150  // Don't do anything if the library has already been loaded.
151 
152  if (fTheLibrary == 0) {
153 
154  if ((!fLibraryBuildJob) && (!fDoNotLoadLibrary)) {
155  std::string LibraryFileWithPath;
156  cet::search_path sp("FW_SEARCH_PATH");
157 
158  if (!sp.find_file(fLibraryFile, LibraryFileWithPath))
159  throw cet::exception("PhotonVisibilityService")
160  << "Unable to find photon library in " << sp.to_string() << "\n";
161 
162  if (!fParameterization) {
163  art::ServiceHandle<geo::Geometry const> geom;
164 
165  mf::LogInfo("PhotonVisibilityService")
166  << "PhotonVisibilityService Loading photon library from file " << LibraryFileWithPath
167  << " for " << GetVoxelDef().GetNVoxels() << " voxels and " << geom->NOpDets()
168  << " optical detectors." << std::endl;
169 
170  if (fHybrid) {
171  fTheLibrary = new PhotonLibraryHybrid(LibraryFileWithPath, GetVoxelDef());
172  }
173  else {
174  PhotonLibrary* lib = new PhotonLibrary;
175  fTheLibrary = lib;
176 
177  size_t NVoxels = GetVoxelDef().GetNVoxels();
178  lib->LoadLibraryFromFile(LibraryFileWithPath,
179  NVoxels,
181  fStoreReflT0,
184 
185  // if the library does not have metadata, we supply some;
186  // otherwise we check that it's compatible with the configured one
187  // (and shrug if it's not); overriding configured metadata
188  // from the one in the library is currently not supported
189  if (!lib->hasVoxelDef())
190  lib->SetVoxelDef(GetVoxelDef());
191  else if (GetVoxelDef() != lib->GetVoxelDef()) {
192  // this might become a fatal error in the future if some protocol
193  // is imposed... it may also be possible to check only the size
194  // rather than the coordinates, which may allow for translations
195  // of the geometry volumes in world space.
196  mf::LogWarning("PhotonVisbilityService")
197  << "Photon library reports the geometry:\n"
198  << lib->GetVoxelDef() << "while PhotonVisbilityService is configured with:\n"
199  << GetVoxelDef();
200  } // if metadata
201  }
202  }
203  }
204  else {
205  art::ServiceHandle<geo::Geometry const> geom;
206 
207  size_t NOpDets = geom->NOpDets();
208  size_t NVoxels = GetVoxelDef().GetNVoxels();
209  if (fLibraryBuildJob) {
210  mf::LogInfo("PhotonVisibilityService")
211  << " Vis service running library build job. Please ensure "
212  << " job contains LightSource, LArG4, SimPhotonCounter";
213  }
214 
215  art::TFileDirectory* pDir = nullptr;
216  try {
217  pDir = art::ServiceHandle<art::TFileService>().get();
218  }
219  catch (art::Exception const& e) {
220  if (e.categoryCode() != art::errors::ServiceNotFound) throw;
221  if (fLibraryBuildJob) {
222  throw art::Exception(e.categoryCode(), "", e)
223  << "PhotonVisibilityService: "
224  "service `TFileService` is required when building a photon library.\n";
225  }
226  }
227 
228  PhotonLibrary* lib = new PhotonLibrary(pDir);
229  fTheLibrary = lib;
230 
232  lib->SetVoxelDef(GetVoxelDef());
233  }
234  }
235  }
236 
237  //--------------------------------------------------------------------
238  void
240  {
241  if (fTheLibrary == 0) LoadLibrary();
242 
243  if (fLibraryBuildJob) {
244 
245  if (fHybrid) {
246  std::cout << "This is would be building a Hybrid Library. Not defined. " << std::endl;
247  }
248  mf::LogInfo("PhotonVisibilityService") << " Vis service "
249  << " Storing Library entries to file..." << std::endl;
250  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
252  }
253  }
254 
255  //--------------------------------------------------------------------
256  void
257  PhotonVisibilityService::reconfigure(fhicl::ParameterSet const& p)
258  {
259 
260  art::ServiceHandle<geo::Geometry const> geom;
261 
262  // Library details
263  fLibraryBuildJob = p.get<bool>("LibraryBuildJob", false);
264  fParameterization = p.get<bool>("DUNE10ktParameterization", false);
265  fHybrid = p.get<bool>("HybridLibrary", false);
266  fLibraryFile = p.get<std::string>("LibraryFile", "");
267  fDoNotLoadLibrary = p.get<bool>("DoNotLoadLibrary");
268  fStoreReflected = p.get<bool>("StoreReflected", false);
269  fStoreReflT0 = p.get<bool>("StoreReflT0", false);
270  // Parametrizations (time and Nhits)
271  fIncludePropTime = p.get<bool>("IncludePropTime", false);
272  fUseNhitsModel = p.get<bool>("UseNhitsModel", false);
273  fApplyVISBorderCorrection = p.get<bool>("ApplyVISBorderCorrection", false);
274  fVISBorderCorrectionType = p.get<std::string>("VIS_BORDER_correction_type", "");
275 
276  // Voxel parameters
277  fUseCryoBoundary = p.get<bool>("UseCryoBoundary", false);
278  fInterpolate = p.get<bool>("Interpolate", false);
279  fReflectOverZeroX = p.get<bool>("ReflectOverZeroX", false);
280 
281  fParPropTime = p.get<bool>("ParametrisedTimePropagation", false);
282  fParPropTime_npar = p.get<size_t>("ParametrisedTimePropagationNParameters", 0);
283  fParPropTime_formula = p.get<std::string>("ParametrisedTimePropagationFittedFormula", "");
284  fParPropTime_MaxRange = p.get<int>("ParametrisedTimePropagationMaxRange", 200);
285 
286  if (!fParPropTime) { fParPropTime_npar = 0; }
287 
288  if (!fUseNhitsModel) {
289 
290  if (fUseCryoBoundary) {
291  double CryoBounds[6];
292  geom->CryostatBoundaries(CryoBounds);
293  fXmin = CryoBounds[0];
294  fXmax = CryoBounds[1];
295  fYmin = CryoBounds[2];
296  fYmax = CryoBounds[3];
297  fZmin = CryoBounds[4];
298  fZmax = CryoBounds[5];
299  }
300  else {
301  fXmin = p.get<double>("XMin");
302  fXmax = p.get<double>("XMax");
303  fYmin = p.get<double>("YMin");
304  fYmax = p.get<double>("YMax");
305  fZmin = p.get<double>("ZMin");
306  fZmax = p.get<double>("ZMax");
307  }
308 
309  fNx = p.get<int>("NX");
310  fNy = p.get<int>("NY");
311  fNz = p.get<int>("NZ");
312 
314  }
315 
316  if (fIncludePropTime) {
317 
318  // load VUV arrival time distribution parametrization (no detector dependent at first order)
319  std::cout << "Loading the VUV time parametrization" << std::endl;
320  fDistances_landau = p.get<std::vector<double>>("Distances_landau");
321  fNorm_over_entries = p.get<std::vector<std::vector<double>>>("Norm_over_entries");
322  fMpv = p.get<std::vector<std::vector<double>>>("Mpv");
323  fWidth = p.get<std::vector<std::vector<double>>>("Width");
324  fDistances_exp = p.get<std::vector<double>>("Distances_exp");
325  fSlope = p.get<std::vector<std::vector<double>>>("Slope");
326  fExpo_over_Landau_norm = p.get<std::vector<std::vector<double>>>("Expo_over_Landau_norm");
327  fstep_size = p.get<double>("step_size");
328  fmax_d = p.get<double>("max_d");
329  fmin_d = p.get<double>("min_d");
330  fvuv_vgroup_mean = p.get<double>("vuv_vgroup_mean");
331  fvuv_vgroup_max = p.get<double>("vuv_vgroup_max");
332  finflexion_point_distance = p.get<double>("inflexion_point_distance");
333  fangle_bin_timing_vuv = p.get<double>("angle_bin_timing_vuv");
334 
335  if (fStoreReflected) {
336 
337  // load VIS arrival time distribution paramterisation
338  std::cout << "Loading the VIS time paramterisation" << std::endl;
339  fDistances_refl = p.get<std::vector<double>>("Distances_refl");
340  fDistances_radial_refl = p.get<std::vector<double>>("Distances_radial_refl");
341  fCut_off = p.get<std::vector<std::vector<std::vector<double>>>>("Cut_off");
342  fTau = p.get<std::vector<std::vector<std::vector<double>>>>("Tau");
343  fvis_vmean = p.get<double>("vis_vmean");
344  fangle_bin_timing_vis = p.get<double>("angle_bin_timing_vis");
345  }
346  }
347 
348  if (fUseNhitsModel) {
349  std::cout << "Loading semi-analytic mode models" << std::endl;
350  // VUV
351  fIsFlatPDCorr = p.get<bool>("FlatPDCorr", false);
352  fIsDomePDCorr = p.get<bool>("DomePDCorr", false);
353  fdelta_angulo_vuv = p.get<double>("delta_angulo_vuv");
354  if(fIsFlatPDCorr) {
355  fGHvuvpars_flat = p.get<std::vector<std::vector<double>>>("GH_PARS_flat");
356  fborder_corr_angulo_flat = p.get<std::vector<double>>("GH_border_angulo_flat");
357  fborder_corr_flat = p.get<std::vector<std::vector<double>>>("GH_border_flat");
358  }
359  if(fIsDomePDCorr) {
360  fGHvuvpars_dome = p.get<std::vector<std::vector<double>>>("GH_PARS_dome");
361  fborder_corr_angulo_dome = p.get<std::vector<double>>("GH_border_angulo_dome");
362  fborder_corr_dome = p.get<std::vector<std::vector<double>>>("GH_border_dome");
363  }
364 
365  if (fStoreReflected) {
366  fdelta_angulo_vis = p.get<double>("delta_angulo_vis");
367  if(fIsFlatPDCorr) {
368  fvis_distances_x_flat = p.get<std::vector<double>>("VIS_distances_x_flat");
369  fvis_distances_r_flat = p.get<std::vector<double>>("VIS_distances_r_flat");
370  fvispars_flat = p.get<std::vector<std::vector<std::vector<double>>>>("VIS_correction_flat");
371  }
372  if(fIsDomePDCorr) {
373  fvis_distances_x_dome = p.get<std::vector<double>>("VIS_distances_x_dome");
374  fvis_distances_r_dome = p.get<std::vector<double>>("VIS_distances_r_dome");
375  fvispars_dome = p.get<std::vector<std::vector<std::vector<double>>>>("VIS_correction_dome");
376  }
377  }
378  // optical detector information
379  fradius = p.get<double>("PMT_radius", 10.16);
380  }
381 
382  return;
383  }
384 
385  //------------------------------------------------------
386 
387  // Eventually we will calculate the light quenching factor here
388  double
390  {
391  // for now, no quenching
392  return 1.0;
393  }
394 
395  //------------------------------------------------------
396 
397  // Get a vector of the relative visibilities of each OpDet
398  // in the event to a point p
399 
400  auto
402  -> MappedCounts_t
403  {
405 
406  // first we fill a container of visibilities in the library index space
407  // (it is directly the values of the library unless interpolation is
408  // requested)
409  if (fInterpolate) {
410  // this is a punch into multithreading face:
411  static std::vector<float> ret;
412  ret.resize(fMapping->libraryMappingSize(p));
413  for (std::size_t libIndex = 0; libIndex < ret.size(); ++libIndex) {
414  ret[libIndex] = doGetVisibilityOfOpLib(p, LibraryIndex_t(libIndex), wantReflected);
415  }
416  data = &ret.front();
417  }
418  else {
419  auto const VoxID = VoxelAt(p);
420  data = GetLibraryEntries(VoxID, wantReflected);
421  }
422  return fMapping->applyOpDetMapping(p, data);
423  }
424 
425  //------------------------------------------------------
426 
427  // Get distance to optical detector OpDet
428  double
430  {
431  art::ServiceHandle<geo::Geometry const> geom;
432  return geom->OpDetGeoFromOpDet(OpDet).DistanceToPoint(p);
433  }
434 
435  //------------------------------------------------------
436 
437  // Get the solid angle reduction factor for planar optical detector OpDet
438  double
440  {
441  art::ServiceHandle<geo::Geometry const> geom;
442  return geom->OpDetGeoFromOpDet(OpDet).CosThetaFromNormal(p);
443  }
444 
445  //------------------------------------------------------
446 
447  float
449  LibraryIndex_t libIndex,
450  bool wantReflected /* = false */) const
451  {
452  if (!fInterpolate) { return GetLibraryEntry(VoxelAt(p), libIndex, wantReflected); }
453 
454  // In case we're outside the bounding box we'll get a empty optional list.
455  auto const neis = GetVoxelDef().GetNeighboringVoxelIDs(LibLocation(p));
456  if (!neis) return 0.0;
457 
458  // Sum up all the weighted neighbours to get interpolation behaviour
459  float vis = 0.0;
460  for (const sim::PhotonVoxelDef::NeiInfo& n : neis.value()) {
461  if (n.id < 0) continue;
462  vis += n.weight * GetLibraryEntry(n.id, libIndex, wantReflected);
463  }
464  return vis;
465  }
466 
467  //------------------------------------------------------
468 
469  bool
471  bool wantReflected /* = false */) const
472  {
473  return HasLibraryEntries(VoxelAt(p), wantReflected);
474  }
475 
476  //------------------------------------------------------
477 
478  float
480  unsigned int OpChannel,
481  bool wantReflected) const
482  {
483  // here we quietly confuse op. det. channel (interface) and op. det. (library)
484  LibraryIndex_t const libIndex = fMapping->opDetToLibraryIndex(p, OpChannel);
485  return doGetVisibilityOfOpLib(p, libIndex, wantReflected);
486  }
487 
488  //------------------------------------------------------
489 
490  void
492  {
493  fCurrentVoxel = VoxID;
494  fCurrentValue = N;
495  mf::LogInfo("PhotonVisibilityService")
496  << " PVS notes production of " << N << " photons at Vox " << VoxID << std::endl;
497  }
498 
499  //------------------------------------------------------
500 
501  void
503  {
504  N = fCurrentValue;
505  VoxID = fCurrentVoxel;
506  }
507 
508  //------------------------------------------------------
509 
510  void
511  PhotonVisibilityService::SetLibraryEntry(int VoxID, int OpChannel, float N, bool wantReflected)
512  {
513  if (fTheLibrary == 0) LoadLibrary();
514 
515  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
516 
517  if (!wantReflected)
518  lib->SetCount(VoxID, OpChannel, N);
519 
520  else
521  lib->SetReflCount(VoxID, OpChannel, N);
522 
523  //std::cout<< " PVS logging " << VoxID << " " << OpChannel<<std::endl;
524  MF_LOG_DEBUG("PhotonVisibilityService")
525  << " PVS logging " << VoxID << " " << OpChannel << std::endl;
526  }
527 
528  //------------------------------------------------------
529 
531  PhotonVisibilityService::GetLibraryEntries(int VoxID, bool wantReflected) const
532  {
533  if (fTheLibrary == 0) LoadLibrary();
534 
535  if (!wantReflected)
536  return fTheLibrary->GetCounts(VoxID);
537  else
538  return fTheLibrary->GetReflCounts(VoxID);
539  }
540 
541  //------------------------------------------------------
542 
543  bool
545  bool /* wantReflected */ /* = false */) const
546  {
547  if (!fTheLibrary) LoadLibrary();
548  return fTheLibrary->isVoxelValid(VoxID);
549  }
550 
551  //------------------------------------------------------
552 
553  float
555  OpDetID_t libOpChannel,
556  bool wantReflected) const
557  {
558  if (fTheLibrary == 0) LoadLibrary();
559 
560  if (!wantReflected)
561  return fTheLibrary->GetCount(VoxID, libOpChannel);
562  else
563  return fTheLibrary->GetReflCount(VoxID, libOpChannel);
564  }
565  // Methodes to handle the extra library parameter refl T0
566  //------------------------------------------------------
567 
568  // Get a vector of the refl <tfirst> of each OpDet
569  // in the event to a point p
570 
571  auto
573  {
574  // both the input and the output go through mapping to apply needed symmetries.
575  int const VoxID = VoxelAt(p);
576  phot::IPhotonLibrary::Counts_t const& data = GetLibraryReflT0Entries(VoxID);
577  return fMapping->applyOpDetMapping(p, data);
578  }
579 
580  //------------------------------------------------------
581 
584  {
585  if (fTheLibrary == 0) LoadLibrary();
586 
587  return fTheLibrary->GetReflT0s(VoxID);
588  }
589 
590  //------------------------------------------------------
591 
592  void
593  PhotonVisibilityService::SetLibraryReflT0Entry(int VoxID, int OpChannel, float T0)
594  {
595  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
596  if (fTheLibrary == 0) LoadLibrary();
597 
598  lib->SetReflT0(VoxID, OpChannel, T0);
599 
600  MF_LOG_DEBUG("PhotonVisibilityService")
601  << " PVS logging " << VoxID << " " << OpChannel << std::endl;
602  }
603 
604  //------------------------------------------------------
605 
606  float
608  {
609  if (fTheLibrary == 0) LoadLibrary();
610 
611  return fTheLibrary->GetReflT0(VoxID, libOpChannel);
612  }
613 
614  //------------------------------------------------------
615 
616  /////////////****////////////
617 
618  auto
620  {
621  int const VoxID = VoxelAt(p);
622  phot::IPhotonLibrary::Params_t const& params = GetLibraryTimingParEntries(VoxID);
623  return fMapping->applyOpDetMapping(p, params);
624  }
625 
626  auto
628  {
629  int const VoxID = VoxelAt(p);
630  phot::IPhotonLibrary::Functions_t const& functions = GetLibraryTimingTF1Entries(VoxID);
631  return fMapping->applyOpDetMapping(p, functions);
632  }
633 
634  //------------------------------------------------------
635 
638  {
639  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
640  if (fTheLibrary == 0) LoadLibrary();
641 
642  return lib->GetTimingPars(VoxID);
643  }
644 
645  //------------------------------------------------------
646 
649  {
650  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
651  if (fTheLibrary == 0) LoadLibrary();
652 
653  return lib->GetTimingTF1s(VoxID);
654  }
655 
656  //------------------------------------------------------
657 
658  void
660  int OpChannel,
661  float par,
662  size_t parnum)
663  {
664  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
665  if (fTheLibrary == 0) LoadLibrary();
666 
667  lib->SetTimingPar(VoxID, OpChannel, par, parnum);
668 
669  MF_LOG_DEBUG("PhotonVisibilityService")
670  << " PVS logging " << VoxID << " " << OpChannel << std::endl;
671  }
672 
673  //------------------------------------------------------
674 
675  void
676  PhotonVisibilityService::SetLibraryTimingTF1Entry(int VoxID, int OpChannel, TF1 const& func)
677  {
678  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
679  if (fTheLibrary == 0) LoadLibrary();
680 
681  lib->SetTimingTF1(VoxID, OpChannel, func);
682 
683  MF_LOG_DEBUG("PhotonVisibilityService")
684  << " PVS logging " << VoxID << " " << OpChannel << std::endl;
685  }
686 
687  //------------------------------------------------------
688 
689  float
691  OpDetID_t libOpChannel,
692  size_t npar) const
693  {
694  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
695  if (fTheLibrary == 0) LoadLibrary();
696 
697  return lib->GetTimingPar(VoxID, libOpChannel, npar);
698  }
699 
700  //------------------------------------------------------
701 
702  size_t
704  {
705  // the last word about the number of channels belongs to the mapping;
706  // this should be also the same answer as `geo::GeometryCore::NOpDets()`.
707  return fMapping->opDetMappingSize();
708  }
709 
710  //------------------------------------------------------
711  void
713  double& d_break,
714  double& d_max,
715  double& tf1_sampling_factor) const
716  {
717  functions[0] = fparslogNorm;
718  functions[1] = fparsMPV;
719  functions[2] = fparsWidth;
720  functions[3] = fparsCte;
721  functions[4] = fparsSlope;
722  functions[5] = fparslogNorm_far;
723  functions[6] = fparsMPV_far;
724  functions[7] = fparsCte_far;
725 
726  d_break = fD_break;
727  d_max = fD_max;
728  tf1_sampling_factor = fTF1_sampling_factor;
729  }
730 
731  //------------------------------------------------------
732  void
734  double& t0_max,
735  double& t0_break_point) const
736  {
737  functions[0] = fparslogNorm_refl;
738  functions[1] = fparsMPV_refl;
739  functions[2] = fparsWidth_refl;
740  functions[3] = fparsCte_refl;
741  functions[4] = fparsSlope_refl;
742 
743  t0_max = fT0_max;
744  t0_break_point = fT0_break_point;
745  }
746 
747  //------------------------------------------------------
748  void
750  double& step_size,
751  double& max_d,
752  double& min_d,
753  double& vuv_vgroup_mean,
754  double& vuv_vgroup_max,
755  double& inflexion_point_distance,
756  double& angle_bin_timing_vuv) const
757  {
758  v[0] = std::vector(1, fDistances_landau);
759  v[1] = fNorm_over_entries;
760  v[2] = fMpv;
761  v[3] = fWidth;
762  v[4] = std::vector(1, fDistances_exp);
763  v[5] = fSlope;
764  v[6] = fExpo_over_Landau_norm;
765 
766  step_size = fstep_size;
767  max_d = fmax_d;
768  min_d = fmin_d;
769  vuv_vgroup_mean = fvuv_vgroup_mean;
770  vuv_vgroup_max = fvuv_vgroup_max;
771  inflexion_point_distance = finflexion_point_distance;
772  angle_bin_timing_vuv = fangle_bin_timing_vuv;
773  }
774 
775  void
776  PhotonVisibilityService::LoadTimingsForVISPar(std::vector<double>& distances,
777  std::vector<double>& radial_distances,
778  std::vector<std::vector<std::vector<double>>>& cut_off,
779  std::vector<std::vector<std::vector<double>>>& tau,
780  double& vis_vmean,
781  double& angle_bin_timing_vis) const
782  {
783  distances = fDistances_refl;
784  radial_distances = fDistances_radial_refl;
785  cut_off = fCut_off;
786  tau = fTau;
787 
788  vis_vmean = fvis_vmean;
789  angle_bin_timing_vis = fangle_bin_timing_vis;
790  }
791 
793  bool &isDomePDCorr,
794  double &delta_angulo_vuv,
795  double &radius) const
796  {
797  isFlatPDCorr = fIsFlatPDCorr;
798  isDomePDCorr = fIsDomePDCorr;
799  delta_angulo_vuv = fdelta_angulo_vuv;
800  radius = fradius;
801  }
802  void PhotonVisibilityService::LoadGHFlat( std::vector<std::vector<double>> &GHvuvpars_flat,
803  std::vector<double> &border_corr_angulo_flat,
804  std::vector<std::vector<double>> &border_corr_flat) const
805  {
806  if (!fIsFlatPDCorr) return;
807  GHvuvpars_flat = fGHvuvpars_flat;
808  border_corr_angulo_flat = fborder_corr_angulo_flat;
809  border_corr_flat = fborder_corr_flat;
810  }
811  void PhotonVisibilityService::LoadGHDome( std::vector<std::vector<double>> &GHvuvpars_dome,
812  std::vector<double> &border_corr_angulo_dome,
813  std::vector<std::vector<double>> &border_corr_dome) const
814  {
815  if (!fIsDomePDCorr) return;
816  GHvuvpars_dome = fGHvuvpars_dome;
817  border_corr_angulo_dome = fborder_corr_angulo_dome;
818  border_corr_dome = fborder_corr_dome;
819  }
821  double &radius) const
822  {
823  delta_angulo_vis = fdelta_angulo_vis;
824  radius = fradius;
825  }
826  void PhotonVisibilityService::LoadVisParsFlat(std::vector<double> &vis_distances_x_flat,
827  std::vector<double> &vis_distances_r_flat,
828  std::vector<std::vector<std::vector<double>>> &vispars_flat) const
829  {
830  if (!fIsFlatPDCorr) return;
831  vis_distances_x_flat = fvis_distances_x_flat;
832  vis_distances_r_flat = fvis_distances_r_flat;
833  vispars_flat = fvispars_flat;
834  }
835  void PhotonVisibilityService::LoadVisParsDome(std::vector<double> &vis_distances_x_dome,
836  std::vector<double> &vis_distances_r_dome,
837  std::vector<std::vector<std::vector<double>>> &vispars_dome) const
838  {
839  if (!fIsDomePDCorr) return;
840  vis_distances_x_dome = fvis_distances_x_dome;
841  vis_distances_r_dome = fvis_distances_r_dome;
842  vispars_dome = fvispars_dome;
843  }
844 
845  //------------------------------------------------------
846  /***
847  * Preform any necessary transformations on the coordinates before trying to access
848  * a voxel ID.
849  **/
852  {
853  return fMapping->detectorToLibrary(p);
854  }
855 
856 } // namespace
std::vector< std::vector< std::vector< double > > > fTau
MappedParams_t doGetTimingPar(geo::Point_t const &p) const
void LoadTimingsForVISPar(std::vector< double > &distances, std::vector< double > &radial_distances, std::vector< std::vector< std::vector< double >>> &cut_off, std::vector< std::vector< std::vector< double >>> &tau, double &vis_vmean, double &angle_bin_timing_vis) const
phot::IPhotonLibrary::Counts_t GetLibraryReflT0Entries(int VoxID) const
createEngine fYmax
MappedCounts_t doGetAllVisibilities(geo::Point_t const &p, bool wantReflected=false) const
phot::IPhotonLibrary::Counts_t GetLibraryEntries(int VoxID, bool wantReflected=false) const
std::unique_ptr< phot::IPhotonMappingTransformations > fMapping
Mapping of detector space into library space.
void RetrieveLightProd(int &VoxID, double &N) const
void SetDirectLightPropFunctions(TF1 const *functions[8], double &d_break, double &d_max, double &tf1_sampling_factor) const
void LoadLibraryFromFile(std::string LibraryFile, size_t NVoxels, bool storeReflected=false, bool storeReflT0=false, size_t storeTiming=0, int maxrange=200)
virtual float GetReflT0(size_t Voxel, size_t OpChannel) const =0
void LoadTimingsForVUVPar(std::vector< std::vector< double >>(&v)[7], double &step_size, double &max_d, double &min_d, double &vuv_vgroup_mean, double &vuv_vgroup_max, double &inflexion_point_distance, double &angle_bin_timing_vuv) const
void CreateEmptyLibrary(size_t NVoxels, size_t NChannels, bool storeReflected=false, bool storeReflT0=false, size_t storeTiming=0)
void SetLibraryTimingParEntry(int VoxID, int OpChannel, float value, size_t parnum)
sim::PhotonVoxelDef const & GetVoxelDef() const
bool hasVoxelDef() const
Returns whether voxel metadata is available.
virtual float GetReflCount(size_t Voxel, size_t OpChannel) const =0
pdgs p
Definition: selectors.fcl:22
phot::IPhotonMappingTransformations::LibraryIndex_t LibraryIndex_t
Type of optical library index.
Representation of a region of space diced into voxels.
float GetLibraryTimingParEntry(int VoxID, OpDetID_t libOpChannel, size_t npar) const
bool HasLibraryEntries(int VoxID, bool wantReflected=false) const
void LoadVisParsFlat(std::vector< double > &vis_distances_x_flat, std::vector< double > &vis_distances_r_flat, std::vector< std::vector< std::vector< double >>> &vispars_flat) const
float doGetVisibilityOfOpLib(geo::Point_t const &p, LibraryIndex_t libIndex, bool wantReflected=false) const
virtual float GetCount(size_t Voxel, size_t OpChannel) const =0
bool doHasVisibility(geo::Point_t const &p, bool wantReflected=false) const
void StoreLibraryToFile(std::string LibraryFile, bool storeReflected=false, bool storeReflT0=false, size_t storeTiming=0) const
std::vector< std::vector< std::vector< double > > > fCut_off
float GetLibraryEntry(int VoxID, OpDetID_t libOpChannel, bool wantReflected=false) const
void LoadGHDome(std::vector< std::vector< double >> &GHvuvpars_dome, std::vector< double > &border_corr_angulo_dome, std::vector< std::vector< double >> &border_corr_dome) const
void SetCount(size_t Voxel, size_t OpChannel, float Count)
void SetTimingTF1(size_t Voxel, size_t OpChannel, TF1 func)
void LoadVUVSemiAnalyticProperties(bool &isFlatPDCorr, bool &isDomePDCorr, double &delta_angulo_vuv, double &radius) const
double GetQuenchingFactor(double dQdx) const
unsigned int GetNVoxels() const
Returns the total number of voxels in the volume.
float GetLibraryReflT0Entry(int VoxID, OpDetID_t libOpChannel) const
void SetReflectedCOLightPropFunctions(TF1 const *functions[5], double &t0_max, double &t0_break_point) const
void StoreLightProd(int VoxID, double N)
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
void SetVoxelDef(sim::PhotonVoxelDef const &voxelDef)
PhotonVisibilityService(fhicl::ParameterSet const &pset)
Definitions of voxel data structures.
createEngine fXmin
void LoadGHFlat(std::vector< std::vector< double >> &GHvuvpars_flat, std::vector< double > &border_corr_angulo_flat, std::vector< std::vector< double >> &border_corr_flat) const
MappedFunctions_t doGetTimingTF1(geo::Point_t const &p) const
geo::Point_t LibLocation(geo::Point_t const &p) const
void SetReflT0(size_t Voxel, size_t OpChannel, float reflT0)
void LoadVisParsDome(std::vector< double > &vis_distances_x_dome, std::vector< double > &vis_distances_r_dome, std::vector< std::vector< std::vector< double >>> &vispars_dome) const
createEngine fXmax
void SetLibraryReflT0Entry(int VoxID, int OpChannel, float value)
virtual T0s_t GetReflT0s(size_t Voxel) const =0
void reconfigure(fhicl::ParameterSet const &p)
const std::vector< float > * GetTimingPars(size_t Voxel) const
phot::IPhotonMappingTransformations::OpDetID_t OpDetID_t
Type of (global) optical detector ID.
createEngine fYmin
std::vector< float > const * Params_t
virtual bool isVoxelValid(size_t Voxel) const
Encapsulate the geometry of an optical detector.
void SetTimingPar(size_t Voxel, size_t OpChannel, float Count, size_t parnum)
virtual Counts_t GetCounts(size_t Voxel) const =0
Returns a pointer to NOpChannels() visibility values, one per channel.
void LoadVisSemiAnalyticProperties(double &delta_angulo_vis, double &radius) const
std::optional< std::array< NeiInfo, 8U > > GetNeighboringVoxelIDs(Point const &v) const
Returns IDs of the eight neighboring voxels around v.
void SetLibraryEntry(int VoxID, OpDetID_t libOpChannel, float N, bool wantReflected=false)
A container for photon visibility mapping data.
process_name largeant stream1 can override from command line with o or output physics producers generator N
createEngine fZmin
float GetTimingPar(size_t Voxel, size_t OpChannel, size_t parnum) const
float doGetVisibility(geo::Point_t const &p, unsigned int OpChannel, bool wantReflected=false) const
do i e
void SetLibraryTimingTF1Entry(int VoxID, int OpChannel, TF1 const &func)
std::vector< std::vector< std::vector< double > > > fvispars_dome
const float * Counts_t
Type for visibility count per optical channel.
static double DistanceToOpDetImpl(geo::Point_t const &p, unsigned int OpDet)
MappedT0s_t doGetReflT0s(geo::Point_t const &p) const
createEngine fZmax
static double SolidAngleFactorImpl(geo::Point_t const &p, unsigned int OpDet)
void SetReflCount(size_t Voxel, size_t OpChannel, float Count)
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:184
phot::IPhotonLibrary::Functions_t GetLibraryTimingTF1Entries(int VoxID) const
phot::IPhotonLibrary::Params_t GetLibraryTimingParEntries(int VoxID) const
std::vector< std::vector< std::vector< double > > > fvispars_flat
virtual Counts_t GetReflCounts(size_t Voxel) const =0
art framework interface to geometry description
BEGIN_PROLOG could also be cout
process_name opdaq physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator physics producers generator T0
Definition: gen_protons.fcl:45