BayesianAnalysisToolkit  0.9.3
Classes | Public Member Functions | Protected Attributes | List of all members
BCMVCombination Class Reference

#include <BCMVCombination.h>

Inheritance diagram for BCMVCombination:
Inheritance graph
[legend]
Collaboration diagram for BCMVCombination:
Collaboration graph
[legend]

Classes

struct  NuisanceParameter
 

Public Member Functions

 BCMVCombination ()
 
 ~BCMVCombination ()
 
void AddObservable (std::string name, double min, double max)
 
void AddUncertainty (std::string name)
 
void AddMeasurement (std::string name, std::string observable, double value, std::vector< double > uncertainties)
 
int GetNObservables ()
 
int GetNUncertainties ()
 
int GetNMeasurements ()
 
int GetNActiveMeasurements ()
 
BCMVCUncertaintyGetUncertainty (int index)
 
BCMVCMeasurementGetMeasurement (int index)
 
TMatrixD GetCovarianceMatrix ()
 
TMatrixD GetBLUEWeights ()
 
TVectorD GetBLUECentralValues ()
 
TVectorD GetBLUEUncertainties ()
 
TVectorD GetBLUEUncertainties (int index)
 
TMatrixD GetBLUECovarianceMatrix ()
 
TMatrixD GetBLUECovarianceMatrix (int index)
 
TMatrixD GetBLUECorrelationMatrix (int index)
 
TMatrixD GetBLUECorrelationMatrix ()
 
std::vector< int > GetVectorObservable ()
 
TVectorD GetVectorMeasurements ()
 
int GetIndexMeasurement (std::string measurement, std::string observable)
 
int GetIndexUncertainty (std::string name)
 
int GetIndexObservable (std::string name)
 
int ReadInput (std::string filename)
 
void CalculateCorrelationMatrix (int index)
 
void CalculateCovarianceMatrix (std::vector< double > nuisance=std::vector< double >(0))
 
void CalculateHelperVectors ()
 
bool PositiveDefinite ()
 
void CalculateBLUE ()
 
void PrepareAnalysis ()
 
void PrintBLUEResults (std::string filename)
 
void PrintMatrix (TMatrixT< double > &matrix, std::string name="matrix")
 
void PrintVector (TVectorD &vector, std::string name="vector")
 
double LogLikelihood (const std::vector< double > &parameters)
 
- Public Member Functions inherited from BCModel
 BCModel (const char *name="model")
 
 BCModel (const BCModel &bcmodel)
 
virtual ~BCModel ()
 
BCModeloperator= (const BCModel &bcmodel)
 
const std::string & GetName () const
 
double GetModelAPrioriProbability () const
 
double GetModelAPosterioriProbability () const
 
BCDataSetGetDataSet () const
 
BCDataPointGetDataPointLowerBoundaries () const
 
BCDataPointGetDataPointUpperBoundaries () const
 
double GetDataPointLowerBoundary (unsigned int index) const
 
double GetDataPointUpperBoundary (unsigned int index) const
 
bool GetFlagBoundaries () const
 
unsigned GetNDataPoints () const
 
BCDataPointGetDataPoint (unsigned int index) const
 
void SetName (const char *name)
 
void SetModelAPrioriProbability (double probability)
 
void SetModelAPosterioriProbability (double probability)
 
virtual int AddParameter (BCParameter *parameter)
 
void SetDataSet (BCDataSet *dataset)
 
void SetSingleDataPoint (BCDataPoint *datapoint)
 
void SetSingleDataPoint (BCDataSet *dataset, unsigned int index)
 
void SetDataBoundaries (unsigned int index, double lowerboundary, double upperboundary, bool fixed=false)
 
void SetDataPointLowerBoundaries (BCDataPoint *datasetlowerboundaries)
 
void SetDataPointUpperBoundaries (BCDataPoint *datasetupperboundaries)
 
void SetDataPointLowerBoundary (int index, double lowerboundary)
 
void SetDataPointUpperBoundary (int index, double upperboundary)
 
int SetPrior (int index, TF1 *f)
 
int SetPrior (const char *name, TF1 *f)
 
int SetPriorDelta (int index, double value)
 
int SetPriorDelta (const char *name, double value)
 
int SetPriorGauss (int index, double mean, double sigma)
 
int SetPriorGauss (const char *name, double mean, double sigma)
 
int SetPriorGauss (int index, double mean, double sigmadown, double sigmaup)
 
int SetPriorGauss (const char *name, double mean, double sigmadown, double sigmaup)
 
int SetPrior (int index, TH1 *h, bool flag=false)
 
int SetPrior (const char *name, TH1 *h, bool flag=false)
 
int SetPriorConstant (int index)
 
int SetPriorConstant (const char *name)
 
int SetPriorConstantAll ()
 
void Copy (const BCModel &bcmodel)
 
double APrioriProbability (const std::vector< double > &parameters)
 
virtual double LogAPrioriProbability (const std::vector< double > &parameters)
 
virtual double Likelihood (const std::vector< double > &params)
 
double ProbabilityNN (const std::vector< double > &params)
 
double LogProbabilityNN (const std::vector< double > &parameters)
 
double Probability (const std::vector< double > &parameter)
 
double LogProbability (const std::vector< double > &parameter)
 
virtual double SamplingFunction (const std::vector< double > &parameters)
 
double Eval (const std::vector< double > &parameters)
 
virtual double LogEval (const std::vector< double > &parameters)
 
virtual void CorrelateDataPointValues (std::vector< double > &x)
 
double GetPvalueFromChi2 (const std::vector< double > &par, int sigma_index)
 
double GetPvalueFromKolmogorov (const std::vector< double > &par, int index)
 
double GetPvalueFromChi2NDoF (std::vector< double > par, int sigma_index)
 
BCH1DCalculatePValue (std::vector< double > par, bool flag_histogram=false)
 
double GetPValue ()
 
double GetPValueNDoF ()
 
double GetChi2NDoF ()
 
std::vector< double > GetChi2Runs (int dataIndex, int sigmaIndex)
 
void SetGoFNIterationsMax (int n)
 
void SetGoFNIterationsRun (int n)
 
void SetGoFNChains (int n)
 
double HessianMatrixElement (const BCParameter *parameter1, const BCParameter *parameter2, std::vector< double > point)
 
void PrintSummary ()
 
void PrintResults (const char *file)
 
void PrintShortFitSummary (int chi2flag=0)
 
void PrintHessianMatrix (std::vector< double > parameters)
 
virtual double CDF (const std::vector< double > &, int, bool)
 
- Public Member Functions inherited from BCIntegrate
 BCIntegrate ()
 
 BCIntegrate (const BCIntegrate &bcintegrate)
 
virtual ~BCIntegrate ()
 
BCIntegrateoperator= (const BCIntegrate &bcintegrate)
 
int ReadMarginalizedFromFile (const char *file)
 
BCH1DGetMarginalized (const BCParameter *parameter)
 
BCH1DGetMarginalized (const char *name)
 
BCH1DGetMarginalized (unsigned index)
 
BCH2DGetMarginalized (const BCParameter *parameter1, const BCParameter *parameter2)
 
BCH2DGetMarginalized (const char *name1, const char *name2)
 
BCH2DGetMarginalized (unsigned index1, unsigned index2)
 
int PrintAllMarginalized1D (const char *filebase)
 
int PrintAllMarginalized2D (const char *filebase)
 
int PrintAllMarginalized (const char *file, std::string options1d="BTsiB3CS1D0pdf0Lmeanmode", std::string options2d="BTfB3CS1meangmode", unsigned int hdiv=1, unsigned int ndiv=1)
 
