Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Attributes | Friends | List of all members
JDETECTOR::JPMTAnalogueSignalProcessor Struct Reference

PMT analogue signal processor. More...

#include <JPMTAnalogueSignalProcessor.hh>

Inheritance diagram for JDETECTOR::JPMTAnalogueSignalProcessor:
JDETECTOR::JPMTSignalProcessorInterface JDETECTOR::JPMTParameters

Public Types

enum  JThresholdDomains { BELOW_THRESHOLDBAND = -1, IN_THRESHOLDBAND = 0, ABOVE_THRESHOLD = 2 }
 Threshold domain specifiers. More...
 

Public Member Functions

 JPMTAnalogueSignalProcessor (const JPMTParameters &parameters=JPMTParameters())
 Constructor. More...
 
void configure ()
 Configure internal parameters. More...
 
double getDecayTime () const
 Get decay time. More...
 
double getT1 () const
 Get time at transition point from Gaussian to exponential. More...
 
double getY1 () const
 Get amplitude at transition point from Gaussian to exponential. More...
 
double getStartOfLinearisation () const
 Get transition point from a model-dependent to linear relation between time-over-threshold and number of photo-electrons. More...
 
double getAmplitude (const double t1_ns) const
 Get amplitude at given time for a one photo-electron pulse. More...
 
double getRiseTime (const double npe, const double th) const
 Get time to pass from threshold to top of analogue pulse. More...
 
double getDecayTime (const double npe, const double th) const
 Get time to pass from top of analogue pulse to threshold. More...
 
double applySaturation (const double T) const
 Get time-over-threshold with saturation. More...
 
double removeSaturation (const double tot_ns) const
 Get time-over-threshold without saturation. More...
 
double getDerivativeOfSaturation (const double T) const
 Get derivative of saturation factor. More...
 
double getGainSpread (int NPE) const
 Get gain spread for given number of photo-electrons. More...
 
virtual double getChargeProbability (const double npe, const int NPE) const
 Get probability density for given charge. More...
 
double getIntegralOfChargeProbability (const double xmin, const double xmax, const int NPE) const
 Get integral of probability. More...
 
double getIntegralOfChargeProbability (JThresholdDomains domain, const int NPE) const
 Get integral of probability in specific threshold domain. More...
 
void setPMTParameters (const JPMTParameters &parameters)
 Set PMT parameters. More...
 
virtual bool applyQE () const
 Apply relative QE. More...
 
virtual double getRandomTime (const double t_ns) const
 Get randomised time according transition time distribution. More...
 
virtual bool compare (const JPhotoElectron &first, const JPhotoElectron &second) const
 Compare (arrival times of) photo-electrons. More...
 
virtual double getRandomCharge (const int NPE) const
 Get randomised charge according to gain and gain spread. More...
 
virtual double getNPE (const double tot_ns, const double eps=1.0e-3) const
 Get number of photo-electrons. More...
 
virtual JThresholdDomains applyThreshold (const double npe) const
 Apply threshold. More...
 
virtual double getRiseTime (const double npe) const
 Get time to reach threshold. More...
 
virtual double getTimeOverThreshold (const double npe) const
 Get time-over-threshold (ToT). More...
 
virtual double getDerivative (const double npe) const
 Get derivative of number of photo-electrons to time-over-threshold. More...
 
double getTimeOverThresholdProbability (const double tot_ns, const int NPE) const
 Get probability of having a pulse with specific time-over-threshold. More...
 
double getIntegralOfTimeOverThresholdProbability (const double Tmin, const double Tmax, const int NPE) const
 Get cumulative probability of time-over-threshold distribution. More...
 
virtual double getSurvivalProbability (const int NPE) const
 Probability that a hit survives the simulation of the PMT. More...
 
void operator() (const JCalibration &calibration, const JPMTData< JPMTSignal > &input, JPMTData< JPMTPulse > &output) const
 Process hits. More...
 
virtual void merge (JPMTData< JPMTPulse > &data) const
 Merging of PMT hits. More...
 
const JPMTParametersgetPMTParameters () const
 Get PMT parameters. More...
 
bool is_valid () const
 Check validity of PMT parameters. More...
 
JProperties getProperties (const JEquationParameters &equation=JPMTParameters::getEquationParameters())
 Get properties of this class. More...
 
JProperties getProperties (const JEquationParameters &equation=JPMTParameters::getEquationParameters()) const
 Get properties of this class. More...
 

Static Public Member Functions

static double getMaximalRiseTime (const double th)
 Get maximal rise time for given threshold. More...
 
static double getTH0 ()
 Get lower threshold for rise time evaluation. More...
 
static double getTH1 ()
 Get upper threshold for rise time evaluation. More...
 
static double getTmin ()
 Get two photo-electron resolution for time-over-threshold. More...
 
static double getQmin ()
 Get minimal width of charge distribution. More...
 
static JEquationParametersgetEquationParameters ()
 Get equation parameters. More...
 
static void setEquationParameters (const JEquationParameters &equation)
 Set equation parameters. More...
 

Public Attributes

double QE
 relative quantum efficiency More...
 
double gain
 gain [unit] More...
 
double gainSpread
 gain spread [unit] More...
 
double riseTime_ns
 rise time of analogue pulse [ns] More...
 
double TTS_ns
 transition time spread [ns] More...
 
double threshold
 threshold [npe] More...
 
double PunderAmplified
 probability of underamplified hit More...
 
double thresholdBand
 threshold-band [npe] More...
 
double mean_ns
 mean time-over-threshold of threshold-band hits [ns] More...
 
double sigma_ns
 time-over-threshold standard deviation of threshold-band hits [ns] More...
 
double slope
 slope [ns/npe] More...
 
double saturation
 saturation [ns] More...
 
bool slewing
 time slewing of analogue signal More...
 

Protected Attributes

double decayTime_ns
 decay time [ns] More...
 
double t1
 time at match point [ns] More...
 
double y1
 amplitude at match point [npe] More...
 
double x1
 Transition point from a logarithmic to a linear relation between time-over-threshold and number of photo-electrons. More...
 
double MAX_CHARGE
 Maximum charge [npe]. More...
 

Friends

std::istream & operator>> (std::istream &in, JPMTAnalogueSignalProcessor &object)
 Read PMT signal from input. More...
 

Detailed Description

PMT analogue signal processor.

