Go to the documentation of this file.    1 #ifndef H_MARKOV_GENERATOR 
    2 #define H_MARKOV_GENERATOR 
   67     JUniformGenerator( 
double xmin, 
double ymin, 
double zmin, 
double xmax, 
double ymax, 
double zmax ) : posmin(xmin,ymin,zmin), posmax(xmax,ymax,zmax) {
 
   80       double x = gRandom->Uniform( posmin.getX(), posmax.getX() ) ;
 
   81       double y = gRandom->Uniform( posmin.getY(), posmax.getY() ) ;
 
   82       double z = gRandom->Uniform( posmin.getZ(), posmax.getZ() ) ;
 
   87       if( pos.
getX() < posmin.getX() ||  pos.
getX() > posmax.getX() ) 
return 0.0 ;
 
   88       if( pos.
getY() < posmin.getY() ||  pos.
getY() > posmax.getY() ) 
return 0.0 ;
 
   89       if( pos.
getZ() < posmin.getZ() ||  pos.
getZ() > posmax.getZ() ) 
return 0.0 ;
 
   96       V = (posmax.getX()-posmin.getX())*(posmax.getY()-posmin.getY())*(posmax.getZ()-posmin.getZ()) ;
 
  121       double rcube = gRandom->Uniform(0,Rcube) ;
 
  122       double r = pow(rcube,1.0/3.0) ;
 
  124       gRandom->Sphere(x,y,z,
r) ;
 
  160       double xi = gRandom->Uniform() ;
 
  161       double r = -lambda*log(1-xi) ;
 
  163       gRandom->Sphere(x,y,z,
r) ;
 
  169       return lambda*exp(-
r/lambda)/(
r*
r*4*M_PI) ;
 
  193       const int nbins = 1000 ;
 
  194       const double xmax = 100 ;
 
  195       const int nsamples = 10000000 ;
 
  196       h = 
new TH1F(
"h",
"",nbins,0,xmax) ;
 
  197       for( 
int i=0 ; i<nsamples ; ++i ) {
 
  203       h->Scale( 1.0/h->Integral(
"width") ) ;
 
  211       double r = h->GetRandom() ;
 
  213       gRandom->Sphere(x,y,z,
r) ;
 
  219       Int_t bin = h->FindBin(
r) ;
 
  220       return h->GetBinContent(bin)/(
r*
r) ;
 
  246       double r = gRandom->Uniform(0,R) ;
 
  248       gRandom->Sphere(x,y,z,
r) ;
 
  286       r1(_r1), r2(_r2), L(_L) 
 
  289       ximin = getIntegrand(r1) ;
 
  290       ximax = getIntegrand(r2) ;
 
  301       gRandom->Sphere(x,y,z,
r) ;
 
  307       double xi = gRandom->Uniform(ximin,ximax) ;
 
  309       return getInvertedIntegrand(xi) ;
 
  313       return exp(-
r/L) / 
C ;
 
  320       return -4*M_PI*L*exp(-
r/L)*(2*L*L+2*L*
r+
r*
r) ;
 
  325       const double precision = 1e-10 ;
 
  328       double vl = getIntegrand(rl)-x ;
 
  330       double vr = getIntegrand(rr)-x ;
 
  334       while( rr-rl > precision ) {
 
  336         vc = getIntegrand(rc)-x ;
 
  337         if( (vc>0) == (vr>0) ) {
 
  371       C = 1.0 / ( 2*M_PI*sqrt(2*M_PI) * sigma * sigma * sigma ) ;
 
  375       double x = gRandom->Gaus(0,sigma) ;
 
  376       double y = gRandom->Gaus(0,sigma) ;
 
  377       double z = gRandom->Gaus(0,sigma) ;
 
  420       if( gRandom->Uniform()<
p1 ) {
 
  421         return g1->getPosition() ;
 
  423         return g2->getPosition() ;
 
  428       double w1 = 
g1->getWeight(pos) ;
 
  429       double w2 = g2->getWeight(pos) ;
 
  430       return p1*w1 + p2*w2 ;
 
  464       gsub(
c1,
g1,c2,g2), g(
c1+c2,&gsub,c3,g3) {}
 
  494       g(_g), shift(_shift) {}
 
  501       return g->getWeight( pos-shift ) ;
 
  528       hx = (TH1*) _hx->Clone( _hx->GetName() ) ;
 
  529       hy = (TH1*) _hy->Clone( _hy->GetName() ) ;
 
  530       hz = (TH1*) _hz->Clone( _hz->GetName() ) ;
 
  532       hx->Scale( 1.0/hx->Integral(
"width") ) ;
 
  533       hy->Scale( 1.0/hy->Integral(
"width") ) ;
 
  534       hz->Scale( 1.0/hz->Integral(
"width") ) ;
 
  544       double x = hx->GetRandom() ;
 
  545       double y = hy->GetRandom() ;
 
  546       double z = hz->GetRandom() ;
 
  554       bin = hx->GetXaxis()->FindBin( pos.
getX() ) ;
 
  555       if( bin==0 || bin==hx->GetNbinsX()+1 ) 
return 0.0 ; 
 
  556       w *= hx->GetBinContent(bin) ;
 
  558       bin = hy->GetXaxis()->FindBin( pos.
getY() ) ;
 
  559       if( bin==0 || bin==hy->GetNbinsX()+1 ) 
return 0.0 ; 
 
  560       w *= hy->GetBinContent(bin) ;
 
  562       bin = hz->GetXaxis()->FindBin( pos.
getZ() ) ;
 
  563       if( bin==0 || bin==hz->GetNbinsX()+1 ) 
return 0.0 ; 
 
  564       w *= hz->GetBinContent(bin) ;
 
  597         h = (TH2*) _h->Clone( _h->GetName() ) ;
 
  599         if( h->GetXaxis()->GetXmin()<-1 || h->GetXaxis()->GetXmax()>1 ) {
 
  600           cerr << 
"FATAL ERROR in JSphereGenerator. Invalid x-axis range." << endl ;
 
  603         if( h->GetYaxis()->GetXmin() < -M_PI || h->GetYaxis()->GetXmax() > M_PI ) {
 
  604           cerr << 
"FATAL ERROR in JSphereGenerator. Invalid y-axis range." << endl ;
 
  608         h->Scale( 1.0/h->Integral(
"width") ) ;
 
  609         h->SetOption(
"colz") ;
 
  618       if( 
r==0 ) 
return x0 ;
 
  620       double costheta, phi ;
 
  621       h->GetRandom2(costheta,phi) ;
 
  626       if( 
r==0 ) 
return 1 ;
 
  629       const double tolerance = 1e-5 ;
 
  635       double ct = dir.
getDZ() ;
 
  636       double phi = dir.
getPhi() ;
 
  638       if( ct < h->GetXaxis()->GetXmin() || ct >=h->GetXaxis()->GetXmax() ) 
return 0.0 ;
 
  639       if( phi< h->GetYaxis()->GetXmin() || phi>=h->GetYaxis()->GetXmax() ) 
return 0.0 ;
 
  641       Int_t bin = h->FindBin(ct,phi) ;
 
  642       return h->GetBinContent(bin) / (
r*
r) ;
 
  667       h = (TH3*) _h->Clone( _h->GetName() ) ;
 
  669       h->Scale( 1.0/h->Integral(
"width") ) ;
 
  678       h->GetRandom3(x,y,z) ;
 
  684       if( pos.
getX()<h->GetXaxis()->GetXmin() || pos.
getX()>=h->GetXaxis()->GetXmax() ) 
return 0.0 ;
 
  685       if( pos.
getY()<h->GetYaxis()->GetXmin() || pos.
getY()>=h->GetYaxis()->GetXmax() ) 
return 0.0 ;
 
  686       if( pos.
getZ()<h->GetZaxis()->GetXmin() || pos.
getZ()>=h->GetZaxis()->GetXmax() ) 
return 0.0 ;
 
  688       Int_t bin = h->FindBin( pos.
getX(), pos.
getY(), pos.
getZ() ) ;
 
  689       return h->GetBinContent(bin) ;
 
 
double getDZ() const
Get z direction.
 
Implementation of the JGenerator interface.
 
JPosition3D getPosition()
Return a randomly generated position.
 
double getLengthSquared() const
Get length squared.
 
Implementation of the JGenerator interface.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
Double_t g1(const Double_t x)
Function.
 
Abstract interface for the generation of points in 3D space.
 
JPosition3D getPosition()
Return a randomly generated position.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
double getZ() const
Get z position.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
JPosition3D getPosition()
Return a randomly generated position.
 
JSphereGenerator(const JPosition3D &_x0, double _r=0, TH2 *_h=NULL)
Constructor.
 
TCanvas * c1
Global variables to handle mouse events.
 
JExponentialGenerator(double _r1, double _r2, double _L)
Constructor.
 
Implementation of the JGenerator interface.
 
double getDistance(const JVector3D &pos) const
Get distance to point.
 
JSingularityGenerator(double _R, JPosition3D _x0)
Constructor.
 
JPosition3D getPosition()
Return a randomly generated position.
 
Data structure for normalised vector in three dimensions.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
double getIntegrand(double r)
 
JMagicalDistribution(unsigned int _N, double _lambda)
 
JPosition3D getPosition()
Return a randomly generated position.
 
JPosition3D getPosition()
Return a randomly generated position.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
 
JGaussianGenerator(double _sigma)
Constructor.
 
JHistGenerator(TH1 *_hx, TH1 *_hy, TH1 *_hz)
Constructor.
 
Implementation of the JGenerator interface.
 
JExpRsqInvGenerator(double _lambda)
Constructor.
 
Implementation of the JGenerator interface.
 
Data structure for vector in three dimensions.
 
JShiftedGenerator(JGenerator *_g, JPosition3D _shift)
Constructor.
 
J3DhistGenerator(TH3 *_h)
Constructor.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
Data structure for position in three dimensions.
 
Implementation of the JGenerator interface.
 
Implementation of the JGenerator interface.
 
Implementation of the JGenerator interface.
 
double getInvertedIntegrand(double x)
return value y such that getIntegrand(y) = x
 
The 'magical distributions' are a class of distributions.
 
JBallGenerator(double _R)
Constructor.
 
JCombinedGenerator(double _c1, JGenerator *_g1, double _c2, JGenerator *_g2)
Constructor.
 
double getWeight(double r)
 
JPosition3D getPosition()
Return a randomly generated position.
 
JPosition3D getPosition()
Return a randomly generated position.
 
Data structure for angles in three dimensions.
 
const JPosition3D & getPosition() const
Get position.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
JPosition3D getPosition(const Vec &v)
Get position.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
Implementation of the JGenerator interface.
 
double getY() const
Get y position.
 
JPosition3D getPosition()
Return a randomly generated position.
 
Implementation of the JGenerator interface.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
Auxiliary classes and methods for 3D geometrical objects and operations.
 
double getWeight(JPosition3D pos)
return the weight (=probability density dP/dV) for the given position.
 
JTripleGenerator(double c1, JGenerator *g1, double c2, JGenerator *g2, double c3, JGenerator *g3)
Constructor.
 
double getX() const
Get x position.
 
JPosition3D getPosition()
Return a randomly generated position.
 
JPosition3D getPosition()
Return a randomly generated position.
 
JPosition3D getPosition()
Return a randomly generated position.
 
double getLength() const
Get length.
 
Implementation of the JGenerator interface.
 
double getPhi() const
Get phi angle.