Abstract base class for implementations of the JMarkovIntegrator interface, where the sample distribution is based on histograms filled from an ensemble of representative paths.  
 More...
#include <JMarkovIntegrator.hh>
 | 
|   | JMarkovEnsembleIntegrator (const vector< JPhotonPath > &_ensemble, JTargetModel *_trg, int _nbinsx, int _nbinsy, int _nbinsz) | 
|   | 
| virtual  | ~JMarkovEnsembleIntegrator () | 
|   | 
| virtual void  | writeHistograms ()=0 | 
|   | write the histograms filled from the ensemble  
  | 
|   | 
| vector< double >  | integrate (int N, int nscat, JSourceModel *src, JScatteringModel *sm, JTargetModel *trg, double lambda_abs) | 
|   | Integrate with N samples.  
  | 
|   | 
| vector< double >  | dummy_integrate (int N, int nscat, JSourceModel *src, JScatteringModel *sm, JTargetModel *trg, double lambda_abs) | 
|   | Integrate a test function with N samples.  
  | 
|   | 
| vector< JPhotonPath >  | get_diagnostic_ensemble (int N, int nscat, JSourceModel *src, JScatteringModel *sm, JTargetModel *trg, double lambda_abs) | 
|   | Return photon paths generated with the generatePath method.  
  | 
|   | 
Abstract base class for implementations of the JMarkovIntegrator interface, where the sample distribution is based on histograms filled from an ensemble of representative paths. 
Definition at line 198 of file JMarkovIntegrator.hh.
 
◆ JMarkovEnsembleIntegrator()
  
  
      
        
          | JMARKOV::JMarkovEnsembleIntegrator::JMarkovEnsembleIntegrator  | 
          ( | 
          const vector< JPhotonPath > & |           _ensemble,  | 
         
        
           | 
           | 
          JTargetModel * |           _trg,  | 
         
        
           | 
           | 
          int |           _nbinsx,  | 
         
        
           | 
           | 
          int |           _nbinsy,  | 
         
        
           | 
           | 
          int |           _nbinsz ) | 
         
       
   | 
  
inline   | 
  
 
Definition at line 202 of file JMarkovIntegrator.hh.
  203      
  205        cerr << "FATAL ERROR in JMarkov3DensembleIntegrator constructor. Ensemble size = 0." << endl ;
  206        exit(1) ;
  207      }
  208 
  209      
  210      for( vector<JPhotonPath>::iterator it=
ensemble.begin(); it!=
ensemble.end(); ++it ) {
 
  211        if( (*it)[0].getDistance(
ensemble.front()[0]) > 0 ) {
 
  212          cerr << "Fatal error in JMarkovEnsembleIntegrator, first vertex positions are not all the same." << endl ;
  213          exit(1) ;
  214        }
  215      }
  216      
  218    }
vector< JPhotonPath > ensemble
 
JPosition3D source_position
 
 
 
 
◆ ~JMarkovEnsembleIntegrator()
  
  
      
        
          | virtual JMARKOV::JMarkovEnsembleIntegrator::~JMarkovEnsembleIntegrator  | 
          ( | 
           | ) | 
           | 
         
       
   | 
  
inlinevirtual   | 
  
 
 
◆ writeHistograms()
  
  
      
        
          | virtual void JMARKOV::JMarkovEnsembleIntegrator::writeHistograms  | 
          ( | 
           | ) | 
           | 
         
       
   | 
  
pure virtual   | 
  
 
 
◆ generatePosition()
  
  
      
        
          | JPosition3D JMARKOV::JMarkovEnsembleIntegrator::generatePosition  | 
          ( | 
          int |           nscat,  | 
         
        
           | 
           | 
          int |           nv,  | 
         
        
           | 
           | 
          double & |           winv ) | 
         
       
   | 
  