This class provides for an implementation of the JDETECTOR::JPMTSignalProcessorInterface using a specific model for the analogue pulse of the PMT.
In this, the leading edge of the analogue pulse from the PMT is assumed to be a Gaussian and the tail an exponential.
The width of the Gaussian is referred to as the rise time and the inverse slope of the exponential to the decay time.
The two functions are matched at a point where the values and first derivatives are identical.
Note that the decay time is related to the rise time via the specification JDETECTOR::TIME_OVER_THRESHOLD_NS.

The charge distribution is assumed to be a Gaussian which is centered at the specified gain and truncated by the specified threshold.

The transition times are generated according the specified spread as follows.

Definition at line 48 of file JPMTAnalogueSignalProcessor.hh.

Member Enumeration Documentation

Threshold domain specifiers.

Enumerator
BELOW_THRESHOLDBAND 

below threshold band

IN_THRESHOLDBAND 

inside threshold band

ABOVE_THRESHOLD 

above threshold band

Definition at line 89 of file JPMTSignalProcessorInterface.hh.

89  {
90  BELOW_THRESHOLDBAND = -1, //!< below threshold band
91  IN_THRESHOLDBAND = 0, //!< inside threshold band
92  ABOVE_THRESHOLD = 2 //!< above threshold band
93  };

Constructor & Destructor Documentation

JDETECTOR::JPMTAnalogueSignalProcessor::JPMTAnalogueSignalProcessor ( const JPMTParameters parameters = JPMTParameters())
inline

Constructor.

Parameters
parametersPMT parameters

Definition at line 57 of file JPMTAnalogueSignalProcessor.hh.

57  :
60  decayTime_ns(0.0),
61  t1(0.0),
62  y1(0.0),
63  x1(std::numeric_limits<double>::max())
64  {
65  configure();
66  }
double x1
Transition point from a logarithmic to a linear relation between time-over-threshold and number of ph...
JPMTSignalProcessorInterface(const double max_charge=100.0)
Default constructor.
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
void configure()
Configure internal parameters.
JPMTParameters()
Default constructor.
double y1
amplitude at match point [npe]

Member Function Documentation

void JDETECTOR::JPMTAnalogueSignalProcessor::configure ( )
inline

Configure internal parameters.

This method provides the implementations for

  • matching of the leading edge of the analogue pulse (Gaussian) and the tail (exponential); and
  • determination of number of photo-electrons above which the time-over-threshold linearly depends on the number of photo-electrons (apart from saturation).

Note that this method will throw an error if the value of the rise time (i.e. width of the Gaussian) is too large with respect to the specification JDETECTOR::TIME_OVER_THRESHOLD_NS.

Definition at line 80 of file JPMTAnalogueSignalProcessor.hh.

81  {
82  static const int N = 100;
83  static const double precision = 1.0e-3;
84 
85  // check thresholdband
86 
87  if (threshold - thresholdBand < getTH0()) {
88  THROW(JValueOutOfRange, "JPMTAnalogueSignalProcessor::configure(): Invalid thresholdband [npe] " << thresholdBand);
89  }
90 
91  // check rise time
92 
94  THROW(JValueOutOfRange, "JPMTAnalogueSignalProcessor::configure(): Invalid rise time [ns] " << riseTime_ns);
95  }
96 
97  // decay time
98 
99  const double y = -log(threshold);
100 
101  const double a = y;
102  const double b = riseTime_ns * sqrt(2.0*y) - TIME_OVER_THRESHOLD_NS;
103  const double c = 0.5*riseTime_ns*riseTime_ns;
104 
105  const double Q = b*b - 4.0*a*c;
106 
107  if (Q > 0.0)
108  decayTime_ns = (-b + sqrt(Q)) / (2.0*a);
109  else
110  decayTime_ns = -b / (2.0*a);
111 
112  // fix matching of Gaussian and exponential
113 
114  const double x = riseTime_ns / decayTime_ns;
115 
116  t1 = riseTime_ns*x;
117  y1 = exp(-0.5*x*x);
118 
119  // determine transition point to linear dependence of time-over-threshold as a function of number of photo-electrons
120 
121  const double v = slope * getDerivativeOfSaturation(getTimeOverThreshold(1.0));
122 
123  x1 = std::numeric_limits<double>::max(); // disable linearisation
124 
125  double xmin = 1.0;
126  double xmax = 1.0 / (getDerivative(1.0) * slope);
127 
128  for (int i = 0; i != N; ++i) {
129 
130  const double x = 0.5 * (xmin + xmax);
131  const double u = getDerivative(x) * v;
132 
133  if (fabs(1.0 - u) < precision) {
134  break;
135  }
136 
137  if (u < 1.0)
138  xmin = x;
139  else
140  xmax = x;
141  }
142 
143  x1 = 0.5 * (xmin + xmax);
144  }
double x1
Transition point from a logarithmic to a linear relation between time-over-threshold and number of ph...
double getDerivativeOfSaturation(const double T) const
Get derivative of saturation factor.
const double TIME_OVER_THRESHOLD_NS
Specification for time-over-threshold corresponding to a one photo-electron pulse.
virtual double getTimeOverThreshold(const double npe) const
Get time-over-threshold (ToT).
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
double threshold
threshold [npe]
static double getTH0()
Get lower threshold for rise time evaluation.
fi JEventTimesliceWriter a
double thresholdBand
threshold-band [npe]
double riseTime_ns
rise time of analogue pulse [ns]
static double getMaximalRiseTime(const double th)
Get maximal rise time for given threshold.
data_type v[N+1][M+1]
Definition: JPolint.hh:707
double u[N+1]
Definition: JPolint.hh:706
virtual double getDerivative(const double npe) const
Get derivative of number of photo-electrons to time-over-threshold.
then usage $script[input file[working directory[option]]] nWhere option can be N
Definition: JMuonPostfit.sh:37
double y1
amplitude at match point [npe]
double slope
slope [ns/npe]
then set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable OUTPUT_FILE histogram.root JHistogram1D -o $WORKDIR/$OUTPUT_FILE -F "$FORMULA" -
double JDETECTOR::JPMTAnalogueSignalProcessor::getDecayTime ( ) const
inline

Get decay time.

Returns
decay time [ns]

Definition at line 152 of file JPMTAnalogueSignalProcessor.hh.

153  {
154  return decayTime_ns;
155  }
double JDETECTOR::JPMTAnalogueSignalProcessor::getT1 ( ) const
inline

Get time at transition point from Gaussian to exponential.

Returns
time [ns]

Definition at line 163 of file JPMTAnalogueSignalProcessor.hh.

164  {
165  return t1;
166  }
double JDETECTOR::JPMTAnalogueSignalProcessor::getY1 ( ) const
inline

Get amplitude at transition point from Gaussian to exponential.

Returns
amplitude [npe]

Definition at line 174 of file JPMTAnalogueSignalProcessor.hh.

175  {
176  return y1;
177  }
double y1
amplitude at match point [npe]
double JDETECTOR::JPMTAnalogueSignalProcessor::getStartOfLinearisation ( ) const
inline

Get transition point from a model-dependent to linear relation between time-over-threshold and number of photo-electrons.

Returns
number of photo-electrons [npe]

Definition at line 185 of file JPMTAnalogueSignalProcessor.hh.

186  {
187  return x1;
188  }
double x1
Transition point from a logarithmic to a linear relation between time-over-threshold and number of ph...
double JDETECTOR::JPMTAnalogueSignalProcessor::getAmplitude ( const double  t1_ns) const
inline

Get amplitude at given time for a one photo-electron pulse.

Parameters
t1_nstime [ns]
Returns
amplitude [npe]

Definition at line 197 of file JPMTAnalogueSignalProcessor.hh.

198  {
199  if (t1_ns < t1) {
200 
201  const double x = t1_ns / riseTime_ns;
202 
203  return exp(-0.5*x*x); // Gaussian
204 
205  } else {
206 
207  const double x = t1_ns / decayTime_ns;
208 
209  return exp(-x) / y1; // exponential
210  }
211  }
double riseTime_ns
rise time of analogue pulse [ns]
double y1
amplitude at match point [npe]
then set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable OUTPUT_FILE histogram.root JHistogram1D -o $WORKDIR/$OUTPUT_FILE -F "$FORMULA" -
double JDETECTOR::JPMTAnalogueSignalProcessor::getRiseTime ( const double  npe,
const double  th 
) const
inline

Get time to pass from threshold to top of analogue pulse.


In this, the leading edge of the analogue pulse is assumed to be Gaussian.

Parameters
npenumber of photo-electrons
ththreshold [npe]
Returns
time [ns]

Definition at line 222 of file JPMTAnalogueSignalProcessor.hh.

223  {
224  return riseTime_ns * sqrt(2.0*log(npe/th)); // Gaussian
225  }
double riseTime_ns
rise time of analogue pulse [ns]
double JDETECTOR::JPMTAnalogueSignalProcessor::getDecayTime ( const double  npe,
const double  th 
) const
inline

Get time to pass from top of analogue pulse to threshold.


In this, the trailing edge of the analogue pulse is assumed to be exponential.

Parameters
npenumber of photo-electrons
ththreshold [npe]
Returns
time [ns]

Definition at line 236 of file JPMTAnalogueSignalProcessor.hh.

237  {
238  if (npe*y1 > th)
239  return decayTime_ns * (log(npe/th) - log(y1)); // exponential
240  else
241  return riseTime_ns * sqrt(2.0*log(npe/th)); // Gaussian
242  }
double riseTime_ns
rise time of analogue pulse [ns]
double y1
amplitude at match point [npe]
static double JDETECTOR::JPMTAnalogueSignalProcessor::getMaximalRiseTime ( const double  th)
inlinestatic

Get maximal rise time for given threshold.

Note that the rise time is entirely constrained by the specification JDETECTOR::TIME_OVER_THRESHOLD_NS.

Parameters
ththreshold [npe]
Returns
rise time [ns]

Definition at line 253 of file JPMTAnalogueSignalProcessor.hh.

254  {
255  if (th > 0.0 && th < 1.0)
256  return 0.5 * TIME_OVER_THRESHOLD_NS / sqrt(-2.0*log(th));
257  else
258  THROW(JValueOutOfRange, "JPMTAnalogueSignalProcessor::getMaximalRiseTime(): Invalid threshold " << th);
259  }
const double TIME_OVER_THRESHOLD_NS
Specification for time-over-threshold corresponding to a one photo-electron pulse.
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
double JDETECTOR::JPMTAnalogueSignalProcessor::applySaturation ( const double  T) const
inline

Get time-over-threshold with saturation.

Parameters
Ttime-over-threshold without saturation
Returns
time-over-threshold with saturation

Definition at line 268 of file JPMTAnalogueSignalProcessor.hh.

269  {
270  return saturation / sqrt(T*T + saturation*saturation) * T;
271  }
do set_variable OUTPUT_DIRECTORY $WORKDIR T
double saturation
saturation [ns]
double JDETECTOR::JPMTAnalogueSignalProcessor::removeSaturation ( const double  tot_ns) const
inline

Get time-over-threshold without saturation.

Parameters
tot_nstime-over-threshold with saturation
Returns
time-over-threshold without saturation

Definition at line 280 of file JPMTAnalogueSignalProcessor.hh.

281  {
282  return saturation / sqrt(saturation*saturation - tot_ns*tot_ns) * tot_ns;
283  }
double saturation
saturation [ns]
double JDETECTOR::JPMTAnalogueSignalProcessor::getDerivativeOfSaturation ( const double  T) const
inline

Get derivative of saturation factor.

Parameters
Ttime-over-threshold without saturation
Returns
derivative of saturation factor

Definition at line 292 of file JPMTAnalogueSignalProcessor.hh.

293  {
295  }
do set_variable OUTPUT_DIRECTORY $WORKDIR T
double saturation
saturation [ns]
double JDETECTOR::JPMTAnalogueSignalProcessor::getGainSpread ( int  NPE) const
inline

Get gain spread for given number of photo-electrons.

Parameters
NPEnumber of photo-electrons
Returns
gain spread

Definition at line 304 of file JPMTAnalogueSignalProcessor.hh.

305  {
306  return sqrt((double) NPE * gain) * gainSpread;
307  }
double gain
gain [unit]
double gainSpread
gain spread [unit]
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getChargeProbability ( const double  npe,
const int  NPE 
) const
inlinevirtual

Get probability density for given charge.

Parameters
npeobserved number of photo-electrons
NPEtrue number of photo-electrons
Returns
probability [npe^-1]

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 317 of file JPMTAnalogueSignalProcessor.hh.

318  {
319  if (applyThreshold(npe) > BELOW_THRESHOLDBAND) {
320 
321  const double f = gainSpread * gainSpread;
322 
323  double norm = 0.0;
324  double prob = 0.0;
325  double weight = (PunderAmplified > 0.0 ? pow(PunderAmplified, NPE) : 1.0);
326 
327  for (int k = (PunderAmplified > 0.0 ? 0 : NPE); k <= NPE; ++k) {
328 
329  const double mu = ((NPE-k) * f * gain) + (k * gain);
330  const double sigma = sqrt((NPE-k) * f + k) * getGainSpread(1);
331 
332  norm += weight * (0.5 * TMath::Erfc(((threshold-thresholdBand) - mu) / sqrt(2.0) / sigma));
333  prob += weight * TMath::Gaus(npe, mu, sigma, kTRUE);
334 
335  weight *= ( (NPE - k) / ((double) (k+1)) *
336  (1.0 - PunderAmplified) / PunderAmplified );
337  }
338 
339  return prob / norm;
340  }
341 
342  return 0.0;
343  }
*fatal Wrong number of arguments esac JCookie sh JRuns D $DETECTOR d sort n k
Definition: JRunrange.sh:16
double gain
gain [unit]
double gainSpread
gain spread [unit]
double threshold
threshold [npe]
double PunderAmplified
probability of underamplified hit
then JPizza f
Definition: JPizza.sh:46
double thresholdBand
threshold-band [npe]
virtual JThresholdDomains applyThreshold(const double npe) const
Apply threshold.
double getGainSpread(int NPE) const
Get gain spread for given number of photo-electrons.
std::vector< double > weight
Definition: JAlgorithm.hh:428
double JDETECTOR::JPMTAnalogueSignalProcessor::getIntegralOfChargeProbability ( const double  xmin,
const double  xmax,
const int  NPE 
) const
inline

Get integral of probability.

Parameters
xminminimum number of photo-electrons
xmaxmaximum number of photo-electrons
NPEtrue number of photo-electrons
Returns
probability

Definition at line 354 of file JPMTAnalogueSignalProcessor.hh.

355  {
356  double zmin = xmin;
357  double zmax = xmax;
358 
359  const double th = threshold - thresholdBand;
360  const double f = gainSpread * gainSpread;
361 
362  if (zmin < th) { zmin = th; }
363  if (zmax < th) { zmax = th; }
364 
365  double norm = 0.0;
366  double cumulP = 0.0;
367  double weight = (PunderAmplified > 0.0 ? pow(PunderAmplified, NPE) : 1.0);
368 
369  for (int k = (PunderAmplified > 0.0 ? 0 : NPE); k <= NPE; ++k) {
370 
371  const double mu = ((NPE-k) * f * gain) + (k * gain);
372  const double sigma = sqrt((NPE-k) * f + k) * getGainSpread(1);
373 
374  norm += weight * (0.5 * TMath::Erfc((th - mu) / sqrt(2.0) / sigma));
375  cumulP += weight * (0.5 * TMath::Erfc((zmin - mu) / sqrt(2.0) / sigma) -
376  0.5 * TMath::Erfc((zmax - mu) / sqrt(2.0) / sigma));
377 
378  weight *= ( (NPE - k) / ((double) (k+1)) *
379  (1.0 - PunderAmplified) / PunderAmplified);
380  }
381 
382  return cumulP / norm;
383  }
*fatal Wrong number of arguments esac JCookie sh JRuns D $DETECTOR d sort n k
Definition: JRunrange.sh:16
double gain
gain [unit]
double gainSpread
gain spread [unit]
double threshold
threshold [npe]
double PunderAmplified
probability of underamplified hit
then JPizza f
Definition: JPizza.sh:46
double thresholdBand
threshold-band [npe]
double getGainSpread(int NPE) const
Get gain spread for given number of photo-electrons.
std::vector< double > weight
Definition: JAlgorithm.hh:428
double JDETECTOR::JPMTAnalogueSignalProcessor::getIntegralOfChargeProbability ( JThresholdDomains  domain,
const int  NPE 
) const
inline

Get integral of probability in specific threshold domain.

Parameters
domainthreshold domain
NPEtrue number of photo-electrons
Returns
probability

Definition at line 393 of file JPMTAnalogueSignalProcessor.hh.

394  {
395  switch (domain) {
396  case ABOVE_THRESHOLD:
398 
399  case IN_THRESHOLDBAND:
401 
402  default:
403  return 0.0;
404  }
405  }
double getIntegralOfChargeProbability(const double xmin, const double xmax, const int NPE) const
Get integral of probability.
double threshold
threshold [npe]
double thresholdBand
threshold-band [npe]
void JDETECTOR::JPMTAnalogueSignalProcessor::setPMTParameters ( const JPMTParameters parameters)
inline

Set PMT parameters.

Parameters
parametersPMT parameters

Definition at line 413 of file JPMTAnalogueSignalProcessor.hh.

414  {
415  static_cast<JPMTParameters&>(*this).setPMTParameters(parameters);
416 
417  configure();
418  }
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
void configure()
Configure internal parameters.
JPMTParameters()
Default constructor.
virtual bool JDETECTOR::JPMTAnalogueSignalProcessor::applyQE ( ) const
inlinevirtual

Apply relative QE.

Returns
true if accepted; false if rejected

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 443 of file JPMTAnalogueSignalProcessor.hh.

444  {
445  if (QE <= 0.0)
446  return false;
447  else if (QE < 1.0)
448  return gRandom->Rndm() < QE;
449  else
450  return true;
451  }
double QE
relative quantum efficiency
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getRandomTime ( const double  t_ns) const
inlinevirtual

Get randomised time according transition time distribution.

Parameters
t_nstime [ns]
Returns
time [ns]

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 460 of file JPMTAnalogueSignalProcessor.hh.

461  {
462  if (TTS_ns < 0.0)
463  return t_ns + getTransitionTime(gRandom->Rndm());
464  else
465  return gRandom->Gaus(t_ns, TTS_ns);
466  }
static const JTransitionTimeGenerator getTransitionTime
Function object to generate transition time.
double TTS_ns
transition time spread [ns]
virtual bool JDETECTOR::JPMTAnalogueSignalProcessor::compare ( const JPhotoElectron first,
const JPhotoElectron second 
) const
inlinevirtual

Compare (arrival times of) photo-electrons.

Parameters
firstfirst photo-electron
secondsecond photo-electron
Returns
true if arrival times of photo-electrons are within rise time of analogue pulses; else false

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 476 of file JPMTAnalogueSignalProcessor.hh.

477  {
478  return second.t_ns < first.t_ns + riseTime_ns;
479  }
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
double riseTime_ns
rise time of analogue pulse [ns]
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getRandomCharge ( const int  NPE) const
inlinevirtual

Get randomised charge according to gain and gain spread.

Parameters
NPEnumber of photo-electrons
Returns
number of photo-electrons

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 488 of file JPMTAnalogueSignalProcessor.hh.

489  {
490  double q;
491  do {
492 
493  // Determine which contribution to sample from
494  int k = NPE;
495  if (PunderAmplified > 0.0) {
496 
497  const double X = gRandom->Uniform();
498  double weight = pow(1.0 - PunderAmplified, NPE);
499 
500  for (double sum_p = 0.0; k > 0; --k) {
501 
502  sum_p += weight;
503  if (sum_p > X) { break; }
504 
505  weight *= ( k / ((double) (NPE - (k-1))) *
506  PunderAmplified / (1.0 - PunderAmplified) );
507  }
508  }
509 
510  // Sample from chosen contribution
511  const double f = gainSpread * gainSpread;
512  const double mu = ((NPE-k) * f * gain) + (k * gain);
513  const double sigma = sqrt((NPE-k) * f + k) * getGainSpread(1);
514 
515  q = gRandom->Gaus(mu,sigma);
516 
517  } while (q < 0.0);
518 
519  return q;
520  }
*fatal Wrong number of arguments esac JCookie sh JRuns D $DETECTOR d sort n k
Definition: JRunrange.sh:16
double gain
gain [unit]
double gainSpread
gain spread [unit]
double PunderAmplified
probability of underamplified hit
then JPizza f
Definition: JPizza.sh:46
double getGainSpread(int NPE) const
Get gain spread for given number of photo-electrons.
std::vector< double > weight
Definition: JAlgorithm.hh:428
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getNPE ( const double  tot_ns,
const double  eps = 1.0e-3 
) const
inlinevirtual

Get number of photo-electrons.

Parameters
tot_nstime over threshold [ns]
epsprecision
Returns
number of photo-electrons

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 530 of file JPMTAnalogueSignalProcessor.hh.

532  {
533  const double t = removeSaturation(tot_ns);
534 
535  double npe;
536  if (t <= 2*getRiseTime(threshold/y1,threshold)) { // Gaussian + Gaussian
537 
538  npe = threshold * exp(t*t/riseTime_ns/riseTime_ns/8.0);
539 
540  } else if (t <= (getRiseTime(getStartOfLinearisation(), threshold) +
542  // Gaussian + Exponential
543 
544  const double A = decayTime_ns;
545  const double B = sqrt(2.0) * riseTime_ns;
546  const double C = -(decayTime_ns*log(y1) + t);
547  const double z = (-B + sqrt(B*B - 4*A*C)) / (2*A);
548 
549  npe = threshold * exp(z*z);
550 
551  } else { // linear
552 
553  const double T = (getRiseTime(getStartOfLinearisation(), threshold) +
555 
556  npe = getStartOfLinearisation() + (t - T) / slope;
557  }
558 
559  return (npe <= MAX_CHARGE ? npe : MAX_CHARGE);
560  }
double threshold
threshold [npe]
double riseTime_ns
rise time of analogue pulse [ns]
do set_variable OUTPUT_DIRECTORY $WORKDIR T
double removeSaturation(const double tot_ns) const
Get time-over-threshold without saturation.
static const double C
Speed of light in vacuum [m/ns].
Definition: JConstants.hh:22
double getStartOfLinearisation() const
Get transition point from a model-dependent to linear relation between time-over-threshold and number...
source $JPP_DIR setenv csh $JPP_DIR eval JShellParser o a A
double y1
amplitude at match point [npe]
double slope
slope [ns/npe]
double getRiseTime(const double npe, const double th) const
Get time to pass from threshold to top of analogue pulse.
then set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable OUTPUT_FILE histogram.root JHistogram1D -o $WORKDIR/$OUTPUT_FILE -F "$FORMULA" -
virtual JThresholdDomains JDETECTOR::JPMTAnalogueSignalProcessor::applyThreshold ( const double  npe) const
inlinevirtual

Apply threshold.

Parameters
npenumber of photo-electrons
Returns
true if pass; else false

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 569 of file JPMTAnalogueSignalProcessor.hh.

570  {
571  if (npe > threshold) {
572 
573  return ABOVE_THRESHOLD;
574 
575  } else if (npe > threshold - thresholdBand) {
576 
577  return IN_THRESHOLDBAND;
578 
579  } else {
580 
581  return BELOW_THRESHOLDBAND;
582  }
583  }
double threshold
threshold [npe]
double thresholdBand
threshold-band [npe]
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getRiseTime ( const double  npe) const
inlinevirtual

Get time to reach threshold.

Note that the rise time is defined to be zero for a one photo-electron signal.

Parameters
npenumber of photo-electrons
Returns
time [ns]

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 594 of file JPMTAnalogueSignalProcessor.hh.

595  {
596  if (slewing) {
597 
598  switch (applyThreshold(npe)) {
599  case IN_THRESHOLDBAND:
600  return ((getRiseTime(npe, getTH0()) - getRiseTime(npe, threshold-thresholdBand)) -
602 
603  case ABOVE_THRESHOLD:
604  return ((getRiseTime(npe, getTH0()) - getRiseTime(npe, threshold)) -
605  (getRiseTime(1.0, getTH0()) - getRiseTime(1.0, threshold)));
606 
607  default:
608  THROW(JValueOutOfRange, "JPMTAnalogueSignalProcessor::getRiseTime: Invalid charge " << npe);
609  }
610 
611  } else {
612 
613  return 0.0;
614  }
615  }
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
double threshold
threshold [npe]
static double getTH0()
Get lower threshold for rise time evaluation.
double thresholdBand
threshold-band [npe]
bool slewing
time slewing of analogue signal
virtual JThresholdDomains applyThreshold(const double npe) const
Apply threshold.
double getRiseTime(const double npe, const double th) const
Get time to pass from threshold to top of analogue pulse.
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getTimeOverThreshold ( const double  npe) const
inlinevirtual

Get time-over-threshold (ToT).

Parameters
npenumber of photo-electrons
Returns
ToT [ns]

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 624 of file JPMTAnalogueSignalProcessor.hh.

625  {
626 
627  switch (applyThreshold(npe)) {
628  case IN_THRESHOLDBAND: {
629 
630  return gRandom->Gaus(mean_ns, sigma_ns);
631  }
632 
633  case ABOVE_THRESHOLD: {
634 
635  double tot = 0.0;
636 
637  if (npe*y1 <= threshold) {
638 
639  tot += getRiseTime(npe, threshold); // Gaussian
640  tot += getRiseTime(npe, threshold); // Gaussian
641 
642  } else if (npe <= getStartOfLinearisation()) {
643 
644  tot += getRiseTime (npe, threshold); // Gaussian
645  tot += getDecayTime(npe, threshold); // exponential
646 
647  } else {
648 
649  tot += getRiseTime (getStartOfLinearisation(), threshold); // Gaussian
650  tot += getDecayTime(getStartOfLinearisation(), threshold); // exponential
651 
652  tot += slope * (npe - getStartOfLinearisation()); // linear
653  }
654 
655  return applySaturation(tot);
656  }
657 
658  default: {
659 
660  THROW(JValueOutOfRange, "JPMTAnalogueSignalProcessor::getTimeOverThreshold: Invalid charge " << npe);
661  }}
662  }
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
double threshold
threshold [npe]
double mean_ns
mean time-over-threshold of threshold-band hits [ns]
double sigma_ns
time-over-threshold standard deviation of threshold-band hits [ns]
virtual JThresholdDomains applyThreshold(const double npe) const
Apply threshold.
double applySaturation(const double T) const
Get time-over-threshold with saturation.
double getStartOfLinearisation() const
Get transition point from a model-dependent to linear relation between time-over-threshold and number...
double y1
amplitude at match point [npe]
double slope
slope [ns/npe]
double getRiseTime(const double npe, const double th) const
Get time to pass from threshold to top of analogue pulse.
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getDerivative ( const double  npe) const
inlinevirtual

Get derivative of number of photo-electrons to time-over-threshold.

Parameters
npenumber of photo-electrons
Returns
dnpe/dToT [ns^-1]

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 671 of file JPMTAnalogueSignalProcessor.hh.

672  {
673  switch(applyThreshold(npe)) {
674  case ABOVE_THRESHOLD: {
675 
676  const double z = riseTime_ns / sqrt(2.0 * log(npe/threshold));
677 
678  double y = 0.0;
679 
680  if (npe*y1 > threshold) {
681 
682  if (npe <= getStartOfLinearisation())
683  y = npe / (z + decayTime_ns); // Gaussian + exponential
684  else
685  y = 1.0 / slope; // linear
686 
687  } else {
688 
689  y = npe / (2.0 * z); // Gaussian + Gaussian
690  }
691 
692  const double tot_ns = getTimeOverThreshold(npe);
693 
694  return y / getDerivativeOfSaturation(removeSaturation(tot_ns));
695 
696  }
697 
698  default:
699  return 0.0;
700  }
701  }
double getDerivativeOfSaturation(const double T) const
Get derivative of saturation factor.
virtual double getTimeOverThreshold(const double npe) const
Get time-over-threshold (ToT).
double threshold
threshold [npe]
double riseTime_ns
rise time of analogue pulse [ns]
virtual JThresholdDomains applyThreshold(const double npe) const
Apply threshold.
double removeSaturation(const double tot_ns) const
Get time-over-threshold without saturation.
double getStartOfLinearisation() const
Get transition point from a model-dependent to linear relation between time-over-threshold and number...
double y1
amplitude at match point [npe]
double slope
slope [ns/npe]
double JDETECTOR::JPMTAnalogueSignalProcessor::getTimeOverThresholdProbability ( const double  tot_ns,
const int  NPE 
) const
inline

Get probability of having a pulse with specific time-over-threshold.

Parameters
tot_nstime-over-threshold with saturation [ns]
NPEtrue number of photo-electrons
Returns
probability [ns^-1]

Definition at line 711 of file JPMTAnalogueSignalProcessor.hh.

712  {
713 
714  const double PthBand = getIntegralOfChargeProbability(IN_THRESHOLDBAND, NPE);
715 
716  const double npe = getNPE(tot_ns);
717  const double y = getChargeProbability(npe, NPE);
718  const double v = getDerivative(npe);
719 
720  const double RthBand = PthBand * TMath::Gaus(tot_ns, mean_ns, sigma_ns, kTRUE);
721  const double RaboveTh = y * v;
722 
723  return RthBand + RaboveTh;
724  }
double getIntegralOfChargeProbability(const double xmin, const double xmax, const int NPE) const
Get integral of probability.
virtual double getChargeProbability(const double npe, const int NPE) const
Get probability density for given charge.
double mean_ns
mean time-over-threshold of threshold-band hits [ns]
double sigma_ns
time-over-threshold standard deviation of threshold-band hits [ns]
virtual double getNPE(const double tot_ns, const double eps=1.0e-3) const
Get number of photo-electrons.
data_type v[N+1][M+1]
Definition: JPolint.hh:707
virtual double getDerivative(const double npe) const
Get derivative of number of photo-electrons to time-over-threshold.
double JDETECTOR::JPMTAnalogueSignalProcessor::getIntegralOfTimeOverThresholdProbability ( const double  Tmin,
const double  Tmax,
const int  NPE 
) const
inline

Get cumulative probability of time-over-threshold distribution.

Parameters
Tminminimum time-over-threshold (with saturation) [ns]
Tmaxmaximum time-over-threshold (with saturation) [ns]
NPEtrue number of photo-electrons
Returns
probability [ns^-1]

Definition at line 735 of file JPMTAnalogueSignalProcessor.hh.

736  {
737 
738  const double PthBand = getIntegralOfChargeProbability(IN_THRESHOLDBAND, NPE);
739 
740  const double IthBand = PthBand * (0.5 * TMath::Erfc((Tmin - mean_ns) / sqrt(2.0) / sigma_ns) -
741  0.5 * TMath::Erfc((Tmax - mean_ns) / sqrt(2.0) / sigma_ns));
742  const double IaboveTh = getIntegralOfChargeProbability(getNPE(Tmin), getNPE(Tmax), NPE);
743 
744  return IthBand + IaboveTh;
745  }
double getIntegralOfChargeProbability(const double xmin, const double xmax, const int NPE) const
Get integral of probability.
double mean_ns
mean time-over-threshold of threshold-band hits [ns]
double sigma_ns
time-over-threshold standard deviation of threshold-band hits [ns]
virtual double getNPE(const double tot_ns, const double eps=1.0e-3) const
Get number of photo-electrons.
virtual double JDETECTOR::JPMTAnalogueSignalProcessor::getSurvivalProbability ( const int  NPE) const
inlinevirtual

Probability that a hit survives the simulation of the PMT.

Parameters
NPEnumber of photo-electrons
Returns
probability

Reimplemented from JDETECTOR::JPMTSignalProcessorInterface.

Definition at line 754 of file JPMTAnalogueSignalProcessor.hh.

755  {
756  if (QE <= 0.0) {
757 
758  return 0.0;
759 
760  } else if (QE <= 1.0) {
761 
762  double P = 0.0;
763 
764  const double f = gainSpread * gainSpread;
765 
766  for (int i = 1; i <= NPE; ++i) {
767 
768  const double p = (TMath::Binomial(NPE, i) *
769  TMath::Power(QE, i) * TMath::Power(1.0 - QE, NPE - i));
770 
771  double Ptotal = 0.0;
772  double Pabove = 0.0;
773  double weight = (PunderAmplified > 0.0 ? pow(PunderAmplified, i) : 1.0);
774 
775  for (int k = (PunderAmplified > 0.0 ? 0 : NPE); k <= NPE; ++k) {
776 
777  const double mu = ((i-k) * f * gain) + (k * gain);
778  const double sigma = sqrt((i-k) * f + k) * getGainSpread(1);
779 
780  Ptotal += weight * (0.5 * TMath::Erfc(( 0.0 - mu) / sqrt(2.0) / sigma));
781  Pabove += weight * (0.5 * TMath::Erfc((threshold - mu) / sqrt(2.0) / sigma));
782 
783  weight *= ( (i - k) / ((double) (k+1)) *
784  (1.0 - PunderAmplified) / PunderAmplified );
785  }
786 
787  P += p*Pabove/Ptotal;
788  }
789 
790  return P;
791 
792  } else {
793 
794  THROW(JValueOutOfRange, "JPMTAnalogueSignalProcessor::getSurvivalProbability: Invalid QE " << QE);
795  }
796  }
*fatal Wrong number of arguments esac JCookie sh JRuns D $DETECTOR d sort n k
Definition: JRunrange.sh:16
double gain
gain [unit]
double gainSpread
gain spread [unit]
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
Definition: JException.hh:670
double threshold
threshold [npe]
double PunderAmplified
probability of underamplified hit
then JPizza f
Definition: JPizza.sh:46
double getGainSpread(int NPE) const
Get gain spread for given number of photo-electrons.
double QE
relative quantum efficiency
std::vector< double > weight
Definition: JAlgorithm.hh:428
then $DIR JPlotNPE PDG P
Definition: JPlotNPE-PDG.sh:60
static double JDETECTOR::JPMTAnalogueSignalProcessor::getTH0 ( )
inlinestatic

Get lower threshold for rise time evaluation.

Returns
threshold [npe]

Definition at line 804 of file JPMTAnalogueSignalProcessor.hh.

805  {
806  return 0.1;
807  }
static double JDETECTOR::JPMTAnalogueSignalProcessor::getTH1 ( )
inlinestatic

Get upper threshold for rise time evaluation.

Returns
threshold [npe]

Definition at line 815 of file JPMTAnalogueSignalProcessor.hh.

816  {
817  return 0.9;
818  }
void JDETECTOR::JPMTSignalProcessorInterface::operator() ( const JCalibration calibration,
const JPMTData< JPMTSignal > &  input,
JPMTData< JPMTPulse > &  output 
) const
inlineinherited

Process hits.

Parameters
calibrationPMT calibration
inputPMT signals
outputPMT hits

Definition at line 103 of file JPMTSignalProcessorInterface.hh.

106  {
107  // apply transition time distribution to each photo-electron.
108 
109  buffer.clear();
110 
111  for (JPMTData<JPMTSignal>::const_iterator hit = input.begin(); hit != input.end(); ++hit) {
112  for (int i = 0; i != hit->npe; ++i) {
113  if (applyQE()) {
114  buffer.push_back(JPhotoElectron(getRandomTime(hit->t_ns)));
115  }
116  }
117  }
118 
119  if (!buffer.empty()) {
120 
122 
123  buffer.sort();
124 
125 
126  // generate PMT hits from time sequence of photo-electrons.
127 
128  for (JPMTData<JPhotoElectron>::const_iterator q = buffer.begin(), p = q++; q != buffer.end(); ++q) {
129 
130  while (compare(*p,*q)) {
131  ++q;
132  }
133 
134  const double npe = getRandomCharge(distance(p,q));
135  if (applyThreshold(npe) > BELOW_THRESHOLDBAND) {
136  output.push_back(JPMTPulse(putTime(p->t_ns + getRiseTime(npe), calibration), getTimeOverThreshold(npe)));
137  }
138 
139  p = q;
140  }
141 
142  // merge overlapping PMT hits.
143 
144  merge(output);
145  }
146  }
Data structure for PMT digital pulse.
virtual double getRandomCharge(const int NPE) const
Get randomised charge in accordance with gain and gain spread.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Data structure for single photo-electron.
virtual void merge(JPMTData< JPMTPulse > &data) const
Merging of PMT hits.
virtual bool applyQE() const
Apply relative QE.
virtual double getRiseTime(const double npe) const
Get time to pass threshold.
double putTime(const T &t1, const JCalibration &cal)
Get de-calibrated time.
virtual JThresholdDomains applyThreshold(const double npe) const
Apply threshold.
virtual bool compare(const JPhotoElectron &first, const JPhotoElectron &second) const
Compare arrival times of photo-electrons.
static JPhotoElectron getEndMarker()
Get end marker.
virtual double getRandomTime(const double t_ns) const
Get randomised time according transition time distribution.
virtual double getTimeOverThreshold(const double npe) const
Get time over threshold (ToT).
std::vector< JElement_t >::const_iterator const_iterator
virtual void JDETECTOR::JPMTSignalProcessorInterface::merge ( JPMTData< JPMTPulse > &  data) const
inlinevirtualinherited

Merging of PMT hits.

Hits with overlapping time-over-threshold signals should -de facto- be combined. In this, the leading edge is maintained and the time-over-threshold is set to the difference between the overall trailing and leading edges. As a result, the number of PMT hits may be reduced.

Parameters
dataPMT hits

Definition at line 294 of file JPMTSignalProcessorInterface.hh.

295  {
296  using namespace std;
297 
298  JPMTData<JPMTPulse>::iterator out = data.begin();
299 
300  for (JPMTData<JPMTPulse>::iterator i = data.begin(); i != data.end(); ) {
301 
302  double t1 = i->t_ns;
303  double t2 = i->t_ns + i->tot_ns;
304 
305  while (++i != data.end() && i->t_ns < t2 + getTmin()) {
306  t2 = max(t2, i->t_ns + i->tot_ns);
307  }
308 
309  out->t_ns = t1;
310  out->tot_ns = t2 - t1;
311 
312  ++out;
313  }
314 
315  data.resize(distance(data.begin(), out));
316  }
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
static double getTmin()
Get two photo-electron resolution for time-over-threshold.
std::vector< JElement_t >::iterator iterator
static double JDETECTOR::JPMTSignalProcessorInterface::getTmin ( )
inlinestaticinherited

Get two photo-electron resolution for time-over-threshold.

Returns
minimal time [ns]

Definition at line 380 of file JPMTSignalProcessorInterface.hh.

381  {
382  return 1.0;
383  }
static double JDETECTOR::JPMTSignalProcessorInterface::getQmin ( )
inlinestaticinherited

Get minimal width of charge distribution.

Returns
minimal charge [npe]

Definition at line 391 of file JPMTSignalProcessorInterface.hh.

392  {
393  return 1.0e-3;
394  }
const JPMTParameters& JDETECTOR::JPMTParameters::getPMTParameters ( ) const
inlineinherited

Get PMT parameters.

Returns
PMT parameters

Definition at line 111 of file JPMTParameters.hh.

112  {
113  return static_cast<const JPMTParameters&>(*this);
114  }
Data structure for PMT parameters.
bool JDETECTOR::JPMTParameters::is_valid ( ) const
inlineinherited

Check validity of PMT parameters.

Returns
true if valid; else false

Definition at line 133 of file JPMTParameters.hh.

134  {
135  if (this->QE < 0.0 ||
136  this->gain < 0.0 ||
137  this->gainSpread < 0.0 ||
138  this->threshold < 0.0 ||
139  this->thresholdBand < 0.0) {
140  return false;
141  }
142 
143  return true;
144  }
double gain
gain [unit]
double gainSpread
gain spread [unit]
double threshold
threshold [npe]
double thresholdBand
threshold-band [npe]
double QE
relative quantum efficiency
static JEquationParameters& JDETECTOR::JPMTParameters::getEquationParameters ( )
inlinestaticinherited

Get equation parameters.

Returns
equation parameters

Definition at line 192 of file JPMTParameters.hh.

193  {
194  static JEquationParameters equation("=", ",", "./", "#");
195 
196  return equation;
197  }
Simple data structure to support I/O of equations (see class JLANG::JEquation).
static void JDETECTOR::JPMTParameters::setEquationParameters ( const JEquationParameters equation)
inlinestaticinherited

Set equation parameters.

Parameters
equationequation parameters

Definition at line 205 of file JPMTParameters.hh.

206  {
207  getEquationParameters() = equation;
208  }
static JEquationParameters & getEquationParameters()
Get equation parameters.
JProperties JDETECTOR::JPMTParameters::getProperties ( const JEquationParameters equation = JPMTParameters::getEquationParameters())
inlineinherited

Get properties of this class.

Parameters
equationequation parameters

Definition at line 216 of file JPMTParameters.hh.

217  {
218  return JPMTParametersHelper(*this, equation);
219  }
JProperties JDETECTOR::JPMTParameters::getProperties ( const JEquationParameters equation = JPMTParameters::getEquationParameters()) const
inlineinherited

Get properties of this class.

Parameters
equationequation parameters

Definition at line 227 of file JPMTParameters.hh.

228  {
229  return JPMTParametersHelper(*this, equation);
230  }

Friends And Related Function Documentation

std::istream& operator>> ( std::istream &  in,
JPMTAnalogueSignalProcessor object 
)
friend

Read PMT signal from input.

Parameters
ininput stream
objectPMT signal
Returns
input stream

Definition at line 428 of file JPMTAnalogueSignalProcessor.hh.

429  {
430  in >> static_cast<JPMTParameters&>(object);
431 
432  object.configure();
433 
434  return in;
435  }
esac print_variable DETECTOR INPUT_FILE OUTPUT_FILE CDF for TYPE in
Definition: JSirene.sh:45
JPMTParameters()
Default constructor.

Member Data Documentation

double JDETECTOR::JPMTAnalogueSignalProcessor::decayTime_ns
protected

decay time [ns]

Definition at line 822 of file JPMTAnalogueSignalProcessor.hh.

double JDETECTOR::JPMTAnalogueSignalProcessor::t1
protected

time at match point [ns]

Definition at line 823 of file JPMTAnalogueSignalProcessor.hh.

double JDETECTOR::JPMTAnalogueSignalProcessor::y1
protected

amplitude at match point [npe]

Definition at line 824 of file JPMTAnalogueSignalProcessor.hh.

double JDETECTOR::JPMTAnalogueSignalProcessor::x1
protected

Transition point from a logarithmic to a linear relation between time-over-threshold and number of photo-electrons.


Measurements by B. Schermer and R. Bruijn at Nikhef.

Definition at line 830 of file JPMTAnalogueSignalProcessor.hh.

double JDETECTOR::JPMTSignalProcessorInterface::MAX_CHARGE
protectedinherited

Maximum charge [npe].

Definition at line 398 of file JPMTSignalProcessorInterface.hh.

double JDETECTOR::JPMTParameters::QE
inherited

relative quantum efficiency

Definition at line 233 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::gain
inherited

gain [unit]

Definition at line 234 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::gainSpread
inherited

gain spread [unit]

Definition at line 235 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::riseTime_ns
inherited

rise time of analogue pulse [ns]

Definition at line 236 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::TTS_ns
inherited

transition time spread [ns]

Definition at line 237 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::threshold
inherited

threshold [npe]

Definition at line 238 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::PunderAmplified
inherited

probability of underamplified hit

Definition at line 239 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::thresholdBand
inherited

threshold-band [npe]

Definition at line 240 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::mean_ns
inherited

mean time-over-threshold of threshold-band hits [ns]

Definition at line 241 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::sigma_ns
inherited

time-over-threshold standard deviation of threshold-band hits [ns]

Definition at line 242 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::slope
inherited

slope [ns/npe]

Definition at line 243 of file JPMTParameters.hh.

double JDETECTOR::JPMTParameters::saturation
inherited

saturation [ns]

Definition at line 244 of file JPMTParameters.hh.

bool JDETECTOR::JPMTParameters::slewing
inherited

time slewing of analogue signal

Definition at line 245 of file JPMTParameters.hh.


The documentation for this struct was generated from the following file: