LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
PhotonVisibilityService_service.cc
Go to the documentation of this file.
1 // ////////////////////////////////////////////////////////////////////////
2 //
3 // \file PhotonVisibilityService_service.cc
4 //
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 // Framework includes
25 
26 // LArSoft includes
34 
37 
38 #include "TF1.h"
39 
40 #include <math.h>
41 
42 namespace phot{
43 
45  {
46  delete fparslogNorm;
47  delete fparslogNorm_far;
48  delete fparsMPV;
49  delete fparsMPV_far;
50  delete fparsWidth;
51  delete fparsCte;
52  delete fparsCte_far;
53  delete fparsSlope;
54  delete fparslogNorm_refl;
55  delete fparsMPV_refl;
56  delete fparsWidth_refl;
57  delete fparsCte_refl;
58  delete fparsSlope_refl;
59  delete fTheLibrary;
60  }
61 
62  //--------------------------------------------------------------------
64 
65  fCurrentVoxel(0),
66  fCurrentValue(0.),
67  fXmin(0.),
68  fXmax(0.),
69  fYmin(0.),
70  fYmax(0.),
71  fZmin(0.),
72  fZmax(0.),
73  fNx(0),
74  fNy(0),
75  fNz(0),
76  fUseCryoBoundary(false),
77  fLibraryBuildJob(false),
78  fDoNotLoadLibrary(false),
79  fParameterization(false),
80  fHybrid(false),
81  fStoreReflected(false),
82  fStoreReflT0(false),
83  fIncludePropTime(false),
84  fParPropTime(false),
88  fInterpolate(false),
89  fReflectOverZeroX(false),
90  fparslogNorm(nullptr),
91  fparslogNorm_far(nullptr),
92  fparsMPV(nullptr),
93  fparsMPV_far(nullptr),
94  fparsWidth(nullptr),
95  fparsCte(nullptr),
96  fparsCte_far(nullptr),
97  fparsSlope(nullptr),
98  fD_break(0.0),
99  fD_max(0.0),
101  fparslogNorm_refl(nullptr),
102  fparsMPV_refl(nullptr),
103  fparsWidth_refl(nullptr),
104  fparsCte_refl(nullptr),
105  fparsSlope_refl(nullptr),
106  fT0_max(0.0),
107  fT0_break_point(0.0),
108  fLibraryFile(),
109  fTheLibrary(nullptr),
110  fVoxelDef()
111  {
112  this->reconfigure(pset);
113  mf::LogInfo("PhotonVisibilityService")<<"PhotonVisbilityService initializing"<<std::endl;
114  }
115 
116  //--------------------------------------------------------------------
118  {
119  // Don't do anything if the library has already been loaded.
120 
121  if(fTheLibrary == 0) {
122 
124  std::string LibraryFileWithPath;
125  cet::search_path sp("FW_SEARCH_PATH");
126 
127  if( !sp.find_file(fLibraryFile, LibraryFileWithPath) )
128  throw cet::exception("PhotonVisibilityService") << "Unable to find photon library in " << sp.to_string() << "\n";
129 
130  if(!fParameterization) {
132 
133  mf::LogInfo("PhotonVisibilityService") << "PhotonVisibilityService Loading photon library from file "
134  << LibraryFileWithPath
135  << " for "
136  << GetVoxelDef().GetNVoxels()
137  << " voxels and "
138  << geom->NOpDets()
139  << " optical detectors."
140  << std::endl;
141 
142  if(fHybrid){
143  fTheLibrary = new PhotonLibraryHybrid(LibraryFileWithPath,
144  GetVoxelDef());
145  }
146  else{
147  PhotonLibrary* lib = new PhotonLibrary;
148  fTheLibrary = lib;
149 
150  size_t NVoxels = GetVoxelDef().GetNVoxels();
152  }
153  }
154  }
155  else {
157 
158  size_t NOpDets = geom->NOpDets();
159  size_t NVoxels = GetVoxelDef().GetNVoxels();
160  mf::LogInfo("PhotonVisibilityService") << " Vis service running library build job. Please ensure "
161  << " job contains LightSource, LArG4, SimPhotonCounter"<<std::endl;
162  PhotonLibrary* lib = new PhotonLibrary;
163  fTheLibrary = lib;
164 
166  }
167 
168  }
169  }
170 
171  //--------------------------------------------------------------------
173  {
174  if(fTheLibrary == 0)
175  LoadLibrary();
176 
177  if(fLibraryBuildJob )
178  {
179 
180  if(fHybrid){
181  std::cout<< "This is would be building a Hybrid Library. Not defined. "<<std::endl;
182  }
183  mf::LogInfo("PhotonVisibilityService") << " Vis service "
184  << " Storing Library entries to file..." <<std::endl;
185  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
187  }
188  }
189 
190 
191  //--------------------------------------------------------------------
193  {
194 
196 
197  // Library details
198  fLibraryBuildJob = p.get< bool >("LibraryBuildJob" );
199  fParameterization = p.get< bool >("DUNE10ktParameterization", false);
200  fHybrid = p.get< bool >("HybridLibrary", false);
201  fLibraryFile = p.get< std::string >("LibraryFile" );
202  fDoNotLoadLibrary = p.get< bool >("DoNotLoadLibrary" );
203  fStoreReflected = p.get< bool >("StoreReflected", false);
204  fStoreReflT0 = p.get< bool >("StoreReflT0", false);
205  fIncludePropTime = p.get< bool >("IncludePropTime", false);
206  // Voxel parameters
207  fUseCryoBoundary = p.get< bool >("UseCryoBoundary" );
208  fInterpolate = p.get< bool >("Interpolate", false);
209  fReflectOverZeroX = p.get< bool >("ReflectOverZeroX", false);
210 
211  fParPropTime = p.get< bool >("ParametrisedTimePropagation", false);
212  fParPropTime_npar = p.get< size_t >("ParametrisedTimePropagationNParameters", 0);
213  fParPropTime_formula = p.get< std::string >("ParametrisedTimePropagationFittedFormula","");
214  fParPropTime_MaxRange = p.get< int >("ParametrisedTimePropagationMaxRange", 200);
215 
216  if (!fParPropTime)
217  {
218  fParPropTime_npar=0;
219  }
220 
221  if(fUseCryoBoundary)
222  {
223  double CryoBounds[6];
224  geom->CryostatBoundaries(CryoBounds);
225  fXmin = CryoBounds[0];
226  fXmax = CryoBounds[1];
227  fYmin = CryoBounds[2];
228  fYmax = CryoBounds[3];
229  fZmin = CryoBounds[4];
230  fZmax = CryoBounds[5];
231  }
232  else
233  {
234  fXmin = p.get< double >("XMin" );
235  fXmax = p.get< double >("XMax" );
236  fYmin = p.get< double >("YMin" );
237  fYmax = p.get< double >("YMax" );
238  fZmin = p.get< double >("ZMin" );
239  fZmax = p.get< double >("ZMax" );
240  }
241 
242  fNx = p.get< int >("NX" );
243  fNy = p.get< int >("NY" );
244  fNz = p.get< int >("NZ" );
245 
247 
248  if(fIncludePropTime)
249  {
250  // Construct parameterized model parameter functions.
251  std::cout<< "Getting direct light parameters from .fcl file"<<std::endl;
252  std::vector<std::string> direct_functions = p.get<std::vector<std::string> >("Direct_functions");
253  //range of distances where the parametrization is valid
254  fD_break = p.get<double>("D_break");
255  fD_max = p.get<double>("D_max");
256 
257  fTF1_sampling_factor = p.get<double>("TF1_sampling_factor");
258 
259  std::vector<double> direct_landauNormpars = p.get<std::vector<double> >("Direct_landauNormpars");
260  fparslogNorm = new TF1("fparslogNorm", direct_functions[0].c_str(), 0., fD_break);
261  for(unsigned int i=0; i<direct_landauNormpars.size(); ++i)
262  fparslogNorm->SetParameter(i, direct_landauNormpars[i]);
263 
264  std::vector<double> direct_landauMPVpars = p.get<std::vector<double> >("Direct_landauMPVpars");
265  fparsMPV = new TF1("fparsMPV", direct_functions[1].c_str(), 0., fD_break);
266  for(unsigned int i=0; i<direct_landauMPVpars.size(); ++i)
267  fparsMPV->SetParameter(i, direct_landauMPVpars[i]);
268 
269  std::vector<double> direct_landauWidthpars = p.get<std::vector<double> >("Direct_landauWidthpars");
270  fparsWidth = new TF1("fparsWidth", direct_functions[2].c_str(), 0., fD_break);
271  for(unsigned int i=0; i<direct_landauWidthpars.size(); ++i)
272  fparsWidth->SetParameter(i, direct_landauWidthpars[i]);
273 
274  std::vector<double> direct_expoCtepars = p.get<std::vector<double> >("Direct_expoCtepars");
275  fparsCte = new TF1("fparsCte", direct_functions[3].c_str(), 0., fD_break);
276  for(unsigned int i=0; i<direct_expoCtepars.size(); ++i)
277  fparsCte->SetParameter(i, direct_expoCtepars[i]);
278 
279  std::vector<double> direct_expoSlopepars = p.get<std::vector<double> >("Direct_expoSlopepars");
280  fparsSlope = new TF1("fparsSlope", direct_functions[4].c_str(), 0., fD_break);
281  for(unsigned int i=0; i<direct_expoSlopepars.size(); ++i)
282  fparsSlope->SetParameter(i, direct_expoSlopepars[i]);
283 
284  std::vector<double> direct_landauNormpars_far = p.get<std::vector<double> >("Direct_landauNormpars_far");
285  fparslogNorm_far = new TF1("fparslogNorm_far", direct_functions[5].c_str(), fD_break, fD_max);
286  for(unsigned int i=0; i<direct_landauNormpars_far.size(); ++i)
287  fparslogNorm_far->SetParameter(i, direct_landauNormpars_far[i]);
288 
289  std::vector<double> direct_landauMPVpars_far = p.get<std::vector<double> >("Direct_landauMPVpars_far");
290  fparsMPV_far = new TF1("fparsMPV_far", direct_functions[6].c_str(), fD_break, fD_max);
291  for(unsigned int i=0; i<direct_landauMPVpars_far.size(); ++i)
292  fparsMPV_far->SetParameter(i, direct_landauMPVpars_far[i]);
293 
294  std::vector<double> direct_expoCtepars_far = p.get<std::vector<double> >("Direct_expoCtepars_far");
295  fparsCte_far = new TF1("fparsCte_far", direct_functions[7].c_str(), fD_break - 50., fD_max);
296  for(unsigned int i=0; i<direct_expoCtepars_far.size(); ++i)
297  fparsCte_far->SetParameter(i, direct_expoCtepars_far[i]);
298 
299  std::vector<std::string> reflected_functions = p.get<std::vector<std::string> >("Reflected_functions");
300  //times where the parametrizations are valid or change
301  fT0_max = p.get<double>("T0_max");
302  fT0_break_point = p.get<double>("T0_break_point");
303 
304  std::vector<double> reflected_landauNormpars = p.get<std::vector<double> >("Reflected_landauNormpars");
305  fparslogNorm_refl = new TF1("fparslogNorm_refl", reflected_functions[0].c_str(), 0., fT0_max);
306  for(unsigned int i=0; i<reflected_landauNormpars.size(); ++i)
307  fparslogNorm_refl->SetParameter(i, reflected_landauNormpars[i]);
308 
309  std::vector<double> reflected_landauMPVpars = p.get<std::vector<double> >("Reflected_landauMPVpars");
310  fparsMPV_refl = new TF1("fparsMPV_refl", reflected_functions[1].c_str(), 0., fT0_max);
311  for(unsigned int i=0; i<reflected_landauMPVpars.size(); ++i)
312  fparsMPV_refl->SetParameter(i, reflected_landauMPVpars[i]);
313 
314  std::vector<double> reflected_landauWidthpars = p.get<std::vector<double> >("Reflected_landauWidthpars");
315  fparsWidth_refl = new TF1("fparsWidth_refl", reflected_functions[2].c_str(), 0., fT0_max);
316  for(unsigned int i=0; i<reflected_landauWidthpars.size(); ++i)
317  fparsWidth_refl->SetParameter(i, reflected_landauWidthpars[i]);
318 
319  std::vector<double> reflected_expoCtepars = p.get<std::vector<double> >("Reflected_expoCtepars");
320  fparsCte_refl = new TF1("fparsCte_refl", reflected_functions[3].c_str(), 0., fT0_max);
321  for(unsigned int i=0; i<reflected_expoCtepars.size(); ++i)
322  fparsCte_refl->SetParameter(i, reflected_expoCtepars[i]);
323 
324  std::vector<double> reflected_expoSlopepars = p.get<std::vector<double> >("Reflected_expoSlopepars");
325  fparsSlope_refl = new TF1("fparsSlope_refl", reflected_functions[4].c_str(), 0., fT0_max);
326  for(unsigned int i=0; i<reflected_expoSlopepars.size(); ++i)
327  fparsSlope_refl->SetParameter(i, reflected_expoSlopepars[i]);
328 
329 
330  }
331 
332 
333  return;
334 
335  }
336 
337 
338 
339  //------------------------------------------------------
340 
341  // Eventually we will calculate the light quenching factor here
342  double PhotonVisibilityService::GetQuenchingFactor(double /* dQdx */) const
343  {
344  // for now, no quenching
345  return 1.0;
346 
347  }
348 
349 
350  //------------------------------------------------------
351 
352  // Get a vector of the relative visibilities of each OpDet
353  // in the event to a point xyz
354 
355  float const* PhotonVisibilityService::GetAllVisibilities(double const* xyz, bool wantReflected) const
356  {
357  if(fInterpolate){
358  static std::vector<float> ret;
359  if(ret.size() != NOpChannels()) ret.resize(NOpChannels());
360  for(unsigned int i = 0; i < NOpChannels(); ++i) ret[i] = GetVisibility(xyz, i, wantReflected);
361  return &ret.front();
362  }
363  else{
364  size_t VoxID = fVoxelDef.GetVoxelID(LibLocation(xyz));
365  return GetLibraryEntries(VoxID, wantReflected);
366  }
367  }
368 
369 
370  //------------------------------------------------------
371 
372  // Get distance to optical detector OpDet
373  double PhotonVisibilityService::DistanceToOpDet( double const* xyz, unsigned int OpDet )
374  {
376  return geom->OpDetGeoFromOpDet(OpDet).DistanceToPoint(xyz);
377 
378  }
379 
380 
381  //------------------------------------------------------
382 
383 
384  // Get the solid angle reduction factor for planar optical detector OpDet
385  double PhotonVisibilityService::SolidAngleFactor( double const* xyz, unsigned int OpDet )
386  {
388  return geom->OpDetGeoFromOpDet(OpDet).CosThetaFromNormal(xyz);
389  }
390 
391  //------------------------------------------------------
392 
393  float PhotonVisibilityService::GetVisibility(double const* xyz, unsigned int OpChannel, bool wantReflected) const
394  {
395  // Static to avoid reallocating this buffer between calls
396  // (GetNeighbouringVoxelIDs makes sure to clear it).
397  static std::vector<sim::PhotonVoxelDef::NeiInfo> neis;
398 
399  if(fInterpolate){
400  // In case we're outside the bounding box we'll get an empty vector here
401  // and return visibility 0, which seems OK.
403  }
404  else{
405  // For no interpolation, use a single entry with weight 1
406  neis.clear();
407  neis.emplace_back(fVoxelDef.GetVoxelID(LibLocation(xyz)), 1);
408  }
409 
410  // Sum up all the weighted neighbours to get interpolation behaviour
411  float vis = 0;
412  for(const sim::PhotonVoxelDef::NeiInfo& n: neis){
413  vis += n.weight * GetLibraryEntry(n.id, OpChannel, wantReflected);
414  }
415  return vis;
416  }
417 
418 
419  //------------------------------------------------------
420 
421  void PhotonVisibilityService::StoreLightProd(int VoxID, double N)
422  {
423  fCurrentVoxel = VoxID;
424  fCurrentValue = N;
425  mf::LogInfo("PhotonVisibilityService") << " PVS notes production of " << N << " photons at Vox " << VoxID<<std::endl;
426  }
427 
428 
429  //------------------------------------------------------
430 
431 
432  void PhotonVisibilityService::RetrieveLightProd(int& VoxID, double& N) const
433  {
434  N = fCurrentValue;
435  VoxID = fCurrentVoxel;
436  }
437 
438  //------------------------------------------------------
439 
440  void PhotonVisibilityService::SetLibraryEntry(int VoxID, int OpChannel, float N, bool wantReflected)
441  {
442  if(fTheLibrary == 0)
443  LoadLibrary();
444 
445  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
446 
447  if(!wantReflected)
448  lib->SetCount(VoxID,OpChannel, N);
449 
450  else
451  lib->SetReflCount(VoxID,OpChannel, N);
452 
453  //std::cout<< " PVS logging " << VoxID << " " << OpChannel<<std::endl;
454  LOG_DEBUG("PhotonVisibilityService") << " PVS logging " << VoxID << " " << OpChannel<<std::endl;
455  }
456 
457  //------------------------------------------------------
458 
459  float const* PhotonVisibilityService::GetLibraryEntries(int VoxID, bool wantReflected) const
460  {
461  if(fTheLibrary == 0)
462  LoadLibrary();
463 
464  if(!wantReflected)
465  return fTheLibrary->GetCounts(VoxID);
466  else
467  return fTheLibrary->GetReflCounts(VoxID);
468  }
469 
470  //------------------------------------------------------
471 
472  float PhotonVisibilityService::GetLibraryEntry(int VoxID, int Channel, bool wantReflected) const
473  {
474  if(fTheLibrary == 0)
475  LoadLibrary();
476 
477  if(!wantReflected)
478  return fTheLibrary->GetCount(VoxID, Channel);
479  else
480  return fTheLibrary->GetReflCount(VoxID, Channel);
481  }
482  // Methodes to handle the extra library parameter refl T0
483  //------------------------------------------------------
484 
485  // Get a vector of the refl <tfirst> of each OpDet
486  // in the event to a point xyz
487 
488  float const* PhotonVisibilityService::GetReflT0s(double const* xyz) const
489  {
490  int VoxID = fVoxelDef.GetVoxelID(LibLocation(xyz));
491  return GetLibraryReflT0Entries(VoxID);
492  }
493 
494  //------------------------------------------------------
495 
497  {
498  if(fTheLibrary == 0)
499  LoadLibrary();
500 
501  return fTheLibrary->GetReflT0s(VoxID);
502  }
503 
504  //------------------------------------------------------
505 
506  void PhotonVisibilityService::SetLibraryReflT0Entry(int VoxID, int OpChannel, float T0)
507  {
508  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
509  if(fTheLibrary == 0)
510  LoadLibrary();
511 
512  lib->SetReflT0(VoxID,OpChannel,T0);
513 
514  LOG_DEBUG("PhotonVisibilityService") << " PVS logging " << VoxID << " " << OpChannel<<std::endl;
515  }
516 
517  //------------------------------------------------------
518 
519  float PhotonVisibilityService::GetLibraryReflT0Entry(int VoxID, int Channel) const
520  {
521  if(fTheLibrary == 0)
522  LoadLibrary();
523 
524  return fTheLibrary->GetReflT0(VoxID, Channel);
525  }
526 
527  //------------------------------------------------------
528 
529 
531 
532  const std::vector<float>* PhotonVisibilityService::GetTimingPar(double const* xyz) const
533  {
534  int VoxID = fVoxelDef.GetVoxelID(LibLocation(xyz));
535  return GetLibraryTimingParEntries(VoxID);
536  }
537 
538  TF1* PhotonVisibilityService::GetTimingTF1(double const* xyz) const
539  {
540  int VoxID = fVoxelDef.GetVoxelID(LibLocation(xyz));
541  return GetLibraryTimingTF1Entries(VoxID);
542  }
543 
544 
545  //------------------------------------------------------
546 
547  const std::vector<float>* PhotonVisibilityService::GetLibraryTimingParEntries(int VoxID) const
548  {
549  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
550  if(fTheLibrary == 0)
551  LoadLibrary();
552 
553  return lib->GetTimingPars(VoxID);
554  }
555 
556  //------------------------------------------------------
557 
559  {
560  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
561  if(fTheLibrary == 0)
562  LoadLibrary();
563 
564  return lib->GetTimingTF1s(VoxID);
565  }
566 
567  //------------------------------------------------------
568 
569  void PhotonVisibilityService::SetLibraryTimingParEntry(int VoxID, int OpChannel, float par, size_t parnum)
570  {
571  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
572  if(fTheLibrary == 0)
573  LoadLibrary();
574 
575  lib->SetTimingPar(VoxID,OpChannel,par, parnum);
576 
577  LOG_DEBUG("PhotonVisibilityService") << " PVS logging " << VoxID << " " << OpChannel<<std::endl;
578  }
579 
580  //------------------------------------------------------
581 
582  void PhotonVisibilityService::SetLibraryTimingTF1Entry(int VoxID, int OpChannel, TF1 func)
583  {
584  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
585  if(fTheLibrary == 0)
586  LoadLibrary();
587 
588  lib->SetTimingTF1(VoxID,OpChannel,func);
589 
590  LOG_DEBUG("PhotonVisibilityService") << " PVS logging " << VoxID << " " << OpChannel<<std::endl;
591  }
592 
593 
594  //------------------------------------------------------
595 
596  float PhotonVisibilityService::GetLibraryTimingParEntry(int VoxID, int Channel, size_t npar) const
597  {
598  PhotonLibrary* lib = dynamic_cast<PhotonLibrary*>(fTheLibrary);
599  if(fTheLibrary == 0)
600  LoadLibrary();
601 
602  return lib->GetTimingPar(VoxID, Channel,npar);
603  }
604 
605  //------------------------------------------------------
606 
608  {
609  if(fTheLibrary == 0)
610  LoadLibrary();
611 
612  return fTheLibrary->NOpChannels();
613  }
614 
615  //------------------------------------------------------
616  void PhotonVisibilityService::SetDirectLightPropFunctions(TF1 const* functions[8], double& d_break, double& d_max, double& tf1_sampling_factor) const
617  {
618  functions[0] = fparslogNorm;
619  functions[1] = fparsMPV;
620  functions[2] = fparsWidth;
621  functions[3] = fparsCte;
622  functions[4] = fparsSlope;
623  functions[5] = fparslogNorm_far;
624  functions[6] = fparsMPV_far;
625  functions[7] = fparsCte_far;
626 
627  d_break = fD_break;
628  d_max = fD_max;
629  tf1_sampling_factor = fTF1_sampling_factor;
630  }
631 
632  //------------------------------------------------------
633  void PhotonVisibilityService::SetReflectedCOLightPropFunctions(TF1 const* functions[5], double& t0_max, double& t0_break_point) const
634  {
635  functions[0] = fparslogNorm_refl;
636  functions[1] = fparsMPV_refl;
637  functions[2] = fparsWidth_refl;
638  functions[3] = fparsCte_refl;
639  functions[4] = fparsSlope_refl;
640 
641  t0_max = fT0_max;
642  t0_break_point = fT0_break_point;
643  }
644 
645 
646  //------------------------------------------------------
647  /***
648  * Preform any necessary transformations on the coordinates before trying to access
649  * a voxel ID.
650  **/
651  const TVector3 PhotonVisibilityService::LibLocation(const double * xyz) const
652  {
653  TVector3 location(xyz);
654 
655  // Always use postive X coordinate if set
656  if (fReflectOverZeroX && location.x() < 0) {
657  location.SetX( fabs(location.x() ) );
658  }
659  return location;
660  }
661 
662 } // namespace
663 
664 namespace phot{
665 
667 
668 } // namespace phot
void GetNeighboringVoxelIDs(const TVector3 &v, std::vector< NeiInfo > &ret) const
OpDetGeo const & OpDetGeoFromOpDet(unsigned int OpDet) const
Number of OpDets in the whole detector.
virtual float GetReflT0(size_t Voxel, size_t OpChannel) const =0
double CosThetaFromNormal(geo::Point_t const &point) const
Get cos(angle) to normal of this detector - used for solid angle calcs.
Definition: OpDetGeo.cxx:121
void RetrieveLightProd(int &VoxID, double &N) const
Encapsulate the construction of a single cyostat.
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)
const std::vector< float > * GetLibraryTimingParEntries(int VoxID) const
virtual const float * GetReflCounts(size_t Voxel) const =0
void CreateEmptyLibrary(size_t NVoxels, size_t NChannels, bool storeReflected=false, bool storeReflT0=false, size_t storeTiming=0)
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
#define DEFINE_ART_SERVICE(svc)
Definition: ServiceMacros.h:93
void SetLibraryTimingParEntry(int VoxID, int OpChannel, float value, size_t parnum)
virtual float GetReflCount(size_t Voxel, size_t OpChannel) const =0
const std::vector< float > * GetTimingPar(double const *xyz) const
float const * GetLibraryReflT0Entries(int VoxID) const
void StoreLibraryToFile(std::string LibraryFile, bool storeReflected=false, bool storeReflT0=false, size_t storeTiming=0) const
void SetCount(size_t Voxel, size_t OpChannel, float Count)
PhotonVisibilityService(fhicl::ParameterSet const &pset, art::ActivityRegistry &reg)
void SetLibraryTimingTF1Entry(int VoxID, int OpChannel, TF1 func)
void CryostatBoundaries(double *boundaries, geo::CryostatID const &cid) const
Returns the boundaries of the specified cryostat.
void SetTimingTF1(size_t Voxel, size_t OpChannel, TF1 func)
int GetNVoxels() const
float const * GetLibraryEntries(int VoxID, bool wantReflected=false) const
void SetReflectedCOLightPropFunctions(TF1 const *functions[5], double &t0_max, double &t0_break_point) const
void SetReflT0(size_t Voxel, size_t OpChannel, float reflT0)
float const * GetAllVisibilities(double const *xyz, bool wantReflected=false) const
void SetLibraryReflT0Entry(int VoxID, int OpChannel, float value)
T get(std::string const &key) const
Definition: ParameterSet.h:231
virtual float GetCount(size_t Voxel, size_t OpChannel) const =0
TF1 * GetTimingTF1(double const *xyz) const
void reconfigure(fhicl::ParameterSet const &p)
virtual const float * GetCounts(size_t Voxel) const =0
Returns a pointer to NOpChannels() visibility values, one per channel.
TF1 * GetTimingTF1s(size_t Voxel) const
const std::vector< float > * GetTimingPars(size_t Voxel) const
void SetLibraryEntry(int VoxID, int OpChannel, float N, bool wantReflected=false)
int GetVoxelID(const TVector3 &) const
unsigned int NOpDets() const
Number of OpDets in the whole detector.
virtual int NOpChannels() const =0
Encapsulate the geometry of an optical detector.
void SetTimingPar(size_t Voxel, size_t OpChannel, float Count, size_t parnum)
const sim::PhotonVoxelDef & GetVoxelDef() const
double DistanceToPoint(geo::Point_t const &point) const
Returns the distance of the specified point from detector center [cm].
Definition: OpDetGeo.cxx:114
General LArSoft Utilities.
virtual const float * GetReflT0s(size_t Voxel) const =0
float GetLibraryTimingParEntry(int VoxID, int Channel, size_t npar) const
static double SolidAngleFactor(double const *xyz, unsigned int OpDet)
float GetTimingPar(size_t Voxel, size_t OpChannel, size_t parnum) const
float GetLibraryReflT0Entry(int VoxID, int Channel) const
#define LOG_DEBUG(id)
Char_t n[5]
const TVector3 LibLocation(const double *xyz) const
void SetReflCount(size_t Voxel, size_t OpChannel, float Count)
art framework interface to geometry description
float GetLibraryEntry(int VoxID, int OpChannel, bool wantReflected=false) const
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
float const * GetReflT0s(double const *xyz) const
float GetVisibility(double const *xyz, unsigned int OpChannel, bool wantReflected=false) const
static double DistanceToOpDet(double const *xyz, unsigned int OpDet)