double GetIntegral () const
 
BCIntegrate::BCOptimizationMethod GetOptimizationMethod () const
 
BCIntegrate::BCIntegrationMethod GetIntegrationMethod () const
 
BCIntegrate::BCMarginalizationMethod GetMarginalizationMethod () const
 
BCIntegrate::BCSASchedule GetSASchedule () const
 
void GetRandomVectorUnitHypercube (std::vector< double > &x) const
 
void GetRandomVectorInParameterSpace (std::vector< double > &x) const
 
double GetRandomPoint (std::vector< double > &x)
 
int GetNIterationsMin () const
 
int GetNIterationsMax () const
 
int GetNIterationsPrecisionCheck () const
 
int GetNIterationsOutput () const
 
int GetNIterations () const
 
double GetRelativePrecision () const
 
double GetAbsolutePrecision () const
 
BCCubaMethod GetCubaIntegrationMethod () const
 
const BCCubaOptions::VegasGetCubaVegasOptions () const
 
const BCCubaOptions::SuaveGetCubaSuaveOptions () const
 
const BCCubaOptions::DivonneGetCubaDivonneOptions () const
 
const BCCubaOptions::CuhreGetCubaCuhreOptions () const
 
BCH1DGetSlice (const BCParameter *parameter, const std::vector< double > parameters=std::vector< double >(0), int bins=0, bool flag_norm=true)
 
BCH1DGetSlice (const char *name, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
BCH2DGetSlice (const BCParameter *parameter1, const BCParameter *parameter2, const std::vector< double > parameters=std::vector< double >(0), int bins=0, bool flag_norm=true)
 
BCH2DGetSlice (const char *name1, const char *name2, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
BCH2DGetSlice (unsigned index1, unsigned index2, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool flag_norm=true)
 
double GetError () const
 
TMinuit * GetMinuit ()
 
int GetMinuitErrorFlag () const
 
double GetSAT0 () const
 
double GetSATmin () const
 
double GetBestFitParameter (unsigned index) const
 
double GetBestFitParameterError (unsigned index) const
 
double GetLogMaximum ()
 
const std::vector< double > & GetBestFitParameters () const
 
const std::vector< double > & GetBestFitParameterErrors () const
 
void SetMinuitArlist (double *arglist)
 
void SetFlagIgnorePrevOptimization (bool flag)
 
void SetOptimizationMethod (BCIntegrate::BCOptimizationMethod method)
 
void SetIntegrationMethod (BCIntegrate::BCIntegrationMethod method)
 
void SetMarginalizationMethod (BCIntegrate::BCMarginalizationMethod method)
 
void SetSASchedule (BCIntegrate::BCSASchedule schedule)
 
void SetNIterationsMin (int niterations)
 
void SetNIterationsMax (int niterations)
 
void SetNIterationsPrecisionCheck (int niterations)
 
void SetNIterationsOutput (int niterations)
 
void SetRelativePrecision (double relprecision)
 
void SetAbsolutePrecision (double absprecision)
 
void SetCubaIntegrationMethod (BCCubaMethod type)
 
void SetCubaOptions (const BCCubaOptions::Vegas &options)
 
void SetCubaOptions (const BCCubaOptions::Suave &options)
 
void SetCubaOptions (const BCCubaOptions::Divonne &options)
 
void SetCubaOptions (const BCCubaOptions::Cuhre &options)
 
void SetSAT0 (double T0)
 
void SetSATmin (double Tmin)
 
void SetFlagWriteSAToFile (bool flag)
 
TTree * GetSATree ()
 
void InitializeSATree ()
 
double Normalize ()
 
double Integrate (BCIntegrationMethod intmethod)
 
double Integrate ()
 
double Integrate (BCIntegrationMethod type, tRandomizer randomizer, tEvaluator evaluator, tIntegralUpdater updater, std::vector< double > &sums)
 
double EvaluatorMC (std::vector< double > &sums, const std::vector< double > &point, bool &accepted)
 
int MarginalizeAll ()
 
int MarginalizeAll (BCMarginalizationMethod margmethod)
 
virtual void MarginalizePreprocess ()
 
virtual void MarginalizePostprocess ()
 
void SAInitialize ()
 
std::vector< double > FindMode (std::vector< double > start=std::vector< double >())
 
std::vector< double > FindMode (BCIntegrate::BCOptimizationMethod optmethod, std::vector< double > start=std::vector< double >())
 
double SATemperature (double t)
 
double SATemperatureBoltzmann (double t)
 
double SATemperatureCauchy (double t)
 
virtual double SATemperatureCustom (double t)
 
std::vector< double > GetProposalPointSA (const std::vector< double > &x, int t)
 
std::vector< double > GetProposalPointSABoltzmann (const std::vector< double > &x, int t)
 
std::vector< double > GetProposalPointSACauchy (const std::vector< double > &x, int t)
 
virtual std::vector< double > GetProposalPointSACustom (const std::vector< double > &x, int t)
 
std::vector< double > SAHelperGetRandomPointOnHypersphere ()
 
double SAHelperGetRadialCauchy ()
 
double SAHelperSinusToNIntegral (int dim, double theta)
 
virtual void ResetResults ()
 
std::string DumpIntegrationMethod (BCIntegrationMethod type)
 
std::string DumpCurrentIntegrationMethod ()
 
std::string DumpUsedIntegrationMethod ()
 
std::string DumpMarginalizationMethod (BCMarginalizationMethod type)
 
std::string DumpCurrentMarginalizationMethod ()
 
std::string DumpUsedMarginalizationMethod ()
 
std::string DumpOptimizationMethod (BCOptimizationMethod type)
 
std::string DumpCurrentOptimizationMethod ()
 
std::string DumpUsedOptimizationMethod ()
 
std::string DumpCubaIntegrationMethod (BCCubaMethod type)
 
std::string DumpCubaIntegrationMethod ()
 
void SetBestFitParameters (const std::vector< double > &x)
 
void SetBestFitParameters (const std::vector< double > &x, const double &new_value, double &old_value)
 
unsigned GetNIntegrationVariables ()
 
double CalculateIntegrationVolume ()
 
bool CheckMarginalizationAvailability (BCMarginalizationMethod type)
 
bool CheckMarginalizationIndices (TH1 *hist, const std::vector< unsigned > &index)
 
- Public Member Functions inherited from BCEngineMCMC
void WriteMarkovChain (bool flag)
 
 BCEngineMCMC ()
 
 BCEngineMCMC (const BCEngineMCMC &enginemcmc)
 
virtual ~BCEngineMCMC ()
 
BCEngineMCMCoperator= (const BCEngineMCMC &engineMCMC)
 
unsigned MCMCGetNChains () const
 
unsigned MCMCGetNLag () const
 
const std::vector< unsigned > & MCMCGetNIterations () const
 
unsigned MCMCGetCurrentIteration () const
 
unsigned MCMCGetCurrentChain () const
 
unsigned MCMCGetNIterationsConvergenceGlobal () const
 
bool MCMCGetFlagConvergenceGlobal () const
 
unsigned MCMCGetNIterationsMax () const
 
unsigned MCMCGetNIterationsRun () const
 
unsigned MCMCGetNIterationsPreRunMin () const
 
unsigned MCMCGetNIterationsUpdate () const
 
unsigned MCMCGetNIterationsUpdateMax () const
 
std::vector< int > MCMCGetNTrialsTrue () const
 
int MCMCGetNTrials () const
 
const std::vector< double > & MCMCGetprobMean () const
 
const std::vector< double > & MCMCGetVariance () const
 
const std::vector< double > & MCMCGetTrialFunctionScaleFactor () const
 
std::vector< double > MCMCGetTrialFunctionScaleFactor (unsigned ichain) const
 
double MCMCGetTrialFunctionScaleFactor (unsigned ichain, unsigned ipar)
 
const std::vector< double > & MCMCGetx () const
 
std::vector< double > MCMCGetx (unsigned ichain)
 
double MCMCGetx (unsigned ichain, unsigned ipar) const
 
const std::vector< double > & MCMCGetLogProbx () const
 
double MCMCGetLogProbx (unsigned ichain)
 
int MCMCGetPhase () const
 
const std::vector< double > & MCMCGetMaximumPoints () const
 
std::vector< double > MCMCGetMaximumPoint (unsigned i) const
 
const std::vector< double > & MCMCGetMaximumLogProb () const
 
int MCMCGetFlagInitialPosition () const
 
double MCMCGetRValueCriterion () const
 
double MCMCGetRValueParametersCriterion () const
 
double MCMCGetRValue () const
 
double MCMCGetRValueParameters (unsigned i)
 
bool MCMCGetRValueStrict () const
 
bool MCMCGetFlagRun () const
 
TTree * MCMCGetMarkovChainTree (unsigned i)
 
BCH1DMCMCGetH1Marginalized (unsigned i)
 
BCH2DMCMCGetH2Marginalized (unsigned i, unsigned j)
 
BCParameterGetParameter (int index) const
 
BCParameterGetParameter (const char *name) const
 
unsigned int GetNParameters () const
 
unsigned int GetNFixedParameters ()
 
unsigned int GetNFreeParameters ()
 
const std::vector< double > & GetBestFitParametersMarginalized () const
 
void MCMCSetTrialFunctionScaleFactor (std::vector< double > scale)
 
void MCMCSetNChains (unsigned n)
 
void MCMCSetNLag (unsigned n)
 
void MCMCSetNIterationsMax (unsigned n)
 
void MCMCSetNIterationsRun (unsigned n)
 
void MCMCSetNIterationsPreRunMin (unsigned n)
 
void MCMCSetNIterationsUpdate (unsigned n)
 
void MCMCSetNIterationsUpdateMax (unsigned n)
 
void MCMCSetMinimumEfficiency (double efficiency)
 
void MCMCSetMaximumEfficiency (double efficiency)
 
void MCMCSetRandomSeed (unsigned seed)
 
void MCMCSetWriteChainToFile (bool flag)
 
void MCMCSetWritePreRunToFile (bool flag)
 
void MCMCSetInitialPositions (const std::vector< double > &x0s)
 
void MCMCSetInitialPositions (std::vector< std::vector< double > > x0s)
 
void MCMCSetFlagInitialPosition (int flag)
 
void MCMCSetFlagOrderParameters (bool flag)
 
void MCMCSetFlagFillHistograms (bool flag)
 
void MCMCSetFlagPreRun (bool flag)
 
void MCMCSetRValueCriterion (double r)
 
void MCMCSetRValueParametersCriterion (double r)
 
void MCMCSetRValueStrict (bool strict=true)
 
void MCMCSetMarkovChainTrees (const std::vector< TTree * > &trees)
 
void MCMCInitializeMarkovChainTrees ()
 
int SetMarginalized (unsigned index, TH1D *h)
 
int SetMarginalized (unsigned index1, unsigned index2, TH2D *h)
 
void MCMCSetValuesDefault ()
 
void MCMCSetValuesQuick ()
 
void MCMCSetValuesDetail ()
 
void MCMCSetPrecision (BCEngineMCMC::Precision precision)
 
void SetNbins (unsigned int nbins)
 
void Copy (const BCEngineMCMC &enginemcmc)
 
virtual int AddParameter (const char *name, double min, double max, const char *latexname="")
 
virtual void MCMCTrialFunction (unsigned ichain, std::vector< double > &x)
 
virtual double MCMCTrialFunctionSingle (unsigned ichain, unsigned ipar)
 
bool MCMCGetProposalPointMetropolis (unsigned chain, std::vector< double > &x)
 
bool MCMCGetProposalPointMetropolis (unsigned chain, unsigned parameter, std::vector< double > &x)
 
bool MCMCGetNewPointMetropolis (unsigned chain=0)
 
bool MCMCGetNewPointMetropolis (unsigned chain, unsigned parameter)
 
void MCMCInChainCheckMaximum ()
 
void MCMCInChainUpdateStatistics ()
 
void MCMCInChainFillHistograms ()
 
void MCMCInChainTestConvergenceAllChains ()
 
void MCMCInChainWriteChains ()
 
int MCMCMetropolis ()
 
int MCMCMetropolisPreRun ()
 
void MCMCResetRunStatistics ()
 
void MCMCInitializeMarkovChains ()
 
int MCMCInitialize ()
 
void ClearParameters (bool hard=false)
 
virtual void MCMCIterationInterface ()
 
virtual void MCMCUserIterationInterface ()
 
virtual void MCMCCurrentPointInterface (std::vector< double > &, int, bool)
 

Protected Attributes

std::vector< BCMVCUncertainty * > fUncertainties
 
std::vector< BCMVCMeasurement * > fMeasurements
 
TMatrixD fCovarianceMatrix
 
TMatrixD fInvCovarianceMatrix
 
double fDetCovariance
 
TVectorD fVectorMeasurements
 
TVectorD fVectorActiveMeasurements
 
std::vector< int > fVectorObservable
 
std::vector< int > fVectorActiveObservable
 
TMatrixD fBLUEWeights
 
TVectorD fBLUECentral
 
TVectorD fBLUEUncertainties
 
std::vector< TVectorD > fBLUEUncertaintiesPerSource
 
TMatrixD fBLUECovarianceMatrix
 
std::vector< TMatrixD > fBLUECovarianceMatrices
 
std::vector< TMatrixD > fBLUECorrelationMatrices
 
TMatrixD fBLUECorrelationMatrix
 
int fNObservables
 
int fNNuisanceCorrelation
 
std::vector< NuisanceParameterfNuisanceCorrelation
 
std::vector< BCMVCObservable * > fObservables
 
- Protected Attributes inherited from BCModel
std::string fName
 
double fModelAPriori
 
double fModelAPosteriori
 
BCDataSetfDataSet
 
BCDataPointfDataPointLowerBoundaries
 
BCDataPointfDataPointUpperBoundaries
 
std::vector< bool > fDataFixedValues
 
double fPValue
 
double fChi2NDoF
 
double fPValueNDoF
 
bool flag_discrete
 
int fGoFNIterationsMax
 
int fGoFNIterationsRun
 
int fGoFNChains
 
std::vector< TNamed * > fPriorContainer
 
bool fPriorConstantAll
 
std::vector< bool > fPriorContainerConstant
 
std::vector< bool > fPriorContainerInterpolate
 
- Protected Attributes inherited from BCIntegrate
int fID
 
TMinuit * fMinuit
 
double fMinuitArglist [2]
 
int fMinuitErrorFlag
 
bool fFlagIgnorePrevOptimization
 
double fSAT0
 
double fSATmin
 
TTree * fTreeSA
 
bool fFlagWriteSAToFile
 
int fSANIterations
 
double fSATemperature
 
double fSALogProb
 
std::vector< double > fSAx
 
std::vector< BCH1D * > fMarginalized1D
 
std::vector< BCH2D * > fMarginalized2D
 
bool fFlagMarginalized
 
- Protected Attributes inherited from BCEngineMCMC
BCParameterSet fParameters
 
unsigned fMCMCNChains
 
unsigned fMCMCNLag
 
std::vector< unsigned > fMCMCNIterations
 
int fMCMCCurrentIteration
 
int fMCMCCurrentChain
 
unsigned fMCMCNIterationsUpdate
 
unsigned fMCMCNIterationsUpdateMax
 
int fMCMCNIterationsConvergenceGlobal
 
bool fMCMCFlagConvergenceGlobal
 
unsigned fMCMCNIterationsMax
 
unsigned fMCMCNIterationsRun
 
unsigned fMCMCNIterationsPreRunMin
 
std::vector< int > fMCMCNTrialsTrue
 
unsigned fMCMCNTrials
 
bool fMCMCFlagWriteChainToFile
 
bool fMCMCFlagWritePreRunToFile
 
std::vector< double > fMCMCTrialFunctionScaleFactor
 
std::vector< double > fMCMCTrialFunctionScaleFactorStart
 
bool fMCMCFlagPreRun
 
bool fMCMCFlagRun
 
std::vector< double > fMCMCInitialPosition
 
double fMCMCEfficiencyMin
 
double fMCMCEfficiencyMax
 
int fMCMCFlagInitialPosition
 
bool fMCMCFlagOrderParameters
 
int fMCMCPhase
 
std::vector< double > fMCMCx
 
std::vector< double > fMCMCxMax
 
std::vector< double > fMCMCxMean
 
std::vector< double > fMCMCxVar
 
std::vector< double > fMCMCprob
 
std::vector< double > fMCMCprobMax
 
std::vector< double > fMCMCprobMean
 
std::vector< double > fMCMCprobVar
 
bool fMCMCRValueUseStrict
 
double fMCMCRValueCriterion
 
double fMCMCRValueParametersCriterion
 
double fMCMCRValue
 
std::vector< double > fMCMCRValueParameters
 
TRandom3 * fRandom
 
std::vector< TH1D * > fMCMCH1Marginalized
 
std::vector< TH2D * > fMCMCH2Marginalized
 
std::vector< TTree * > fMCMCTrees
 
std::vector< double > fMCMCBestFitParameters
 
double fMCMCLogMaximum
 
std::vector< double > fMarginalModes
 

Additional Inherited Members

- Public Types inherited from BCIntegrate
enum  BCOptimizationMethod {
  kOptEmpty, kOptSimAnn, kOptMetropolis, kOptMinuit,
  kOptDefault, NOptMethods
}
 
enum  BCIntegrationMethod {
  kIntEmpty, kIntMonteCarlo, kIntCuba, kIntGrid,
  kIntDefault, NIntMethods
}
 
enum  BCMarginalizationMethod {
  kMargEmpty, kMargMetropolis, kMargMonteCarlo, kMargGrid,
  kMargDefault, NMargMethods
}
 
enum  BCSASchedule { kSACauchy, kSABoltzmann, kSACustom, NSAMethods }
 
enum  BCCubaMethod {
  kCubaVegas, kCubaSuave, kCubaDivonne, kCubaCuhre,
  NCubaMethods
}
 
typedef void(BCIntegrate::* tRandomizer )(std::vector< double > &) const
 
typedef double(BCIntegrate::* tEvaluator )(std::vector< double > &, const std::vector< double > &, bool &)
 
typedef void(* tIntegralUpdater )(const std::vector< double > &, const int &, double &, double &)
 
- Static Public Member Functions inherited from BCIntegrate
static void IntegralUpdaterMC (const std::vector< double > &sums, const int &nIterations, double &integral, double &absprecision)
 
static int CubaIntegrand (const int *ndim, const double xx[], const int *ncomp, double ff[], void *userdata)
 
static void FCNLikelihood (int &npar, double *grad, double &fval, double *par, int flag)
 
- Protected Member Functions inherited from BCIntegrate
unsigned IntegrationOutputFrequency () const
 
void LogOutputAtStartOfIntegration (BCIntegrationMethod type, BCCubaMethod cubatype)
 
void LogOutputAtIntegrationStatusUpdate (BCIntegrationMethod type, double integral, double absprecision, int nIterations)
 
void LogOutputAtEndOfIntegration (double integral, double absprecision, double relprecision, int nIterations)
 
void Copy (const BCIntegrate &bcintegrate)
 

Detailed Description

Definition at line 26 of file BCMVCombination.h.

Constructor & Destructor Documentation

BCMVCombination::BCMVCombination ( )

Definition at line 26 of file BCMVCombination.cxx.

27  : BCModel("BCMVCombination")
28  , fDetCovariance(0)
29  , fNObservables(0)
31 {
32 }
BCMVCombination::~BCMVCombination ( )

Definition at line 35 of file BCMVCombination.cxx.

36 {
37  int nuncertainties = GetNUncertainties();
38  int nmeasurements = GetNMeasurements();
39 
40  for (int i = 0; i < nuncertainties; ++i) {
42  delete u;
43  }
44  fUncertainties.clear();
45 
46  for (int i = 0; i < nmeasurements; ++i) {
48  delete m;
49  }
50  fMeasurements.clear();
51 }

Member Function Documentation

void BCMVCombination::AddMeasurement ( std::string  name,
std::string  observable,
double  value,
std::vector< double >  uncertainties 
)

Definition at line 83 of file BCMVCombination.cxx.

84 {
85  // get index of the corresponding observable
86  int index = GetIndexObservable(observable);
87 
88  // check if observable exists
89  if (index < 0) {
90  BCLog::OutWarning(Form("BCMVCombination::AddMeasurement. Observable \"%s\" does not exist. Measurement was not added.", observable.c_str()));
91  return;
92  }
93 
94  BCMVCMeasurement* m = new BCMVCMeasurement(name);
95  m->SetObservable(index);
96  m->SetCentralValue(value);
97  m->SetUncertainties(uncertainties);
98 
99  fMeasurements.push_back(m);
100 
101  fVectorObservable.push_back(index);
102 
103  int n = GetNMeasurements();
104  fVectorMeasurements.ResizeTo(n);
105  fVectorMeasurements[n-1]=value;
106 }
void BCMVCombination::AddObservable ( std::string  name,
double  min,
double  max 
)

Definition at line 54 of file BCMVCombination.cxx.

55 {
56  // check if observable exists already
57  int index = GetIndexObservable(name);
58 
59  if (index >= 0)
60  return;
61 
62  BCMVCObservable* obs = new BCMVCObservable();
63  obs->SetName(name);
64  obs->SetMinMax(min, max);
65  fObservables.push_back(obs);
66 
67  fNObservables++;
68 
69  AddParameter(name.c_str(), min, max);
70 
71  SetPriorConstant(name.c_str());
72 }
void BCMVCombination::AddUncertainty ( std::string  name)

Definition at line 75 of file BCMVCombination.cxx.

76 {
77  BCMVCUncertainty* u = new BCMVCUncertainty(name);
78 
79  fUncertainties.push_back(u);
80 }
void BCMVCombination::CalculateBLUE ( )

Definition at line 492 of file BCMVCombination.cxx.

493 {
494  int nmeas = GetNMeasurements();
495  int nactivemeas = GetNActiveMeasurements();
496  int nobs = GetNObservables();
497  int nunc = GetNUncertainties();
498 
499  // calculate U matrix
500  // TMatrixD u(nmeas, nobs);
501  TMatrixD u(nactivemeas, nobs);
502 
503  int counter = 0;
504  for (int i = 0; i < nmeas; ++i) {
506 
507  // if measurement is active fill matrix u
508  if (m->GetFlagActive()) {
509  for (int j = 0; j < nobs; ++j) {
510  if (m->GetObservable() == j)
511  u[counter][j] = 1;
512  else
513  u[counter][j] = 0;
514  }
515  counter++;
516  }
517  }
518 
519  // calculate weight matrix
520  TMatrixD ut = u;
521  ut.Transpose(ut);
522 
523  TMatrixD m1 = ut * fInvCovarianceMatrix;
524  TMatrixD m2 = m1 * u;
525  m2.Invert();
526 
527  fBLUEWeights.Clear();
528  // fBLUEWeights.ResizeTo(nobs, nmeas);
529  fBLUEWeights.ResizeTo(nobs, nactivemeas);
530 
531  fBLUEWeights = m2*m1;
532 
533  // calculate central values
534  fBLUECentral.Clear();
535  fBLUECentral.ResizeTo(nobs);
536 
537  // fBLUECentral = fBLUEWeights * fVectorMeasurements;
539 
540  // calculate covariance matrix
541  fBLUECovarianceMatrix.Clear();
542  fBLUECovarianceMatrix.ResizeTo(nobs, nobs);
543 
544  TMatrixD weightt = fBLUEWeights;
545  weightt.Transpose(weightt);
546 
548 
549  // calculate uncertainties, covariance and correlation matrices for each uncertainty
550  for (int i = 0; i < nunc; ++i) {
551 
552  // calculate covariance matrix
554  if (!u->GetFlagActive())
555  continue;
556 
557  TMatrixD cov = u->GetCovarianceMatrix();
558  TMatrixD mat = fBLUEWeights * cov * weightt;
559  fBLUECovarianceMatrices.push_back(mat);
560 
561  // calculate uncertainties
562  TVectorD vec(nobs);
563  for (int j = 0; j < nobs; ++j) {
564  double sigma_j = sqrt(mat[j][j]);
565  vec[j] = sigma_j;
566  }
567  fBLUEUncertaintiesPerSource.push_back(vec);
568 
569  TMatrixD mat2(nobs, nobs);
570  for (int j = 0; j < nobs; ++j)
571  for (int k = 0; k < nobs; ++k) {
572  mat2[j][k] = mat[j][k]/vec[j]/vec[k];
573  }
574  fBLUECorrelationMatrices.push_back(mat2);
575  }
576 
577  // calculate uncertainties
578  fBLUEUncertainties.Clear();
579  fBLUEUncertainties.ResizeTo(nobs);
580 
581  for (int i = 0; i < nobs; ++i)
582  fBLUEUncertainties[i] = sqrt(fBLUECovarianceMatrix[i][i]);
583 
584  // calculate correlation matrix
585  fBLUECorrelationMatrix.Clear();
586  fBLUECorrelationMatrix.ResizeTo(nobs, nobs);
587 
588  for (int i = 0; i < nobs; ++i)
589  for (int j = 0; j < nobs; ++j)
591 }
void BCMVCombination::CalculateCorrelationMatrix ( int  index)

Definition at line 354 of file BCMVCombination.cxx.

355 {
356  BCMVCUncertainty* u = fUncertainties.at(index);
357 
358  int n = GetNMeasurements();
359  int nactive = GetNActiveMeasurements();
360 
361  TMatrixD cov(nactive, nactive);
362 
363  TMatrixD corr = u->GetCorrelationMatrix();
364 
365  int counteri = 0;
366  for (int i = 0; i < n; ++i) {
368  double sigma_i = mi->GetUncertainty(index);
369 
370  // skip line if not active
371  if (!mi->GetFlagActive())
372  continue;
373 
374  int counterj = 0;
375  for (int j = 0; j < n; ++j) {
377  double sigma_j = mj->GetUncertainty(index);
378 
379  if (mj->GetFlagActive()) {
380  cov[counteri][counterj] = corr[i][j]*sigma_i*sigma_j;
381  counterj++;
382  }
383  }
384  counteri++;
385  }
386 
387  u->SetCovarianceMatrix(cov);
388 }
void BCMVCombination::CalculateCovarianceMatrix ( std::vector< double >  nuisance = std::vector<double>(0))

Definition at line 391 of file BCMVCombination.cxx.

392 {
393  int n = GetNUncertainties();
394  int nmeasurements = GetNMeasurements();
395  int nnuisance = fNNuisanceCorrelation;
396 
397  fCovarianceMatrix.Clear();
399 
400  // loop over all uncertainties
401  for (int unc_i = 0; unc_i < n; ++unc_i) {
402  BCMVCUncertainty* u = fUncertainties.at(unc_i);
403 
404  // shrink covariance matrix such that it fits only active measurements
405  TMatrixD mat_small = u->GetCovarianceMatrix();
406 
407  // loop over all measurements (i)
408  int counteri = 0;
409  for (int meas_i = 0; meas_i < nmeasurements; ++meas_i) {
410  BCMVCMeasurement* mi = GetMeasurement(meas_i);
411 
412  // skip line if not active
413  if (!mi->GetFlagActive())
414  continue;
415 
416  // loop over all measurements (j)
417  int counterj = 0;
418  for (int meas_j = 0; meas_j < nmeasurements; ++meas_j) {
419  BCMVCMeasurement* mj = GetMeasurement(meas_j);
420 
421  // skip line if not active
422  if (!mj->GetFlagActive())
423  continue;
424 
425  // modify matrix if nuisance parameter present
426  if (parameters.size() > 0) {
427 
428  // loop over all nuisance parameters
429  for (int nuis_k = 0; nuis_k < nnuisance; ++nuis_k) {
430 
431  // get nuisance parameter
432  NuisanceParameter p = fNuisanceCorrelation.at(nuis_k);
433 
434  // compare
435  if (p.index_uncertainty == unc_i) {
436  double sigma_i = GetMeasurement(p.index_measurement1)->GetUncertainty(unc_i);
437  double sigma_j = GetMeasurement(p.index_measurement2)->GetUncertainty(unc_i);
438  double pre = p.pre;
439 
440  // check if indices agree
441  if (meas_i == p.index_measurement1 && meas_j == p.index_measurement2) {
442  mat_small[counteri][counterj] = pre * parameters.at(p.index_rhoparameter) * sigma_i*sigma_j;
443  mat_small[counterj][counteri] = mat_small[counteri][counterj];
444  }
445  }
446  } // end loop: nuisance parameters
447  }
448 
449  // increase counter of active measurements
450  counterj++;
451  } // end loop: measurements j
452 
453  // increase counter of active measurements
454  counteri++;
455  } // end loop: measurements i
456 
457  // add matrix if active
458  if (u->GetFlagActive())
459  fCovarianceMatrix += mat_small;
460  }
461  fInvCovarianceMatrix.Clear();
464  fInvCovarianceMatrix.Invert();
465 
466  fDetCovariance = fCovarianceMatrix.Determinant();
467 }
void BCMVCombination::CalculateHelperVectors ( )

Definition at line 321 of file BCMVCombination.cxx.

322 {
323  int nmeasurements = GetNMeasurements();
324 
325  fVectorMeasurements.Clear();
326  fVectorMeasurements.ResizeTo(nmeasurements);
327  fVectorObservable.clear();
328 
329  for (int i = 0; i < nmeasurements; ++i) {
332  fVectorObservable.push_back(m->GetObservable());
333  }
334 
335  int nactive = GetNActiveMeasurements();
336 
338  fVectorActiveMeasurements.ResizeTo(nactive);
339  fVectorActiveObservable.clear();
340 
341  int counter = 0;
342  for (int i = 0; i < nmeasurements; ++i) {
344  if (m->GetFlagActive()) {
346  fVectorActiveObservable.push_back(m->GetObservable());
347  counter++;
348  }
349  }
350 
351 }
TVectorD BCMVCombination::GetBLUECentralValues ( )
inline

Definition at line 88 of file BCMVCombination.h.

89  { return fBLUECentral; };
TMatrixD BCMVCombination::GetBLUECorrelationMatrix ( int  index)
inline

Definition at line 114 of file BCMVCombination.h.

115  { return fBLUECorrelationMatrices.at(index); };
TMatrixD BCMVCombination::GetBLUECorrelationMatrix ( )
inline

Definition at line 118 of file BCMVCombination.h.

119  { return fBLUECorrelationMatrix; };
TMatrixD BCMVCombination::GetBLUECovarianceMatrix ( )
inline

Definition at line 102 of file BCMVCombination.h.

103  { return fBLUECovarianceMatrix; };
TMatrixD BCMVCombination::GetBLUECovarianceMatrix ( int  index)
inline

Definition at line 108 of file BCMVCombination.h.

109  { return fBLUECovarianceMatrices.at(index); };
TVectorD BCMVCombination::GetBLUEUncertainties ( )
inline

Definition at line 92 of file BCMVCombination.h.

93  { return fBLUEUncertainties; };
TVectorD BCMVCombination::GetBLUEUncertainties ( int  index)
inline

Definition at line 98 of file BCMVCombination.h.

99  { return fBLUEUncertaintiesPerSource.at(index); };
TMatrixD BCMVCombination::GetBLUEWeights ( )
inline

Definition at line 84 of file BCMVCombination.h.

85  { return fBLUEWeights; };
TMatrixD BCMVCombination::GetCovarianceMatrix ( )
inline

Definition at line 80 of file BCMVCombination.h.

81  { return fCovarianceMatrix; };
int BCMVCombination::GetIndexMeasurement ( std::string  measurement,
std::string  observable 
)

Definition at line 802 of file BCMVCombination.cxx.

803 {
804  int index_observable = GetIndexObservable(observable);
805 
806  int nmeasurements = GetNMeasurements();
807 
808  for (int i = 0; i < nmeasurements; ++i) {
809  if ((measurement == GetMeasurement(i)->GetName()) && (index_observable == GetMeasurement(i)->GetObservable()))
810  return i;
811  }
812 
813  return -1;
814 }
int BCMVCombination::GetIndexObservable ( std::string  name)

Definition at line 832 of file BCMVCombination.cxx.

833 {
834  int n = GetNObservables();
835 
836  // go through the list of parameters and compare strings
837  for (int i = 0; i < n; ++i) {
838  // if (name == std::string(GetParameter(i)->GetName()))
839  if (name == std::string(fObservables.at(i)->GetName()))
840  return i;
841  }
842 
843  // return -1 if not in the list
844  return -1;
845 }
int BCMVCombination::GetIndexUncertainty ( std::string  name)

Definition at line 817 of file BCMVCombination.cxx.

818 {
819  int nuncertainties = GetNUncertainties();
820 
821  int index = -1;
822 
823  for (int i = 0; i < nuncertainties; ++i) {
824  if (name == GetUncertainty(i)->GetName())
825  index = i;
826  }
827 
828  return index;
829 }
BCMVCMeasurement* BCMVCombination::GetMeasurement ( int  index)
inline

Definition at line 76 of file BCMVCombination.h.

77  { return fMeasurements.at(index); }
int BCMVCombination::GetNActiveMeasurements ( )

Definition at line 306 of file BCMVCombination.cxx.

307 {
308  int n = GetNMeasurements();
309 
310  int counter = 0;
311 
312  for (int i = 0; i < n; ++i) {
314  if (m->GetFlagActive())
315  counter++;
316  }
317  return counter;
318 }
int BCMVCombination::GetNMeasurements ( )
inline

Definition at line 65 of file BCMVCombination.h.

66  { return int(fMeasurements.size()); };
int BCMVCombination::GetNObservables ( )
inline

Definition at line 57 of file BCMVCombination.h.

58  { return fNObservables; };
int BCMVCombination::GetNUncertainties ( )
inline

Definition at line 61 of file BCMVCombination.h.

62  { return int(fUncertainties.size()); };
BCMVCUncertainty* BCMVCombination::GetUncertainty ( int  index)
inline

Definition at line 72 of file BCMVCombination.h.

73  { return fUncertainties.at(index); }
TVectorD BCMVCombination::GetVectorMeasurements ( )
inline

Definition at line 126 of file BCMVCombination.h.

127  { return fVectorMeasurements; };
std::vector<int> BCMVCombination::GetVectorObservable ( )
inline

Definition at line 122 of file BCMVCombination.h.

123  { return fVectorObservable; };
double BCMVCombination::LogLikelihood ( const std::vector< double > &  params)
virtual

Calculates natural logarithm of the likelihood. Method needs to be overloaded by the user.

Parameters
paramsA set of parameter values
Returns
Natural logarithm of the likelihood

Implements BCModel.

Reimplemented in BCMVCPhysicsModel.

Definition at line 109 of file BCMVCombination.cxx.

110 {
111  double logprob = 0.;
112 
113  if (fNNuisanceCorrelation > 0) {
114  CalculateCovarianceMatrix(parameters);
115  if (!PositiveDefinite())
116  return -1e90;
117  }
118 
119  int nmeas = GetNActiveMeasurements();
120 
121  // copy parameters into a vector
122  TVectorD observables(nmeas);
123 
124  for (int i = 0; i < nmeas; ++i) {
125  observables[i] = parameters[fVectorActiveObservable[i]];
126  }
127 
128  TVectorD prod1 = observables - fVectorActiveMeasurements;
129  TVectorD prod2 = fInvCovarianceMatrix * prod1;
130  double prod = prod1 * prod2;
131 
132  logprob = -0.5 * prod - log(TMath::Power(2*TMath::Pi(), nmeas/2.) * sqrt(fabs(fDetCovariance)));
133 
134  return logprob;
135 }
bool BCMVCombination::PositiveDefinite ( )

Definition at line 470 of file BCMVCombination.cxx.

471 {
472  TMatrixDEigen m(fCovarianceMatrix);
473 
474  // calculate real part of all eigenvalues
475  TVectorD eigen_re = m.GetEigenValuesRe();
476 
477  int n_eigen = eigen_re.GetNoElements();
478 
479  bool flag_ispositive = true;
480 
481  // check if eigenvalues are positive
482  for (int i = 0; i < n_eigen; ++i) {
483  if (eigen_re[i] < 0)
484  flag_ispositive = false;
485  }
486 
487  // true if all eigenvalues are positive
488  return flag_ispositive;
489 }
void BCMVCombination::PrepareAnalysis ( )

Definition at line 287 of file BCMVCombination.cxx.

288 {
289  int nuncertainties = GetNUncertainties();
290 
291  // prepare analysis
292  for (int i = 0; i < nuncertainties; ++i)
294 
296 
297  if (!PositiveDefinite()) {
298  BCLog::OutWarning("BCMVCombination::PrepareAnalysis. Covariance matrix is not positive definite.");
299  }
300 
302 
303 }
void BCMVCombination::PrintBLUEResults ( std::string  filename)

Definition at line 594 of file BCMVCombination.cxx.

595 {
596  // open file
597  std::ofstream ofi(filename.c_str());
598 
599  // check if file is open
600  if (!ofi.is_open()) {
601  std::cerr << "Couldn't open file " << filename << std::endl;
602  return;
603  }
604 
605  int nobs = GetNObservables();
606  int nmeas = GetNMeasurements();
607  int nunc = GetNUncertainties();
608 
609  ofi << std::endl;
610  ofi << "Summary of the combination" << std::endl;
611  ofi << "==========================" << std::endl << std::endl;
612 
613  ofi << "* Observables:" << std::endl;
614  ofi << " Observable (range): " << std::endl;
615  for (int i = 0; i < nobs; ++i)
616  ofi << " " << std::setiosflags(std::ios::left) << GetParameter(i)->GetName()
617  << " (" << GetParameter(i)->GetLowerLimit() << " - " << GetParameter(i)->GetUpperLimit() << ")" << std::endl;
618  ofi << std::endl;
619 
620  ofi << "* Measurements:" << std::endl;
621  ofi << " Measurement (observable): central value +- total uncertainty" << std::endl;
622  for (int i = 0; i < nmeas; ++i) {
624  if (m->GetFlagActive()) {
625  double total2 = 0;
626  for (int j = 0; j < nunc; ++j) {
627  if (GetUncertainty(j)->GetFlagActive())
628  total2+=m->GetUncertainty(j)*m->GetUncertainty(j);
629  }
630  ofi << " " << std::setiosflags(std::ios::left) << std::setw(20) << m->GetName()
631  << std::setiosflags(std::ios::left) << " (" << GetParameter(m->GetObservable())->GetName() << ")"
632  << ": " << std::setiosflags(std::ios::left) << std::setw(7) << std::setprecision(4) << m->GetCentralValue()
633  << " +- " << std::setiosflags(std::ios::left) << std::setw(7) << std::setprecision(4) << sqrt(total2) << std::endl;
634  }
635  }
636  ofi << std::endl;
637 
638  ofi << "* Uncertainties:" << std::endl;
639  ofi << " Measurement (observable): Uncertainty (";
640  for (int j = 0; j < nunc-1; ++j )
641  if (GetUncertainty(j)->GetFlagActive())
642  ofi << GetUncertainty(j)->GetName() << ", ";
643  if (GetUncertainty(nunc-1)->GetFlagActive())
644  ofi << GetUncertainty(nunc-1)->GetName() << ")" << std::endl;
645  else
646  ofi << ")" << std::endl;
647 
648  for (int i = 0; i < nmeas; ++i) {
650  if (m->GetFlagActive()) {
651  ofi << " " << std::setiosflags(std::ios::left) << std::setw(20) << m->GetName()
652  << std::setiosflags(std::ios::left) << " (" << GetParameter(m->GetObservable())->GetName() << "): ";
653  for (int j = 0; j < nunc; ++j )
654  if (GetUncertainty(j)->GetFlagActive())
655  ofi << std::setiosflags(std::ios::left) << std::setw(7) << m->GetUncertainty(j);
656  ofi << std::endl;
657  }
658  }
659  ofi << std::endl;
660 
661  for (int i = 0; i < nunc; ++i) {
663 
664  if (!u->GetFlagActive())
665  continue;
666 
667  ofi << " " << u->GetName() << " " << "(correlation matrix)" << std::endl;
668  TMatrixD mat = u->GetCorrelationMatrix();
669 
670  int counterk = 0;
671  for (int k = 0; k < nmeas; ++k) {
673 
674  // skip line if not active
675  if (!mk->GetFlagActive())
676  continue;
677 
678  ofi << " ";
679 
680  int counterj = 0;
681  for (int j = 0; j < nmeas; ++j) {
683 
684  if (mj->GetFlagActive()) {
685  ofi << std::setw(7) << std::showpos << mat[k][j] << " ";
686  counterj++;
687  }
688  }
689  ofi << std::noshowpos << std::endl;
690  counterk++;
691  }
692  ofi << std::endl;
693  }
694 
695  ofi << "* BLUE results: " << std::endl;
696  ofi << " Observable: estimate +- total uncertainty" << std::endl;
697  for (int i = 0; i < nobs; ++i) {
698  if (i < fBLUECentral.GetNoElements())
699  ofi << " " << GetParameter(i)->GetName() << ": " << fBLUECentral[i] << " +- " << std::setprecision(4) << fBLUEUncertainties[i] << std::endl;
700  }
701  ofi << std::endl;
702 
703  ofi << " Observable: Uncertainty (";
704  for (int j = 0; j < nunc-1; ++j )
705  if (GetUncertainty(j)->GetFlagActive())
706  ofi << GetUncertainty(j)->GetName() << ", ";
707  if (GetUncertainty(nunc-1)->GetFlagActive())
708  ofi << GetUncertainty(nunc-1)->GetName() << ")" << std::endl;
709  else
710  ofi << ")" << std::endl;
711 
712  for (int i = 0; i < nobs; ++i) {
713  ofi << " " << std::setiosflags(std::ios::left) << GetParameter(i)->GetName()<< ":";
714  int counterj = 0;
715  for (int j = 0; j < nunc; ++j )
716  if (GetUncertainty(j)->GetFlagActive()) {
717  ofi << " " << std::setiosflags(std::ios::left) << std::setw(7) << std::setprecision(4) << GetBLUEUncertainties(counterj)[i];
718  counterj++;
719  }
720  ofi << std::endl;
721  }
722  ofi << std::endl;
723 
724  if (nobs > 1) {
725  ofi << " Individual correlation matrices " << std::endl;
726  for (int i = 0; i < nunc; ++i) {
727  if (GetUncertainty(i)->GetFlagActive()) {
728  TMatrixD mat = GetBLUECorrelationMatrix(i);
729  ofi << " " << GetUncertainty(i)->GetName() << std::endl;
730  for (int j = 0; j < nobs; ++j) {
731  ofi << " ";
732  for (int k = 0; k < nobs; ++k) {
733  ofi << std::setw(7) << std::setprecision(4) << std::showpos << mat[j][k] << " ";
734  }
735  ofi << std::noshowpos << std::endl;
736  }
737  ofi << std::endl;
738  }
739  }
740  }
741 
742  if (nobs > 1) {
743  ofi << " Overall correlation matrix" << std::endl;
744  TMatrixD mat = fBLUECorrelationMatrix;
745  for (int j = 0; j < nobs; ++j) {
746  ofi << " ";
747  for (int k = 0; k < nobs; ++k)
748  ofi << std::setw(7) << std::setprecision(4) << std::showpos << mat[j][k] << " ";
749  ofi << std::noshowpos << std::endl;
750  }
751  ofi << std::endl;
752  }
753 
754  ofi << " Weights [%]:" <<std::endl;
755  int counter = 0;
756  for (int j = 0; j < nmeas; ++j) {
757  if (GetMeasurement(j)->GetFlagActive()) {
758  ofi << " " << std::setw(20) << GetMeasurement(j)->GetName() << " : ";
759  for (int k = 0; k < nobs; ++k)
760  ofi << std::setw(7) << std::setprecision(4) << std::showpos << fBLUEWeights[k][counter]*100. << " ";
761  ofi << std::endl;
762  counter++;
763  }
764  }
765  ofi << std::endl;
766 
767  // close file
768  ofi.close();
769 }
void BCMVCombination::PrintMatrix ( TMatrixT< double > &  matrix,
std::string  name = "matrix" 
)

Definition at line 772 of file BCMVCombination.cxx.

773 {
774  int nrows = matrix.GetNrows();
775  int ncols = matrix.GetNcols();
776 
777  std::cout << std::endl;
778  std::cout << name.c_str() << " (" << nrows << "x" << ncols << "):" << std::endl;
779 
780  for (int i = 0; i < nrows; ++i) {
781  for (int j = 0; j < ncols; ++j)
782  std::cout << std::setprecision(3) << std::setw(7) << matrix[i][j] << " ";
783  std::cout << std::endl;
784  }
785 }
void BCMVCombination::PrintVector ( TVectorD &  vector,
std::string  name = "vector" 
)

Definition at line 788 of file BCMVCombination.cxx.

789 {
790  int nrows = vector.GetNoElements();
791 
792  std::cout << std::endl;
793  std::cout << name.c_str() << " (" << nrows << "):" << std::endl;
794 
795  for (int i = 0; i < nrows; ++i) {
796  std::cout << std::setprecision(3) << std::setw(7) << vector[i] << " ";
797  std::cout << std::endl;
798  }
799 }
int BCMVCombination::ReadInput ( std::string  filename)

Definition at line 138 of file BCMVCombination.cxx.

139 {
140  // open input file
141  ifstream infile;
142  infile.open(filename.c_str(), std::ifstream::in);
143 
144  // check if file is open
145  if (!infile.is_open()) {
146  BCLog::OutWarning(Form("BCMVCombination::ReadInput. Could not open input file %s.", filename.c_str()));
147  return 0;
148  }
149 
150  int nobservables;
151  int nmeasurements;
152  int nuncertainties;
153  int nnuisance;
154 
155  infile >> nobservables >> nmeasurements >> nuncertainties >> nnuisance;
156 
157  std::vector<std::string> observable_names;
158 
159  for (int i = 0; i < nobservables; ++i) {
160  std::string name;
161  double min;
162  double max;
163  infile >> name >> min >> max;
164 
165  // add observable
166  AddObservable(name.c_str(), min, max);
167  }
168 
169  for (int i = 0; i < nuncertainties; ++i) {
170  std::string name;
171  infile >> name;
172 
173  // add uncertainty
174  AddUncertainty(name);
175  }
176 
177  for (int i = 0; i < nmeasurements; ++i) {
178  std::string name;
179  std::string observable;
180  double central;
181  std::vector<double> uncertainties(0);
182 
183  infile >> name;
184  infile >> observable;
185  infile >> central;
186 
187  for (int j = 0; j < nuncertainties; ++j) {
188  double uncertainty;
189  infile >> uncertainty;
190  uncertainties.push_back(uncertainty);
191  }
192 
193  // add measurement
194  AddMeasurement(name, observable, central, uncertainties);
195  }
196 
197  for (int i = 0; i < nuncertainties; ++i) {
198  TMatrixD mat(nmeasurements, nmeasurements);
199 
200  for (int j = 0; j < nmeasurements; ++j)
201  for (int k = 0; k < nmeasurements; ++k) {
202  double corr;
203  infile >> corr;
204  mat[j][k] = corr;
205  }
206 
207  // set correlation matrix
209  }
210 
211  for (int i = 0; i < nnuisance; ++i) {
212  std::string uncertainty;
213  std::string measurement1;
214  std::string observable1;
215  std::string measurement2;
216  std::string observable2;
217  std::string parname;
218  double min;
219  double max;
220  double pre;
221  std::string priorshape;
222 
223  infile >> uncertainty >> measurement1 >> observable1 >> measurement2 >> observable2 >> parname;
224 
225  // check if parameter exists already
226  int index = -1;
227 
228  for (unsigned int i = 0; i < GetNParameters(); i++)
229  if (parname.c_str() == GetParameter(i)->GetName())
230  index = i;
231 
232  if (index >= 0)
233  infile >> pre;
234 
235  else if (index < 0) {
236  // read properties of parameter
237  infile >> min >> max >> priorshape;
238 
239  // add nuisance parameter
240  AddParameter(parname.c_str(), min, max);
241 
242  // set pre-factor
243  pre = 1;
244 
245  // set index
246  index = GetNParameters()-1;
247 
248  if (priorshape == "flat") {
249  SetPriorConstant(parname.c_str());
250  }
251  else if (priorshape == "gauss") {
252  double mean;
253  double std;
254 
255  infile >> mean >> std;
256 
257  SetPriorGauss(parname.c_str(), mean, std);
258  }
259  else {
260  BCLog::OutWarning(Form("BCMVCombination::ReadInput. Unknown prior shape %s.", priorshape.c_str()));
261  }
262  }
263 
264  // increase counter of nuisance parametera
266 
267  NuisanceParameter p;
268  p.index_uncertainty = GetIndexUncertainty(uncertainty);
269  p.index_measurement1 = GetIndexMeasurement(measurement1, observable1);
270  p.index_measurement2 = GetIndexMeasurement(measurement2, observable2);
271  p.index_rhoparameter = index;
272  p.pre = pre;
273 
274  fNuisanceCorrelation.push_back(p);
275  }
276 
277  // close input file
278  infile.close();
279 
280  PrepareAnalysis();
281 
282  // no error
283  return 1;
284 }

Member Data Documentation

TVectorD BCMVCombination::fBLUECentral
protected

Definition at line 215 of file BCMVCombination.h.

std::vector<TMatrixD> BCMVCombination::fBLUECorrelationMatrices
protected

Definition at line 230 of file BCMVCombination.h.

TMatrixD BCMVCombination::fBLUECorrelationMatrix
protected

Definition at line 233 of file BCMVCombination.h.

std::vector<TMatrixD> BCMVCombination::fBLUECovarianceMatrices
protected

Definition at line 227 of file BCMVCombination.h.

TMatrixD BCMVCombination::fBLUECovarianceMatrix
protected

Definition at line 224 of file BCMVCombination.h.

TVectorD BCMVCombination::fBLUEUncertainties
protected

Definition at line 218 of file BCMVCombination.h.

std::vector<TVectorD> BCMVCombination::fBLUEUncertaintiesPerSource
protected

Definition at line 221 of file BCMVCombination.h.

TMatrixD BCMVCombination::fBLUEWeights
protected

Definition at line 212 of file BCMVCombination.h.

TMatrixD BCMVCombination::fCovarianceMatrix
protected

Definition at line 191 of file BCMVCombination.h.

double BCMVCombination::fDetCovariance
protected

Definition at line 197 of file BCMVCombination.h.

TMatrixD BCMVCombination::fInvCovarianceMatrix
protected

Definition at line 194 of file BCMVCombination.h.

std::vector<BCMVCMeasurement*> BCMVCombination::fMeasurements
protected

Definition at line 188 of file BCMVCombination.h.

int BCMVCombination::fNNuisanceCorrelation
protected

Definition at line 239 of file BCMVCombination.h.

int BCMVCombination::fNObservables
protected

Definition at line 236 of file BCMVCombination.h.

std::vector<NuisanceParameter> BCMVCombination::fNuisanceCorrelation
protected

Definition at line 242 of file BCMVCombination.h.

std::vector<BCMVCObservable*> BCMVCombination::fObservables
protected

Definition at line 245 of file BCMVCombination.h.

std::vector<BCMVCUncertainty*> BCMVCombination::fUncertainties
protected

Definition at line 185 of file BCMVCombination.h.

TVectorD BCMVCombination::fVectorActiveMeasurements
protected

Definition at line 203 of file BCMVCombination.h.

std::vector<int> BCMVCombination::fVectorActiveObservable
protected

Definition at line 209 of file BCMVCombination.h.

TVectorD BCMVCombination::fVectorMeasurements
protected

Definition at line 200 of file BCMVCombination.h.

std::vector<int> BCMVCombination::fVectorObservable
protected

Definition at line 206 of file BCMVCombination.h.


The documentation for this class was generated from the following files: