All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
flashmatch::PhotonLibHypothesis Class Reference

#include <PhotonLibHypothesis.h>

Inheritance diagram for flashmatch::PhotonLibHypothesis:
flashmatch::BaseFlashHypothesis flashmatch::BaseAlgorithm flashmatch::LoggerFeature

Public Member Functions

 PhotonLibHypothesis (const std::string name="PhotonLibHypothesis")
 Default constructor. More...
 
virtual ~PhotonLibHypothesis ()
 Default destructor. More...
 
void FillEstimate (const QCluster_t &, Flash_t &) const
 Method to simply fill provided reference of flashmatch::Flash_t. More...
 
- Public Member Functions inherited from flashmatch::BaseFlashHypothesis
 BaseFlashHypothesis (const std::string name="noname")
 Default constructor. More...
 
 ~BaseFlashHypothesis ()
 Default destructor. More...
 
Flash_t GetEstimate (const QCluster_t &) const
 Method to create flashmatch::Flash_t object and return. More...
 
void SetChannelMask (std::vector< int > ch_mask)
 Sets the channels to use. More...
 
void SetUncoatedPMTs (std::vector< int > ch_uncoated)
 Sets the channels sensitive to visible light. More...
 
void SetSemiAnalyticalModel (std::unique_ptr< SemiAnalyticalModel > model)
 Sets the semi analytical model. More...
 
- Public Member Functions inherited from flashmatch::BaseAlgorithm
 BaseAlgorithm (const Algorithm_t type, const std::string name)
 Default constructor. More...
 
 ~BaseAlgorithm ()
 Default destructor. More...
 
void Configure (const Config_t &pset)
 Function to accept configuration. More...
 
Algorithm_t AlgorithmType () const
 Algorithm type. More...
 
const std::string & AlgorithmName () const
 Algorithm name. More...
 
- Public Member Functions inherited from flashmatch::LoggerFeature
 LoggerFeature (const std::string logger_name="LoggerFeature")
 Default constructor. More...
 
 LoggerFeature (const LoggerFeature &original)
 Default copy constructor. More...
 
virtual ~LoggerFeature ()
 Default destructor. More...
 
const flashmatch::loggerlogger () const
 Logger getter. More...
 
void set_verbosity (::flashmatch::msg::Level_t level)
 Verbosity level. More...
 
const std::string & name () const
 Name getter, defined in a logger instance attribute. More...
 

Protected Member Functions

void _Configure_ (const Config_t &pset)
 

Protected Attributes

double _global_qe
 Global QE for direct light. More...
 
double _global_qe_refl
 Global QE for reflected light. More...
 
double _sigma_qe
 Sigma for Gaussian centered on Global QE. More...
 
std::vector< double > _qe_v
 PMT-wise relative QE. More...
 
bool _use_semi_analytical
 If the semi-analytical approach should be used. More...
 
phot::PhotonVisibilityService
const *const 
_vis
 
- Protected Attributes inherited from flashmatch::BaseFlashHypothesis
std::vector< int > _channel_mask
 The list of channels to use. More...
 
std::vector< int > _uncoated_pmt_list
 A list of opdet sensitive to visible (reflected) light. More...
 
std::unique_ptr
< SemiAnalyticalModel
_semi_model
 

Private Member Functions

void FillEstimateSemiAnalytical (const QCluster_t &, Flash_t &) const
 Fills the estimate using the semi analytical approach (SBND) More...
 
void FillEstimateLibrary (const QCluster_t &, Flash_t &) const
 Fills the estimate using the photon library (ICARUS, SBND) More...
 

Detailed Description

User custom analysis class made by SHELL_USER_NAME

Definition at line 43 of file PhotonLibHypothesis.h.

Constructor & Destructor Documentation

flashmatch::PhotonLibHypothesis::PhotonLibHypothesis ( const std::string  name = "PhotonLibHypothesis")

Default constructor.

Definition at line 18 of file PhotonLibHypothesis.cxx.

20  , _vis(art::ServiceHandle<phot::PhotonVisibilityService const>().get())
21  {}
BaseFlashHypothesis(const std::string name="noname")
Default constructor.
phot::PhotonVisibilityService const *const _vis
const std::string & name() const
Name getter, defined in a logger instance attribute.
Definition: LoggerFeature.h:51
virtual flashmatch::PhotonLibHypothesis::~PhotonLibHypothesis ( )
inlinevirtual

Default destructor.

Definition at line 51 of file PhotonLibHypothesis.h.

51 {}

Member Function Documentation

void flashmatch::PhotonLibHypothesis::_Configure_ ( const Config_t pset)
protectedvirtual

Implements flashmatch::BaseAlgorithm.

Definition at line 23 of file PhotonLibHypothesis.cxx.

