Classes | Public Member Functions | Protected Attributes | List of all members
phot::PhotonLibraryHybrid Class Reference

TODO doc. More...

#include <PhotonLibraryHybrid.h>

Inheritance diagram for phot::PhotonLibraryHybrid:
phot::IPhotonLibrary

Classes

struct  Exception
 
struct  FitFunc
 
struct  OpDetRecord
 

Public Member Functions

 PhotonLibraryHybrid (const std::string &fname, const sim::PhotonVoxelDef &voxdef)
 
virtual ~PhotonLibraryHybrid ()
 
virtual float GetCount (size_t Voxel, size_t OpChannel) const override
 
virtual const float * GetCounts (size_t Voxel) const override
 Returns a pointer to NOpChannels() visibility values, one per channel. More...
 
virtual bool hasReflected () const override
 Don't implement reflected light. More...
 
virtual const float * GetReflCounts (size_t Voxel) const override
 
virtual float GetReflCount (size_t Voxel, size_t OpChannel) const override
 
virtual bool hasReflectedT0 () const override
 Don't implement reflected light timing. More...
 
virtual const float * GetReflT0s (size_t Voxel) const override
 
virtual float GetReflT0 (size_t Voxel, size_t OpChannel) const override
 
virtual int NOpChannels () const override
 
virtual int NVoxels () const override
 
- Public Member Functions inherited from phot::IPhotonLibrary
virtual ~IPhotonLibrary ()=default
 
virtual bool isVoxelValid (size_t Voxel) const
 
size_t LibrarySize () const
 Returns the number of elements in the library. More...
 

Protected Attributes

const sim::PhotonVoxelDeffVoxDef
 
std::vector< OpDetRecordfRecords
 

Additional Inherited Members

- Public Types inherited from phot::IPhotonLibrary
using Counts_t = const float *
 Type for visibility count per optical channel. More...
 
using T0s_t = const float *
 Type for time of arrival per optical channel. More...
 
using Params_t = std::vector< float > const *
 
using Functions_t = TF1 *
 

Detailed Description

TODO doc.

Definition at line 18 of file PhotonLibraryHybrid.h.

Constructor & Destructor Documentation

phot::PhotonLibraryHybrid::PhotonLibraryHybrid ( const std::string fname,
const sim::PhotonVoxelDef voxdef 
)

Definition at line 27 of file PhotonLibraryHybrid.cxx.

29  : fVoxDef(voxdef)
30  {
31  TFile f(fname.c_str());
32  !f.IsZombie() || fatal("Could not open PhotonLibrary "+fname);
33 
34  for(int opdetIdx = 0; true; ++opdetIdx){
35  const std::string dirname = TString::Format("opdet_%d", opdetIdx).Data();
36  TDirectory* dir = (TDirectory*)f.Get(dirname.c_str());
37  if(!dir) break; // Ran out of opdets
38 
39  OpDetRecord rec;
40 
41  TVectorD* fit = (TVectorD*)dir->Get("fit");
42  fit || fatal("Didn't find "+dirname+"/fit in "+fname);
43  rec.fit = FitFunc((*fit)[0], (*fit)[1]);
44 
45  TTree* tr = (TTree*)dir->Get("tr");
46  tr || fatal("Didn't find "+dirname+"/tr in "+fname);
47  int vox;
48  float vis;
49  tr->SetBranchAddress("vox", &vox);
50  tr->SetBranchAddress("vis", &vis);
51 
52  rec.exceptions.reserve(tr->GetEntries());
53  for(int i = 0; i < tr->GetEntries(); ++i){
54  tr->GetEntry(i);
55  vox < NVoxels() || fatal("Voxel out of range");
56  rec.exceptions.emplace_back(vox, vis);
57  }
58  rec.exceptions.shrink_to_fit(); // In case we don't trust reserve()
59 
60  // In case they weren't sorted in the file
61  std::sort(rec.exceptions.begin(), rec.exceptions.end());
62 
63  fRecords.push_back(rec);
64  } // end for opdetIdx
65 
66  !fRecords.empty() || fatal("No opdet_*/ directories in "+fname);
67 
69  geom->NOpDets() == fRecords.size() || fatal("Number of opdets mismatch");
70 
71  fRecords.shrink_to_fit(); // save memory
72  }
rec
Definition: tracks.py:88
Format
Definition: utils.h:7
std::string string
Definition: nybbler.cc:12
string dir
virtual int NVoxels() const override
void fatal(const char *msg,...)
Definition: qglobal.cpp:465
const sim::PhotonVoxelDef & fVoxDef
unsigned int NOpDets() const
Number of OpDets in the whole detector.
std::vector< OpDetRecord > fRecords
phot::PhotonLibraryHybrid::~PhotonLibraryHybrid ( )
virtual

Definition at line 75 of file PhotonLibraryHybrid.cxx.

76  {
77  }

Member Function Documentation

float phot::PhotonLibraryHybrid::GetCount ( size_t  Voxel,
size_t  OpChannel 
) const
overridevirtual

Implements phot::IPhotonLibrary.

Definition at line 98 of file PhotonLibraryHybrid.cxx.

99  {
100  int(vox) < NVoxels() || fatal("GetCount(): Voxel out of range");
101  int(opchan) < NOpChannels() || fatal("GetCount(): OpChan out of range");
102 
103  const OpDetRecord& rec = fRecords[opchan];
104 
105  auto it2 = std::lower_bound(rec.exceptions.begin(),
106  rec.exceptions.end(),
107  vox);
108  if(it2->vox == vox) return it2->vis;
109 
110  // Otherwise, we use the interpolation, which requires a distance
111 
113  const geo::OpDetGeo& opdet = geom->OpDetGeoFromOpDet(opchan);
114 
115  const auto voxvec = fVoxDef.GetPhotonVoxel(vox).GetCenter();
116 
117  const double dist = opdet.DistanceToPoint(voxvec);
118 
119  return rec.fit.Eval(dist);
120  }
rec
Definition: tracks.py:88
OpDetGeo const & OpDetGeoFromOpDet(unsigned int OpDet) const
Returns the geo::OpDetGeo object for the given detector number.
virtual int NVoxels() const override
void fatal(const char *msg,...)
Definition: qglobal.cpp:465
const sim::PhotonVoxelDef & fVoxDef
virtual int NOpChannels() const override
double DistanceToPoint(geo::Point_t const &point) const
Returns the distance of the specified point from detector center [cm].
Definition: OpDetGeo.cxx:123
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
Point GetCenter() const
Returns the center of the voxel (type Point).
Definition: PhotonVoxels.h:199
std::vector< OpDetRecord > fRecords
PhotonVoxel GetPhotonVoxel(int ID) const
const float * phot::PhotonLibraryHybrid::GetCounts ( size_t  Voxel) const
overridevirtual

Returns a pointer to NOpChannels() visibility values, one per channel.

Implements phot::IPhotonLibrary.

Definition at line 86 of file PhotonLibraryHybrid.cxx.

87  {
88  // The concept of GetCounts() conflicts fairly badly with how this library
89  // works. This is probably the best we can do.
90 
91  static float* counts = 0;
92  if(!counts) counts = new float[NOpChannels()];
93  for(int od = 0; od < NOpChannels(); ++od) counts[od] = GetCount(vox, od);
94  return counts;
95  }
counts_as<> counts
Number of ADC counts, represented by signed short int.
Definition: electronics.h:116
virtual int NOpChannels() const override
virtual float GetCount(size_t Voxel, size_t OpChannel) const override
virtual float phot::PhotonLibraryHybrid::GetReflCount ( size_t  Voxel,
size_t  OpChannel 
) const
inlineoverridevirtual

Implements phot::IPhotonLibrary.

Definition at line 33 of file PhotonLibraryHybrid.h.

33 {return 0;}
virtual const float* phot::PhotonLibraryHybrid::GetReflCounts ( size_t  Voxel) const
inlineoverridevirtual

Implements phot::IPhotonLibrary.

Definition at line 32 of file PhotonLibraryHybrid.h.

32 {return 0;}
virtual float phot::PhotonLibraryHybrid::GetReflT0 ( size_t  Voxel,
size_t  OpChannel 
) const
inlineoverridevirtual

Implements phot::IPhotonLibrary.

Definition at line 38 of file PhotonLibraryHybrid.h.

38 {return 0;}
virtual const float* phot::PhotonLibraryHybrid::GetReflT0s ( size_t  Voxel) const
inlineoverridevirtual

Implements phot::IPhotonLibrary.

Definition at line 37 of file PhotonLibraryHybrid.h.

37 {return 0;}
virtual bool phot::PhotonLibraryHybrid::hasReflected ( ) const
inlineoverridevirtual

Don't implement reflected light.

Implements phot::IPhotonLibrary.

Definition at line 31 of file PhotonLibraryHybrid.h.

31 {return false;}
virtual bool phot::PhotonLibraryHybrid::hasReflectedT0 ( ) const
inlineoverridevirtual

Don't implement reflected light timing.

Implements phot::IPhotonLibrary.

Definition at line 36 of file PhotonLibraryHybrid.h.

36 {return false;}
virtual int phot::PhotonLibraryHybrid::NOpChannels ( ) const
inlineoverridevirtual

Implements phot::IPhotonLibrary.

Definition at line 40 of file PhotonLibraryHybrid.h.

40 {return fRecords.size();}
std::vector< OpDetRecord > fRecords
int phot::PhotonLibraryHybrid::NVoxels ( ) const
overridevirtual

Implements phot::IPhotonLibrary.

Definition at line 80 of file PhotonLibraryHybrid.cxx.

81  {
82  return fVoxDef.GetNVoxels();
83  }
const sim::PhotonVoxelDef & fVoxDef
unsigned int GetNVoxels() const
Returns the total number of voxels in the volume.

Member Data Documentation

std::vector<OpDetRecord> phot::PhotonLibraryHybrid::fRecords
protected

Definition at line 74 of file PhotonLibraryHybrid.h.

const sim::PhotonVoxelDef& phot::PhotonLibraryHybrid::fVoxDef
protected

Definition at line 44 of file PhotonLibraryHybrid.h.


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