inlineprotectedvirtual   | 
  
 
generate a position for vertex nv for the given number of scatterings, winv will be set to 1/weight 
Implements JMARKOV::JMarkovIntegrator.
Definition at line 228 of file JMarkovIntegrator.hh.
  228                                                                    {
  229      
  230      if( nv==0 ) {
  231        winv = 1.0 ;
  233      }
  234      
  235      if( nv == nscat+1 ) {
  238        }
  239        JPosition3D pos = 
target_gens[nscat]->getPosition() ;
 
  241        return pos ;
  242      } else {
  243        
  245      }
  246    }
virtual JPosition3D generateScatteringVertexPosition(int nscat, int nv, double &winv)=0
 
void initgenerators(int nscat)
initialize the position generators for a given number of scatterings (i.e. create histograms and fill...
 
vector< JSphereGenerator * > target_gens
 
 
 
 
◆ generateScatteringVertexPosition()
  
  
      
        
          | virtual JPosition3D JMARKOV::JMarkovEnsembleIntegrator::generateScatteringVertexPosition  | 
          ( | 
          int |           nscat,  | 
         
        
           | 
           | 
          int |           nv,  | 
         
        
           | 
           | 
          double & |           winv ) | 
         
       
   | 
  
protectedpure virtual   | 
  
 
 
◆ getIndex()
  
  
      
        
          | int JMARKOV::JMarkovEnsembleIntegrator::getIndex  | 
          ( | 
          int |           nv,  | 
         
        
           | 
           | 
          int |           nscat ) | 
         
       
   | 
  
inlineprotected   | 
  
 
return the internal index for a given number of scatterings and vertex number 
Definition at line 251 of file JMarkovIntegrator.hh.
  251                                      {
  252      
  253      
  254      
  255      
  256      
  257      
  258      
  259      return (nscat*(nscat-1))/2 + nv - 1 ;
  260    }
 
 
 
◆ initgenerators()
  
  
      
        
          | void JMARKOV::JMarkovEnsembleIntegrator::initgenerators  | 
          ( | 
          int |           nscat | ) | 
           | 
         
       
   | 
  
inlineprotected   | 
  
 
initialize the position generators for a given number of scatterings (i.e. create histograms and fill them from the ensemble) 
Definition at line 263 of file JMarkovIntegrator.hh.
  263                                     {
  264      
  266 
  267      
  269    }
virtual void init_scattering_vertex_generators(int nscat)=0
 
void init_target_generator(int nscat)
initialize the generator for a target vertex for a given number of scatterings
 
 
 
 
◆ init_target_generator()
  
  
      
        
          | void JMARKOV::JMarkovEnsembleIntegrator::init_target_generator  | 
          ( | 
          int |           nscat | ) | 
           | 
         
       
   | 
  
inlineprotected   | 
  
 
initialize the generator for a target vertex for a given number of scatterings 
Definition at line 272 of file JMarkovIntegrator.hh.
  272                                            {
  274 
  276        
  277        char hname[200] ;
  278        sprintf( hname, "htarget_nscat%i", nscat ) ;
  279        TH2D* ht = new TH2D(hname,";cos(#theta);#phi",100,-1,1,100,-M_PI,M_PI) ;
  280        for( vector<JPhotonPath>::iterator it=
ensemble.begin() ; it!=
ensemble.end() ; ++it ) {
 
  281          if( it->n-2 != nscat ) continue ;
  283          ht->Fill( dir.getDZ(), dir.getPhi() ) ;
  284        }
  286        delete ht ;
  287      } else {
  288        
  290      }
  291    }
const JPosition3D & getPosition() const
 
 
 
 
◆ init_scattering_vertex_generators()
  
  
      
        
          | virtual void JMARKOV::JMarkovEnsembleIntegrator::init_scattering_vertex_generators  | 
          ( | 
          int |           nscat | ) | 
           | 
         
       
   | 
  
protectedpure virtual   | 
  
 
 
◆ free_target_gens()
  
  
      
        
          | void JMARKOV::JMarkovEnsembleIntegrator::free_target_gens  | 
          ( | 
           | ) | 
           | 
         
       
   | 
  
inlineprotected   | 
  
 
 
