ATLAS Offline Software
Loading...
Searching...
No Matches
TripleGaussCollFit.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5/*** Collisions monitoring algorithm by Tudor Costin & Peter Onyisi
6
7Revision history:
8
9v1.1a (TC):
10
11-multLB_writeTags() will now query the bins to get lower and upper good LBs
12instead of assuming a 1-1 correspondence between bins and LBs
13-added a couple of #ifdef statements:
14 TC_WORKBENCH - when set, alg. uses a different heade in local mode
15 This is used for testing purposes ONLY and should never
16 be used in the live environment.
17
18 PARAM_IS_AREA - when set, the three-Gaussian formula is defined so that the
19 parameters multiplying the Gaussian exponents represent the
20 integrated area under the Gaussian; otherwise, they
21 represent the height of the peaks.
22
23v1.1 (TC):
24
25-improved code to pick starting height of Gaussians. Will now query the 1D
26histogram for given ranges of ns (was: bin ranges)
27-starting zeros are added to the output LB_(#S)_(#F) to improve formatting
28-Mean_zero_gaus parameter will now translate the entire function, not just middle peak
29-Changed names of some configurable parameters to reflect this.
30-All parameters except the heights/areas of gaussians are now fixed by default
31-Reasonable(?) limits imposed in case the user overrides the above fixing.
32 (cannot be overriden at present)
33
34v1.0a,b(TC)
35
36-various bug fixes.
37
38v1.0 (TC):
39-first release with multiple LB functionality
40-This should actually be labelled v2.0, but too late now
41
42
43*/
44
45#ifndef TC_WORKBENCH
47#include <dqm_core/AlgorithmManager.h>
48#include <dqm_core/exceptions.h>
50#include "ers/ers.h"
51
52namespace {
53 dqm_algorithms::TripleGaussCollFit three_gauss_fit( "default" );
54}
55
56#else
57#include "TripleGaussCollFit_standalone.h"
58#endif
59
60#include "TH2.h"
61#include "TH1.h"
62#include "TF1.h"
63#include "TClass.h" //shocking that this doesn't fall in by now
64#include <iostream>
65#include <cmath>
66
67
68
69//for some reason, FindBin() is not a const method, so I have to discard the const part. Sigh.
70Double_t getStartValue(TH1* thehist, Double_t minx, Double_t maxx) {
71
72 Double_t result = 0.0;
73 int minidx = thehist->FindBin(minx), maxidx = thehist->FindBin(maxx);
74 for (int idx = minidx; idx <= maxidx; idx++) {
75 Double_t val = thehist->GetBinContent(idx);
76 if (val > result)
77 result = val;
78 }
79
80 return result;
81}
82
83
84const char* const fmt = "%0*d";
85
86void multLB_writeTags(std::vector<int>* goodLB, dqm_core::Result* result, const TH1* thehist) {
87
88 /* This method assumes that the entries in goodLB are in ascending order!! */
89
90 char buf[10];
91 std::string lb_str = "LB_";
92 int firstLB = (*goodLB)[0], offset = 1;
93 unsigned int c = 1, sz = goodLB->size(), psz = static_cast<unsigned int>( std::log10(sz))+1;
94 Double_t lo, hi;
95
96 while (c < sz) {
97 if ((*goodLB)[c] == firstLB + offset)
98 offset++;
99
100 else {
101 std::string tagname = lb_str;
102 //find the low LB
103
104 lo = thehist->GetXaxis()->GetBinLowEdge(firstLB);
105
106 if (lo == std::floor(lo))
107 sprintf(buf, fmt, psz, static_cast<Int_t>(lo));
108 else
109 sprintf(buf, fmt, psz, static_cast<Int_t>(std::floor(lo))+1);
110
111 tagname += buf;
112 tagname += "_";
113
114 hi = thehist->GetXaxis()->GetBinUpEdge(firstLB+offset-1);
115
116 if (hi == std::floor(hi))
117 sprintf(buf, fmt, psz, static_cast<Int_t>(std::floor(hi))-1);
118 else
119 sprintf(buf, fmt, psz, static_cast<Int_t>(std::floor(hi)));
120
121 tagname += buf;
122 result->tags_ [tagname] = 1.0;
123 firstLB = (*goodLB)[c];
124 offset = 1;
125 }
126
127 c++; //I've always wanted to write this
128 }
129
130 lo = thehist->GetXaxis()->GetBinLowEdge(firstLB);
131
132 if (lo == std::floor(lo))
133 sprintf(buf, fmt, psz, static_cast<Int_t>(lo));
134 else
135 sprintf(buf, fmt, psz, static_cast<Int_t>(std::floor(lo))+1);
136
137 lb_str += buf;
138 lb_str += "_";
139 hi = thehist->GetXaxis()->GetBinUpEdge(firstLB+offset-1);
140 if (hi == std::floor(hi))
141 sprintf(buf, fmt, psz, static_cast<Int_t>(std::floor(hi))-1);
142 else
143 sprintf(buf, fmt, psz, static_cast<Int_t>(std::floor(hi)));
144 lb_str += buf;
145 result->tags_ [lb_str] = 1.0;
146
147}
148
149/****
150
151The following parameters are held fixed by default (whether user specified or using built-in values):
152 -Mean_zero_gaus (default is 0.0)
153 -Mean_offzero_gaus (default is 23.2, obtained by empirical fits. Should be overriden for certain detectors!)
154 -Width_all (default is 4.3, also obtained empirically)
155
156If they are allowed to float, the following limits are imposed. This is currently not configurable.
157 -Mean_zero_gaus: specified or default value, plus or minus 5
158 -Mean_offzero_gaus: specified or default value, plus or minus 10
159 -Width_all: between 1.0 and 10.0
160**/
161
162typedef std::map<std::string, double> param_map;
163typedef std::map<std::string, double> err_map;
164
165#define PARAM_IS_AREA //if this is not in, the constants multiplying the Gaussians are defined as max_ht, not integr. area
166
167namespace dqm_algorithms {
168
169 TripleGaussCollFit::TripleGaussCollFit(const std::string& name) : m_name(name) {
170
171 //Register our instance with the Alg. Manager
172 m_gaus3_fn = 0;
173 dqm_core::AlgorithmManager::instance().registerAlgorithm( "Gaus3_" + m_name +"_Fit", this );
174
175 }
176
178
179 if (m_gaus3_fn != 0) {
180 delete m_gaus3_fn;
181 m_gaus3_fn = 0;
182 }
183
184 }
185
189
190// This method probably need not be implemented, but I should have it anyway...
191 void TripleGaussCollFit::printDescription(std::ostream& out) {
192
193 out << "Gaus3_"+ m_name +"_Fit: specialized fit for the Tile timining difference plot." << std::endl;
194 out << "Current parameter values: " << std::endl;
195
196 if (m_gaus3_fn != 0)
197 for (int i=0; i < m_gaus3_fn->GetNpar(); i++)
198 out << " Param. " << i << ": " << m_gaus3_fn->GetParName(i) << " , current value: "
199 << m_gaus3_fn->GetParameter(i) << std::endl;
200
201 else
202 out << " execute() has not been called" << std::endl;
203
204 out << "Parameter list: " << std::endl
205 << " allowFitMean: 1 allows fitting for a shift in the entire function around zero by up to 5ns" << std::endl
206 << " allowFitMeanNonzeroPeaks: 1 allows fitting of the distance between center and off-center Gaussians" << std::endl
207 << " allowFitResolution: 1 allows fitting the (shared) width of the Gaussians" << std::endl
208 << " reportSignificance: 1 reports the central peak significance for each LB as a tag" << std::endl
209 << " meanZeroPeak: starting or fixed mean of the Gaussian centered at zero" << std::endl
210 << " meanNonZeroPeaks: starting or fixed mean of the other Gaussians" << std::endl
211 << " resolutionAll: starting or fixed value of the width of the Gaussians" << std::endl
212 << " minSignificance: number of events required for checks to be valid, undefined if not set." << std::endl
213 << " centralPeakMinSignificance: significance of central peak required for a LB to pass threshold; default is 3.0" << std::endl
214 << " minEvents: number of events required to attempt a fit; default is 5." << std::endl ;
215
216 out << "Currently, no thresholds are defined." << std::endl;
217
218 }
219
220 dqm_core::Result* TripleGaussCollFit::execute (const std::string& /* s */, const TObject& object, const dqm_core::AlgorithmConfig& cfg) {
221
222 const TH1* histogram;
223 if(object.IsA()->InheritsFrom("TH1")) //or this can be done via dynamic casting
224 {
225 histogram = static_cast<const TH1*>(&object);
226 if (histogram->GetDimension() > 2 ){
227 throw dqm_core::BadConfig( ERS_HERE, m_name, "called with histogram of dimension > 2" );
228 }
229 }
230 else
231 throw dqm_core::BadConfig( ERS_HERE, m_name, "called with object that does not inherit from TH1" );
232
233 //Setup our TF1: find the binwidth
234 Double_t binwidth = 0;
235 if (histogram->GetDimension() == 1) //running in 1D mode
236 if (histogram->GetNbinsX() >= 1) {
237 binwidth = histogram->GetBinWidth(1);
238 } else {
239 throw dqm_core::BadConfig( ERS_HERE, m_name, "called with histogram with no bins");
240 }
241
242 else { //running in 2D mode {
243 const TH1* proj = static_cast<const TH2*>(histogram)->ProjectionY("", 0, 0);
244 if (proj->GetNbinsX() >= 1) {
245 binwidth = proj ->GetBinWidth(1);
246 } else {
247 throw dqm_core::BadConfig( ERS_HERE, m_name, "called with histogram with no bins");
248 }
249 }
250 #ifdef PARAM_IS_AREA
251 m_gaus3_fn = new TF1("gaus3", Form("[3]/[2]/2.5066*%f*exp(-0.5*((x-[0]+[1])/[2])**2) + [4]/[2]/2.5066*%f*exp(-0.5*((x-[0])/[2])**2) + [5]/[2]/2.5066*%f*exp(-0.5*((x-[0]-[1])/[2])**2)", binwidth, binwidth, binwidth), -100, 100);
252 #else
253 m_gaus3_fn = new TF1("gaus3", "[3]*exp(-0.5*((x-[0]+[1])/[2])**2) + [4]*exp(-0.5*((x-[0])/[2])**2) + [5]*exp(-0.5*((x-[0]-[1])/[2])**2)", -100, 100);
254 #endif
255 m_gaus3_fn->SetParName(0, "Mean_zero_gaus"); //mean of gaussian centered at zero, should be close to zero
256 m_gaus3_fn->SetParName(1, "Mean_offzero_gaus"); //mean of other two gaussians, should be close to 20
257 m_gaus3_fn->SetParName(2, "Width_all"); //shared width, should be about 5
258 m_gaus3_fn->SetParName(3, "N_negative_gaus"); //normalization for negative-centered gaussian
259 m_gaus3_fn->SetParName(4, "N_zero_gaus"); //normalization for zero-centered gaussian
260 m_gaus3_fn->SetParName(5, "N_positive_gaus"); //normalization for positive-centered gaussian
261
262 Double_t resolutionAll, meanZeroPeak, meanNonzeroPeaks, minEvents;
263 param_map mymap = cfg.getParameters();
264
265 meanZeroPeak = dqm_algorithms::tools::GetFirstFromMap("meanZeroPeak", mymap, 0.0);
266 meanNonzeroPeaks = dqm_algorithms::tools::GetFirstFromMap("meanNonzeroPeaks", mymap, 23.2);
267 resolutionAll = dqm_algorithms::tools::GetFirstFromMap("resolutionAll", mymap, 4.57);
268 minEvents = dqm_algorithms::tools::GetFirstFromMap("minEvents", mymap, 5.0);
269
270 //we DO NOT allow fits by default. If the parameter is not found, the retval is 0
271 if ( dqm_algorithms::tools::GetFirstFromMap("allowFitMean", mymap, 0) == 1 ) {
272
273 //according to Peter, this can behave strangely, so we don't allow too much horizontal translation
274 m_gaus3_fn->SetParameter(0, meanZeroPeak);
275 m_gaus3_fn->SetParLimits(0, meanZeroPeak-5, meanZeroPeak+5);
276 }
277 else
278 m_gaus3_fn->FixParameter(0, meanZeroPeak);
279
280 if ( dqm_algorithms::tools::GetFirstFromMap("allowFitMeanNonzeroPeaks", mymap, 0) == 1 ) {
281 m_gaus3_fn->SetParameter(1, meanNonzeroPeaks);
282 m_gaus3_fn->SetParLimits(1, meanNonzeroPeaks-10, meanNonzeroPeaks+10);
283 }
284 else
285 m_gaus3_fn->FixParameter(1, meanNonzeroPeaks);
286
287 if ( dqm_algorithms::tools::GetFirstFromMap("allowFitResolution", mymap, 0) == 1 ) {
288 m_gaus3_fn->SetParameter(2, resolutionAll);
289 m_gaus3_fn->SetParLimits(2, 1, 10);
290 }
291 else {
292 m_gaus3_fn->FixParameter(2, resolutionAll);
293
294 }
295 // get some starting values for the histogram
296
297 dqm_core::Result* result = new dqm_core::Result();
298
299 if (histogram->GetDimension() == 1) {
300 fitSingle(const_cast<TH1*>(histogram));
301
302 try {
303
306 for (err_map::const_iterator peit = paramErrors.begin(); peit != paramErrors.end(); ++peit) {
307 result->tags_[(*peit).first + " Error"] = (*peit).second;
308 }
309
310 }
311 catch ( dqm_core::Exception & ex) {
312 throw dqm_core::BadConfig( ERS_HERE, m_name, ex.what(), ex);
313 }
314 }
315 else { //working with a 2D histogram
316
317 int nxbin = histogram->GetNbinsX();
318 Double_t minsig = dqm_algorithms::tools::GetFirstFromMap("centralPeakMinSignificance", mymap, 3.0);
319
320 bool writeSig = (dqm_algorithms::tools::GetFirstFromMap("reportSignificance", mymap, 0) == 1);
321
322
323 unsigned int psz = (unsigned int) std::log10(nxbin)+1;
324
325 std::vector<int> goodLB;
326 char buf[10];
327 const std::string sig_str = "SIG_LB_";
328
329 for (int lb = 0; lb < nxbin; lb++) {
330
331 sprintf(buf, fmt, psz, lb); //to avoid possible leaks (?) - should see how DQ framework behaves
332 TH1D* proj = static_cast<const TH2*>(histogram)->ProjectionY(buf, lb, lb);
333 Double_t good_entries = proj->GetEntries() - proj->GetBinContent(0) - proj->GetBinContent(proj->GetNbinsX()+1);
334 if (good_entries < minEvents) {
335 if (writeSig)
336 result->tags_[sig_str+buf] = 0;
337 continue; //nothing to fit here, move along
338 }
339
340 // I probably need to reset parameters 0-2 if they are not fixed
341 if ( dqm_algorithms::tools::GetFirstFromMap("allowFitMeanZeroPeak", mymap, 0) == 1 )
342 m_gaus3_fn->SetParameter(0, meanZeroPeak);
343
344 if ( dqm_algorithms::tools::GetFirstFromMap("allowFitMeanNonzeroPeaks", mymap, 0) == 1 )
345 m_gaus3_fn->SetParameter(1, meanNonzeroPeaks);
346
347 if ( dqm_algorithms::tools::GetFirstFromMap("allowFitResolution", mymap, 0) == 1 )
348 m_gaus3_fn->SetParameter(2, resolutionAll);
349
350 fitSingle(proj);
351
352#ifdef TC_WORKBENCH
353 result->histos_.push_back((TH1D*)proj->Clone());
354#endif
355
356 if (m_gaus3_fn->GetParameter(4) == 0) {
357 if (writeSig)
358 result->tags_[sig_str+buf] = 0;
359 continue;
360 }
361 Double_t sig = std::abs(m_gaus3_fn->GetParameter(4) / m_gaus3_fn->GetParError(4));
362 if (writeSig)
363 result->tags_[sig_str+buf] = sig;
364
365 if (sig > minsig)
366 goodLB.push_back(lb); //this is grossly non-optimal
367
368
369 } //end loop over lumi blocks
370
371 if (goodLB.size() != 0)
373
374 }
375 return result;
376 }
377
378
379#ifdef PARAM_IS_AREA
380#define ADJUST_VAL startVal = startVal * hist1D->GetBinWidth(1) / m_gaus3_fn->GetParameter(2) / 2.5066
381#else
382#define ADJUST_VAL
383#endif
384
386
387 Double_t startVal;
388 const Double_t window = 5.0;
389
390 Double_t offset_all = m_gaus3_fn->GetParameter(0);
391 Double_t offset_nz = m_gaus3_fn->GetParameter(1);
392
393 startVal = getStartValue(hist1D, offset_all-offset_nz-window, offset_all-offset_nz+window); ADJUST_VAL;
394 if (startVal == 0.0)
395 m_gaus3_fn->FixParameter(3, 0.0);
396 else {
397 m_gaus3_fn->SetParameter(3, startVal);
398 // DO NOT do something like this if we are using areas!!
399 // m_gaus3_fn->SetParLimits(3, 0, 5 * startVal);
400 }
401
402 startVal = getStartValue(hist1D, offset_all-window, offset_all+window); ADJUST_VAL;
403 if (startVal == 0.0)
404 m_gaus3_fn->FixParameter(4, 0.0);
405 else {
406 m_gaus3_fn->SetParameter(4, startVal);
407 // m_gaus3_fn->SetParLimits(4, 0, 100 * startVal);
408 }
409 startVal = getStartValue(hist1D, offset_all+offset_nz-window , offset_all+offset_nz+window); ADJUST_VAL;
410 if (startVal == 0.0)
411 m_gaus3_fn->FixParameter(5, 0.0);
412 else {
413 m_gaus3_fn->SetParameter(5, startVal);
414 }
415
416
417
418 hist1D->Fit(m_gaus3_fn, "QN");
419
420 }
421
422} //close namespace
static Double_t sz
Double_t getStartValue(TH1 *thehist, Double_t minx, Double_t maxx)
#define ADJUST_VAL
std::map< std::string, double > err_map
void multLB_writeTags(std::vector< int > *goodLB, dqm_core::Result *result, const TH1 *thehist)
std::map< std::string, double > param_map
const char *const fmt
std::string histogram
Definition chains.cxx:52
dqm_core::Result * execute(const std::string &s, const TObject &object, const dqm_core::AlgorithmConfig &cfg)
TripleGaussCollFit(const std::string &_name)
int lb
Definition globals.cxx:23
bool binwidth
Definition listroot.cxx:58
double GetFirstFromMap(const std::string &paramName, const std::map< std::string, double > &params)
std::map< std::string, double > GetFitParamErrors(const TF1 *func)
std::map< std::string, double > GetFitParams(const TF1 *func)
#define IsA
Declare the TObject style functions.