Jpp  18.2.0-rc.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Public Attributes | Protected Attributes | List of all members
JPHYSICS::JLED_C Class Referenceabstract

Probability Density Functions of the time response of a PMT (C-like interface) More...

#include <JLED.hh>

Inheritance diagram for JPHYSICS::JLED_C:
JPHYSICS::JLED JPHYSICS::JDispersion JPHYSICS::JDispersionInterface JPHYSICS::JAbstractPMT JPHYSICS::JAbstractLED JPHYSICS::JAbstractMedium JPHYSICS::JDispersionInterface

Public Member Functions

 JLED_C (const double Area, const JAbstractLED *LED, double(*QE)(const double), double(*Pmt)(const double), const double L_abs, const double L_s, double(*Km3)(const double), const double P_atm, const double wavelength, const double Tmin_ns, const double Tmax_ns, const JQuadrature &engine=JCotangent(20), const int numberOfPoints=20, const double epsilon=1e-12)
 Constructor. More...
 
virtual double getLightFromLED (const double ct, const double phi, const double dt) const override
 Light yield from LED (number of p.e. More...
 
virtual double getPhotocathodeArea () const override
 Photo-cathode area of PMT. More...
 
virtual double getQE (const double lambda) const override
 Quantum efficiency of PMT (incl. More...
 
virtual double getAngularAcceptance (const double ct) const override
 Angular acceptance of PMT. More...
 
virtual double getAbsorptionLength (const double lambda) const override
 Absorption length. More...
 
virtual double getScatteringLength (const double lambda) const override
 Scattering length. More...
 
virtual double getScatteringProbability (const double ct) const override
 Model specific function to describe light scattering in water (integral over full solid angle normalised to one). More...
 
double getDirectLightFromLED (const double D_m, const double cd, const double theta, const double phi, const double t_ns) const
 Probability density function for direct light from LED. More...
 
double getScatteredLightFromLED (const double D_m, const double cd, const double theta, const double phi, const double t_ns) const
 Probability density function for scattered light from LED. More...
 
virtual double getIndexOfRefractionPhase (const double lambda) const =0
 Index of refraction for phase velocity. More...
 
virtual double getDispersionPhase (const double lambda) const =0
 Dispersion of light for phase velocity. More...
 
virtual double getIndexOfRefractionGroup (const double lambda) const
 Index of refraction for group velocity. More...
 
virtual double getDispersionGroup (const double lambda) const =0
 Dispersion of light for group velocity. More...
 
double getKappa (const double lambda) const
 Get effective index of refraction for muon light. More...
 
double getKmin (const double lambda) const
 Get smallest index of refraction for Bremsstrahlung light (i.e. point at which dt/dz = 0). More...
 
virtual double getIndexOfRefractionPhase (const double lambda) const
 Index of refraction (phase velocity). More...
 
virtual double getDispersionPhase (const double lambda) const
 Dispersion of light for phase velocity. More...
 
virtual double getDispersionGroup (const double lambda) const
 Dispersion of light for group velocity. More...
 

Public Attributes

const double P
 Dispersion parameters (x = 1/lambda) More...
 
const double a0
 offset More...
 
const double a1
 dn/dP More...
 
const double a2
 d^1n/(dx)^1 More...
 
const double a3
 d^2n/(dx)^2 More...
 
const double a4
 d^3n/(dx)^3 More...
 

Protected Attributes

const double A
 photo-cathode area [m2] More...
 
const JAbstractLEDled
 Pointer to interface for emission profile from LED. More...
 
double(* qe )(const double lambda)
 Quantum efficiency of PMT (incl. More...
 
double l_abs
 Absorption length. More...
 
double ls
 Scattering length. More...
 
double(* pmt )(const double ct)
 Angular acceptance of PMT. More...
 
double(* km3 )(const double ct)
 Model specific function to describe light scattering in water. More...
 
double wavelength
 
double tmin
 
double tmax
 
JQuadrature main_engine
 
JQuadrature beta_engine
 
std::vector< JElement3D_tphi_engine
 

Detailed Description

Probability Density Functions of the time response of a PMT (C-like interface)

Definition at line 275 of file JLED.hh.

Constructor & Destructor Documentation

JPHYSICS::JLED_C::JLED_C ( const double  Area,
const JAbstractLED LED,
double(*)(const double)  QE,
double(*)(const double)  Pmt,
const double  L_abs,
const double  L_s,
double(*)(const double)  Km3,
const double  P_atm,
const double  wavelength,
const double  Tmin_ns,
const double  Tmax_ns,
const JQuadrature engine = JCotangent(20),
const int  numberOfPoints = 20,
const double  epsilon = 1e-12 
)
inline

Constructor.

Parameters
Areaphoto-cathode area [m^2]
LEDpointer to interface for emission profile from LED
QEpointer to function for quantum efficiency of PMT
Pmtpointer to function for angular acceptance of PMT
L_absabsorption length [m]
L_sscattering length [m]
Km3pointer to model specific function to describe light scattering in water
P_atmambient pressure [atm]
wavelengthwavelength of light [nm]
Tmin_nsminimal time of emmision [ns]
Tmax_nsminimal time of emmision [ns]
enginescattering angle integrator
numberOfPointsnumber of points for integration
epsilonprecision of points for integration

Definition at line 298 of file JLED.hh.

317  :
318  JLED(wavelength, Tmin_ns, Tmax_ns, engine, numberOfPoints, epsilon),
319  JDispersion(P_atm),
320  A (Area ),
321  led (LED ),
322  qe (QE ),
323  l_abs(L_abs),
324  ls (L_s ),
325  pmt (Pmt ),
326  km3 (Km3 )
327  {}
const JAbstractLED * led
Pointer to interface for emission profile from LED.
Definition: JLED.hh:426
double l_abs
Absorption length.
Definition: JLED.hh:439
double(* km3)(const double ct)
Model specific function to describe light scattering in water.
Definition: JLED.hh:460
JDispersion(const double P_atm)
Constructor.
Definition: JDispersion.hh:35
double ls
Scattering length.
Definition: JLED.hh:444
double(* pmt)(const double ct)
Angular acceptance of PMT.
Definition: JLED.hh:452
const double A
photo-cathode area [m2]
Definition: JLED.hh:421
int numberOfPoints
Definition: JResultPDF.cc:22
double(* qe)(const double lambda)
Quantum efficiency of PMT (incl.
Definition: JLED.hh:434
const double epsilon
Definition: JQuadrature.cc:21
JLED(const double lambda, const double Tmin_ns, const double Tmax_ns, const JQuadrature &engine=JCotangent(20), const int numberOfPoints=20, const double epsilon=1e-12)
Constructor.
Definition: JLED.hh:78
double wavelength
Definition: JLED.hh:261

Member Function Documentation

virtual double JPHYSICS::JLED_C::getLightFromLED ( const double  ct,
const double  phi,
const double  dt 
) const
inlineoverridevirtual

Light yield from LED (number of p.e.

per unit solid angle per unit time).

Parameters
ctzenith angle of emission
phiazimuth angle of emission
dttime of emission [ns]
Returns
d^2P / dOmega dt [npe/ns/sr]

Implements JPHYSICS::JAbstractLED.

Definition at line 338 of file JLED.hh.

341  {
342  return led->getLightFromLED(ct, phi, dt);
343  }
const JAbstractLED * led
Pointer to interface for emission profile from LED.
Definition: JLED.hh:426
virtual double getLightFromLED(const double ct, const double phi, const double dt) const =0
Light yield from LED (number of p.e.
virtual double JPHYSICS::JLED_C::getPhotocathodeArea ( ) const
inlineoverridevirtual

Photo-cathode area of PMT.

Returns
photo-cathode area [m^2]

Implements JPHYSICS::JAbstractPMT.

Definition at line 351 of file JLED.hh.

352  {
353  return A;
354  }
const double A
photo-cathode area [m2]
Definition: JLED.hh:421
virtual double JPHYSICS::JLED_C::getQE ( const double  lambda) const
inlineoverridevirtual

Quantum efficiency of PMT (incl.

absorption in glass, gel, etc.).

Parameters
lambdawavelenth [nm]
Returns
QE

Implements JPHYSICS::JAbstractPMT.

Definition at line 363 of file JLED.hh.

364  {
365  return qe(lambda);
366  }
double(* qe)(const double lambda)
Quantum efficiency of PMT (incl.
Definition: JLED.hh:434
virtual double JPHYSICS::JLED_C::getAngularAcceptance ( const double  ct) const
inlineoverridevirtual

Angular acceptance of PMT.

Parameters
ctcosine angle of incidence
Returns
relative efficiency

Implements JPHYSICS::JAbstractPMT.

Definition at line 375 of file JLED.hh.

376  {
377  return pmt(ct);
378  }
double(* pmt)(const double ct)
Angular acceptance of PMT.
Definition: JLED.hh:452
virtual double JPHYSICS::JLED_C::getAbsorptionLength ( const double  lambda) const
inlineoverridevirtual

Absorption length.

Parameters
lambdawavelenth [nm]
Returns
absorption length [m]

Implements JPHYSICS::JAbstractMedium.

Definition at line 387 of file JLED.hh.

388  {
389  return l_abs;
390  }
double l_abs
Absorption length.
Definition: JLED.hh:439
virtual double JPHYSICS::JLED_C::getScatteringLength ( const double  lambda) const
inlineoverridevirtual

Scattering length.

Parameters
lambdawavelenth [nm]
Returns
scattering length [m]

Implements JPHYSICS::JAbstractMedium.

Definition at line 399 of file JLED.hh.

400  {
401  return ls;
402  }
double ls
Scattering length.
Definition: JLED.hh:444
virtual double JPHYSICS::JLED_C::getScatteringProbability ( const double  ct) const
inlineoverridevirtual

Model specific function to describe light scattering in water (integral over full solid angle normalised to one).

Parameters
ctcosine scattering angle
Returns
probability

Implements JPHYSICS::JAbstractMedium.

Definition at line 412 of file JLED.hh.

413  {
414  return km3(ct);
415  }
double(* km3)(const double ct)
Model specific function to describe light scattering in water.
Definition: JLED.hh:460
double JPHYSICS::JLED::getDirectLightFromLED ( const double  D_m,
const double  cd,
const double  theta,
const double  phi,
const double  t_ns 
) const
inlineinherited

Probability density function for direct light from LED.

Parameters
D_mdistance between LED and PMT [m]
cdcosine angle LED orientation and LED - PMT position
thetazenith angle orientation PMT
phiazimuth angle orientation PMT
t_nstime difference relative to direct light [ns]
Returns
dP/dt [npe/ns]

Definition at line 108 of file JLED.hh.

113  {
114  const double sd = sqrt((1.0 + cd)*(1.0 - cd));
115  const double d = D_m; // photon path [m]
116 
117  const double A = getPhotocathodeArea();
118 
119  const double px = sin(theta)*cos(phi);
120  //const double py = sin(theta)*sin(phi);
121  const double pz = cos(theta);
122 
123  const double ct = sd*px + cd*pz; // cosine angle of incidence on PMT
124 
125  const double qe = getQE(wavelength);
126  const double l_abs = getAbsorptionLength(wavelength);
127  const double ls = getScatteringLength(wavelength);
128 
129  const double U = getAngularAcceptance(ct);
130  const double V = exp(-d/l_abs) * exp(-d/ls); // absorption & scattering
131  const double W = A/(d*d); // solid angle
132 
133  return qe * getLightFromLED(cd, 0.0, t_ns) * U * V * W;
134  }
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
then usage set_variable ACOUSTICS_WORKDIR $WORKDIR set_variable FORMULA sin([0]+2 *$PI *([1]+[2]*x)*x)" set_variable DY 1.0e-8 mkdir $WORKDIR for DETECTOR in $DETECTORS[*]
virtual double getScatteringLength(const double lambda) const =0
Scattering length.
virtual double getAbsorptionLength(const double lambda) const =0
Absorption length.
virtual double getPhotocathodeArea() const =0
Photo-cathode area of PMT.
then JMuonMCEvt f $INPUT_FILE o $INTERMEDIATE_FILE d
Definition: JMuonPath.sh:47
virtual double getAngularAcceptance(const double ct) const =0
Angular acceptence of PMT.
virtual double getLightFromLED(const double ct, const double phi, const double dt) const =0
Light yield from LED (number of p.e.
virtual double getQE(const double lambda) const =0
Quantum efficiency of PMT (incl.
source $JPP_DIR setenv csh $JPP_DIR &dev null eval JShellParser o a A
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
double wavelength
Definition: JLED.hh:261
double JPHYSICS::JLED::getScatteredLightFromLED ( const double  D_m,
const double  cd,
const double  theta,
const double  phi,
const double  t_ns 
) const
inlineinherited

Probability density function for scattered light from LED.

Parameters
D_mdistance between LED and PMT [m]
cdcosine angle LED orientation and LED - PMT position
thetazenith angle orientation PMT
phiazimuth angle orientation PMT
t_nstime difference relative to direct light [ns]
Returns
dP/dt [npe/ns]

Definition at line 147 of file JLED.hh.

152  {
153  using namespace std;
154 
155  double value = 0;
156 
157  if (t_ns >= tmin) {
158 
159  const double ng = getIndexOfRefractionGroup(wavelength);
160 
161  const double sd = sqrt((1.0 + cd)*(1.0 - cd));
162  const double l = D_m; // distance [m]
163 
164  const double A = getPhotocathodeArea();
165  //const double sr = sqrt(A/PI) / D_m;
166  //const double cr = sqrt((1.0 + sr)*(1.0 - sr));
167 
168  const double px = sin(theta)*cos(phi);
169  const double py = sin(theta)*sin(phi);
170  const double pz = cos(theta);
171 
172  const double qx = cd*px + 0 - sd*pz;
173  const double qy = 1*py;
174  const double qz = sd*px + 0 + cd*pz;
175 
176  const double qe = getQE(wavelength);
177  const double l_abs = getAbsorptionLength(wavelength);
178  const double ls = getScatteringLength(wavelength);
179 
180  const double Jc = 1.0 / ls; // dN/dx
181 
182  const double xmin = tmin;
183  const double xmax = std::min(tmax, t_ns);
184 
185  JQuadrature buffer;
186 
187  if (xmax > xmin) {
188 
189  for (JQuadrature::const_iterator i = main_engine.begin(); i != main_engine.end(); ++i) {
190 
191  const double t = 0.5 * (xmax + xmin) + i->getX() * 0.5 * (xmax - xmin);
192  const double dt = i->getY() * 0.5 * (xmax - xmin);
193 
194  buffer[t] = dt;
195  }
196 
197  } else {
198 
199  buffer[tmin] = 1.0;
200  }
201 
202  for (JQuadrature::const_iterator i = buffer.begin(); i != buffer.end(); ++i) {
203 
204  const double t = i->getX();
205  const double dt = i->getY();
206 
207  const double d = l + C*(t_ns - t)/ng; // photon path
208  const double V = exp(-d/l_abs); // absorption
209 
210  for (JQuadrature::const_iterator j = beta_engine.begin(); j != beta_engine.end(); ++j) {
211 
212  const double cb = j->getX();
213  const double dc = j->getY();
214 
215  const double sb = sqrt((1.0 + cb)*(1.0 - cb));
216 
217  const double v = 0.5 * (d + l) * (d - l) / (d - l*cb);
218  const double u = d - v;
219 
220  if (u <= 0) continue;
221  if (v <= 0) continue;
222 
223  const double cts = (l*cb - v) / u; // cosine scattering angle
224 
225  const double W = min(A/(v*v), 2.0*PI); // solid angle
226  const double Ja = getScatteringProbability(cts); // P(cos(),phi)
227  const double Jd = ng * (1.0 - cts) / C; // dt/du
228 
229  const double ca = (l - v*cb) / u;
230  const double sa = v*sb / u;
231 
232  for (std::vector<JElement3D_t>::const_iterator k = phi_engine.begin(); k != phi_engine.end(); ++k) {
233 
234  const double cp = k->getX();
235  const double sp = k->getY();
236  const double dp = k->getZ();
237 
238  const double dom = dc * dp * v*v / (u*u);
239 
240  const double ct0 = cd*ca - sd*sa*cp;
241  const double ph0 = atan2(sa*sp, cd*sa*cp + sd*ca);
242 
243  const double vx = -sb*cp * qx;
244  const double vy = -sb*sp * qy;
245  const double vz = cb * qz;
246 
247  const double U =
248  getAngularAcceptance(vx + vy + vz) * getLightFromLED(ct0, +ph0, t) +
249  getAngularAcceptance(vx - vy + vz) * getLightFromLED(ct0, -ph0, t);
250 
251  value += qe * dt * dom * Ja * Jc * U * V * W / Jd;
252  }
253  }
254  }
255  }
256 
257  return value;
258  }
const double xmax
Definition: JQuadrature.cc:24
double tmin
Definition: JLED.hh:263
then fatal No hydrophone data file $HYDROPHONE_TXT fi sort gr k
JQuadrature main_engine
Definition: JLED.hh:266
double tmax
Definition: JLED.hh:264
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
then usage set_variable ACOUSTICS_WORKDIR $WORKDIR set_variable FORMULA sin([0]+2 *$PI *([1]+[2]*x)*x)" set_variable DY 1.0e-8 mkdir $WORKDIR for DETECTOR in $DETECTORS[*]
static const double C
Physics constants.
std::vector< JElement3D_t > phi_engine
Definition: JLED.hh:268
virtual double getScatteringLength(const double lambda) const =0
Scattering length.
JQuadrature beta_engine
Definition: JLED.hh:267
virtual double getAbsorptionLength(const double lambda) const =0
Absorption length.
static const double PI
Mathematical constants.
virtual double getPhotocathodeArea() const =0
Photo-cathode area of PMT.
const double xmin
Definition: JQuadrature.cc:23
virtual double getIndexOfRefractionGroup(const double lambda) const
Index of refraction for group velocity.
then JMuonMCEvt f $INPUT_FILE o $INTERMEDIATE_FILE d
Definition: JMuonPath.sh:47
virtual double getAngularAcceptance(const double ct) const =0
Angular acceptence of PMT.
virtual double getLightFromLED(const double ct, const double phi, const double dt) const =0
Light yield from LED (number of p.e.
then cp
container_type::const_iterator const_iterator
Definition: JCollection.hh:91
int j
Definition: JPolint.hh:703
virtual double getQE(const double lambda) const =0
Quantum efficiency of PMT (incl.
data_type v[N+1][M+1]
Definition: JPolint.hh:777
double u[N+1]
Definition: JPolint.hh:776
source $JPP_DIR setenv csh $JPP_DIR &dev null eval JShellParser o a A
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
virtual double getScatteringProbability(const double ct) const =0
Model specific function to describe light scattering in water (integral over full solid angle normali...
double wavelength
Definition: JLED.hh:261
virtual double JPHYSICS::JDispersionInterface::getIndexOfRefractionPhase ( const double  lambda) const
pure virtualinherited

Index of refraction for phase velocity.

Parameters
lambdawavelenth [nm]
Returns
index of refraction

Implemented in JPHYSICS::JDispersion.

virtual double JPHYSICS::JDispersionInterface::getDispersionPhase ( const double  lambda) const
pure virtualinherited

Dispersion of light for phase velocity.

Parameters
lambdawavelength of light [nm]
Returns
dn/dlambda

Implemented in JPHYSICS::JDispersion.

virtual double JPHYSICS::JDispersionInterface::getIndexOfRefractionGroup ( const double  lambda) const
inlinevirtualinherited

Index of refraction for group velocity.

Parameters
lambdawavelenth [nm]
Returns
index of refraction

Definition at line 52 of file JDispersionInterface.hh.

53  {
54  const double n = getIndexOfRefractionPhase(lambda);
55  const double y = getDispersionPhase(lambda);
56 
57  return n / (1.0 + y*lambda/n);
58  }
const int n
Definition: JPolint.hh:697
virtual double getIndexOfRefractionPhase(const double lambda) const =0
Index of refraction for phase velocity.
virtual double getDispersionPhase(const double lambda) const =0
Dispersion of light for phase velocity.
virtual double JPHYSICS::JDispersionInterface::getDispersionGroup ( const double  lambda) const
pure virtualinherited

Dispersion of light for group velocity.

Parameters
lambdawavelength of light [nm]
Returns
dn/dlambda

Implemented in JPHYSICS::JDispersion.

double JPHYSICS::JDispersionInterface::getKappa ( const double  lambda) const
inlineinherited

Get effective index of refraction for muon light.

Parameters
lambdawavelength of light [nm]
Returns
index of refraction

Definition at line 76 of file JDispersionInterface.hh.

77  {
78  const double n = getIndexOfRefractionPhase(lambda);
79  const double ng = getIndexOfRefractionGroup(lambda);
80 
81  return (ng * n - 1.0) / sqrt(n*n - 1.0);
82  }
const int n
Definition: JPolint.hh:697
virtual double getIndexOfRefractionPhase(const double lambda) const =0
Index of refraction for phase velocity.
virtual double getIndexOfRefractionGroup(const double lambda) const
Index of refraction for group velocity.
double JPHYSICS::JDispersionInterface::getKmin ( const double  lambda) const
inlineinherited

Get smallest index of refraction for Bremsstrahlung light (i.e. point at which dt/dz = 0).

Parameters
lambdawavelength of light [nm]
Returns
index of refraction

Definition at line 91 of file JDispersionInterface.hh.

92  {
93  const double ng = getIndexOfRefractionGroup(lambda);
94 
95  return sqrt(ng*ng - 1.0);
96  }
virtual double getIndexOfRefractionGroup(const double lambda) const
Index of refraction for group velocity.
virtual double JPHYSICS::JDispersion::getIndexOfRefractionPhase ( const double  lambda) const
inlinevirtualinherited

Index of refraction (phase velocity).

Parameters
lambdawavelenth [nm]
Returns
index of refraction

Implements JPHYSICS::JDispersionInterface.

Definition at line 51 of file JDispersion.hh.

52  {
53  const double x = 1.0 / lambda;
54 
55  return a0 + a1*P + x*(a2 + x*(a3 + x*a4));
56  }
const double a4
d^3n/(dx)^3
Definition: JDispersion.hh:100
const double a1
dn/dP
Definition: JDispersion.hh:97
const double P
Dispersion parameters (x = 1/lambda)
Definition: JDispersion.hh:95
const double a0
offset
Definition: JDispersion.hh:96
const double a2
d^1n/(dx)^1
Definition: JDispersion.hh:98
const double a3
d^2n/(dx)^2
Definition: JDispersion.hh:99
virtual double JPHYSICS::JDispersion::getDispersionPhase ( const double  lambda) const
inlinevirtualinherited

Dispersion of light for phase velocity.

Parameters
lambdawavelength of light [nm]
Returns
dn/dlambda

Implements JPHYSICS::JDispersionInterface.

Definition at line 65 of file JDispersion.hh.

66  {
67  const double x = 1.0 / lambda;
68 
69  return -x*x*(a2 + x*(2.0*a3 + x*3.0*a4));
70  }
const double a4
d^3n/(dx)^3
Definition: JDispersion.hh:100
const double a2
d^1n/(dx)^1
Definition: JDispersion.hh:98
const double a3
d^2n/(dx)^2
Definition: JDispersion.hh:99
virtual double JPHYSICS::JDispersion::getDispersionGroup ( const double  lambda) const
inlinevirtualinherited

Dispersion of light for group velocity.

Parameters
lambdawavelength of light [nm]
Returns
dn/dlambda

Implements JPHYSICS::JDispersionInterface.

Definition at line 79 of file JDispersion.hh.

80  {
81  const double x = 1.0 / lambda;
82 
83  const double n = getIndexOfRefractionPhase(lambda);
84  const double np = getDispersionPhase(lambda);
85  const double npp = x*x*x*(2.0*a2 + x*(6.0*a3 + x*12.0*a4));
86  const double ng = n / (1.0 + np*lambda/n);
87 
88  return ng*ng * (2*np*np - n*npp) * lambda / (n*n*n);
89  }
const double a4
d^3n/(dx)^3
Definition: JDispersion.hh:100
const int n
Definition: JPolint.hh:697
virtual double getIndexOfRefractionPhase(const double lambda) const
Index of refraction (phase velocity).
Definition: JDispersion.hh:51
virtual double getDispersionPhase(const double lambda) const
Dispersion of light for phase velocity.
Definition: JDispersion.hh:65
const double a2
d^1n/(dx)^1
Definition: JDispersion.hh:98
const double a3
d^2n/(dx)^2
Definition: JDispersion.hh:99

Member Data Documentation

const double JPHYSICS::JLED_C::A
protected

photo-cathode area [m2]

Definition at line 421 of file JLED.hh.

const JAbstractLED* JPHYSICS::JLED_C::led
protected

Pointer to interface for emission profile from LED.

Definition at line 426 of file JLED.hh.

double(* JPHYSICS::JLED_C::qe)(const double lambda)
protected

Quantum efficiency of PMT (incl.

absorption in glass, gel, etc.)

Parameters
lambdawavelenth [nm]
Returns
QE

Definition at line 434 of file JLED.hh.

double JPHYSICS::JLED_C::l_abs
protected

Absorption length.

Definition at line 439 of file JLED.hh.

double JPHYSICS::JLED_C::ls
protected

Scattering length.

Definition at line 444 of file JLED.hh.

double(* JPHYSICS::JLED_C::pmt)(const double ct)
protected

Angular acceptance of PMT.

Parameters
ctcosine angle of incidence
Returns
relative efficiency

Definition at line 452 of file JLED.hh.

double(* JPHYSICS::JLED_C::km3)(const double ct)
protected

Model specific function to describe light scattering in water.

Parameters
ctcosine scattering angle
Returns
probability

Definition at line 460 of file JLED.hh.

double JPHYSICS::JLED::wavelength
protectedinherited

Definition at line 261 of file JLED.hh.

double JPHYSICS::JLED::tmin
protectedinherited

Definition at line 263 of file JLED.hh.

double JPHYSICS::JLED::tmax
protectedinherited

Definition at line 264 of file JLED.hh.

JQuadrature JPHYSICS::JLED::main_engine
protectedinherited

Definition at line 266 of file JLED.hh.

JQuadrature JPHYSICS::JLED::beta_engine
protectedinherited

Definition at line 267 of file JLED.hh.

std::vector<JElement3D_t> JPHYSICS::JLED::phi_engine
protectedinherited

Definition at line 268 of file JLED.hh.

const double JPHYSICS::JDispersion::P
inherited

Dispersion parameters (x = 1/lambda)

ambient pressure [atm]

Definition at line 95 of file JDispersion.hh.

const double JPHYSICS::JDispersion::a0
inherited

offset

Definition at line 96 of file JDispersion.hh.

const double JPHYSICS::JDispersion::a1
inherited

dn/dP

Definition at line 97 of file JDispersion.hh.

const double JPHYSICS::JDispersion::a2
inherited

d^1n/(dx)^1

Definition at line 98 of file JDispersion.hh.

const double JPHYSICS::JDispersion::a3
inherited

d^2n/(dx)^2

Definition at line 99 of file JDispersion.hh.

const double JPHYSICS::JDispersion::a4
inherited

d^3n/(dx)^3

Definition at line 100 of file JDispersion.hh.


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