ATLAS Offline Software
Loading...
Searching...
No Matches
AlgorithmHelper.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
8
9#ifndef DQM_ALGORITHMS_TOOLS_ALGORITHMHELPER_H
10#define DQM_ALGORITHMS_TOOLS_ALGORITHMHELPER_H
11
12#include <dqm_core/Result.h>
13#include <dqm_core/AlgorithmConfig.h>
14#include <dqm_core/exceptions.h>
15
16#include <map>
17#include <vector>
18#include <deque>
19#include <list>
20#include <utility>
21
22#include <TH1.h>
23#include <TObject.h>
24#include <TClass.h>
25#include <boost/thread/once.hpp>
26#include <boost/thread/mutex.hpp>
27#include <boost/thread/locks.hpp>
28#include <boost/thread.hpp>
29#include <boost/utility.hpp>
30#include <iostream>
31
32class TF1;
33class TAxis;
34
35namespace dqm_algorithms
36{
37 namespace tools
38 {
40
42
43 //Structure to hold the essence of a 2D histogram bin, and be sortable:
44 struct binContainer {
45 double value;
46 double error;
47 int test;
48 int ix;
49 int iy;
50 double x;
51 double y;
52 //Comparison function, for sorting
53 inline static bool comp(const binContainer &lhs, const binContainer &rhs) { return fabs(lhs.value) < fabs(rhs.value); }
54 };
55
56 struct binCluster {
57 double value;
58 double error;
59 double x;
60 double y;
61 double radius;
62 int n;
63 int ixmin;
64 int ixmax;
65 int iymin;
66 int iymax;
67 inline static bool comp(const binCluster &lhs, const binCluster &rhs) { return fabs(lhs.value) < fabs(rhs.value); }
68 };
69
70 std::map<std::string, double > GetFitParams(const TF1 * func);
71
72 std::map<std::string, double > GetFitParamErrors(const TF1 * func);
73
74 dqm_core::Result * MakeComparisons( const std::map<std::string,double> & algparams,
75 const std::map<std::string,double> & gthreshold,
76 const std::map<std::string,double> & rthreshold );
77
78 dqm_core::Result * CompareWithErrors( const std::map<std::string,double> & algparams,
79 const std::map<std::string,double> & paramErrors,
80 const std::map<std::string,double> & gthreshold,
81 const std::map<std::string,double> & rthreshold, double minSig);
82
83 dqm_core::Result * GetFitResult (const TF1 * func, const dqm_core::AlgorithmConfig & config, double minSig = 0) ;
84
85 double GetFirstFromMap(const std::string &paramName, const std::map<std::string, double > &params);
86 // mandatory: throws an exception if the parameter is not found
87
88 double GetFirstFromMap(const std::string &paramName, const std::map<std::string, double > &params, double defaultValue);
89 // optional: returns defaultValue if the parameter is not found
90
91 //string overloads
92 const std::string& GetFirstFromMap(const std::string &paramName, const std::map<std::string, std::string > &params);
93 // mandatory: throws an exception if the parameter is not found
94
95 const std::string& GetFirstFromMap(const std::string &paramName, const std::map<std::string, std::string > &params, const std::string& defaultValue);
96 // optional: returns defaultValue if the parameter is not found
97
98 std::vector<int> GetBinRange(const TH1* histogram, const std::map<std::string, double > & params);
99
100 void PublishBin(const TH1 * histogram, int xbin, int ybin, double content, dqm_core::Result *result);
101
102 TH1* DivideByHistogram(const TH1* hNumerator, const TH1* hDenominator);
103
104 void ModifyHistogram(TH1 * histogram, const dqm_core::AlgorithmConfig & config);
105
106 std::string ExtractAlgorithmName(const dqm_core::AlgorithmConfig& config);
107
108 dqm_core::Result * ExecuteNamedAlgorithm(const std::string & name, const TObject & object,
109 const dqm_core::AlgorithmConfig & config);
110
111 TH1* BookHistogramByExample(const TH1* histogram, const std::string& name, const std::string& title, AxisType axisType);
112
113 template <class T>
114 const T & GetFromMap( const std::string & pname, const std::map<std::string,T> & params )
115 {
116 typename std::map<std::string,T>::const_iterator it = params.find( pname );
117 if ( it != params.end() ){
118 return it->second;
119 }else {
120 throw dqm_core::BadConfig( ERS_HERE, "None", pname );
121 }
122 }
123
137 void handleReference( const TObject& inputReference , const TObject*& firstReference , TObject*& secondReference);
138
139 // Function to find outliers in input; iterates over values nIteration times, recalculating mean each time and
140 // removing values that are beyond threshold * scale, where:
141 //
142 // scale = ( sum_in[ abs( value - mean )^ exponent ] / (Nin - 1 - SBCF * Nout) ) ^ ( 1 / exponent ).
143 //
144 // If all bins are in, and the exponent is two, this is just and unbiased estimator of the standard variance.
145 // SBCF, or the Scale Bias Correction Factor, is an empirical quantity intended to correct for the bias induced from the
146 // bin exclusion process (in a sense, the act of excluding bins could be thought of as decreasing the number of degrees
147 // of freedom). In practice, the SBCF serves to impose an upper bound on the fraction of bins that can be excluded.
148
149
150 void findOutliers( std::vector<binContainer>& input, double& mean, double& scale, int& nIn, int nIterations,
151 double exponent, double threshold, double SBCF = 1., double nStop = 8. );
152
153 void findOutliersUsingErrors( std::vector<binContainer>& input, double& mean, double& meanError, int& nIn,
154 double mindiff = 0, int minNin = 4);
155
156 // Method for building a cluster:
157 binCluster buildCluster( binContainer& seed, const std::vector<std::vector<binContainer*> >& binMap,
158 const std::vector<double>& xValues, const std::vector<double>& yValues,
159 double threhold, int topology = CylinderX);
160
161 // Method for mapping binContainer object by their relative positions:
162 std::vector<std::vector<binContainer*> >
163 makeBinMap(std::vector<dqm_algorithms::tools::binContainer>& bins, int ixmax, int iymax, int topology = CylinderX);
164
165 dqm_core::Result::Status
166 WorstCaseAddStatus(dqm_core::Result::Status baseStatus, dqm_core::Result::Status addedStatus, float weight = 1.0);
167
168 dqm_core::Result::Status
169 BestCaseAddStatus(dqm_core::Result::Status baseStatus, dqm_core::Result::Status addedStatus, float weight = 1.0);
170
171 std::pair<double,double> CalcBinsProbChisq(const std::vector<double>& inputval,const std::vector<double>& inputerr,
172 double x0, double x0_err);
173 std::pair<double,double> CalcBinsProbChisq(const std::vector<double>& inputval,const std::vector<double>& inputerr,
174 const std::vector<double>& x0,const std::vector<double>& x0_err);
175
176 void MergePastMinStat(std::vector<std::vector<tools::binContainer> >& strips, int minStat);
177
178 void MakeBinTag( const binContainer& bin, std::string & tag );
179
180 void FormatToSize( double value, int size, std::string & str, bool showSign = true );
181
182 }
183}
184
185
186#endif // #ifndef DQM_ALGORITHMS_TOOLS_ALGORITHMHELPER_H
static const std::vector< std::string > bins
std::string histogram
Definition chains.cxx:52
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
void FormatToSize(double value, int size, std::string &str, bool showSign=true)
std::vector< int > GetBinRange(const TH1 *histogram, const std::map< std::string, double > &params)
dqm_core::Result::Status WorstCaseAddStatus(dqm_core::Result::Status baseStatus, dqm_core::Result::Status addedStatus, float weight=1.0)
void findOutliersUsingErrors(std::vector< binContainer > &input, double &mean, double &meanError, int &nIn, double mindiff=0, int minNin=4)
std::pair< double, double > CalcBinsProbChisq(const std::vector< double > &inputval, const std::vector< double > &inputerr, double x0, double x0_err)
double GetFirstFromMap(const std::string &paramName, const std::map< std::string, double > &params)
TH1 * DivideByHistogram(const TH1 *hNumerator, const TH1 *hDenominator)
std::vector< std::vector< binContainer * > > makeBinMap(std::vector< dqm_algorithms::tools::binContainer > &bins, int ixmax, int iymax, int topology=CylinderX)
dqm_core::Result * CompareWithErrors(const std::map< std::string, double > &algparams, const std::map< std::string, double > &paramErrors, const std::map< std::string, double > &gthreshold, const std::map< std::string, double > &rthreshold, double minSig)
TH1 * BookHistogramByExample(const TH1 *histogram, const std::string &name, const std::string &title, AxisType axisType)
const T & GetFromMap(const std::string &pname, const std::map< std::string, T > &params)
void MakeBinTag(const binContainer &bin, std::string &tag)
void MergePastMinStat(std::vector< std::vector< tools::binContainer > > &strips, int minStat)
dqm_core::Result * GetFitResult(const TF1 *func, const dqm_core::AlgorithmConfig &config, double minSig=0)
std::string ExtractAlgorithmName(const dqm_core::AlgorithmConfig &config)
void findOutliers(std::vector< binContainer > &input, double &mean, double &scale, int &nIn, int nIterations, double exponent, double threshold, double SBCF=1., double nStop=8.)
std::map< std::string, double > GetFitParamErrors(const TF1 *func)
dqm_core::Result * MakeComparisons(const std::map< std::string, double > &algparams, const std::map< std::string, double > &gthreshold, const std::map< std::string, double > &rthreshold)
dqm_core::Result::Status BestCaseAddStatus(dqm_core::Result::Status baseStatus, dqm_core::Result::Status addedStatus, float weight=1.0)
dqm_core::Result * ExecuteNamedAlgorithm(const std::string &name, const TObject &object, const dqm_core::AlgorithmConfig &config)
void ModifyHistogram(TH1 *histogram, const dqm_core::AlgorithmConfig &config)
std::map< std::string, double > GetFitParams(const TF1 *func)
binCluster buildCluster(binContainer &seed, const std::vector< std::vector< binContainer * > > &binMap, const std::vector< double > &xValues, const std::vector< double > &yValues, double threhold, int topology=CylinderX)
void handleReference(const TObject &inputReference, const TObject *&firstReference, TObject *&secondReference)
Helper function used to handle complex reference histograms This function gets as input a reference o...
void PublishBin(const TH1 *histogram, int xbin, int ybin, double content, dqm_core::Result *result)
static bool comp(const binCluster &lhs, const binCluster &rhs)
static bool comp(const binContainer &lhs, const binContainer &rhs)