◆ integrate()
Integrate with N samples. 
Returns a vector with the contribution to the integral of each sample. The mean of those values is the estimate of the result of the integral, while the distribution itself can be used to estimate the stability of the result. In this distribution, you want to avoid
- long tails (because they make the result unstable)
 
- small contributions (because it means that parts of the parameter space are being oversampled, so it is less efficient) This can be achieved by tuning the sample distribution to the problem at hand. 
 
Definition at line 130 of file JMarkovIntegrator.hh.
  130                                                                                                                                               {
  131    vector<double> contributions(N,-1) ;
  132 
  133    for( int i=0 ; i<N ; ++i ) {
  134      double winv ;
  136 
  138      contributions[i] = rho * winv ; 
  139    }
  140    return contributions ;
  141  }
virtual JPhotonPath generatePath(int nscat, double &winv)
Generate a random photon path with a given number of scatterings.
 
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.
 
 
 
 
◆ dummy_integrate()
Integrate a test function with N samples. 
This can be used as a sanity check for derived classes of JMarkovIntegrator.
The integral should yield 1 when the complete relevant part of the volume is taken into account. If it does not, it may be a sign that the implementation is not correct.
Returns a vector with the contribution to the integral of each sample. 
Definition at line 143 of file JMarkovIntegrator.hh.
  143                                                                                                                                                     {
  144    vector<double> contributions(N,-1) ;
  145 
  146    for( int i=0 ; i<N ; ++i ) {
  147      const double r = 10 ;
  148      double winv ;
  150      double rho = 1.0/(4.0/3.0*M_PI*r*r*r) ;
  151      if( p[1].getLength()>r ) rho = 0 ;
  152      contributions[i] = rho * winv ; 
  153    }
  154    return contributions ;
  155  }
 
 
 
◆ get_diagnostic_ensemble()
Return photon paths generated with the generatePath method. 
This can be used to identify the parts of parameter space that are over- or undersampled in a given problem so that the integrator may be optimized to handle those better. 
Definition at line 157 of file JMarkovIntegrator.hh.
  157                                                                                                                                                                  {
  158    vector<JPhotonPath> paths ;
  159 
  160    for( int i=0 ; i<N ; ++i ) {
  161      double winv ;
  163      paths.push_back(p) ;
  164    }
  165    return paths ;
  166  }
 
 
 
◆ generatePath()
  
  
      
        
          | virtual JPhotonPath JMARKOV::JMarkovIntegrator::generatePath  | 
          ( | 
          int |           nscat,  | 
         
        
           | 
           | 
          double & |           winv ) | 
         
       
   | 
  
inlineprotectedvirtualinherited   | 
  
 
Generate a random photon path with a given number of scatterings. 
winv must be set to the inverted probability density to generate this particular path. 
Definition at line 94 of file JMarkovIntegrator.hh.
   94                                                                {
   95      
   96      
   97      JPhotonPath p(nscat) ;
   98      double _winv = 1 ;
   99      for( int nv=0 ; nv<nscat+2 ; ++nv ) {
  100        double part_winv ;
  102        _winv *= part_winv ;
  103      }
  104      winv = _winv ;
  105      return p ;
  106    }
virtual JPosition3D generatePosition(int nscat, int nv, double &winv)=0
Generate a random position for vertex nv.
 
 
 
 
◆ target_gens
◆ ensemble
◆ nbinsx
  
  
      
        
          | int JMARKOV::JMarkovEnsembleIntegrator::nbinsx | 
         
       
   | 
  
protected   | 
  
 
 
◆ nbinsy
  
  
      
        
          | int JMARKOV::JMarkovEnsembleIntegrator::nbinsy | 
         
       
   | 
  
protected   | 
  
 
 
◆ nbinsz
  
  
      
        
          | int JMARKOV::JMarkovEnsembleIntegrator::nbinsz | 
         
       
   | 
  
protected   | 
  
 
 
◆ trg
◆ source_position
  
  
      
        
          | JPosition3D JMARKOV::JMarkovEnsembleIntegrator::source_position | 
         
       
   | 
  
protected   | 
  
 
 
The documentation for this class was generated from the following file: