class Teddy
Use MC techniques to apply a function \(f\) to a distribution $x$ with (non-diagonal) covariance matrix \(\mathbf{V}\).
Procedure:
- Diagonalise the covariance matrix
 
- Generate a multi-dimensional Gaussian event \(\delta\) distributed according to the (sqrt of the) eigenvalues \(\epsilon_i\) 
\[ \delta_i \sim \mathcal{N}(0,\sqrt{\epsilon_i}) \]
 
- Transform back to the original base (where the cov matrix is non-diagonal)
 
- Add rotated vector to the input distribution an apply the function \(f\): 
\[ \mathbf{y}_k = \]
 
- Go back to 2. until a satisfying statistics is reached, with 
Teddy::play() (the criterion of satisfaction and the loop are to be defined outside of the Teddy instance) 
- Get output histogram and covariance matrix with 
Teddy::buy()  
 
#include <Teddy.h>
 | 
|   | Teddy (const Eigen::VectorXd &, const Eigen::MatrixXd &, std::function< Eigen::VectorXd(const Eigen::VectorXd &)>=identity, int=42) | 
|   | 
|   | Teddy (const std::unique_ptr< TH1D > &, const std::unique_ptr< TH2D > &, std::function< Eigen::VectorXd(const Eigen::VectorXd &)>=identity, int=42) | 
|   | 
| void  | play () | 
|   | 
| std::pair< Eigen::VectorXd, Eigen::MatrixXd >  | buy () const | 
|   | 
| template<typename... Args>  | 
| std::pair< std::unique_ptr< TH1D >, std::unique_ptr< TH2D > >  | buy (TString name, Args... args) const | 
|   | 
 | 
| static Eigen::VectorXd  | H2V (const std::unique_ptr< TH1D > &) | 
|   | 
| static Eigen::MatrixXd  | H2M (const std::unique_ptr< TH2D > &) | 
|   | 
| static int  | GuessOutDim (int, std::function< Eigen::VectorXd(const Eigen::VectorXd &)>) | 
|   | 
| static void  | SanityCheck (const Eigen::VectorXd &, Eigen::VectorXd &, double=0) | 
|   | 
◆ Teddy() [1/2]
      
        
          | Teddy  | 
          ( | 
          const Eigen::VectorXd &  | 
          ,  | 
        
        
           | 
           | 
          const Eigen::MatrixXd &  | 
          ,  | 
        
        
           | 
           | 
          std::function< Eigen::VectorXd(const Eigen::VectorXd &)>  | 
           = identity,  | 
        
        
           | 
           | 
          int  | 
           = 42  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
The rotation and the eigenvalues are determined from the input covariance matrix. Default function is the identity (useful for testing). 
 
 
◆ Teddy() [2/2]
      
        
          | Teddy  | 
          ( | 
          const std::unique_ptr< TH1D > &  | 
          ,  | 
        
        
           | 
           | 
          const std::unique_ptr< TH2D > &  | 
          ,  | 
        
        
           | 
           | 
          std::function< Eigen::VectorXd(const Eigen::VectorXd &)>  | 
           = identity,  | 
        
        
           | 
           | 
          int  | 
           = 42  | 
        
        
           | 
          ) | 
           |  | 
        
      
 
The rotation and the eigenvalues are determined from the input covariance matrix. Default function is the identity (useful for testing). 
 
 
◆ buy() [1/2]
      
        
          | pair< VectorXd, MatrixXd > buy  | 
          ( | 
           | ) | 
           const | 
        
      
 
Get output distributions, can be called several times. 
  155     VectorXd outputDist = 
sum*(1./
N);
 
  157     MatrixXd avSq = 
sum2*(1./
N);
 
  158     MatrixXd outputCov = avSq - outputDist * outputDist.transpose();
 
  160     return {outputDist, outputCov};
 
 
 
 
◆ buy() [2/2]
  
  
      
        
          | std::pair<std::unique_ptr<TH1D>,std::unique_ptr<TH2D> > buy  | 
          ( | 
          TString  | 
          name,  | 
         
        
           | 
           | 
          Args...  | 
          args  | 
         
        
           | 
          ) | 
           |  const | 
         
       
   | 
  
inline   | 
  
 
Get output distributions, can be called several times. Input histrogams are defined from scratch and filled. The statistical uncertainty of the 1D histogram correspondings to the square-root of the diagonal elements of the covariance matrix 
 
  108         using namespace Eigen;
 
  110         TAxis * axis = 
new TAxis(
args...);
 
  111         int ndiv = axis->GetNbins();
 
  112         double * edges = 
new double[ndiv];
 
  113         axis->GetLowEdge(edges);
 
  114         edges[ndiv-1] = axis->GetBinUpEdge(ndiv);
 
  115         auto h   = make_unique<TH1D>(
name         , 
name, ndiv-1, edges);
 
  116         auto cov = make_unique<TH2D>(
name + 
"_cov", 
name, ndiv-1, edges, ndiv-1, edges);
 
  123             double content = x(i),
 
  124                    error = sqrt(err(i,i));
 
  126             h->SetBinContent(i+1, content);
 
  127             h->SetBinError  (i+1, error  );
 
  130                 cov->SetBinContent(i+1,j+1, err(i,j));
 
  133         return make_pair<unique_ptr<TH1D>,unique_ptr<TH2D>>(std::move(h),std::move(
cov));
 
 
 
 
◆ GuessOutDim()
  
  
      
        
          | int GuessOutDim  | 
          ( | 
          int  | 
          ,  | 
         
        
           | 
           | 
          std::function< Eigen::VectorXd(const Eigen::VectorXd &)>  | 
            | 
         
        
           | 
          ) | 
           |  | 
         
       
   | 
  
staticprivate   | 
  
 
Guess size of the output distribution by calling once the function with dummy input. 
   48     if (
verbose) cout << __FILE__ << 
':' << __LINE__ << 
'\t' << __func__ << endl;
 
 
 
 
◆ H2M()
  
  
      
        
          | MatrixXd H2M  | 
          ( | 
          const std::unique_ptr< TH2D > &  | 
          h | ) | 
           | 
         
       
   | 
  
staticprivate   | 
  
 
conversion from (2D) histogram to matrix 
   36     if (
verbose) cout << __FILE__ << 
':' << __LINE__ << 
'\t' << __func__ << endl;
 
   37     int sizx = h->GetNbinsX(),
 
   38         sizy = h->GetNbinsY();
 
   39     MatrixXd mat(sizx,sizy);
 
   40     for(
int j = 1; j <= sizx; ++j)
 
   41     for(
int i = 1; i <= sizy; ++i)
 
   42         mat(i-1,j-1) = h->GetBinContent(j,i);
 
 
 
 
◆ H2V()
  
  
      
        
          | VectorXd H2V  | 
          ( | 
          const std::unique_ptr< TH1D > &  | 
          h | ) | 
           | 
         
       
   | 
  
staticprivate   | 
  
 
conversion from (1D) histogram to vector 
   27     if (
verbose) cout << __FILE__ << 
':' << __LINE__ << 
'\t' << __func__ << endl;
 
   28     VectorXd v(h->GetNbinsX());
 
   29     for(
int i = 1; i <= h->GetNbinsX(); ++i)
 
   30         v(i-1) = h->GetBinContent(i);
 
 
 
 
◆ play()
Needs to be called from a loop, but requires no argument. 
  139     VectorXd x = 
vec + delta;
 
  142     VectorXd y = 
func(x); 
 
  146     MatrixXd square = y * y.transpose(); 
 
 
 
 
◆ SanityCheck()
  
  
      
        
          | void SanityCheck  | 
          ( | 
          const Eigen::VectorXd &  | 
          ,  | 
         
        
           | 
           | 
          Eigen::VectorXd &  | 
          ,  | 
         
        
           | 
           | 
          double  | 
           = 0  | 
         
        
           | 
          ) | 
           |  | 
         
       
   | 
  
staticprivate   | 
  
 
Check the number of negative eigenvalues and fail if any is found. 
   56     int nEigv = 0, nNegEigVals = 0, nEmptyBins = 0;
 
   57     for (
long k = 0; k < 
eigVals.size(); ++k) {
 
   64             cout << __FILE__ << 
':' << __LINE__ << 
'\t' << k << 
' ' << 
eigVals(k) << 
'\n';
 
   73         BOOST_THROW_EXCEPTION( invalid_argument(Form(
"%d negative eigenvalues (%d empty bins in input dist)",
 
   74                                                     nNegEigVals,                nEmptyBins)) );
 
 
 
 
◆ cov
  
  
      
        
          | const Eigen::MatrixXd cov | 
         
       
   | 
  
private   | 
  
 
 
◆ eigVals
eigenvalues of input covariance matrix 
 
 
◆ func
  
  
      
        
          | std::function<Eigen::VectorXd(const Eigen::VectorXd &)> func | 
         
       
   | 
  
private   | 
  
 
transformation defining output as a function of the input 
 
 
◆ nNewBins
number of bins of output distributions 
 
 
◆ nOldBins
number of bins of input distributions 
 
 
◆ random
◆ rotation
rotation from diagonal basis to original basis 
 
 
◆ sigma
Gaussian widths in diagonal base. 
 
 
◆ sum
sum of events, increased at each call and from which output distribution is estimated 
 
 
◆ sum2
sum of tensor products of each event, increased at each call and from which output covariance is estimated 
 
 
◆ vec
  
  
      
        
          | const Eigen::VectorXd vec | 
         
       
   | 
  
private   | 
  
 
 
◆ verbose
The documentation for this class was generated from the following files:
- /builds/cms-analysis/general/DasAnalysisSystem/Core/Installer/Core/Uncertainties/interface/Teddy.h
 
- /builds/cms-analysis/general/DasAnalysisSystem/Core/Installer/Core/Uncertainties/src/Teddy.cc