Classes | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
gar::detinfo::DetectorPropertiesStandard Class Reference

#include <DetectorPropertiesStandard.h>

Inheritance diagram for gar::detinfo::DetectorPropertiesStandard:
gar::detinfo::DetectorProperties

Classes

struct  Configuration_t
 Structure for configuration parameters. More...
 
struct  SternheimerParameters_t
 Parameters for Sternheimer density effect corrections. More...
 

Public Types

using providers_type = gar::ProviderPack< geo::GeometryCore, detinfo::GArProperties, detinfo::ECALProperties, detinfo::DetectorClocks >
 List of service providers we depend on. More...
 

Public Member Functions

 DetectorPropertiesStandard ()
 
 DetectorPropertiesStandard (fhicl::ParameterSet const &pset, const geo::GeometryCore *geo, const detinfo::GArProperties *gp, const detinfo::ECALProperties *ecalp, const detinfo::DetectorClocks *c, std::set< std::string > const &ignore_params={})
 
 DetectorPropertiesStandard (fhicl::ParameterSet const &pset, providers_type providers, std::set< std::string > const &ignore_params={})
 Constructs the provider and sets up the dependencies. More...
 
 DetectorPropertiesStandard (DetectorPropertiesStandard const &)=delete
 
virtual ~DetectorPropertiesStandard ()=default
 
void ValidateAndConfigure (fhicl::ParameterSet const &p, std::set< std::string > const &ignore_params={})
 Configures the provider, first validating the configuration. More...
 
void Configure (Configuration_t const &config)
 Extracts the relevant configuration from the specified object. More...
 
Configuration_t ValidateConfiguration (fhicl::ParameterSet const &p, std::set< std::string > const &ignore_params={})
 Validates the specified configuration. More...
 
bool Update (uint64_t ts)
 
bool UpdateClocks (const detinfo::DetectorClocks *clks)
 
void Setup (providers_type providers)
 Sets all the providers at once. More...
 
void SetGeometry (const geo::GeometryCore *g)
 
void SetGArProperties (const detinfo::GArProperties *gp)
 
void SetECALProperties (const detinfo::ECALProperties *ecalp)
 
void SetDetectorClocks (const detinfo::DetectorClocks *clks)
 
void SetNumberTimeSamples (unsigned int nsamp)
 
virtual double Efield (unsigned int planegap=0) const override
 kV/cm More...
 
virtual double DriftVelocity (double efield=0., double temperature=0., bool cmPerns=true) const override
 cm/ns if true, otherwise cm/us More...
 
virtual double ElectronLifetime () const override
 dQ/dX in electrons/cm, returns dE/dX in MeV/cm. More...
 
virtual double Density (double temperature) const override
 Returns argon density at a given temperature. More...
 
virtual double Density () const override
 Returns argon density at the temperature from Temperature() More...
 
virtual double Temperature () const override
 In kelvin. More...
 
virtual double Eloss (double mom, double mass, double tcut) const override
 Restricted mean energy loss (dE/dx) More...
 
virtual double ElossVar (double mom, double mass) const override
 Energy loss fluctuation ( $ \sigma_{E}^2 / x $) More...
 
virtual double SamplingRate () const override
 
virtual double ElectronsToADC () const override
 
virtual unsigned int NumberTimeSamples () const override
 
virtual int TriggerOffset () const override
 
virtual double ConvertXToTicks (double X) const override
 
virtual double ConvertTicksToX (double ticks) const override
 
virtual double ConvertTDCToTicks (double tdc) const override
 
virtual double ConvertTicksToTDC (double ticks) const override
 
virtual double EffectivePixel () const override
 
virtual double LightYield () const override
 
virtual double SiPMGain () const override
 
virtual double IntercalibrationFactor () const override
 
virtual double ADCSaturation () const override
 
virtual double TimeResolution () const override
 
virtual double MeVtoMIP () const override
 
virtual double NoisePx () const override
 
void CheckIfConfigured () const
 
- Public Member Functions inherited from gar::detinfo::DetectorProperties
 DetectorProperties (const DetectorProperties &)=delete
 
 DetectorProperties (DetectorProperties &&)=delete
 
DetectorPropertiesoperator= (const DetectorProperties &)=delete
 
DetectorPropertiesoperator= (DetectorProperties &&)=delete
 
virtual ~DetectorProperties ()=default
 

Protected Member Functions

void CalculateXTicksParams ()
 
- Protected Member Functions inherited from gar::detinfo::DetectorProperties
 DetectorProperties ()=default
 

Protected Attributes

const detinfo::GArPropertiesfGP
 
const detinfo::ECALPropertiesfECALP
 
const detinfo::DetectorClocksfClocks
 
const geo::GeometryCorefGeo
 
std::vector< double > fEfield
 kV/cm (per inter-plane volume) More...
 
double fElectronlifetime
 microseconds More...
 
double fTemperature
 kelvin More...
 
double fDriftVelocity
 centimeters / microsecond More...
 
double fSamplingRate
 in ns More...
 
double fElectronsToADC
 conversion factor for # of ionization electrons to 1 ADC count More...
 
unsigned int fNumberTimeSamples
 number of clock ticks per event (= readout window) More...
 
SternheimerParameters_t fSternheimerParameters
 Sternheimer parameters. More...
 
double fXTicksCoefficient
 Parameters for x<–>ticks. More...
 
detinfo::ElecClock fTPCClock
 TPC electronics clock. More...
 

Detailed Description

Definition at line 35 of file DetectorPropertiesStandard.h.

Member Typedef Documentation

List of service providers we depend on.

Definition at line 41 of file DetectorPropertiesStandard.h.

Constructor & Destructor Documentation

gar::detinfo::DetectorPropertiesStandard::DetectorPropertiesStandard ( )

Definition at line 35 of file DetectorPropertiesStandard.cxx.

36  : fGP(0)
37  , fECALP(0)
38  , fClocks(0)
39  , fGeo(0)
40  {
41 
42  }
gar::detinfo::DetectorPropertiesStandard::DetectorPropertiesStandard ( fhicl::ParameterSet const &  pset,
const geo::GeometryCore geo,
const detinfo::GArProperties gp,
const detinfo::ECALProperties ecalp,
const detinfo::DetectorClocks c,
std::set< std::string > const &  ignore_params = {} 
)

Definition at line 45 of file DetectorPropertiesStandard.cxx.

51  : fGP(gp)
52  , fECALP(ecalp)
53  , fClocks(c)
54  , fGeo(geo)
55  {
56  {
57  mf::LogInfo debug("setupProvider<DetectorPropertiesStandard>");
58 
59  debug << "Asked to ignore " << ignore_params.size() << " keys:";
60  for (auto const& key: ignore_params) debug << " '" << key << "'";
61  }
62 
63  ValidateAndConfigure(pset, ignore_params);
64 
66 
67  MF_LOG_WARNING("DetectorPropertiesStandard")
68  << "DetectorPropertiesStandard Warning!"
69  << "\t The following functions need to be verified for gaseous argon:"
70  << "\t\t BirksCorrection"
71  << "\t\t Density"
72  << "\t\t DriftVelocity"
73  << "\t\t EField";
74 
75 
76  }
void ValidateAndConfigure(fhicl::ParameterSet const &p, std::set< std::string > const &ignore_params={})
Configures the provider, first validating the configuration.
def key(type, name=None)
Definition: graph.py:13
virtual const ElecClock & TPCClock() const =0
Borrow a const TPC clock with time set to Trigger time [ns].
#define MF_LOG_WARNING(category)
detinfo::ElecClock fTPCClock
TPC electronics clock.
gar::detinfo::DetectorPropertiesStandard::DetectorPropertiesStandard ( fhicl::ParameterSet const &  pset,
providers_type  providers,
std::set< std::string > const &  ignore_params = {} 
)

Constructs the provider and sets up the dependencies.

Parameters
psetFHiCL parameter set for provider configuration
providerspack of providers DetectorPropertiesStandard depends on
See also
Setup()

Definition at line 79 of file DetectorPropertiesStandard.cxx.

83  providers.get<geo::GeometryCore>(),
84  providers.get<detinfo::GArProperties>(),
85  providers.get<detinfo::ECALProperties>(),
86  providers.get<detinfo::DetectorClocks>(),
87  ignore_params)
88  {}
Description of geometry of one entire detector.
Class used for the conversion of times between different formats and references.
gar::detinfo::DetectorPropertiesStandard::DetectorPropertiesStandard ( DetectorPropertiesStandard const &  )
delete
virtual gar::detinfo::DetectorPropertiesStandard::~DetectorPropertiesStandard ( )
virtualdefault

Member Function Documentation

virtual double gar::detinfo::DetectorPropertiesStandard::ADCSaturation ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 244 of file DetectorPropertiesStandard.h.

244 { return fECALP->ADCSaturation(); }
virtual double ADCSaturation() const =0
void gar::detinfo::DetectorPropertiesStandard::CalculateXTicksParams ( )
protected

Definition at line 361 of file DetectorPropertiesStandard.cxx.

362  {
364 
365  double samplingRate = SamplingRate();
366  double efield = Efield();
367  double temperature = Temperature();
368  double driftVelocity = DriftVelocity(efield, temperature);
369 
370  fXTicksCoefficient = 0.001 * driftVelocity * samplingRate;
371 
372  return;
373  }
virtual double DriftVelocity(double efield=0., double temperature=0., bool cmPerns=true) const override
cm/ns if true, otherwise cm/us
virtual double Temperature() const override
In kelvin.
double fXTicksCoefficient
Parameters for x<–>ticks.
virtual double Efield(unsigned int planegap=0) const override
kV/cm
void gar::detinfo::DetectorPropertiesStandard::CheckIfConfigured ( ) const

Verifies that the provider is in a fully configured status

Exceptions
cet::exception(category DetectorPropertiesStandard) if not ok

Definition at line 351 of file DetectorPropertiesStandard.cxx.

352  {
353  if (!fGeo) throw cet::exception(__FUNCTION__) << "Geometry is uninitialized!";
354  if (!fGP) throw cet::exception(__FUNCTION__) << "GArPropertiesStandard is uninitialized!";
355  if (!fECALP) throw cet::exception(__FUNCTION__) << "ECALPropertiesStandard is uninitialized!";
356  if (!fClocks) throw cet::exception(__FUNCTION__) << "DetectorClocks is uninitialized!";
357  }
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
void gar::detinfo::DetectorPropertiesStandard::Configure ( Configuration_t const &  config)

Extracts the relevant configuration from the specified object.

Definition at line 122 of file DetectorPropertiesStandard.cxx.

122  {
123 
124  fEfield = config.Efield();
125  fElectronlifetime = config.Electronlifetime();
126  fTemperature = config.Temperature();
127  fDriftVelocity = config.DriftVelocity();
128  fElectronsToADC = config.ElectronsToADC();
129  fNumberTimeSamples = config.NumberTimeSamples();
130 
131  fSternheimerParameters.a = config.SternheimerA();
132  fSternheimerParameters.k = config.SternheimerK();
133  fSternheimerParameters.x0 = config.SternheimerX0();
134  fSternheimerParameters.x1 = config.SternheimerX1();
135  fSternheimerParameters.cbar = config.SternheimerCbar();
136 
138 
139  } // DetectorPropertiesStandard::Configure()
SternheimerParameters_t fSternheimerParameters
Sternheimer parameters.
static Config * config
Definition: config.cpp:1054
double fElectronsToADC
conversion factor for # of ionization electrons to 1 ADC count
std::vector< double > fEfield
kV/cm (per inter-plane volume)
unsigned int fNumberTimeSamples
number of clock ticks per event (= readout window)
double fDriftVelocity
centimeters / microsecond
double gar::detinfo::DetectorPropertiesStandard::ConvertTDCToTicks ( double  tdc) const
overridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 109 of file DetectorPropertiesStandard.cxx.

110  {
111  return fClocks->TPCTDC2Tick(tdc);
112  }
virtual double TPCTDC2Tick(double tdc) const =0
Given electronics clock count [tdc] returns TPC time-tick.
double gar::detinfo::DetectorPropertiesStandard::ConvertTicksToTDC ( double  ticks) const
overridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 115 of file DetectorPropertiesStandard.cxx.

116  {
117  return fClocks->TPCTick2TDC(ticks);
118  }
tick ticks
Alias for common language habits.
Definition: electronics.h:78
virtual double TPCTick2TDC(double tick) const =0
Given TPC time-tick (waveform index), returns electronics clock count [tdc].
double gar::detinfo::DetectorPropertiesStandard::ConvertTicksToX ( double  ticks) const
overridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 344 of file DetectorPropertiesStandard.cxx.

345  {
346  return ticks * fXTicksCoefficient;
347  }
tick ticks
Alias for common language habits.
Definition: electronics.h:78
double fXTicksCoefficient
Parameters for x<–>ticks.
double gar::detinfo::DetectorPropertiesStandard::ConvertXToTicks ( double  X) const
overridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 336 of file DetectorPropertiesStandard.cxx.

337  {
338  return (X / fXTicksCoefficient);
339  }
double fXTicksCoefficient
Parameters for x<–>ticks.
double gar::detinfo::DetectorPropertiesStandard::Density ( double  temperature) const
overridevirtual

Returns argon density at a given temperature.

Parameters
temperaturethe temperature in kelvin
Returns
argon density in g/cm^3

Density is nearly a linear function of temperature. See the NIST tables for details Slope is between -6.2 and -6.1, intercept is 1928 kg/m^3. This parameterization will be good to better than 0.5%.g/cm^3

Implements gar::detinfo::DetectorProperties.

Definition at line 189 of file DetectorPropertiesStandard.cxx.

190  {
191  // Default temperature use internal value.
192  if(temperature == 0.)
193  temperature = Temperature();
194 
195  double density = -0.00615*temperature + 1.928;
196 
197  return density;
198  } // DetectorPropertiesStandard::Density()
virtual double Temperature() const override
In kelvin.
virtual double gar::detinfo::DetectorPropertiesStandard::Density ( ) const
inlineoverridevirtual

Returns argon density at the temperature from Temperature()

Reimplemented from gar::detinfo::DetectorProperties.

Definition at line 196 of file DetectorPropertiesStandard.h.

196 { return Density(Temperature()); }
virtual double Temperature() const override
In kelvin.
virtual double Density() const override
Returns argon density at the temperature from Temperature()
double gar::detinfo::DetectorPropertiesStandard::DriftVelocity ( double  efield = 0.,
double  temperature = 0.,
bool  cmPerns = true 
) const
overridevirtual

cm/ns if true, otherwise cm/us

Implements gar::detinfo::DetectorProperties.

Definition at line 292 of file DetectorPropertiesStandard.cxx.

295  {
296 
297  // Efield should have units of kV/cm
298  // Temperature should have units of Kelvin
299 
300  // Default Efield, use internal value.
301  if(efield == 0.)
302  efield = Efield();
303 
304  //
305  //if(efield > 4.0)
306  // MF_LOG_WARNING("DetectorPropertiesStandard")
307  // << "DriftVelocity Warning! : E-field value of "
308  // << efield
309  // << " kV/cm is outside of range covered by drift"
310  // << " velocity parameterization. Returned value"
311  //<< " may not be correct";
312 
313  // Default temperature use internal value.
314  if(temperature == 0.)
315  temperature = Temperature();
316 
317  // read in from fcl parameter
318 
319  double vd = fDriftVelocity; // cm/us. For now just take it out of a fcl parameter. Calcualted with magboltz and it's a strong function of gas composition
320 
321  if(cmPerns) return vd * 1.e-3; // cm/ns
322 
323  return vd; // in cm/us
324  }
virtual double Temperature() const override
In kelvin.
virtual double Efield(unsigned int planegap=0) const override
kV/cm
double fDriftVelocity
centimeters / microsecond
virtual double gar::detinfo::DetectorPropertiesStandard::EffectivePixel ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 240 of file DetectorPropertiesStandard.h.

240 { return fECALP->EffectivePixel(); }
virtual double EffectivePixel() const =0
double gar::detinfo::DetectorPropertiesStandard::Efield ( unsigned int  planegap = 0) const
overridevirtual

kV/cm

Implements gar::detinfo::DetectorProperties.

Definition at line 178 of file DetectorPropertiesStandard.cxx.

179  {
180  if(planegap >= fEfield.size())
181  throw cet::exception("DetectorPropertiesStandard")
182  << "requesting Electric field in a plane gap that is not defined\n";
183 
184  return fEfield[planegap];
185  }
std::vector< double > fEfield
kV/cm (per inter-plane volume)
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
virtual double gar::detinfo::DetectorPropertiesStandard::ElectronLifetime ( ) const
inlineoverridevirtual

dQ/dX in electrons/cm, returns dE/dX in MeV/cm.

Implements gar::detinfo::DetectorProperties.

Definition at line 181 of file DetectorPropertiesStandard.h.

181 { return fElectronlifetime; } //< microseconds
virtual double gar::detinfo::DetectorPropertiesStandard::ElectronsToADC ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 227 of file DetectorPropertiesStandard.h.

227 { return fElectronsToADC; }
double fElectronsToADC
conversion factor for # of ionization electrons to 1 ADC count
double gar::detinfo::DetectorPropertiesStandard::Eloss ( double  mom,
double  mass,
double  tcut 
) const
overridevirtual

Restricted mean energy loss (dE/dx)

Parameters
mommomentum of incident particle [GeV/c]
massmass of incident particle [GeV/c^2]
tcutmaximum kinetic energy of delta rays [MeV]; 0 for unlimited
Returns
the restricted mean energy loss (dE/dx) in units of MeV/cm

Returned value is always positive. For unrestricted mean energy loss, set tcut = 0 (special case), or tcut large.

Based on Bethe-Bloch formula as contained in particle data book. Material parameters are from the configuration.

Implements gar::detinfo::DetectorProperties.

Definition at line 218 of file DetectorPropertiesStandard.cxx.

221  {
222  // Some constants.
223 
224  double K = 0.307075; // 4 pi N_A r_e^2 m_e c^2 (MeV cm^2/mol).
225  double me = 0.510998918; // Electron mass (MeV/c^2).
226 
227  // Calculate kinematic quantities.
228 
229  double bg = mom / mass; // beta*gamma.
230  double gamma = sqrt(1. + bg*bg); // gamma.
231  double beta = bg / gamma; // beta (velocity).
232  double mer = 0.001 * me / mass; // electron mass / mass of incident particle.
233  double tmax = 2.*me* bg*bg / (1. + 2.*gamma*mer + mer*mer); // Maximum delta ray energy (MeV).
234 
235  // Make sure tcut does not exceed tmax.
236 
237  if(tcut == 0. || tcut > tmax)
238  tcut = tmax;
239 
240  // Calculate density effect correction (delta).
241 
242  double x = std::log10(bg);
243  double delta = 0.;
244  if(x >= fSternheimerParameters.x0) {
245  delta = 2. * std::log(10.) * x - fSternheimerParameters.cbar;
246  if(x < fSternheimerParameters.x1)
248  }
249 
250  // Calculate stopping number.
251 
252  double B = 0.5 * std::log(2.*me*bg*bg*tcut / (1.e-12 * sqr(fGP->ExcitationEnergy())))
253  - 0.5 * beta*beta * (1. + tcut / tmax) - 0.5 * delta;
254 
255  // Don't let the stopping number become negative.
256 
257  if(B < 1.)
258  B = 1.;
259 
260  // Calculate dE/dx.
261 
262  double dedx = Density() * K*fGP->AtomicNumber()*B / (fGP->AtomicMass() * beta*beta);
263 
264  // Done.
265 
266  return dedx;
267  } // DetectorPropertiesStandard::Eloss()
SternheimerParameters_t fSternheimerParameters
Sternheimer parameters.
double beta(double KE, const simb::MCParticle *part)
virtual double ExcitationEnergy() const =0
Mean excitation energy of the gas (eV)
constexpr T pow(T x)
Definition: pow.h:72
virtual double AtomicNumber() const =0
Atomic number of the gas.
virtual double AtomicMass() const =0
Atomic mass of the gas (g/mol)
T sqr(T v)
const double e
virtual double Density() const override
Returns argon density at the temperature from Temperature()
double gamma(double KE, const simb::MCParticle *part)
list x
Definition: train.py:276
double gar::detinfo::DetectorPropertiesStandard::ElossVar ( double  mom,
double  mass 
) const
overridevirtual

Energy loss fluctuation ( $ \sigma_{E}^2 / x $)

Parameters
mommomentum of incident particle in [GeV/c]
Returns
energy loss fluctuation in MeV^2/cm

Based on Bichsel formula referred to but not given in pdg.

Implements gar::detinfo::DetectorProperties.

Definition at line 270 of file DetectorPropertiesStandard.cxx.

272  {
273  // Some constants.
274 
275  double K = 0.307075; // 4 pi N_A r_e^2 m_e c^2 (MeV cm^2/mol).
276  double me = 0.510998918; // Electron mass (MeV/c^2).
277 
278  // Calculate kinematic quantities.
279 
280  double bg = mom / mass; // beta*gamma.
281  double gamma2 = 1. + bg*bg; // gamma^2.
282  double beta2 = bg*bg / gamma2; // beta^2.
283 
284  // Calculate final result.
285 
286  double result = gamma2 * (1. - 0.5 * beta2) * me * (fGP->AtomicNumber() / fGP->AtomicMass()) * K * Density();
287  return result;
288  } // DetectorPropertiesStandard::ElossVar()
static QCString result
virtual double AtomicNumber() const =0
Atomic number of the gas.
virtual double AtomicMass() const =0
Atomic mass of the gas (g/mol)
virtual double Density() const override
Returns argon density at the temperature from Temperature()
virtual double gar::detinfo::DetectorPropertiesStandard::IntercalibrationFactor ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 243 of file DetectorPropertiesStandard.h.

243 { return fECALP->IntercalibrationFactor(); }
virtual double IntercalibrationFactor() const =0
virtual double gar::detinfo::DetectorPropertiesStandard::LightYield ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 241 of file DetectorPropertiesStandard.h.

241 { return fECALP->LightYield(); }
virtual double LightYield() const =0
virtual double gar::detinfo::DetectorPropertiesStandard::MeVtoMIP ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 246 of file DetectorPropertiesStandard.h.

246 { return fECALP->MeVtoMIP(); }
virtual double MeVtoMIP() const =0
virtual double gar::detinfo::DetectorPropertiesStandard::NoisePx ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 247 of file DetectorPropertiesStandard.h.

247 { return fECALP->NoisePx(); }
virtual double NoisePx() const =0
virtual unsigned int gar::detinfo::DetectorPropertiesStandard::NumberTimeSamples ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 228 of file DetectorPropertiesStandard.h.

228 { return fNumberTimeSamples; }
unsigned int fNumberTimeSamples
number of clock ticks per event (= readout window)
virtual double gar::detinfo::DetectorPropertiesStandard::SamplingRate ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 226 of file DetectorPropertiesStandard.h.

226 { return fTPCClock.TickPeriod() * 1.e3; }
double TickPeriod() const
A single tick period in nano-second, frequency is in MHz.
Definition: ElecClock.h:125
detinfo::ElecClock fTPCClock
TPC electronics clock.
void gar::detinfo::DetectorPropertiesStandard::SetDetectorClocks ( const detinfo::DetectorClocks clks)
inline

Definition at line 168 of file DetectorPropertiesStandard.h.

168 { fClocks = clks; }
void gar::detinfo::DetectorPropertiesStandard::SetECALProperties ( const detinfo::ECALProperties ecalp)
inline

Definition at line 167 of file DetectorPropertiesStandard.h.

167 { fECALP = ecalp; }
void gar::detinfo::DetectorPropertiesStandard::SetGArProperties ( const detinfo::GArProperties gp)
inline

Definition at line 166 of file DetectorPropertiesStandard.h.

166 { fGP = gp; }
void gar::detinfo::DetectorPropertiesStandard::SetGeometry ( const geo::GeometryCore g)
inline

Definition at line 165 of file DetectorPropertiesStandard.h.

165 { fGeo = g; }
static constexpr double g
Definition: Units.h:144
void gar::detinfo::DetectorPropertiesStandard::SetNumberTimeSamples ( unsigned int  nsamp)
inline

Definition at line 170 of file DetectorPropertiesStandard.h.

170 { fNumberTimeSamples=nsamp;}
unsigned int fNumberTimeSamples
number of clock ticks per event (= readout window)
void gar::detinfo::DetectorPropertiesStandard::Setup ( providers_type  providers)

Sets all the providers at once.

Parameters
providersthe pack of service providers we depend on

Example:

gar::DetectorPropertiesStandard::providers_type providers;
providers.set(gar::providerFrom<geo::GeometryGAr>());
providers.set(gar::providerFrom<detinfo::GArPropertiesService>());
providers.set(gar::providerFrom<detinfo::DetectorClocksService>());
detprop->Setup(providers);

Definition at line 167 of file DetectorPropertiesStandard.cxx.

167  {
168 
169  SetGeometry(providers.get<geo::GeometryCore>());
170  SetGArProperties(providers.get<detinfo::GArProperties>());
171  SetECALProperties(providers.get<detinfo::ECALProperties>());
172  SetDetectorClocks(providers.get<detinfo::DetectorClocks>());
173 
174  } // DetectorPropertiesStandard::Setup()
void SetDetectorClocks(const detinfo::DetectorClocks *clks)
void SetGArProperties(const detinfo::GArProperties *gp)
Description of geometry of one entire detector.
Class used for the conversion of times between different formats and references.
void SetECALProperties(const detinfo::ECALProperties *ecalp)
virtual double gar::detinfo::DetectorPropertiesStandard::SiPMGain ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 242 of file DetectorPropertiesStandard.h.

242 { return fECALP->SiPMGain(); }
virtual double SiPMGain() const =0
virtual double gar::detinfo::DetectorPropertiesStandard::Temperature ( ) const
inlineoverridevirtual

In kelvin.

Implements gar::detinfo::DetectorProperties.

Definition at line 199 of file DetectorPropertiesStandard.h.

virtual double gar::detinfo::DetectorPropertiesStandard::TimeResolution ( ) const
inlineoverridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 245 of file DetectorPropertiesStandard.h.

245 { return fECALP->TimeResolution(); }
virtual double TimeResolution() const =0
int gar::detinfo::DetectorPropertiesStandard::TriggerOffset ( ) const
overridevirtual

Implements gar::detinfo::DetectorProperties.

Definition at line 327 of file DetectorPropertiesStandard.cxx.

328  {
329  return fTPCClock.Ticks(fClocks->TriggerOffsetTPC() * -1.);
330  }
virtual double TriggerOffsetTPC() const =0
int Ticks() const
of Ticks
Definition: ElecClock.h:95
detinfo::ElecClock fTPCClock
TPC electronics clock.
bool gar::detinfo::DetectorPropertiesStandard::Update ( uint64_t  ts)

Definition at line 91 of file DetectorPropertiesStandard.cxx.

92  {
93 
95  return true;
96  }
bool gar::detinfo::DetectorPropertiesStandard::UpdateClocks ( const detinfo::DetectorClocks clks)

Definition at line 99 of file DetectorPropertiesStandard.cxx.

100  {
101  fClocks = clks;
102 
105  return true;
106  }
virtual const ElecClock & TPCClock() const =0
Borrow a const TPC clock with time set to Trigger time [ns].
detinfo::ElecClock fTPCClock
TPC electronics clock.
void gar::detinfo::DetectorPropertiesStandard::ValidateAndConfigure ( fhicl::ParameterSet const &  p,
std::set< std::string > const &  ignore_params = {} 
)

Configures the provider, first validating the configuration.

Parameters
pconfiguration parameter set
ignore_paramsparameters to be ignored (optional)

This method will validate the parameter set (except for the parameters it's explicitly told to ignore) and extract the useful information out of it.

Definition at line 159 of file DetectorPropertiesStandard.cxx.

161  {
162  Configure(ValidateConfiguration(p, ignore_params));
163  } // ValidateAndConfigure()
p
Definition: test.py:223
Configuration_t ValidateConfiguration(fhicl::ParameterSet const &p, std::set< std::string > const &ignore_params={})
Validates the specified configuration.
void Configure(Configuration_t const &config)
Extracts the relevant configuration from the specified object.
DetectorPropertiesStandard::Configuration_t gar::detinfo::DetectorPropertiesStandard::ValidateConfiguration ( fhicl::ParameterSet const &  p,
std::set< std::string > const &  ignore_params = {} 
)

Validates the specified configuration.

Parameters
pconfiguration parameter set
ignore_paramsparameters to be ignored (optional)
Returns
a parsed configuration object
See also
ValidateAndConfigure(), Configure()

This method will validate the parameter set (except for the parameters it's explicitly told to ignore) and it returns an object ready to be used with Configure().

Definition at line 143 of file DetectorPropertiesStandard.cxx.

145  {
146  std::set<std::string> ignorable_keys = gar::IgnorableProviderConfigKeys();
147  ignorable_keys.insert(ignore_params.begin(), ignore_params.end());
148 
149  gar::debug::printBacktrace(mf::LogDebug("DetectorPropertiesStandard"), 15);
150 
151  // parses and validates the parameter set:
152  fhicl::Table<Configuration_t> config_table { p, ignorable_keys };
153 
154  return std::move(config_table());
155 
156  } // DetectorPropertiesStandard::ValidateConfiguration()
def move(depos, offset)
Definition: depos.py:107
std::set< std::string > const & IgnorableProviderConfigKeys()
Returns a list of configuration keys that providers should ignore.
Definition: ProviderUtil.h:35
p
Definition: test.py:223
void printBacktrace(Stream &&out, unsigned int maxLines=5, std::string indent=" ", CallInfoPrinter::opt const *options=nullptr)
Prints the full backtrace into a stream.
Definition: DebugUtils.h:249

Member Data Documentation

const detinfo::DetectorClocks* gar::detinfo::DetectorPropertiesStandard::fClocks
protected

Definition at line 271 of file DetectorPropertiesStandard.h.

double gar::detinfo::DetectorPropertiesStandard::fDriftVelocity
protected

centimeters / microsecond

Definition at line 277 of file DetectorPropertiesStandard.h.

const detinfo::ECALProperties* gar::detinfo::DetectorPropertiesStandard::fECALP
protected

Definition at line 270 of file DetectorPropertiesStandard.h.

std::vector< double > gar::detinfo::DetectorPropertiesStandard::fEfield
protected

kV/cm (per inter-plane volume)

Definition at line 274 of file DetectorPropertiesStandard.h.

double gar::detinfo::DetectorPropertiesStandard::fElectronlifetime
protected

microseconds

Definition at line 275 of file DetectorPropertiesStandard.h.

double gar::detinfo::DetectorPropertiesStandard::fElectronsToADC
protected

conversion factor for # of ionization electrons to 1 ADC count

Definition at line 279 of file DetectorPropertiesStandard.h.

const geo::GeometryCore* gar::detinfo::DetectorPropertiesStandard::fGeo
protected

Definition at line 272 of file DetectorPropertiesStandard.h.

const detinfo::GArProperties* gar::detinfo::DetectorPropertiesStandard::fGP
protected

Definition at line 269 of file DetectorPropertiesStandard.h.

unsigned int gar::detinfo::DetectorPropertiesStandard::fNumberTimeSamples
protected

number of clock ticks per event (= readout window)

Definition at line 280 of file DetectorPropertiesStandard.h.

double gar::detinfo::DetectorPropertiesStandard::fSamplingRate
protected

in ns

Definition at line 278 of file DetectorPropertiesStandard.h.

SternheimerParameters_t gar::detinfo::DetectorPropertiesStandard::fSternheimerParameters
protected

Sternheimer parameters.

Definition at line 282 of file DetectorPropertiesStandard.h.

double gar::detinfo::DetectorPropertiesStandard::fTemperature
protected

kelvin

Definition at line 276 of file DetectorPropertiesStandard.h.

detinfo::ElecClock gar::detinfo::DetectorPropertiesStandard::fTPCClock
protected

TPC electronics clock.

Definition at line 286 of file DetectorPropertiesStandard.h.

double gar::detinfo::DetectorPropertiesStandard::fXTicksCoefficient
protected

Parameters for x<–>ticks.

Definition at line 284 of file DetectorPropertiesStandard.h.


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