ATLAS Offline Software
FitMonitor.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 
9 #include "LArSamplesMon/Data.h"
10 
11 #include "TTree.h"
12 #include "TFile.h"
13 #include "TMath.h"
14 
15 #include <iostream>
16 using std::cout;
17 using std::endl;
18 
19 using namespace LArSamples;
20 
21 
22 bool FitMonitor::makeSummary(const char* fileName) const
23 {
24  TFile* f = TFile::Open(fileName, "RECREATE");
25  TTree* summary = new TTree("summary", "");
26 
27  int tree_hash, tree_calo, tree_layer, tree_iEta, tree_iPhi, tree_feb, tree_channel, tree_index;
28  float tree_chi2, tree_refitChi2;
29  float tree_time, tree_energy, tree_adcMax, tree_refitTime, tree_refitADCMax;
30 
31  summary->Branch("hash", &tree_hash, "hash/I");
32  summary->Branch("calo", &tree_calo, "calo/I");
33  summary->Branch("layer", &tree_layer, "layer/I");
34  summary->Branch("iEta", &tree_iEta, "iEta/I");
35  summary->Branch("iPhi", &tree_iPhi, "iPhi/I");
36  summary->Branch("feb", &tree_feb, "feb/I");
37  summary->Branch("channel", &tree_channel, "channel/I");
38  summary->Branch("index", &tree_index, "index/I");
39  summary->Branch("time", &tree_time, "time/F");
40  summary->Branch("energy", &tree_energy, "energy/F");
41  summary->Branch("adcMax", &tree_adcMax, "adcMax/F");
42  summary->Branch("refitTime", &tree_refitTime, "refitTime/F");
43  summary->Branch("chi2", &tree_chi2, "chi2/F");
44  summary->Branch("refitChi2", &tree_refitChi2, "refitChi2/F");
45  summary->Branch("refitADCMax", &tree_refitADCMax, "refitADCMax/F");
46 
47  for (unsigned int i = 0; i < nChannels(); i++) {
48  const History* history = cellHistory(i);
49  if (!history || !history->isValid()) continue;
50 
51  if (i%1000 == 0) cout << "-> History hash = " << i << endl;
52  History* refitHistory = history->refit();
53 
54  tree_hash = i;
55  tree_calo = history->cellInfo()->calo();
56  tree_layer = history->cellInfo()->layer();
57  tree_iEta = history->cellInfo()->iEta();
58  tree_iPhi = history->cellInfo()->iPhi();
59  tree_feb = history->cellInfo()->feb();
60  tree_channel = history->cellInfo()->channel();
61 
62  for (unsigned int j = 0; j < history->nData(); j++) {
63  tree_index = j;
64  tree_chi2 = history->chi2(j);
65  tree_time = history->data(j)->ofcTime();
66  tree_adcMax = history->data(j)->adcMax();
67  tree_energy = history->data(j)->energy();
68  tree_refitChi2 = refitHistory->chi2(j);
69  tree_refitTime = refitHistory->data(j)->ofcTime();
70  tree_refitADCMax = refitHistory->data(j)->adcMax();
71  summary->Fill();
72  }
73  delete refitHistory;
74  }
75 
76  summary->Write();
77  cout << "Closing file..." << endl;
78 
79  f->Close();
80  return true;
81 }
82 
LArSamples::CellInfo::feb
short feb() const
Definition: CellInfo.cxx:102
python.CaloRecoConfig.f
f
Definition: CaloRecoConfig.py:127
LArSamples::CellInfo::calo
CaloId calo() const
Definition: CellInfo.h:50
LArSamples::History::cellInfo
const CellInfo * cellInfo() const
Definition: History.h:61
LArSamples::AbsLArCells::cellHistory
virtual const History * cellHistory(unsigned int i) const
Definition: AbsLArCells.cxx:57
LArSamples::Data::adcMax
double adcMax() const
Definition: Data.h:129
LArSamples::History
Definition: History.h:40
LArSamples
Definition: AbsShape.h:24
LArSamples::Data::ofcTime
double ofcTime() const
Definition: Data.h:116
LArSamples::CellInfo::iEta
short iEta() const
Definition: CellInfo.h:56
LArSamples::MonitorBase::nChannels
unsigned int nChannels() const
Definition: LArCalorimeter/LArSamplesMon/src/MonitorBase.cxx:45
LArSamples::History::data
const Data * data(unsigned int i) const
Definition: History.cxx:88
FortranAlgorithmOptions.fileName
fileName
Definition: FortranAlgorithmOptions.py:13
LArSamples::History::refit
History * refit(Chi2Params pars=DefaultChi2) const
Definition: History.cxx:232
FitMonitor.h
lumiFormat.i
int i
Definition: lumiFormat.py:92
LArSamples::History::nData
unsigned int nData() const
Definition: History.h:56
LArSamples::History::chi2
double chi2(int i, int lwb=-1, int upb=-1, int chi2Params=DefaultChi2, ShapeErrorType shapeErrorType=BestShapeError, unsigned int *nDof=0) const
Definition: History.cxx:161
LArSamples::FitMonitor::makeSummary
bool makeSummary(const char *fileName) const
Definition: FitMonitor.cxx:22
LArSamples::CellInfo::iPhi
short iPhi() const
Definition: CellInfo.h:62
LArSamples::CellInfo::channel
short channel() const
Definition: CellInfo.h:76
Data.h
History.h
LArSamples::History::isValid
bool isValid() const
Definition: History.cxx:149
Interface.h
LArSamples::CellInfo::layer
short layer() const
Definition: CellInfo.h:53
LArSamples::Data::energy
double energy() const
Definition: Data.h:113
SCT_Monitoring::summary
@ summary
Definition: SCT_MonitoringNumbers.h:65