1 #ifndef H_J_SCATTERING_MODEL 
    2 #define H_J_SCATTERING_MODEL 
  127       return 1.0/(4.0*M_PI) ;
 
  132       gRandom->Sphere(dx,dy,dz,1) ;
 
  155       source_dir(_source_dir), ctmin(cos(0.5*_alpha)), alpha(_alpha), cutSurface(_cutSurface), surface_normal(_surface_normal) {
 
  162       delta = 0.5*M_PI - acos(surface_normal.getDot(source_dir)) ;
 
  163       omega = 2.0*M_PI*(1.0-ctmin) ; 
 
  166         if( delta>0.0 && delta>0.5*_alpha )  {
 
  169         } 
else if( delta<=0 && fabs(delta)>=0.5*alpha ) {
 
  173           cerr << 
"FATAL ERROR in JDirectedSource: beam points completely into the surface." << endl ;
 
  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 ;
 
  193           omega_absorbed *= __dct ;
 
  194           omega -= omega_absorbed ;
 
  198             omega -= 2.0 * M_PI * (1.0-cos(delta)) ;
 
  211       double ct = dir.
getDot(source_dir) ;
 
  213         if( cutSurface && dir.
getDot(surface_normal)<0 ) 
return 0.0 ;
 
  221         double phi = gRandom->Uniform(0,2*M_PI) ;
 
  222         double ct  = gRandom->Uniform(ctmin,1) ;
 
  226         if( !cutSurface ) 
return dir ;
 
  227         if( dir.
getDot(surface_normal)>0 ) 
return dir ;
 
  244         if( theta<delta ) 
return M_PI ;
 
  245         if( theta+delta > M_PI ) 
return 0.0 ;
 
  248         if( theta<fabs(delta) ) 
return 0.0 ;
 
  249         if( theta+fabs(delta) > M_PI ) 
return M_PI ;
 
  251       double phimax = acos(-tan(delta)/tan(theta)) ;
 
  336       ctmin = cos(0.5*_alpha) ;
 
  343       const double ct = axis.getDirection().getDot(dir) ;
 
  352       return 2.0*
r*
r*M_PI*(1-ctmin) ; 
 
  387       target_dir(_target_dir) {}
 
  390       return max( dir.
getDot(target_dir), 0.0 ) ;
 
  415       lambda_scat  = _lambda_scat ;
 
  421       double den = 1 + g*g - 2*g*ct ;
 
  422       return (1-g*g)/(4*M_PI*den*sqrt(den)) ;
 
  426       double phi = gRandom->Uniform(0,2*M_PI) ;
 
  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) ) ) ;
 
  432         ct = gRandom->Uniform(-1,1) ;
 
  455       lambda_scat  = _lambda_scat ;
 
  459       return 1.0/(4.0*M_PI) ;
 
  463       double phi = gRandom->Uniform(0,2*M_PI) ;
 
  464       double ct  = gRandom->Uniform(-1,1) ;
 
  483       lambda_scat  = _lambda_scat ;
 
  486         cerr << 
"Fatal error in initialization of JRayleighScattering: a<0 !" << endl ;
 
  493       return  (1+
a*ct*ct)/(4.0*M_PI*(1+
a/3.0)) ;
 
  497       double phi = gRandom->Uniform(0,2*M_PI) ;
 
  500         double xi  = gRandom->Uniform(0,1) ;
 
  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) ;
 
  506         ct = gRandom->Uniform(-1,1) ;
 
  535       double l1 = sm1->getScatteringLength() ;
 
  536       double l2 = sm2->getScatteringLength() ;
 
  537       return 1.0 / ( 1.0/l1 + 1.0/l2 ) ;
 
  541       double val1 = sm1->getScatteringProbability(ct) ;
 
  542       double val2 = sm2->getScatteringProbability(ct) ;
 
  543       double l1 = sm1->getScatteringLength() ;
 
  544       double l2 = sm2->getScatteringLength() ;
 
  547       return lambda * ( val1/l1 + val2/l2 ) ;
 
  551       double l1 = sm1->getScatteringLength() ;
 
  552       double l2 = sm2->getScatteringLength() ;
 
  554       double P = l2/(l1+l2) ; 
 
  555       if( gRandom->Uniform(0,1)<P ) {
 
  556         return sm1->generateDirection() ;
 
  558         return sm2->generateDirection() ;
 
  588     double lambda_abs_inv = 0 ;
 
  590     if( lambda_abs == lambda_abs ) lambda_abs_inv = 1.0/lambda_abs ;
 
  593     const int nvert = p.
n   ; 
 
  594     const int nscat = p.
n-2 ; 
 
  603     for( 
int i=0 ; i<nvert-1 ; ++i ) {
 
  604       lengths[i] = (p[i+1]-p[i]).getLength() ;
 
  605       Ltotal += lengths[i] ;
 
  615     for( 
int i=0 ; i<nscat ; ++i ) {
 
  618       cts[i] = seg1.
getDot(seg2) ;
 
  623     rho *= exp(-Ltotal*lambda_abs_inv ) ;
 
  648         pshort.resize( pshort.size()-1 ) ;
 
  664     for( 
int i=0 ; i<nscat ; ++i ) {
 
  676     for( 
int i=0; i<nscat; ++i ) {
 
  678       rho *= lambda_scat_inv / (lengths[i]*lengths[i]) * exp(-lengths[i]*lambda_scat_inv) ;
 
  686     rho *= exp(-lengths[nscat]*lambda_scat_inv) / ( lengths[nscat]*lengths[nscat] ) ;
 
  692       cerr << 
"NaN for probability density" << endl ;
 
  693       cerr << 
"rho = " << rho << endl ;
 
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.
 
Data structure for direction in three dimensions.
 
JDirection3D & rotate_back(const JRotation3D &R)
Rotate back.
 
double getDot(const JAngle3D &angle) const
Get dot product.
 
Data structure for position in three dimensions.
 
const JPosition3D & getPosition() const
Get position.
 
Data structure for normalised vector in three dimensions.
 
double getDot(const JVersor3D &versor) const
Get dot product.
 
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...
 
bool hitsSphere(const JPosition3D &pos, double r)
Returns whether the photon path intersects a sphere of radius r at position pos.
 
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 ~JScatteringModel()
 
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.
 
Auxiliary classes and methods for 3D geometrical objects and operations.
 
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 .
 
static const double C
Physics constants.