24  {
25  #if USING_LARSOFT == 0
26  omp_set_num_threads(NUM_THREADS);
27  #endif
28 
29  _global_qe = pset.get<double>("GlobalQE");
30  _global_qe_refl = pset.get<double>("GlobalQERefl", 0);
31  _use_semi_analytical = pset.get<bool>("UseSemiAnalytical", 0);
32 
33  _qe_v.clear();
34  _qe_v = pset.get<std::vector<double> >("CCVCorrection",_qe_v);
35  if(_qe_v.empty()) _qe_v.resize(DetectorSpecs::GetME().NOpDets(),1.0);
36  if(_qe_v.size() != DetectorSpecs::GetME().NOpDets()) {
37  FLASH_CRITICAL() << "CCVCorrection factor array has size " << _qe_v.size()
38  << " != number of opdet (" << DetectorSpecs::GetME().NOpDets() << ")!" << std::endl;
39  throw OpT0FinderException();
40  }
41 
42  // By default, add all opdets to the channel mask
43  // Note that this may be overridden by the manager
44  // via the SetChannelMask() method.
45  _channel_mask.clear();
47  for (size_t i = 0; i < DetectorSpecs::GetME().NOpDets(); i++) {
48  _channel_mask[i] = i;
49  }
50  }
double _global_qe
Global QE for direct light.
#define NUM_THREADS
double _global_qe_refl
Global QE for reflected light.
size_t NOpDets() const
of PMTs
Definition: FMWKInterface.h:60
static DetectorSpecs & GetME(std::string filename="detector_specs.cfg")
Definition: FMWKInterface.h:44
#define FLASH_CRITICAL()
Compiler macro for CRITICAL message.
bool _use_semi_analytical
If the semi-analytical approach should be used.
std::vector< double > _qe_v
PMT-wise relative QE.
std::vector< int > _channel_mask
The list of channels to use.
void flashmatch::PhotonLibHypothesis::FillEstimate ( const QCluster_t ,
Flash_t  
) const
virtual

Method to simply fill provided reference of flashmatch::Flash_t.

Implements flashmatch::BaseFlashHypothesis.

Definition at line 57 of file PhotonLibHypothesis.cxx.

58  {
59 
60  for ( auto& v : flash.pe_v ) v = 0;
61 
63  FillEstimateSemiAnalytical(trk, flash);
64  } else {
65  FillEstimateLibrary(trk, flash);
66  }
67  }
void FillEstimateLibrary(const QCluster_t &, Flash_t &) const
Fills the estimate using the photon library (ICARUS, SBND)
void FillEstimateSemiAnalytical(const QCluster_t &, Flash_t &) const
Fills the estimate using the semi analytical approach (SBND)
bool _use_semi_analytical
If the semi-analytical approach should be used.
void flashmatch::PhotonLibHypothesis::FillEstimateLibrary ( const QCluster_t trk,
Flash_t flash 
) const
private

Fills the estimate using the photon library (ICARUS, SBND)

Definition at line 134 of file PhotonLibHypothesis.cxx.

135  {
136 
137  static double xyz[3] = {0.};
138 
139  size_t n_pmt = DetectorSpecs::GetME().NOpDets();
140 
141  for ( size_t ipmt = 0; ipmt < n_pmt; ++ipmt) {
142 
143  if (std::find(_channel_mask.begin(), _channel_mask.end(), ipmt) == _channel_mask.end()) {
144  continue;
145  }
146 
147  bool is_uncoated = false;
148  if (std::find(_uncoated_pmt_list.begin(), _uncoated_pmt_list.end(), ipmt) != _uncoated_pmt_list.end()) {
149  is_uncoated = true;
150  }
151 
152  for (size_t ipt = 0; ipt < trk.size(); ++ipt) {
153  auto const& pt = trk[ipt];
154 
155  double q = pt.q;
156 
157  xyz[0] = pt.x;
158  xyz[1] = pt.y;
159  xyz[2] = pt.z;
160 
161  double q_direct = 0, q_refl = 0;
162 
163  // Direct light
164  if (is_uncoated) {
165  q_direct = 0;
166  } else {
167  q_direct = q * _vis->GetVisibility(xyz, ipmt) * _global_qe / _qe_v[ipmt];
168  }
169  // std::cout << "PMT : " << ipmt << " [x,y,z] -> [q] : [" << pt.x << ", " << pt.y << ", " << pt.z << "] -> [" << q << "]" << std::endl;
170 
171  // Reflected light
172  q_refl = q * _vis->GetVisibility(xyz, ipmt, true) * _global_qe_refl / _qe_v[ipmt];
173 
174  flash.pe_v[ipmt] += q_direct;
175  flash.pe_v[ipmt] += q_refl;
176  }
177  }
178  return;
179  }
double _global_qe
Global QE for direct light.
double _global_qe_refl
Global QE for reflected light.
size_t NOpDets() const
of PMTs
Definition: FMWKInterface.h:60
static DetectorSpecs & GetME(std::string filename="detector_specs.cfg")
Definition: FMWKInterface.h:44
std::vector< int > _uncoated_pmt_list
A list of opdet sensitive to visible (reflected) light.
std::vector< double > _qe_v
PMT-wise relative QE.
phot::PhotonVisibilityService const *const _vis
float GetVisibility(Point const &p, unsigned int OpChannel, bool wantReflected=false) const
std::vector< int > _channel_mask
The list of channels to use.
void flashmatch::PhotonLibHypothesis::FillEstimateSemiAnalytical ( const QCluster_t trk,
Flash_t flash 
) const
private

Fills the estimate using the semi analytical approach (SBND)

Get the 3D point in space from where photons should be propagated

Definition at line 70 of file PhotonLibHypothesis.cxx.

71  {
72 
73  for ( size_t ipt = 0; ipt < trk.size(); ++ipt) {
74 
75  /// Get the 3D point in space from where photons should be propagated
76  auto const& pt = trk[ipt];
77 
78  // Get the number of photons produced in such point
79  double n_original_photons = pt.q;
80 
81  geo::Point_t const xyz = {pt.x, pt.y, pt.z};
82 
83  std::vector<double> direct_visibilities;
84  _semi_model->detectedDirectVisibilities(direct_visibilities, xyz);
85 
86  std::vector<double> reflected_visibilities;
87  _semi_model->detectedReflectedVisibilities(reflected_visibilities, xyz);
88 
89  //
90  // Fill Estimate with Direct light
91  //
92  for (size_t op_det=0; op_det<direct_visibilities.size(); ++op_det) {
93 
94  const double visibility = direct_visibilities[op_det];
95 
96  double q = n_original_photons * visibility * _global_qe / _qe_v[op_det];
97 
98  // Coated PMTs don't see direct photons
99  if (std::find(_uncoated_pmt_list.begin(), _uncoated_pmt_list.end(), op_det) != _uncoated_pmt_list.end()) {
100  q = 0;
101  }
102 
103  // std::cout << "OpDet: " << op_det << " [x,y,z] -> [q] : [" << pt.x << ", " << pt.y << ", " << pt.z << "] -> [" << q << "]" << std::endl;
104 
105  if (std::find(_channel_mask.begin(), _channel_mask.end(), op_det) != _channel_mask.end()) {
106  flash.pe_v[op_det] += q;
107  } else {
108  flash.pe_v[op_det] = 0;
109  }
110  }
111 
112  //
113  // Fill Estimate with Reflected light
114  //
115  for (size_t op_det=0; op_det<reflected_visibilities.size(); ++op_det) {
116 
117  const double visibility = reflected_visibilities[op_det];
118 
119  double q = n_original_photons * visibility * _global_qe_refl / _qe_v[op_det];
120 
121  if (std::find(_channel_mask.begin(), _channel_mask.end(), op_det) != _channel_mask.end()) {
122  flash.pe_v[op_det] += q;
123  } else {
124  flash.pe_v[op_det] = 0;
125  }
126 
127  }
128 
129  }
130  }
double _global_qe
Global QE for direct light.
double _global_qe_refl
Global QE for reflected light.
std::vector< int > _uncoated_pmt_list
A list of opdet sensitive to visible (reflected) light.
std::vector< double > _qe_v
PMT-wise relative QE.
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
std::unique_ptr< SemiAnalyticalModel > _semi_model
std::vector< int > _channel_mask
The list of channels to use.

Member Data Documentation

double flashmatch::PhotonLibHypothesis::_global_qe
protected

Global QE for direct light.

Definition at line 66 of file PhotonLibHypothesis.h.

double flashmatch::PhotonLibHypothesis::_global_qe_refl
protected

Global QE for reflected light.

Definition at line 67 of file PhotonLibHypothesis.h.

std::vector<double> flashmatch::PhotonLibHypothesis::_qe_v
protected

PMT-wise relative QE.

Definition at line 69 of file PhotonLibHypothesis.h.

double flashmatch::PhotonLibHypothesis::_sigma_qe
protected

Sigma for Gaussian centered on Global QE.

Definition at line 68 of file PhotonLibHypothesis.h.

bool flashmatch::PhotonLibHypothesis::_use_semi_analytical
protected

If the semi-analytical approach should be used.

Definition at line 70 of file PhotonLibHypothesis.h.

phot::PhotonVisibilityService const* const flashmatch::PhotonLibHypothesis::_vis
protected

Definition at line 72 of file PhotonLibHypothesis.h.


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