Jpp  debug
the software that should make you happy
JScatteringModel.hh
Go to the documentation of this file.
1 #ifndef H_J_SCATTERING_MODEL
2 #define H_J_SCATTERING_MODEL
3 
4 /**
5  * \author mjongen
6  */
7 
8 // c++ standard library includes
9 #include<iostream>
10 #include<iomanip>
11 #include<cmath>
12 #include<cstdlib>
13 
14 // root includes
15 #include "TObject.h"
16 #include "TH1F.h"
17 #include "TH2F.h"
18 #include "TRandom3.h"
19 
20 // JPP includes
21 // these include files are hidden from rootcint
22 #include "JROOT/JRoot.hh"
23 #include "JMarkov/JPhotonPath.hh"
24 #include "JGeometry3D/JAxis3D.hh"
27 
28 #ifndef __CINT__
29  #include "JPhysics/KM3NeT.hh"
30 #endif
31 
32 namespace JMARKOV {
33 
34  // for rootcint, all we need is this class declaration
35  #ifdef __CINT__
36  class JPhotonPath ;
37  #endif
38 }
39 
40 namespace JMARKOV {
41  using namespace JGEOMETRY3D ;
42 
43  /**
44  Virtual base class for a scattering model.
45  **/
47  public :
48 
49  // constructor
51 
52  // destructor
53  virtual ~JScatteringModel() {}
54 
55  /**
56  Return the probability density as a function of cos(theta)
57 
58  dP / dOmega = dP / dcosTheta dPhi
59 
60  to scatter in a given direction. Theta=0 is forward scattering.
61  **/
62  virtual double getScatteringProbability( double ct ) = 0 ;
63 
64  /**
65  Return a randomly generated direction according to the scattering
66  probability distribution.
67 
68  This uses gRandom.
69  **/
70  virtual JVersor3D generateDirection() = 0 ;
71 
72  // return the scattering length in [m]
73  virtual double getScatteringLength() { return lambda_scat ; }
74 
75  protected :
76 
77  double lambda_scat ; // scattering length in [m]
78  } ;
79 
80  /**
81  Virtual base class for a light source.
82  **/
83  class JSourceModel {
84  public:
85 
87 
88  virtual ~JSourceModel() {}
89 
90  /**
91  Return the probability density
92 
93  dP / dOmega = dP / dCosTheta dPhi
94 
95  that a photon from this source is emitted in a given direction,
96  given that a photon is emitted.
97  **/
98  virtual double getEmissionProbability( JVersor3D dir ) = 0 ;
99 
100  /**
101  Return a randomly generated direction according to the emission
102  distribution.
103 
104  This uses gRandom.
105  **/
106  virtual JVersor3D generateDirection() = 0 ;
107 
108  void setPosition( JPosition3D& _pos ) { pos = _pos ; }
109 
110  const JPosition3D& getPosition() const { return pos ; }
111 
112  protected:
113 
115  } ;
116 
117  /**
118  Implementation of the JSourceModel class that represents an
119  isotropic source.
120  **/
122  public :
123 
125 
127  return 1.0/(4.0*M_PI) ;
128  }
129 
131  double dx, dy, dz ;
132  gRandom->Sphere(dx,dy,dz,1) ;
133  return JVersor3D(dx,dy,dz) ;
134  }
135  } ;
136 
137  /**
138  Implementation of the JSourceModel class that represents a
139  directed source with a flat intensity distribution.
140 
141  When _cutSurface = true, photons that would be emitted into
142  the surface (defined by the surface normal) will not be emitted.
143  **/
144  class JDirectedSource : public JSourceModel {
145  public :
146 
147  /**
148  Constructor.
149  _source_dir is the source direction
150  _alpha is the opening angle in radians
151 
152  note that _alpha should be > -1
153  **/
154  JDirectedSource( JVersor3D _source_dir, double _alpha, bool _cutSurface=false, const JVersor3D _surface_normal=JVersor3D() ) :
155  source_dir(_source_dir), ctmin(cos(0.5*_alpha)), alpha(_alpha), cutSurface(_cutSurface), surface_normal(_surface_normal) {
156  // we have to find out which fraction of photons would be emitted into the surface, given
157  // the direction and opening angle
158  // as far as I can tell, this is a non-trivial integral, so I solve it numerically
159 
160  // delta is the angle of the beam direction over the surface
161  // (negative means the beam points into the surface)
162  delta = 0.5*M_PI - acos(surface_normal.getDot(source_dir)) ;
163  omega = 2.0*M_PI*(1.0-ctmin) ; // space angle of beam ("spherical cap")
164 
165  if( cutSurface ) {
166  if( delta>0.0 && delta>0.5*_alpha ) {
167  // no part of the beam points into the surface
168  //fraction_into_surface = 0.0 ;
169  } else if( delta<=0 && fabs(delta)>=0.5*alpha ) {
170  // the beam points completely into the surface, throw an error
171  //fraction_into_surface = 1.0 ;
172  omega = 0.0 ;
173  cerr << "FATAL ERROR in JDirectedSource: beam points completely into the surface." << endl ;
174  exit(1) ;
175  } else {
176  // part of the beam points into the surface, which effectively reduces
177  // the space angle omega
178  // here we numerically compute omega
179 
180  // only at theta=delta to alpha/2 is part of the beam absorbed
181  // we calculate the space angle of the absorbed beam part
182  double omega_absorbed = 0.0 ;
183  const int __nct = 100000 ;
184  const double __ctmin = cos(0.5*alpha) ;
185  const double __ctmax = cos(delta) ;
186  const double __dct = (__ctmax - __ctmin)/__nct ;
187  for( int i=0; i<__nct; ++i ) {
188  double __ct = __ctmin + (i+0.5) * __dct ;
189  double __theta = acos(__ct) ;
190  double __phirange = 2.0*(M_PI-getPhiMax(__theta)) ;
191  omega_absorbed += __phirange ;
192  }
193  omega_absorbed *= __dct ;
194  omega -= omega_absorbed ;
195 
196  // if the beam direction is into the surface, this part gets absorbed completely
197  if( delta<0 ) {
198  omega -= 2.0 * M_PI * (1.0-cos(delta)) ;
199  }
200  }
201  }
202  }
203 
204  /// return the opening angle in radians
205  double getOpeningAngle() const { return alpha ; }
206 
207  /// return the source direction
208  JVersor3D getDirection() const { return source_dir ; }
209 
211  double ct = dir.getDot(source_dir) ;
212  if( ct>ctmin ) {
213  if( cutSurface && dir.getDot(surface_normal)<0 ) return 0.0 ;
214  return 1.0/omega ;
215  }
216  return 0.0 ;
217  }
218 
220  while(true) {
221  double phi = gRandom->Uniform(0,2*M_PI) ;
222  double ct = gRandom->Uniform(ctmin,1) ;
223  JDirection3D _dir( JAngle3D(acos(ct),phi) ) ;
224  JRotation3D R( source_dir ) ;
225  JVersor3D dir(_dir.rotate_back(R)) ;
226  if( !cutSurface ) return dir ;
227  if( dir.getDot(surface_normal)>0 ) return dir ;
228  }
229  }
230 
231  protected :
232 
233  /**
234  For a given angle theta with respect to the beam direction, and a given surface
235  normal, we define phi as the rotation angle around the beam, where phi=0
236  corresponds to the direction sticking out over the surface a much as possible.
237 
238  This function returns 'phi max', which is the largest value of phi for which the
239  direction does not point into the surface.
240 
241  **/
242  double getPhiMax( double theta ) {
243  if( delta>0 ) {
244  if( theta<delta ) return M_PI ;
245  if( theta+delta > M_PI ) return 0.0 ;
246  }
247  if( delta<0 ) {
248  if( theta<fabs(delta) ) return 0.0 ;
249  if( theta+fabs(delta) > M_PI ) return M_PI ;
250  }
251  double phimax = acos(-tan(delta)/tan(theta)) ;
252  return phimax ;
253  }
254 
255 
257  double ctmin ;
258  double alpha ;
259  bool cutSurface ;
261  double delta ;
262  //double fraction_into_surface ; // the fraction of photons that would be emitted into the surface
263  double omega ; // integrated phase angle of photon emission directions
264  } ;
265 
266 
267  /**
268  Virtual base class for a light detector ("photon target").
269  **/
270  class JTargetModel {
271  public:
272 
274  r = 0 ; // infinitely small target by default
275  }
276 
277  virtual ~JTargetModel() {}
278 
279  /**
280  Return the efficiency, which is defined as the probability that
281  a photon with the given direction hitting the target will be
282  registered.
283 
284  Note that we assume by convention that the direction is the
285  PHOTON direction, NOT the direction that you would see the
286  photon coming from!
287 
288  By convention the highest efficiency is 1.
289  **/
290  virtual double getEfficiency( JVersor3D dir ) const = 0 ;
291 
292  /**
293  Return the effective area, i.e. the integral over the whole
294  surface of the target, weighted by the efficiency.
295  **/
296  virtual double getEffectiveArea() { return 0.0 ; }
297 
298  // set the target radius in [m]
299  void setRadius( double _r ) { r = _r ; }
300 
301  // get the target radius in [m]
302  double getRadius() const { return r ; }
303 
304  void setPosition( JPosition3D& _pos ) { axis = JAxis3D(_pos,axis.getDirection()) ; }
305  const JPosition3D& getPosition() const { return axis.getPosition() ; }
306 
307  void setDirection( JVersor3D& _dir ) { axis = JAxis3D(axis.getPosition(),_dir) ; }
308  const JVersor3D& getDirection() const { return axis.getDirection() ; }
309 
310  protected:
311 
312  double r ;
314  } ;
315 
316  /**
317  Implementation of the JTargetModel class that represents a
318  single PMT on a DOM
319  **/
320  class JPMTTarget : public JTargetModel {
321  public :
322 
323  /**
324  Constructor
325 
326  _r is DOM radius, default = 0.2159 [m] (=17 inch diameter glass sphere)
327 
328  _alpha is the opening angle in radians
329  **/
330  JPMTTarget( JAxis3D _axis, double _alpha, double _r=0.2159 ) {
331  setRadius(_r) ;
332  axis = _axis ;
333  alpha = _alpha ;
334 
335  // minimal dot product of PMT direction vs hit direction on DOM
336  ctmin = cos(0.5*_alpha) ;
337  }
338 
339  /// return the PMT opening angle in radians
340  double getOpeningAngle() const { return alpha ; }
341 
342  double getEfficiency( JVersor3D dir ) const {
343  const double ct = axis.getDirection().getDot(dir) ;
344  if( ct > ctmin ) {
345  return 1.0 ;
346  } else {
347  return 0.0 ;
348  }
349  }
350 
351  double getEffectiveArea() const {
352  return 2.0*r*r*M_PI*(1-ctmin) ;
353  }
354 
355  double ctmin ;
356 
357  protected :
358 
359  double alpha ;
360  } ;
361 
362  /**
363  Implementation of the JTargetModel class that represents a
364  spherically symmetric target
365  **/
367  public :
368 
370 
371  double getEfficiency( JVersor3D dir ) const { return 1.0 ; }
372  } ;
373 
374  /**
375  Implementation of the JTargetModel class that represents a
376  directed target with a cos(theta) acceptance.
377  **/
378  class JCosineTarget : public JTargetModel {
379  public :
380 
381  /**
382  Constructor.
383  _target_dir is the direction where the photon has
384  the highest chance to be detected.
385  **/
386  JCosineTarget( JVersor3D _target_dir ) :
387  target_dir(_target_dir) {}
388 
389  double getEfficiency( JVersor3D dir ) const {
390  return max( dir.getDot(target_dir), 0.0 ) ;
391  }
392 
393  protected :
394 
396  } ;
397 
398 
399  /**
400  Implementation of the JScatteringModel interface with scattering
401  according to the Henyey-Greenstein function.
402  **/
404 
405  public:
406 
407  /**
408  Constructor.
409 
410  g is a parameter that determines how forward the scattering is.
411  g = 0 corresponds to isotropic scattering
412  g = 1 corresponds to fully forward scattering
413  **/
414  JHenyeyGreensteinScattering( double _lambda_scat, double _g ) {
415  lambda_scat = _lambda_scat ;
416  g = _g ;
417  }
418 
419  double getScatteringProbability( double ct ) {
420  //double ct = cos(dir.getTheta()) ;
421  double den = 1 + g*g - 2*g*ct ;
422  return (1-g*g)/(4*M_PI*den*sqrt(den)) ;
423  }
424 
426  double phi = gRandom->Uniform(0,2*M_PI) ;
427  double ct ;
428  if( g>0 ) {
429  double xi = gRandom->Uniform(0,1) ;
430  ct = 1.0/(2*g) * ( 1 + g*g - (1-g*g)*(1-g*g)/( (1-g+2*g*xi)*(1-g+2*g*xi) ) ) ;
431  } else {
432  ct = gRandom->Uniform(-1,1) ;
433  }
434  return JDirection3D( JAngle3D( acos(ct),phi) ) ;
435  }
436 
437  protected :
438 
439  double g ;
440 
441  } ;
442 
443  /**
444  Implementation of the JScatteringModel interface with
445  isotropic scattering
446  **/
448 
449  public:
450 
451  /**
452  Constructor.
453  **/
454  JIsotropicScattering( double _lambda_scat ) {
455  lambda_scat = _lambda_scat ;
456  }
457 
458  double getScatteringProbability( double ct ) {
459  return 1.0/(4.0*M_PI) ;
460  }
461 
463  double phi = gRandom->Uniform(0,2*M_PI) ;
464  double ct = gRandom->Uniform(-1,1) ;
465  return JDirection3D( JAngle3D( acos(ct),phi) ) ;
466  }
467  } ;
468 
469  /**
470  Implementation of the JScatteringModel interface with
471  Rayleigh scattering
472  **/
474 
475  public:
476 
477  /**
478  Constructor.
479 
480  note that one should never set a<0
481  **/
482  JRayleighScattering( double _lambda_scat, double _a ) {
483  lambda_scat = _lambda_scat ;
484  a = _a ;
485  if( fabs(a)<0 ) {
486  cerr << "Fatal error in initialization of JRayleighScattering: a<0 !" << endl ;
487  exit(1) ;
488  }
489  }
490 
491  double getScatteringProbability( double ct ) {
492  //double ct = cos(dir.getTheta()) ;
493  return (1+a*ct*ct)/(4.0*M_PI*(1+a/3.0)) ;
494  }
495 
497  double phi = gRandom->Uniform(0,2*M_PI) ;
498  double ct ;
499  if( a>0 ) {
500  double xi = gRandom->Uniform(0,1) ;
501  double d0 = -9/a ;
502  double d1 = 27.0*(1-2*xi)*(3+a)/a ;
503  double C = pow( 0.5*(d1 + sqrt(d1*d1-4*d0*d0*d0) ), 1.0/3.0 ) ;
504  ct = -1.0/(3.0)*(C+d0/C) ;
505  } else {
506  ct = gRandom->Uniform(-1,1) ;
507  }
508  return JDirection3D( JAngle3D( acos(ct),phi) ) ;
509  }
510 
511  protected:
512 
513  double a ;
514  } ;
515 
516  /**
517  Implementation of the JScatteringModel interface with
518  that combines two scattering models into one effective model.
519  **/
521 
522  public:
523 
524  /**
525  Constructor.
526 
527  Note that this only copies only the pointers to the JScatteringModel
528  instances. So if _sm1 or _sm2 are deleted or changed, it will also
529  affect the behaviour of this instance.
530  **/
531  JCombinedScattering( JScatteringModel* _sm1, JScatteringModel* _sm2 ) : sm1(_sm1), sm2(_sm2) { }
532 
533  virtual double getScatteringLength() {
534  // calculate effective scattering length
535  double l1 = sm1->getScatteringLength() ;
536  double l2 = sm2->getScatteringLength() ;
537  return 1.0 / ( 1.0/l1 + 1.0/l2 ) ;
538  }
539 
540  double getScatteringProbability( double ct ) {
541  double val1 = sm1->getScatteringProbability(ct) ;
542  double val2 = sm2->getScatteringProbability(ct) ;
543  double l1 = sm1->getScatteringLength() ;
544  double l2 = sm2->getScatteringLength() ;
545  double lambda = getScatteringLength() ;
546  // effective scattering probability
547  return lambda * ( val1/l1 + val2/l2 ) ;
548  }
549 
551  double l1 = sm1->getScatteringLength() ;
552  double l2 = sm2->getScatteringLength() ;
553  // probability to scatter with process 1
554  double P = l2/(l1+l2) ;
555  if( gRandom->Uniform(0,1)<P ) {
556  return sm1->generateDirection() ;
557  } else {
558  return sm2->generateDirection() ;
559  }
560  }
561 
562  protected:
563 
564  double a ;
567  } ;
568 
569  /**
570  Return the probability density for a photon path with the given
571  ingredients. These are
572  - a model of the source direction distribution
573  - a scattering model
574  - a target model
575  - the absorption length (note that this is allowed to be +INFTY)
576 
577  Multiply by dV1, dV2, etc. and the target cross-section sigma
578  to get a probability.
579 
580  The first (last) vertex of the photon path is used as the source
581  (target) position.
582  **/
583  double getPhotonPathProbabilityDensity( JPhotonPath& p, JSourceModel* src, JScatteringModel* sm, JTargetModel* trg, double lambda_abs ) {
584  //const bool verbose = true ;
585 
586  // calculate inverse scattering and absorption lengths
587  const double lambda_scat_inv = 1.0/sm->getScatteringLength() ;
588  double lambda_abs_inv = 0 ;
589  // if the absorption length is not infinite
590  if( lambda_abs == lambda_abs ) lambda_abs_inv = 1.0/lambda_abs ;
591 
592  double rho = 1 ;
593  const int nvert = p.n ; // number of vertices in the path
594  const int nscat = p.n-2 ; // number of scattering vertices in the path
595 
596  //if( verbose ) cout << "$ rho = " << rho << " (before calculating path segment lengths)" << endl ;
597 
598  // calculate path segment lengths
599  // index i corresponds to the path segment connecting
600  // vertex i and vertex i+1
601  double Ltotal = 0 ;
602  vector<double> lengths(nvert-1) ;
603  for( int i=0 ; i<nvert-1 ; ++i ) {
604  lengths[i] = (p[i+1]-p[i]).getLength() ;
605  Ltotal += lengths[i] ;
606  }
607 
608  //if( verbose ) cout << "$ Total length = " << Ltotal << " m" << endl ;
609 
610  // cosine of scattering angles (angle between previous segment and new segment)
611  // index i is the cosine of the angle of segment i+1 w.r.t segment i
612  vector<double> cts(nscat) ;
613 
614  // scattering angles
615  for( int i=0 ; i<nscat ; ++i ) {
616  JDirection3D seg1(p[i+1]-p[i]) ;
617  JDirection3D seg2(p[i+2]-p[i+1]) ;
618  cts[i] = seg1.getDot(seg2) ;
619  }
620 
621  // term for absorption length
622  // (=probability to not be absorbed within Ltotal)
623  rho *= exp(-Ltotal*lambda_abs_inv ) ;
624 
625  //if( verbose ) cout << "$ rho = " << rho << " (after absorption length)" << endl ;
626 
627  // term for emission profile
628  rho *= src->getEmissionProbability( JDirection3D(p[1]-p[0]) ) ;
629 
630  //if( verbose ) cout << "$ rho = " << rho << " (after emission probability)" << endl ;
631 
632  // term for target efficiency
633  if( trg->getRadius()>0 ) {
634  // finite target (spherical), the efficiency is a measure of how good the
635  // target is at a given spot
636  rho *= trg->getEfficiency( JDirection3D(p[nscat+1]-trg->getPosition()) ) ;
637  // this factor accounts for the angle of incidence
638  double ct = JDirection3D(p[nscat]-p[nscat+1]).getDot( JDirection3D(p[nscat+1]-trg->getPosition()) ) ;
639  rho *= max(0.0,ct) ;
640 
641  //if( verbose ) cout << "$ rho = " << rho << " (after target efficiency and angle of incidence)" << endl ;
642 
643  // check whether the photon path hits the sphere
644  // we do not need to consider the last vertex, because if that intersects the sphere
645  // the angle of incidence factor above would be 0.
646  if( p.n>2 ) {
647  JPhotonPath pshort(p) ;
648  pshort.resize( pshort.size()-1 ) ;
649  --pshort.n ;
650  if( pshort.hitsSphere(trg->getPosition(),trg->getRadius()) ) {
651  /*cout << "Path hits sphere at "
652  << JPosition3D(pshort.getSphereHitPosition(trg->getPosition(),trg->getRadius())-trg->getPosition())
653  << endl ;*/
654  rho *= 0 ;
655  }
656  }
657  } else {
658  rho *= trg->getEfficiency( JDirection3D(p[nscat+1]-p[nscat]) ) ;
659  }
660 
661  //if( verbose ) cout << "$ rho = " << rho << " (after target efficiency)" << endl ;
662 
663  // terms for scattering directions
664  for( int i=0 ; i<nscat ; ++i ) {
665  rho *= sm->getScatteringProbability(cts[i]) ;
666  }
667 
668  //if( verbose ) cout << "$ rho = " << rho << " (after scattering directions)" << endl ;
669 
670  // phase space terms for scattering to each of the scattering
671  // vertices
672  // this is the probability/dV to scatter into a volume element
673  // dV at a given distance
674  // (although it is missing the direction-dependent term, which is
675  // the emission probability or scattering probability)
676  for( int i=0; i<nscat; ++i ) {
677  // reminder: lengths[i] is the length from vertex i to vertex i+1
678  rho *= lambda_scat_inv / (lengths[i]*lengths[i]) * exp(-lengths[i]*lambda_scat_inv) ;
679  }
680 
681  //if( verbose ) cout << "$ rho = " << rho << " (after scattering phase space term)" << endl ;
682 
683  // phase space term for scattering to the target
684  // this is the probability/dsigma to hit a target at a given
685  // distance, without scattering along the way
686  rho *= exp(-lengths[nscat]*lambda_scat_inv) / ( lengths[nscat]*lengths[nscat] ) ;
687 
688  //if( verbose ) cout << "$ rho = " << rho << " (after target phase space term)" << endl ;
689 
690  // check for NaN
691  if( rho!=rho ) {
692  cerr << "NaN for probability density" << endl ;
693  cerr << "rho = " << rho << endl ;
694  exit(1) ;
695  }
696 
697  //if( verbose ) cout << "$ rho = " << rho << endl ;
698 
699  return rho ;
700  }
701 
702 }
703 
704 #endif
This include file serves the purpose of hiding ROOT dependencies and circumphere namespace problems w...
Properties of KM3NeT PMT and deep-sea water.
Data structure for angles in three dimensions.
Definition: JAngle3D.hh:35
Axis object.
Definition: JAxis3D.hh:41
Data structure for direction in three dimensions.
Definition: JDirection3D.hh:35
JDirection3D & rotate_back(const JRotation3D &R)
Rotate back.
double getDot(const JAngle3D &angle) const
Get dot product.
Data structure for position in three dimensions.
Definition: JPosition3D.hh:38
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:130
Rotation matrix.
Definition: JRotation3D.hh:114
Data structure for normalised vector in three dimensions.
Definition: JVersor3D.hh:28
double getDot(const JVersor3D &versor) const
Get dot product.
Definition: JVersor3D.hh:156
Implementation of the JScatteringModel interface with that combines two scattering models into one ef...
virtual double getScatteringLength()
JVersor3D generateDirection()
Return a randomly generated direction according to the scattering probability distribution.
double getScatteringProbability(double ct)
Return the probability density as a function of cos(theta)
JCombinedScattering(JScatteringModel *_sm1, JScatteringModel *_sm2)
Constructor.
Implementation of the JTargetModel class that represents a directed target with a cos(theta) acceptan...
JCosineTarget(JVersor3D _target_dir)
Constructor.
double getEfficiency(JVersor3D dir) const
Return the efficiency, which is defined as the probability that a photon with the given direction hit...
Implementation of the JSourceModel class that represents a directed source with a flat intensity dist...
JVersor3D getDirection() const
return the source direction
JVersor3D generateDirection()
Return a randomly generated direction according to the emission distribution.
JDirectedSource(JVersor3D _source_dir, double _alpha, bool _cutSurface=false, const JVersor3D _surface_normal=JVersor3D())
Constructor.
double getEmissionProbability(JVersor3D dir)
Return the probability density.
double getPhiMax(double theta)
For a given angle theta with respect to the beam direction, and a given surface normal,...
double getOpeningAngle() const
return the opening angle in radians
Implementation of the JScatteringModel interface with scattering according to the Henyey-Greenstein f...
JHenyeyGreensteinScattering(double _lambda_scat, double _g)
Constructor.
double getScatteringProbability(double ct)
Return the probability density as a function of cos(theta)
JVersor3D generateDirection()
Return a randomly generated direction according to the scattering probability distribution.
Implementation of the JScatteringModel interface with isotropic scattering.
JVersor3D generateDirection()
Return a randomly generated direction according to the scattering probability distribution.
double getScatteringProbability(double ct)
Return the probability density as a function of cos(theta)
JIsotropicScattering(double _lambda_scat)
Constructor.
Implementation of the JSourceModel class that represents an isotropic source.
JVersor3D generateDirection()
Return a randomly generated direction according to the emission distribution.
double getEmissionProbability(JVersor3D dir)
Return the probability density.
Implementation of the JTargetModel class that represents a spherically symmetric target.
double getEfficiency(JVersor3D dir) const
Return the efficiency, which is defined as the probability that a photon with the given direction hit...
Implementation of the JTargetModel class that represents a single PMT on a DOM.
JPMTTarget(JAxis3D _axis, double _alpha, double _r=0.2159)
Constructor.
double getOpeningAngle() const
return the PMT opening angle in radians
double getEffectiveArea() const
double getEfficiency(JVersor3D dir) const
Return the efficiency, which is defined as the probability that a photon with the given direction hit...
A photon path.
Definition: JPhotonPath.hh:38
bool hitsSphere(const JPosition3D &pos, double r)
Returns whether the photon path intersects a sphere of radius r at position pos.
Definition: JPhotonPath.hh:176
Implementation of the JScatteringModel interface with Rayleigh scattering.
JVersor3D generateDirection()
Return a randomly generated direction according to the scattering probability distribution.
double getScatteringProbability(double ct)
Return the probability density as a function of cos(theta)
JRayleighScattering(double _lambda_scat, double _a)
Constructor.
Virtual base class for a scattering model.
virtual double getScatteringLength()
virtual JVersor3D generateDirection()=0
Return a randomly generated direction according to the scattering probability distribution.
virtual double getScatteringProbability(double ct)=0
Return the probability density as a function of cos(theta)
Virtual base class for a light source.
void setPosition(JPosition3D &_pos)
virtual JVersor3D generateDirection()=0
Return a randomly generated direction according to the emission distribution.
virtual double getEmissionProbability(JVersor3D dir)=0
Return the probability density.
const JPosition3D & getPosition() const
Virtual base class for a light detector ("photon target").
virtual double getEffectiveArea()
Return the effective area, i.e.
const JVersor3D & getDirection() const
const JPosition3D & getPosition() const
virtual double getEfficiency(JVersor3D dir) const =0
Return the efficiency, which is defined as the probability that a photon with the given direction hit...
void setPosition(JPosition3D &_pos)
void setRadius(double _r)
void setDirection(JVersor3D &_dir)
double getScatteringLength(const double lambda)
Get scattering length.
Definition: Antares.hh:148
const double a
Definition: JQuadrature.cc:42
Auxiliary classes and methods for 3D geometrical objects and operations.
Definition: JAngle3D.hh:19
double getPhotonPathProbabilityDensity(JPhotonPath &p, JSourceModel *src, JScatteringModel *sm, JTargetModel *trg, double lambda_abs)
Return the probability density for a photon path with the given ingredients.
T pow(const T &x, const double y)
Power .
Definition: JMath.hh:97
static const double C
Physics constants.
data_type r[M+1]
Definition: JPolint.hh:868