ATLAS Offline Software
Loading...
Searching...
No Matches
CaloHadDMCoeffFit.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//-----------------------------------------------------------------------
6// File and Version Information:
7// $Id: CaloHadDMCoeffFit.cxx,v 1.1 2009/03/03 17:30:23 pospelov Exp $
8//
9// Description: see CaloHadDMCoeffFit.h
10//
11// Environment:
12// Software developed for the ATLAS Detector at CERN LHC
13//
14// Author List:
15// Gennady Pospelov
16//
17//-----------------------------------------------------------------------
19#include "AthenaKernel/Units.h"
26#include <algorithm>
27#include <cmath>
28#include <fstream>
29#include <iostream>
30#include <sstream>
31
32#include <CLHEP/Vector/LorentzVector.h>
33
34#include "TROOT.h"
35#include "TStyle.h"
36#include "TError.h"
37#include "TFile.h"
38#include "TChain.h"
39#include "TH1.h"
40#include "TH1F.h"
41#include "TF1.h"
42#include "TProfile.h"
43#include "TH2F.h"
44#include "TCanvas.h"
45#include "TPad.h"
46#include "TText.h"
47#include "TLatex.h"
48#include "TGraphErrors.h"
49#include "TProfile2D.h"
50
51
52using CLHEP::HepLorentzVector;
53using Athena::Units::MeV;
54using Athena::Units::GeV;
55
56
57
68
69
75
76
77CaloLocalHadCoeff * CaloHadDMCoeffFit::process(CaloHadDMCoeffData *myData, CaloLocalHadCoeff *myHadDMCoeff, bool isSingleParticle, bool tbflag)
78{
79 m_isTestbeam = tbflag;
80
81 std::cout << std::endl;
82 std::cout << std::endl;
83 std::cout << "--- CaloHadDMCoeffFit::process() --- " << std::endl;
84
85 if(m_isTestbeam) std::cout << "Processing TESTBEAM data" << std::endl;
86
87 if ( m_NormalizationType == "Lin" ) {
88 std::cout << "Using weighting proportional to E_calib" << std::endl;
90 } else if ( m_NormalizationType == "Log" ) {
91 std::cout << "Using weighting proportional to log(E_calib)" << std::endl;
93 } else if ( m_NormalizationType == "NClus" ) {
94 std::cout << "Using weighting proportional to 1/N_Clus_E_calib>0" << std::endl;
96 } else {
97 std::cout << "Using constant weighting" << std::endl;
99 }
100
101 m_data = myData;
102 m_HadDMCoeff = myHadDMCoeff;
103
104 // loading only necessary branches
105 m_data->fChain->SetBranchStatus("*",0);
106 m_data->fChain->SetBranchStatus("mc_ener",1);
107 m_data->fChain->SetBranchStatus("mc_eta",1);
108 m_data->fChain->SetBranchStatus("mc_phi",1);
109 m_data->fChain->SetBranchStatus("mc_pdg",1);
110 m_data->fChain->SetBranchStatus("ncls",1);
111 m_data->fChain->SetBranchStatus("cls_eta",1);
112 m_data->fChain->SetBranchStatus("cls_phi",1);
113 m_data->fChain->SetBranchStatus("cls_lambda",1);
114 m_data->fChain->SetBranchStatus("cls_calib_emfrac",1);
115 m_data->fChain->SetBranchStatus("cls_engcalib",1);
116 m_data->fChain->SetBranchStatus("engClusSumCalib",1);
117 m_data->fChain->SetBranchStatus("cls_ener_unw",1);
118 m_data->fChain->SetBranchStatus("narea",1);
119 m_data->fChain->SetBranchStatus("cls_eprep",1);
120 m_data->fChain->SetBranchStatus("cls_dmener",1);
121 //m_data->fChain->SetBranchStatus("cls_smpener_unw",1);
122 //m_data->fChain->SetBranchStatus("cls_oocener",1);
123 //m_data->fChain->SetBranchStatus("cls_engcalibpres",1);
124
125 if( !m_data->GetEntries() ) {
126 std::cout << "CaloHadDMCoeffFit::process -> Error! No entries in DeadMaterialTree." << std::endl;
127 return nullptr;
128 }
129
130 m_data->GetEntry(0);
131 if(m_data->m_narea != m_HadDMCoeff->getSizeAreaSet()) {
132 std::cout << "CaloHadDMCoeffFit::process -> Error! Different numbers of areas for DM definition" << std::endl;
133 std::cout << "m_data->m_narea:" << m_data->m_narea << " m_HadDMCoeff->getSizeAreaSet():" << m_HadDMCoeff->getSizeAreaSet() << std::endl;
134 return nullptr;
135 }
136
137 clear();
138 for(int i_size=0; i_size<m_HadDMCoeff->getSizeCoeffSet(); i_size++){
139 m_engClus.push_back(new PrepData());
140 m_engPrep.push_back(new PrepData());
141 m_engDm.push_back(new PrepData());
142 m_engDmOverClus.push_back(new PrepData());
143 }
144
145
146 // --------------------------------------------------------------------------
147 // first run to calculate average and rms for histogram limits
148 // --------------------------------------------------------------------------
149 std::cout << "CaloHadDMCoeffFit::process() -> Info. Getting averages..." << std::endl;
150 for(int i_ev=0; m_data->GetEntry(i_ev)>0;i_ev++) {
151 if(i_ev%20000==0) std::cout << " i_ev: " << i_ev << " '" << static_cast<TChain *>(m_data->fChain)->GetFile()->GetName() << "'" << std::endl;
152
153 // checking event quality
154 if(isSingleParticle) {
155 bool GoodClusterFound(false);
156 if( m_data->m_ncls ) {
157 HepLorentzVector hlv_pion(1,0,0,1);
158 hlv_pion.setREtaPhi(1./cosh(m_data->m_mc_eta), m_data->m_mc_eta, m_data->m_mc_phi);
159 for(int i_cls=0; i_cls<m_data->m_ncls; i_cls++){ // loop over clusters
160 HepLorentzVector hlv_cls(1,0,0,1);
161 hlv_cls.setREtaPhi(1./cosh( (*m_data->m_cls_eta)[i_cls] ), (*m_data->m_cls_eta)[i_cls], (*m_data->m_cls_phi)[i_cls]);
162 double r = hlv_pion.angle(hlv_cls.vect());
164 && (*m_data->m_cls_engcalib)[i_cls] > 20.0*MeV
165 //&& (*m_data->m_cls_ener)[i_cls] > 0.01*m_data->m_mc_ener
166 ) {
167 GoodClusterFound = true;
168 break;
169 }
170 } // i_cls
171 }
172 if(!GoodClusterFound) continue;
173 }
174
175 if(m_data->m_engClusSumCalib <=0.0) continue;
176 for(int i_cls=0; i_cls<m_data->m_ncls; i_cls++){ // loop over clusters
177 if( (*m_data->m_cls_ener_unw)[i_cls] < m_energyMin) continue;
178
179 double clus_weight = 1.0;
181 clus_weight = (*m_data->m_cls_engcalib)[i_cls]/m_data->m_engClusSumCalib;
182 }
183 if(clus_weight <= 0) continue;
184 for(int i_dma=0; i_dma<m_data->m_narea; i_dma++){ // loop over DM areas
185 std::vector<float> vars;
186 m_data->PackClusterVars(i_cls, vars);
187 int iBin = m_HadDMCoeff->getBin(i_dma, vars );
188
189 if(iBin < 0) continue;
190 if(iBin >= m_HadDMCoeff->getSizeCoeffSet()) {
191 std::cout << "Panic'! iBin " << iBin << " " << m_HadDMCoeff->getSizeCoeffSet() << std::endl;
192 continue;
193 }
194 m_engClus[iBin]->add( (*m_data->m_cls_ener_unw)[i_cls] );
195 m_engPrep[iBin]->add( (*m_data->m_cls_eprep)[i_cls][i_dma] );
196 m_engDm[iBin]->add( (*m_data->m_cls_dmener)[i_cls][i_dma] );
197 if((*m_data->m_cls_ener_unw)[i_cls] <=0.0) {
198 std::cout << "Panic! Zero cluster energy" << std::endl;
199 continue;
200 }
201 double w = (*m_data->m_cls_dmener)[i_cls][i_dma]/(*m_data->m_cls_ener_unw)[i_cls];
202 if(w>=0 && w<m_weightMax) {
203 m_engDmOverClus[iBin]->add(w, clus_weight);
204 }
205 } // i_dma
206 } // i_cls
207 } // i_ev
208 std::cout << "done" << std::endl;
209 std::vector<double > engDmOverClusLim;
210 engDmOverClusLim.resize(m_HadDMCoeff->getSizeCoeffSet(), 0.0);
211 for(int i_size=0; i_size<m_HadDMCoeff->getSizeCoeffSet(); i_size++){
212 engDmOverClusLim[i_size] = m_engDmOverClus[i_size]->m_aver + 3.0*sqrt(m_engDmOverClus[i_size]->m_rms);
213 }
214
215 // --------------------------------------------------------------------------
216 // creation of histograms
217 // --------------------------------------------------------------------------
218 std::cout << "CaloHadDMCoeffFit::process() -> Info. Creation of histograms..." << std::endl;
219 m_h2_DmVsPrep.resize(m_HadDMCoeff->getSizeCoeffSet(), nullptr);
220 m_hp_DmVsPrep.resize(m_HadDMCoeff->getSizeCoeffSet(), nullptr);
221 m_h1_engDmOverClus.resize(m_HadDMCoeff->getSizeCoeffSet(), nullptr);
222 m_hp2_DmWeight.resize(m_HadDMCoeff->getSizeCoeffSet(), nullptr);
223
224 char hname[1024];
225 int nch_dms=40;
226 // creation of profiles
227 for(int i_size=0; i_size<m_HadDMCoeff->getSizeCoeffSet(); i_size++){
228 const CaloLocalHadCoeff::LocalHadArea *area = m_HadDMCoeff->getAreaFromBin(i_size);
229 if( area->getType() == CaloLocalHadDefs::AREA_DMFIT) {
230 // eprep .vs. edmtrue
231 sprintf(hname,"hp_DmVsPrep_%d",i_size);
232// float x_lim_edm = m_engDm[i_size]->m_aver + 3.0*sqrt(m_engDm[i_size]->m_rms);
233// float x_lim_eprep = m_engPrep[i_size]->m_aver + 3.0*sqrt(m_engPrep[i_size]->m_rms);
234 float x_lim_edm = 2.0*m_engDm[i_size]->m_aver + 3.0*sqrt(m_engDm[i_size]->m_rms);
235 float x_lim_eprep = 2.0*m_engPrep[i_size]->m_aver + 3.0*sqrt(m_engPrep[i_size]->m_rms);
236 if(x_lim_edm <= 0.0) x_lim_edm = 1.0;
237 if(x_lim_eprep <= 0.0) x_lim_eprep = 1.0;
238 m_hp_DmVsPrep[i_size] = new TProfile(hname,hname,nch_dms,0., x_lim_edm);
239 m_hp_DmVsPrep[i_size]->Sumw2();
240 sprintf(hname,"h2_DmVsPrep_%d",i_size);
241 m_h2_DmVsPrep[i_size] = new TH2F(hname,hname,nch_dms/2,0.0, x_lim_edm, nch_dms/2,0., x_lim_eprep);
242 }
243 if( area->getType() == CaloLocalHadDefs::AREA_DMLOOKUP) {
244 sprintf(hname,"h1_engDmOverClus_%d",i_size);
245// float xlim = 2.0*m_engDmOverClus[i_size]->m_aver + 5.0*sqrt(m_engDmOverClus[i_size]->m_rms);
246// if(xlim < 1.5) xlim = 1.5;
247 //m_h1_engDmOverClus[i_size] = new TH1F(hname,hname, nch_dms*2, -0.5, xlim );
248 m_h1_engDmOverClus[i_size] = new TH1F(hname,hname, 150, -0.5, m_weightMax );
249
250 // creation of histograms for inverted dm weight .vs. log10(ener), log10(lambda)
251 int refbin = getFirstEnerLambdaBin(i_size);
252 TProfile2D *hp2=nullptr;
253 if(refbin == i_size) {
254 std::cout << " creating histos:" << hname << " refbin:" << refbin << std::endl;
257 sprintf(hname,"m_hp2_DmWeight_%d",i_size);
258 hp2 = new TProfile2D(hname, hname, dimEner->getNbins(), dimEner->getXmin(), dimEner->getXmax(), dimLambda->getNbins(), dimLambda->getXmin(), dimLambda->getXmax(), 0.0, m_weightMax);
259 }else if (refbin >= 0) {
260 hp2 = m_hp2_DmWeight[refbin];
261 }
262 m_hp2_DmWeight[i_size] = hp2;
263 }
264 }
265
266 // --------------------------------------------------------------------------
267 // Filling of histograms
268 // --------------------------------------------------------------------------
269 std::cout << "CaloHadDMCoeffFit::process() -> Info. Filling histograms..." << std::endl;
270 for(int i_ev=0; m_data->GetEntry(i_ev)>0;i_ev++) {
271 if(i_ev%20000==0) std::cout << " i_ev: " << i_ev << " '" << static_cast<TChain *>(m_data->fChain)->GetFile()->GetName() << "'" << std::endl;
272
273 if(isSingleParticle) {
274 // checking event quality
275 bool GoodClusterFound(false);
276 if( m_data->m_ncls ) {
277 HepLorentzVector hlv_pion(1,0,0,1);
278 hlv_pion.setREtaPhi(1./cosh(m_data->m_mc_eta), m_data->m_mc_eta, m_data->m_mc_phi);
279 for(int i_cls=0; i_cls<m_data->m_ncls; i_cls++){ // loop over clusters
280 HepLorentzVector hlv_cls(1,0,0,1);
281 hlv_cls.setREtaPhi(1./cosh( (*m_data->m_cls_eta)[i_cls] ), (*m_data->m_cls_eta)[i_cls], (*m_data->m_cls_phi)[i_cls]);
282 double r = hlv_pion.angle(hlv_cls.vect());
284 && (*m_data->m_cls_engcalib)[i_cls] > 20.0*MeV
285 //&& (*m_data->m_cls_ener)[i_cls] > 0.01*m_data->m_mc_ener
286 ) {
287 GoodClusterFound = true;
288 break;
289 }
290 } // i_cls
291 }
292 if(!GoodClusterFound) continue;
293 }
294
295 //if(m_data->m_engClusSumCalib <=0.0) continue;
296 for(int i_cls=0; i_cls<m_data->m_ncls; i_cls++){ // loop over clusters
297 if( (*m_data->m_cls_ener_unw)[i_cls] < m_energyMin) continue;
298
299 double clus_weight = 1.0;
301 clus_weight = (*m_data->m_cls_engcalib)[i_cls]/m_data->m_engClusSumCalib;
302 }
303 if(clus_weight <= 0.0) continue;
304 for(int i_dma=0; i_dma<m_data->m_narea; i_dma++){ // loop over DM areas
305 std::vector<float> vars;
306 m_data->PackClusterVars(i_cls, vars);
307 int iBin = m_HadDMCoeff->getBin(i_dma, vars );
308
309 if(iBin < 0) continue;
310 if( (*m_data->m_cls_eprep)[i_cls][i_dma] > 0.0 && m_hp_DmVsPrep[iBin]){
311 m_hp_DmVsPrep[iBin]->Fill( (*m_data->m_cls_dmener)[i_cls][i_dma], (*m_data->m_cls_eprep)[i_cls][i_dma], clus_weight );
312 m_h2_DmVsPrep[iBin]->Fill( (*m_data->m_cls_dmener)[i_cls][i_dma], (*m_data->m_cls_eprep)[i_cls][i_dma], clus_weight );
313 }
314 double w = (*m_data->m_cls_dmener)[i_cls][i_dma]/(*m_data->m_cls_ener_unw)[i_cls];
315 if(m_h1_engDmOverClus[iBin] && w>=0 && w<m_weightMax) {
316 m_h1_engDmOverClus[iBin]->Fill(w, clus_weight);
317 }
318 if( m_hp2_DmWeight[iBin] ) {
319 //double isol = (*m_data->m_cls_isol)[i_cls];
321 }
322
323 } // i_dma
324 } // i_cls
325 } // i_ev
326
327 // --------------------------------------------------------------------------
328 // Fitting histograms
329 // --------------------------------------------------------------------------
330 std::cout << "CaloHadDMCoeffFit::process() -> Info. Fitting histograms..." << std::endl;
331 TF1 *fitFun = new TF1("fitFun","[0]+[1]*pow(x,[2])",1.,200000.);
332 double start_pars[3] = {100.0, 1.0, 1.0};
333 double start_errs[3] = {1.0, 1.0, 0.01};
334 // fitting profiles in silent mode
335 for(int i_size=0; i_size<m_HadDMCoeff->getSizeCoeffSet(); i_size++){
336 const CaloLocalHadCoeff::LocalHadArea *area = m_HadDMCoeff->getAreaFromBin(i_size);
337 if( area->getType() != CaloLocalHadDefs::AREA_DMFIT) continue;
338 TProfile *hp= m_hp_DmVsPrep[i_size];
339 if(hp->GetEntries() <25) continue;
340 for(int ip=0; ip<fitFun->GetNpar(); ip++) fitFun->SetParameter(ip,start_pars[ip]);
341 fitFun->SetParErrors(start_errs);
342 fitFun->FixParameter(2,1.0);
343 std::cout << "----------------------" << std::endl;
344 double fitlim1 = hp->GetBinLowEdge(1);
345 double fitlim2 = ProfileRefiner(hp);
346 //std::cout << "area->getTitle():" << i_size << " " << area->getTitle() << " hlim:"<< hp->GetBinCenter(hp->GetNbinsX()) + hp->GetBinWidth(hp->GetNbinsX())/2. << std::endl;
347 if(area->getTitle()=="ENG_CALIB_DEAD_TILEG3") {
348 //std::cout << " xxx ENG_CALIB_DEAD_TILEG3 " << std::endl;
349 fitlim2 = ProfileRefiner(hp,0.85);
350 }else{
351 fitlim2 = ProfileRefiner(hp, 0.98);
352 }
353 //std::cout << " xxx fitlim2:" << fitlim2<< std::endl;
354 hp->Fit(fitFun,"0Q","",fitlim1,fitlim2); // silent fit (drawing option for fit results is switched off)
355 TF1 *f=hp->GetFunction(fitFun->GetName());
356 if(f) {
357 f->ResetBit(TF1::kNotDraw);
358 f->SetLineWidth(1);
359 f->SetLineColor(2);
360 }
361 } // i_size
362
363 // --------------------------------------------------------------------------
364 // Making output coefficients
365 // --------------------------------------------------------------------------
366 std::cout << "CaloHadDMCoeffFit::process() -> Info. Making new coefficients ... " << std::endl;
368
369 m_FitData.clear();
370 m_FitData.resize(m_HadDMCoeff->getSizeCoeffSet(), nullptr);
371
372 for(int i_size=0; i_size<(int)new_data->getSizeCoeffSet(); i_size++){
373 const CaloLocalHadCoeff::LocalHadArea *dmArea = m_HadDMCoeff->getAreaFromBin(i_size);
374 int n_area;
375 for(n_area=0; n_area<m_HadDMCoeff->getSizeAreaSet(); n_area++){
376 if(dmArea == m_HadDMCoeff->getArea(n_area)) {
377 break;
378 }
379 }
380
381 const CaloLocalHadCoeff::LocalHadCoeff *old_pars = m_HadDMCoeff->getCoeff(i_size);
382 CaloLocalHadCoeff::LocalHadCoeff pars = *old_pars;
383
384 std::vector<int > indexes; // DIM_EMFRAC, DIM_SIDE, DIM_ETA, DIM_PHI, DIM_ENER, DIM_LAMBDA
385 m_HadDMCoeff->bin2indexes(i_size, indexes);
388 float ener = dimEner->getXmin() + dimEner->getDx()*(indexes[CaloLocalHadCoeffHelper::DIM_ENER]);
389 float eta = dimEta->getXmin() + dimEta->getDx()*(indexes[CaloLocalHadCoeffHelper::DIM_ETA]+0.5);
390
391 if( dmArea->getType() == CaloLocalHadDefs::AREA_DMFIT){
392 TF1 *fitFun = m_hp_DmVsPrep[i_size]->GetFunction("fitFun");
393 if(fitFun){
394 float p0 = fitFun->GetParameter(0);
395 float s0 = fitFun->GetParError(0);
396 float p1 = fitFun->GetParameter(1);
397 float s1 = fitFun->GetParError(1);
398 m_FitData[i_size] = new FitData(p0, s0, p1, s1);
399 // checking fit quality
400 if(fitFun->GetNDF() >= 2 && fitFun->GetChisquare()/fitFun->GetNDF() < 100. && p1 !=0.0 && p1>0.2) {
401 //if(fitFun->GetNDF() > 2 && fitFun->GetChisquare()/fitFun->GetNDF() < 15. && p1 !=0.0 && s1<fabs(p1) && p1>0.05) {
402 m_FitData[i_size]->isOK = true;
403 m_FitData[i_size]->descr = std::string("OK");
404 pars[0] = -1.0*p0/p1;
405 pars[1] = 1./p1;
406 pars[2] = 1.0;
407 std::cout << "i_size ok " << i_size << " " << m_FitData[i_size]->descr << std::endl;
408 }else{
409 m_FitData[i_size]->isOK = false;
410 m_FitData[i_size]->descr = std::string("failed");
411 if(fitFun->GetNDF() < 2 || fitFun->GetChisquare()/fitFun->GetNDF() > 100.) m_FitData[i_size]->descr += std::string(" NDFCHI");
412 if(p1==0 || p1<0.2) m_FitData[i_size]->descr += std::string(" p1");
413 //if(s1<fabs(p1)) m_FitData[i_size]->descr += std::string(" s1");
414 std::cout << "i_size failed " << i_size << " " << m_FitData[i_size]->descr << std::endl;
415 }
416 }else{
417 m_FitData[i_size] = new FitData();
418 m_FitData[i_size]->isOK = false;
419 m_FitData[i_size]->descr = "failed noFitFun";
420 std::cout << "i_size nofitfun " << i_size << " " << m_FitData[i_size]->descr << std::endl;
421 }
422
423// if(m_isTestbeam) {
424// // setting coefficients for neutral pions in between emec and hec equal to coefficients of charged pions
425// if(dmArea->getTitle() == "ENG_CALIB_DEAD_HEC0" && indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC]==1) {
426// std::vector<int > m_tmp = indexes;
427// m_tmp[CaloLocalHadCoeffHelper::DIM_EMFRAC]=0;
428// int n_area;
429// for(n_area=0; n_area<m_HadDMCoeff->getSizeAreaSet(); n_area++){
430// if(dmArea == m_HadDMCoeff->getArea(n_area)) {
431// break;
432// }
433// }
434// int iBin = m_HadDMCoeff->getBin( n_area, m_tmp);
435// const CaloLocalHadCoeff::LocalHadCoeff *had_pars = new_data->getCoeff(iBin);
436// pars[0] = (*had_pars)[0];
437// pars[1] = (*had_pars)[1];
438// pars[2] = (*had_pars)[2];
439// }
440// } // isTestbeam
441
442 // neutral pions doesn't reach the material between emec-hec and emb-tile, so normally the fit is screwed there
443 // let's set fit coefficients for neutral pions in this DM area equal to the charged pion
444 if( indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC]==1 && (dmArea->getTitle() == "ENG_CALIB_DEAD_HEC0" || dmArea->getTitle() == "ENG_CALIB_DEAD_TILE0" ) ) {
445 std::vector<int > tmp = std::move(indexes);
447 int iBin = m_HadDMCoeff->getBin( n_area, tmp);
448 const CaloLocalHadCoeff::LocalHadCoeff *had_pars = new_data->getCoeff(iBin);
449 pars[0] = (*had_pars)[0];
450 pars[1] = (*had_pars)[1];
451 pars[2] = (*had_pars)[2];
452 char str[128];
453 sprintf(str, "Afrom %d",iBin);
454 m_FitData[i_size]->descr += std::string(str);
455 //std::cout << "xxx A i_size:" << i_size << " iBin:" << iBin << " " << m_FitData[i_size]->descr << std::endl;
456 }else if( indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC]==1 && (dmArea->getTitle() == "ENG_CALIB_DEAD_TILEG3") ) {
457 // For material before tile scintillator things are weared, in some eta regions neutral pions have good correlation
458 // between signal in scintillator and dead material energy, in some eta regions correlation is absent. So we select here
459 // following approach: if fit failed, take coefficients from charged pion data and use them for neutral
460 if(!m_FitData[i_size]->isOK) {
461 std::vector<int > tmp = std::move(indexes);
463 int iBin = m_HadDMCoeff->getBin( n_area, tmp);
464 const CaloLocalHadCoeff::LocalHadCoeff *had_pars = new_data->getCoeff(iBin);
465 pars[0] = (*had_pars)[0];
466 pars[1] = (*had_pars)[1];
467 pars[2] = (*had_pars)[2];
468 char str[128];
469 sprintf(str, "Afrom %d",iBin);
470 m_FitData[i_size]->descr += std::string(str);
471 std::cout << "xxx A i_size:" << i_size << " iBin:" << iBin << " " << m_FitData[i_size]->descr << std::endl;
472 }
473 }else{
474 // if fit went wrong, lets take fit results from neighboting eta bin
475 if(!m_FitData[i_size]->isOK) {
476 std::vector<int > tmp = std::move(indexes);
479 }else if(tmp[CaloLocalHadCoeffHelper::DIM_ETA]==0){
481 }
482 int iBin = m_HadDMCoeff->getBin( n_area, tmp);
483 const CaloLocalHadCoeff::LocalHadCoeff *had_pars = new_data->getCoeff(iBin);
484 pars[0] = (*had_pars)[0];
485 pars[1] = (*had_pars)[1];
486 pars[2] = (*had_pars)[2];
487 char str[128];
488 sprintf(str, "Bfrom %d",iBin);
489 m_FitData[i_size]->descr += std::string(str);
490 //std::cout << "xxx B i_size:" << i_size << " iBin:" << iBin << " " << m_FitData[i_size]->descr << std::endl;
491 }
492 }
493
494
495 }else if(dmArea->getType() == CaloLocalHadDefs::AREA_DMLOOKUP){
496 //if(dmArea->getTitle() == "ENG_CALIB_DEAD_LEAKAGE") {
497 // right tail cut
498 if( m_h1_engDmOverClus[i_size]->GetEntries() > 0) {
499 double rcut = 0.92;
500 if(m_isTestbeam){ // H6 testbeam
501 rcut = 0.93;
502 if(indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC] == 1) rcut = 0.75;
503
504 if(indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC]==0 && ener < 4.69 && fabs(eta)>2.8) {
505 rcut = 0.85;
506 }
507 if(indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC]==1 && ener < 4.69 && fabs(eta)>3.0) {
508 rcut = 0.85;
509 }
510 }else{ // atlas
511 rcut = 0.95;
512 if(indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC] == 1) {
513 rcut = 0.80;
514 if(fabs(eta) > 1.44 && fabs(eta) < 1.56) {
515 rcut = 0.60;
516 }
517 }
518 if(fabs(eta) > 4.0) rcut = rcut*0.7;
519 }
520 double w = GetAverageWithoutRightTail(m_h1_engDmOverClus[i_size], rcut);
521 pars[0] = 1.0+w;
522 }else{
523 pars[0] = 1.0;
524 }
525 pars[1] = float(m_h1_engDmOverClus[i_size]->GetEntries());
526 pars[2] = float(m_engDm[i_size]->m_aver);
527
528// }else if (dmArea->getTitle() == "ENG_CALIB_DEAD_UNCLASS") {
529// if(m_isTestbeam) {
530// // right tail cut
531// if( m_h1_engDmOverClus[i_size]->GetEntries() > 0) {
532// //double rcut = 0.99;
533// double rcut = 0.93;
534// if(indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC]==0 && ener < 4.69 ) {
535// //rcut = 0.88;
536// rcut = 0.85;
537// }
538// if(indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC] == 1) {
539// //rcut = 0.79; //clsw=lin
540// rcut = 0.85; //clsw=1.0
541// }
542// double w = GetAverageWithoutRightTail(m_h1_engDmOverClus[i_size], rcut);
543// pars[0] = 1.0+w;
544// }else{
545// pars[0] = 1.0;
546// }
547// pars[1] = float(m_h1_engDmOverClus[i_size]->GetEntries());
548// pars[2] = float(m_engDm[i_size]->m_aver);
549//
550// }else{ // Atlas
551// // simple average with weightMax rejection
552// if( m_hp2_DmWeight[i_size]->GetEntries() > 0) {
553// int binx = indexes[CaloLocalHadCoeffHelper::DIM_ENER] + 1;
554// int biny = indexes[CaloLocalHadCoeffHelper::DIM_LAMBDA] + 1;
555// double w = m_hp2_DmWeight[i_size]->GetBinContent(binx, biny);
556// // manual hack to decrease the level of correction
557// if(indexes[CaloLocalHadCoeffHelper::DIM_EMFRAC]==1 && fabs(eta)>1.3 && fabs(eta)<1.7) {
558// if(ener < 4.5){
559// w = w*0.78;
560// }else{
561// w = w*0.92;
562// }
563// }
564// pars[0] = w + 1.0;
565// //pars[0] = w;
566// pars[1] = float(m_hp2_DmWeight[i_size]->GetBinEntries(m_hp2_DmWeight[i_size]->GetBin(binx, biny)));
567// }else{
568// pars[0] = 1.0;
569// }
570// pars[2] = float(m_engDm[i_size]->m_aver);
571// }
572//
573// }else{
574// std::cout << "Warning! Unknown lookup area " << dmArea->getTitle() << std::endl;
575// }
576 }
577 new_data->setCoeff(i_size, pars);
578 }
579
580 return new_data;
581}
582
583
584
585/* ****************************************************************************
586Making nice postscript report
587**************************************************************************** */
588void CaloHadDMCoeffFit::make_report(std::string &sreport)
589{
590 float etaRegions[][2] = { {0.0, 0.1}, {1.4, 1.7}, {2.0, 2.1}, {3.3, 3.4}, {4.0, 4.1} };
591 std::cout << "CaloHadDMCoeffFit::make_report() -> Info. Making report..." << std::endl;
592 gStyle->SetCanvasColor(10);
593 gStyle->SetCanvasBorderMode(0);
594 gStyle->SetPadColor(10);
595 gStyle->SetPadBorderMode(0);
596 gStyle->SetPalette(1);
597 gStyle->SetTitleBorderSize(1);
598 gStyle->SetTitleFillColor(10);
599 int cc_xx = 768, cc_yy = 1024;
600 char str[1024], hname[1024];
601 TText tex;
602 gROOT->SetBatch(kTRUE);
603 gErrorIgnoreLevel=3; // root global variables to supress text output in ->Print() methods
604 std::string sfname = sreport;
605 TCanvas *ctmp = new TCanvas("ctmp","ctmp", cc_xx, cc_yy);
606 sfname += "[";
607 ctmp->Print(sfname.c_str());
608 sfname = sreport;
609
610 TCanvas *c1ps = new TCanvas("c1ps","c1ps", cc_xx, cc_yy);
611 int npage=1;
612 for(int i_dms=0; i_dms<(int)m_HadDMCoeff->getSizeAreaSet(); i_dms++) {
613 const CaloLocalHadCoeff::LocalHadArea *dmArea = m_HadDMCoeff->getArea(i_dms);
614 if( dmArea->getType() != CaloLocalHadDefs::AREA_DMFIT) continue;
615
622
623 std::vector<int > v_indx;
624 v_indx.resize(CaloLocalHadCoeffHelper::DIM_UNKNOWN, 0);
625 for(int i_frac=0; i_frac<dimFrac->getNbins(); i_frac++) {
626 for(int i_ener=0; i_ener<dimEner->getNbins(); i_ener++) {
627 for(int i_lambda=0; i_lambda<dimLambda->getNbins(); i_lambda++) {
628 for(int i_side=0; i_side<dimSide->getNbins(); i_side++) {
629 for(int i_phi=0; i_phi<dimPhi->getNbins(); i_phi++) {
631 v_indx[CaloLocalHadCoeffHelper::DIM_ENER] = i_ener;
632 v_indx[CaloLocalHadCoeffHelper::DIM_LAMBDA] = i_lambda;
633 v_indx[CaloLocalHadCoeffHelper::DIM_PHI] = i_phi;
634 v_indx[CaloLocalHadCoeffHelper::DIM_SIDE] = i_side;
635 for(int k=0; k<2; k++){ // loop over TH2F and TProfile
636 c1ps->Clear();
637 TPad *pad1 = new TPad("p1ps","p1ps",0.0, 0.95, 1.0, 1.0); pad1->Draw();
638 TPad *pad2 = new TPad("p2ps","p2ps",0.0, 0.0, 1.0, 0.95); pad2->Draw();
639 // nice header with energy (in GeV) range on it
640 pad1->cd();
641 float ener1 = pow(10,dimEner->getXmin()+dimEner->getDx()*i_ener)/GeV;
642 float ener2 = pow(10,dimEner->getXmin()+dimEner->getDx()*(i_ener+1))/GeV;
643 sprintf(str,"%s em:%d ener:%d> %5.1f-%6.1f phi:%d s:%d",dmArea->getTitle().c_str(),i_frac,i_ener, ener1, ener2, i_phi, i_side);
644 tex.SetTextSize(0.4);
645 tex.SetTextColor(kBlue);
646 tex.DrawTextNDC(0.05,0.1,str);
647 // number of pad's divisions to display all eta-histograms on one page
648 int ndiv = dimEta->getNbins();
649 if(ndiv <=6){
650 pad2->Divide(2,3);
651 }else if(ndiv <=9){
652 pad2->Divide(3,3);
653 }else if(ndiv <=12){
654 pad2->Divide(3,4);
655 }else if(ndiv <=20){
656 pad2->Divide(4,5);
657 }else if(ndiv <=30){
658 pad2->Divide(5,6);
659 } else {
660 pad2->Divide(6,6);
661 }
662 for(int i_eta=0; i_eta<dimEta->getNbins(); i_eta++){
663 v_indx[CaloLocalHadCoeffHelper::DIM_ETA] = i_eta;
664 int iBin = m_HadDMCoeff->getBin(i_dms, v_indx);
665 if(iBin == -1) {
666 std::cout << "Panic! iBin == -1, i_frac:" << i_frac << " i_ener:" << i_ener << " i_lambda:" << i_lambda << " i_eta:" << i_eta << std::endl;
667 return;
668 }
669 TH1 *hh = nullptr;
670 if(k==0) {
671 hh = m_h2_DmVsPrep[iBin];
672 }else{
673 hh = m_hp_DmVsPrep[iBin];
674 }
675 float eta1 = dimEta->getXmin() + dimEta->getDx()*i_eta;
676 float eta2 = dimEta->getXmin() + dimEta->getDx()*(i_eta+1);
677 pad2->cd(1+i_eta); gPad->SetGrid();
678 hh->GetXaxis()->SetNdivisions(508);
679 hh->GetYaxis()->SetNdivisions(508);
680 hh->GetXaxis()->SetTitle("edmtrue");
681 hh->GetYaxis()->SetTitle("eprep");
682 if(k==0) { // TH2F
683 gStyle->SetOptStat(111111);
684 m_h2_DmVsPrep[iBin]->Draw("box");
685 }else { // TProfile
686 gStyle->SetOptStat(1);
687 gStyle->SetOptFit();
688 //m_hp_DmVsPrep[iBin]->SetMaximum(m_engPrep[iBin]->m_aver + 3.0*sqrt(m_engPrep[iBin]->m_rms));
689 m_hp_DmVsPrep[iBin]->SetMaximum( m_h2_DmVsPrep[iBin]->GetYaxis()->GetXmax() );
690 m_hp_DmVsPrep[iBin]->Draw();
691 }
692 tex.SetTextSize(0.095);
693 tex.SetTextColor(kBlue);
694 sprintf(str,"%4.2f-%4.2f",eta1,eta2);
695 tex.DrawTextNDC(0.2,0.8,str);
696 sprintf(str,"%d",iBin);
697 tex.DrawTextNDC(0.5,0.28,str);
698 // get parameters and draw values nicely
699 if(m_FitData[iBin]) {
700 tex.SetTextColor(kBlue);
701 float p0inv, s0inv, p1inv, s1inv;
702 m_FitData[iBin]->getInverted(p0inv, s0inv, p1inv, s1inv);
703 tex.SetTextColor(kBlue); tex.SetTextSize(0.07);
704 sprintf(str,"p :%6.3f %6.3f",m_FitData[iBin]->p0, m_FitData[iBin]->p1);
705 tex.DrawTextNDC(0.2,0.7,str);
706 sprintf(str,"inv:%6.3f %6.3f",p0inv, p1inv);
707 tex.DrawTextNDC(0.2,0.7-0.1,str);
708 sprintf(str,"err:%6.3f %6.3f",s0inv, s1inv);
709 tex.DrawTextNDC(0.2,0.7-0.2,str);
710 if(!m_FitData[iBin]->isOK) {
711 tex.SetTextColor(kRed);
712 }
713 tex.DrawTextNDC(0.18,0.7-0.3,m_FitData[iBin]->descr.c_str());
714 }else{
715 sprintf(str,"Oops!");
716 tex.SetTextColor(kMagenta); tex.SetTextSize(0.095);
717 tex.DrawTextNDC(0.2,0.7,str);
718 }
719 } // i_eta
720 c1ps->Print(sfname.c_str());
721 printf("page:%d dmArea->m_title:%s frac:%d ener:%d phi:%d side:%d\n",
722 npage++, dmArea->getTitle().c_str(), i_frac, i_ener, i_phi, i_side);
723 } // k
724 } // i_phi
725 } // i_side
726 } // i_lambda
727 } // i_ener
728
729 /* **************************************************
730 as a function of energy
731 ************************************************** */
732 for(int i_side=0; i_side<dimSide->getNbins(); i_side++){
733 v_indx[CaloLocalHadCoeffHelper::DIM_SIDE] = i_side;
734 for(int i_phi=0; i_phi<dimPhi->getNbins(); i_phi++){
735 v_indx[CaloLocalHadCoeffHelper::DIM_PHI] = i_phi;
736 for(int i_par=0; i_par<2; i_par++){
737 c1ps->Clear();
738 TPad *pad1 = new TPad("p1ps","p1ps",0.0, 0.95, 1.0, 1.0); pad1->Draw();
739 TPad *pad2 = new TPad("p2ps","p2ps",0.0, 0.0, 1.0, 0.95); pad2->Draw();
740 // nice header with energy (in GeV) range on it
741 pad1->cd();
742 sprintf(str,"Summary: %s par:%d frac:%d phi:%d side:%d",dmArea->getTitle().c_str(),i_par,i_frac, i_phi, i_side);
743 tex.SetTextSize(0.5); tex.SetTextColor(kBlue);
744 tex.DrawTextNDC(0.05,0.1,str);
745 // number of pad's divisions to display all eta-histograms on one page
746 int ndiv = dimEta->getNbins();
747 if(ndiv <=6){
748 pad2->Divide(2,3);
749 }else if(ndiv <=9){
750 pad2->Divide(3,3);
751 }else if(ndiv <=12){
752 pad2->Divide(3,4);
753 }else if(ndiv <=20){
754 pad2->Divide(4,5);
755 }else if(ndiv <=30){
756 pad2->Divide(5,6);
757 } else {
758 pad2->Divide(6,6);
759 }
760 for(int i_eta=0; i_eta<dimEta->getNbins(); i_eta++){
761 v_indx[CaloLocalHadCoeffHelper::DIM_ETA] = i_eta;
762 TGraphErrors *gr = new TGraphErrors(dimEner->getNbins());
763 for(int i_ener=0; i_ener<dimEner->getNbins(); i_ener++) {
764 v_indx[CaloLocalHadCoeffHelper::DIM_ENER] = i_ener;
765 int iBin = m_HadDMCoeff->getBin(i_dms, v_indx);
766 float y(0), ye(0);
767 if(iBin >= 0 && m_FitData[iBin]) {
768 float p0inv, s0inv, p1inv, s1inv;
769 m_FitData[iBin]->getInverted(p0inv, s0inv, p1inv, s1inv);
770 if(i_par==0) {
771 y = p0inv;
772 ye = s0inv;
773 }else{
774 y = p1inv;
775 ye = s1inv;
776 }
777 }
778 gr->SetPoint(i_ener, dimEner->getXmin()+i_ener*dimEner->getDx(), y);
779 gr->SetPointError(i_ener, 0.0, ye);
780 } // i_ener
781 pad2->cd(1+i_eta); gPad->SetGrid();
782 gr->Draw("apl");
783 float eta1 = dimEta->getXmin() + dimEta->getDx()*i_eta;
784 float eta2 = dimEta->getXmin() + dimEta->getDx()*(i_eta+1);
785 tex.SetTextSize(0.095);
786 tex.SetTextColor(kBlue);
787 sprintf(str,"%4.2f-%4.2f phibin:%d",eta1,eta2,i_phi);
788 tex.DrawTextNDC(0.2,0.8,str);
789 } // i_eta
790 c1ps->Print(sfname.c_str());
791 printf("page:%d dmArea->m_title:%s frac:%d Energy summary\n",
792 npage++, dmArea->getTitle().c_str(), i_frac);
793 } // i_par
794 } // i_phi
795 } // i_side
796 } // i_frac
797 } // i_dms
798
799 /* *********************************************
800 now making TH2F histograms which represents our look-up data
801 ********************************************** */
802 for(int i_dms=0; i_dms<int(m_HadDMCoeff->getSizeAreaSet()); i_dms++) {
803 const CaloLocalHadCoeff::LocalHadArea *dmArea = m_HadDMCoeff->getArea(i_dms);
804 if( dmArea->getType() != CaloLocalHadDefs::AREA_DMLOOKUP) continue;
805
812
813 for(int i_frac=0; i_frac<dimFrac->getNbins(); i_frac++){
814 for(int i_side=0; i_side<dimSide->getNbins(); i_side++){
815 for(int i_phi=0; i_phi<dimPhi->getNbins(); i_phi++){
816 for(int i_eta=0; i_eta<dimEta->getNbins(); i_eta++){
817 // check if eta bin within region of interest
818 bool ShowThisEta = false;
819 for(int i=0; i<int(sizeof(etaRegions)/sizeof(float)/2); i++){
820 float xeta = dimEta->getXmin()+i_eta*dimEta->getDx();
821 if(xeta>=etaRegions[i][0] && xeta <= etaRegions[i][1]) {
822 ShowThisEta = true;
823 break;
824 }
825 }
826 if(!ShowThisEta) continue;
827 std::vector<int > v_indx;
828 v_indx.resize(CaloLocalHadCoeffHelper::DIM_UNKNOWN, 0);
830 v_indx[CaloLocalHadCoeffHelper::DIM_SIDE] = i_side;
831 v_indx[CaloLocalHadCoeffHelper::DIM_ETA] = i_eta;
832 v_indx[CaloLocalHadCoeffHelper::DIM_PHI] = i_phi;
835 // lets find bins for axes
836 Double_t xx[100];
837 bzero(xx,100*sizeof(Double_t));
838 for(int i_ener=0; i_ener<dimEner->getNbins()+1; i_ener++){
839 xx[i_ener] = dimEner->getXmin() + dimEner->getDx()*i_ener;
840 }
841 Double_t yy[100];
842 bzero(yy,100*sizeof(Double_t));
843 for(int i_lambda=0; i_lambda<dimLambda->getNbins()+1; i_lambda++){
844 yy[i_lambda] = dimLambda->getXmin() + dimLambda->getDx()*i_lambda;
845 }
846 int ibin_min = m_HadDMCoeff->getBin(i_dms, v_indx);
847 v_indx[CaloLocalHadCoeffHelper::DIM_ENER] = dimEner->getNbins() - 1;
848 v_indx[CaloLocalHadCoeffHelper::DIM_LAMBDA] = dimLambda->getNbins() - 1;
849 int ibin_max = m_HadDMCoeff->getBin(i_dms, v_indx);
850 const int n_prof= 2;
851 TH2F *hp2[n_prof];
852 sprintf(hname,"%s dm:%d frac:%d eta:%d phi:%d indx:%d-%d<ecls>/<edm>",dmArea->getTitle().c_str(),i_dms, i_frac, i_eta, i_phi, ibin_min, ibin_max);
853 hp2[0] = new TH2F(hname, hname, dimEner->getNbins(), xx, dimLambda->getNbins(), yy);
854 sprintf(hname,"%s dm:%d frac:%d eta:%d phi:%d indx:%d-%d nev",dmArea->getTitle().c_str(),i_dms, i_frac, i_eta, i_phi, ibin_min, ibin_max);
855 hp2[1] = new TH2F(hname, hname, dimEner->getNbins(), xx, dimLambda->getNbins(), yy);
856 for(int i=0; i<n_prof; i++){
857 hp2[i]->GetXaxis()->SetTitle("log10(E_{cls}(MeV))");
858 hp2[i]->GetYaxis()->SetTitle("log10(#lambda)");
859 hp2[i]->GetZaxis()->SetTitle("<E_{cls}>/<E_{dm}>");
860 }
861 // filling histograms
862 float hp_aver[n_prof];
863 bzero(hp_aver,n_prof*sizeof(float));
864 std::vector<std::pair<int, int> > v_occupancy;
865 for(int i_ener=0; i_ener<dimEner->getNbins(); i_ener++) {
866 for(int i_lambda=0; i_lambda<dimLambda->getNbins(); i_lambda++) {
867 v_indx[CaloLocalHadCoeffHelper::DIM_ENER] = i_ener;
868 v_indx[CaloLocalHadCoeffHelper::DIM_LAMBDA] = i_lambda;
869 int iBin = m_HadDMCoeff->getBin(i_dms, v_indx);
870 if(iBin == -1) {
871 std::cout << "Panic in pp3, Wrong iBin=-1" << std::endl;
872 exit(EXIT_FAILURE);
873 }
874 float x = m_engDmOverClus[iBin]->m_aver;
875 hp2[0]->Fill(dimEner->getXmin() + dimEner->getDx()*(i_ener+0.5), dimLambda->getXmin() + dimLambda->getDx()*(i_lambda+0.5), x);
876 hp_aver[0] += x;
877 x = float(m_engDm[iBin]->size());
878 hp2[1]->Fill(dimEner->getXmin() + dimEner->getDx()*(i_ener+0.5), dimLambda->getXmin() + dimLambda->getDx()*(i_lambda+0.5), x );
879 hp_aver[1] += x;
880 std::pair<int, int> pp;
881 pp.first = int(m_engDmOverClus[iBin]->size());
882 pp.second = iBin;
883 v_occupancy.push_back(pp);
884 }
885 }
886 TCanvas *cp2 = new TCanvas("cp2","cp2",cc_xx,cc_yy);
887 cp2->cd();
888 TPad *pad1 = new TPad("cp2_p1ps","cp2_p1ps",0.0, 0.95, 1.0, 1.0); pad1->Draw();
889 TPad *pad2 = new TPad("cp2_p2ps","cp2_p2ps",0.0, 0.0, 1.0, 0.95); pad2->Draw();
890 // nice header with energy (in GeV) range on it
891 pad1->cd();
892 sprintf(str,"Summary: %s frac:%d phi:%d side:%d",dmArea->getTitle().c_str(),i_frac, i_phi, i_side);
893 tex.SetTextSize(0.5); tex.SetTextColor(kBlue);
894 tex.DrawTextNDC(0.05,0.1,str);
895 pad2->Divide(2,3);
896 pad2->cd(1);
897 gPad->SetLogz();
898 gPad->SetRightMargin(0.2);
899 hp2[1]->Draw("colz");
900 pad2->cd(2);
901 if (ibin_min >= 0) {
902 m_hp2_DmWeight[ibin_min]->SetStats(0);
903 m_hp2_DmWeight[ibin_min]->SetMinimum(0.0);
904 m_hp2_DmWeight[ibin_min]->SetMaximum(1.0);
905 m_hp2_DmWeight[ibin_min]->Draw("colz");
906 }
907 // drawing bin number for histograms
908 for(int i_ener=0; i_ener<dimEner->getNbins(); i_ener++) {
909 v_indx[CaloLocalHadCoeffHelper::DIM_ENER] = i_ener;
910 int iBin = m_HadDMCoeff->getBin(i_dms, v_indx);
911 sprintf(str,"%d",iBin);
912 TLatex tex;
913 tex.SetTextSize(0.03);
914 tex.SetTextAngle(90.);
915 tex.DrawLatex(dimEner->getXmin() + dimEner->getDx()*(i_ener+0.5), dimLambda->getXmax()-2.*dimLambda->getDx(), str);
916 }
917
918 // drawing spectras in selected cells
919 std::sort(v_occupancy.begin(), v_occupancy.end());
920 std::reverse(v_occupancy.begin(), v_occupancy.end());
921 for(unsigned int i_spec=0; i_spec<v_occupancy.size(); i_spec++){
922 int iBin = v_occupancy[i_spec].second;
923 if(iBin == -1) {
924 std::cout << "Panic in pp3, Wrong iBin=-1" << std::endl;
925 exit(EXIT_FAILURE);
926 }
927 if(!m_h1_engDmOverClus[iBin]) {
928 std::cout << "Undefined h1 for " << iBin << std::endl;
929 continue;
930 }
931 std::vector<int > indexes;
932 m_HadDMCoeff->bin2indexes(iBin, indexes);
933 sprintf(str, "ibin:%d frac:%d ener:%d lambda:%d eta:%d phi:%d",iBin,
937 pad2->cd(3+i_spec);
938 m_h1_engDmOverClus[iBin]->SetStats(1);
939 m_h1_engDmOverClus[iBin]->SetTitle(str);
940 gStyle->SetOptStat(111111);
941 m_h1_engDmOverClus[iBin]->Draw();
943 sprintf(str,"rmscu:%8.5f evcu:%6.3f",m_engDmOverClus[iBin]->m_aver, av);
944 tex.SetTextSize(0.05);
945 tex.DrawTextNDC(0.2,0.7,str);
946 if(i_spec >= 3) break;
947 }
948 cp2->Print(sfname.c_str());
949 printf("page:%d dmArea->m_title:%s i_frac:%d eta:%d phi:%d side:%d\n",
950 npage, dmArea->getTitle().c_str(), i_frac, i_eta, i_phi, i_side);
951 npage++;
952 }// i_eta
953 }// i_phi
954 } // i_side
955 } // i_frac
956 } // i_dms
957
958 sfname = sreport;
959 sfname += "]";
960 ctmp->Print(sfname.c_str());
961}
962
963
964/* ****************************************************************************
965
966**************************************************************************** */
968{
969 for (PrepData* h : m_engClus) {
970 delete h;
971 }
972 m_engClus.clear();
973
974 for (PrepData* h : m_engPrep) {
975 delete h;
976 }
977 m_engPrep.clear();
978
979 for (PrepData* h : m_engDm) {
980 delete h;
981 }
982 m_engDm.clear();
983
984 for (PrepData* h : m_engDmOverClus) {
985 delete h;
986 }
987 m_engDmOverClus.clear();
988
989 for (TProfile* h : m_hp_DmVsPrep) {
990 delete h;
991 }
992 m_hp_DmVsPrep.clear();
993
994 for (TH2F* h : m_h2_DmVsPrep) {
995 delete h;
996 }
997 m_h2_DmVsPrep.clear();
998
999 for (TH1F* h : m_h1_engDmOverClus) {
1000 delete h;
1001 }
1002 m_h1_engDmOverClus.clear();
1003
1004 int i_size=0;
1005 for (TProfile2D* h : m_hp2_DmWeight) {
1006 if( i_size==getFirstEnerLambdaBin(i_size) ) delete h;
1007 i_size++;
1008 }
1009 m_hp2_DmWeight.clear();
1010
1011}
1012
1013
1014
1015/* ****************************************************************************
1016- reset bin content if number of entries <2 (i.e. with zero errors)
1017- calculates interval on x-axis [0,xLimRight] which contains 92% of all entries
1018**************************************************************************** */
1019double CaloHadDMCoeffFit::ProfileRefiner(TProfile *pH, double ev_ratio)
1020{
1021 double xLimRight(0);
1022// double e(0),cont(0),nent(0);
1023 double nevtot(0);
1024 for(int i=1; i<=pH->GetNbinsX(); i++){
1025 nevtot += pH->GetBinEntries(i);
1026 }
1027 //std::cout << "yyy entries:" << pH->GetEntries() << "bin0:" << pH->GetBinEntries(0) << " bin1:" << pH->GetBinEntries(pH->GetNbinsX()+1)
1028 //<< " nevtot: " << nevtot
1029 //<< " ev_ratio:" << ev_ratio
1030 //<< std::endl;
1031 double nev_in(0);
1032 double inv_nevtot = nevtot ? 1. / nevtot : 1;
1033 for(int i=1; i<=pH->GetNbinsX(); i++){
1034 nev_in += pH->GetBinEntries(i);
1035 if(pH->GetBinEffectiveEntries(i)<2.0){
1036 pH->SetBinEntries(i,0.0);
1037 pH->SetBinContent(i,0.0);
1038 }
1039 //if(xLimRight==0.0 && nev_in/nevtot > ev_ratio) xLimRight=pH->GetBinCenter(i) + pH->GetBinWidth(i)/2.;
1040 //std::cout << "yyy i:" << i << " bin_content:" << pH->GetBinContent(i) << " entries:" << pH->GetBinEntries(i)
1041 //<< " nev_in:" << nev_in << " nevtot:" << nevtot << " nev_in/nevtot:" << nev_in/nevtot
1042 //<< std::endl;
1043 if(nev_in*inv_nevtot > ev_ratio) {
1044 xLimRight=pH->GetBinCenter(i) + pH->GetBinWidth(i)/2.;
1045 //std::cout << "yyy i:" << i << " xLimRight:" << xLimRight << std::endl;
1046 break;
1047 }
1048 }
1049 if(xLimRight==0) xLimRight=pH->GetBinCenter(pH->GetNbinsX()) + pH->GetBinWidth(pH->GetNbinsX())/2.;
1050 //std::cout << "yyy results xLimRight:" << xLimRight << std::endl;
1051 return xLimRight;
1052}
1053
1054
1055
1056/* ****************************************************************************
1057- calculates average of TH1F histogram after rejection of events in the right tail
1058**************************************************************************** */
1059double CaloHadDMCoeffFit::GetAverageWithoutRightTail(TH1F *pH, double ev_ratio)
1060{
1061 double xLimRight=pH->GetBinCenter(pH->GetNbinsX()) + pH->GetBinWidth(pH->GetNbinsX())/2.;
1062 int nbinx = (int)pH->GetNbinsX();
1063// double nevtot = pH->Integral() - pH->GetBinContent(0) - pH->GetBinContent(nbinx+1);
1064 double nevtot = pH->Integral();
1065 double nev = 0;
1066 if(nevtot != 0.0) {
1067 const double inv_nevtot = 1. / nevtot;
1068 for(int i_bin = 0; i_bin <=nbinx; i_bin++){
1069 nev+= pH->GetBinContent(i_bin);
1070
1071 if( nev*inv_nevtot >= ev_ratio){
1072 xLimRight = pH->GetBinCenter(i_bin) + pH->GetBinWidth(i_bin)/2.;
1073 break;
1074 }
1075 }
1076 }
1077 double new_aver = 0.0;
1078 double sum = 0.0;
1079 for(int i_bin = 1; i_bin <=nbinx; i_bin++){
1080 if(pH->GetBinCenter(i_bin) <= xLimRight) {
1081 new_aver += pH->GetBinContent(i_bin)*pH->GetBinCenter(i_bin);
1082 sum += pH->GetBinContent(i_bin);
1083 }
1084 }
1085 if(sum != 0.0) {
1086 new_aver /= sum;
1087 if(new_aver < 0.0) new_aver = 0.0;
1088 }else{
1089 new_aver = pH->GetMean();
1090 }
1091 return new_aver;
1092}
1093
1094
1096{
1097 // DIM_EMFRAC, DIM_SIDE, DIM_ETA, DIM_PHI, DIM_ENER, DIM_LAMBDA
1098 const CaloLocalHadCoeff::LocalHadArea *dmArea = m_HadDMCoeff->getAreaFromBin(ibin);
1099 std::vector<int > indexes;
1100 m_HadDMCoeff->bin2indexes(ibin, indexes);
1103
1104 int n_area;
1105 for(n_area=0; n_area<m_HadDMCoeff->getSizeAreaSet(); n_area++){
1106 if(dmArea == m_HadDMCoeff->getArea(n_area)) {
1107 break;
1108 }
1109 }
1110 int refbin = m_HadDMCoeff->getBin(n_area, indexes);
1111 return refbin;
1112}
1113
1114
Scalar eta() const
pseudorapidity method
double area(double R)
#define gr
static Double_t s0
TGraphErrors * GetEntries(TH2F *histo)
Wrapper to avoid constant divisions when using units.
#define y
#define x
constexpr int pow(int base, int exp) noexcept
Header file for AthHistogramAlgorithm.
Data to read from special DeadMaterialTree.
CaloLocalHadCoeffHelper * m_HadDMHelper
std::vector< PrepData * > m_engClus
std::string m_NormalizationType
CaloHadDMCoeffData * m_data
double ProfileRefiner(TProfile *pH, double ev_ratio=0.92)
CaloLocalHadCoeff * m_HadDMCoeff
std::vector< PrepData * > m_engPrep
std::vector< FitData * > m_FitData
std::vector< PrepData * > m_engDmOverClus
double GetAverageWithoutRightTail(TH1F *pH, double ev_ratio=0.92)
std::vector< PrepData * > m_engDm
std::vector< TProfile * > m_hp_DmVsPrep
std::vector< TProfile2D * > m_hp2_DmWeight
void make_report(std::string &sfname)
std::vector< TH1F * > m_h1_engDmOverClus
std::vector< TH2F * > m_h2_DmVsPrep
CaloLocalHadCoeff * process(CaloHadDMCoeffData *myData, CaloLocalHadCoeff *myHadDMCoeff, bool isSingleParticle=true, bool tbflag=false)
int getFirstEnerLambdaBin(int ibin)
Definition of correction area.
const CaloLocalHadCoeff::LocalHadDimension * getDimension(int n_dim) const
to get dimension
const std::string & getTitle() const
return name
unsigned int getType() const
return area type
Class defines binning for user dimension.
float getDx() const
return size of bin
int getNbins() const
return number of bins
float getXmax() const
return maximum value for the last bin
float getXmin() const
return minimum value for the first bin
Hold binned correction data for local hadronic calibration procedure.
const LocalHadCoeff * getCoeff(const int &iBin) const
get data for given general bin number
void setCoeff(const int iBin, const LocalHadCoeff &theCoeff)
set new data
int getSizeCoeffSet() const
return total number of coefficient sets
std::vector< float > LocalHadCoeff
Correction parameters for one general bin.
static double angle_mollier_factor(double x)
int r
Definition globals.cxx:22
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
void reverse(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of reverse for DataVector/List.