-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathORHistWriter.cc
141 lines (125 loc) · 4.48 KB
/
ORHistWriter.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
// ORHistWriter.cc
#include "ORHistWriter.hh"
#include "TH2.h"
#include "TH3.h"
#include "ORLogger.hh"
#include "TROOT.h"
using namespace std;
ORHistWriter::ORHistWriter(ORVHistDecoder* histDecoder) :
ORDataProcessor(histDecoder)
{
fHistDecoder = histDecoder;
}
ORDataProcessor::EReturnCode ORHistWriter::StartProcessing()
{
if(fHistDecoder == NULL) {
ORLog(kWarning) << "StartProcessing(): fHistDecoder was NULL, can't proceed";
KillProcessor();
return kFailure;
}
return kSuccess;
}
ORDataProcessor::EReturnCode ORHistWriter::StartRun()
{
map<int, TH1*>::iterator iHist;
for (iHist = fHists.begin(); iHist != fHists.end(); iHist++) {
if (iHist->second != NULL) iHist->second->Reset();
}
return kSuccess;
}
ORDataProcessor::EReturnCode ORHistWriter::ProcessMyDataRecord(UInt_t* record)
{
int iHist = fHistDecoder->GetHistIndex(record);
// look up in map only once for a typical record
TH1* hist = fHists[iHist];
if(hist == NULL) {
switch(fHistDecoder->GetNDim()) {
// fHists are owned by parent root file (or gROOT);
// they will be deleted upon fFile->Close() (or end of program)
case 1:
fHists[iHist] = new TH1D(
fHistDecoder->GetHistName(iHist).c_str(), fHistDecoder->GetHistTitle(iHist).c_str(),
fHistDecoder->GetNbinsX(), fHistDecoder->GetXLo(), fHistDecoder->GetXHi()
);
fHists[iHist]->SetXTitle(fHistDecoder->GetXTitle().c_str());
break;
case 2:
fHists[iHist] = new TH2D(
fHistDecoder->GetHistName(iHist).c_str(), fHistDecoder->GetHistTitle(iHist).c_str(),
fHistDecoder->GetNbinsX(), fHistDecoder->GetXLo(), fHistDecoder->GetXHi(),
fHistDecoder->GetNbinsY(), fHistDecoder->GetYLo(), fHistDecoder->GetYHi()
);
fHists[iHist]->SetXTitle(fHistDecoder->GetXTitle().c_str());
fHists[iHist]->SetYTitle(fHistDecoder->GetYTitle().c_str());
break;
case 3:
fHists[iHist] = new TH3D(
fHistDecoder->GetHistName(iHist).c_str(), fHistDecoder->GetHistTitle(iHist).c_str(),
fHistDecoder->GetNbinsX(), fHistDecoder->GetXLo(), fHistDecoder->GetXHi(),
fHistDecoder->GetNbinsY(), fHistDecoder->GetYLo(), fHistDecoder->GetYHi(),
fHistDecoder->GetNbinsZ(), fHistDecoder->GetZLo(), fHistDecoder->GetZHi()
);
fHists[iHist]->SetXTitle(fHistDecoder->GetXTitle().c_str());
fHists[iHist]->SetYTitle(fHistDecoder->GetYTitle().c_str());
fHists[iHist]->SetZTitle(fHistDecoder->GetZTitle().c_str());
break;
default:
ORLog(kWarning) << "ProcessMyDataRecord(): Can't handle more than "
<< "3 dimensions; using 1..." << endl;
fHists[iHist] = new TH1D(
fHistDecoder->GetHistName(iHist).c_str(), fHistDecoder->GetHistTitle(iHist).c_str(),
fHistDecoder->GetNbinsX(), fHistDecoder->GetXLo(), fHistDecoder->GetXHi()
);
}
hist = fHists[iHist];
}
for(size_t iEntry=0; iEntry < fHistDecoder->GetNEntries(record); iEntry++) {
switch(fHistDecoder->GetNDim()) {
case 1:
hist->Fill(fHistDecoder->GetX(record, iEntry),
fHistDecoder->GetWeight(record, iEntry));
break;
case 2:
((TH2D*) hist)->Fill(fHistDecoder->GetX(record, iEntry),
fHistDecoder->GetY(record, iEntry),
fHistDecoder->GetWeight(record, iEntry));
break;
case 3:
((TH3D*) hist)->Fill(fHistDecoder->GetX(record, iEntry),
fHistDecoder->GetY(record, iEntry),
fHistDecoder->GetZ(record, iEntry),
fHistDecoder->GetWeight(record, iEntry));
break;
default:
hist->Fill(fHistDecoder->GetX(record, iEntry),
fHistDecoder->GetWeight(record, iEntry));
break;
}
}
return kSuccess;
}
ORDataProcessor::EReturnCode ORHistWriter::EndRun()
{
map<int, TH1*>::iterator iHist;
for (iHist = fHists.begin(); iHist != fHists.end(); iHist++) {
if (iHist->second != NULL) {
iHist->second->Write();
iHist->second->Reset();
}
}
// if a TFile is open, the histos will be deleted by the file. If not,
// gDirectory will be named "Rint", and the histos will remain in memory.
if(string(gDirectory->GetName()) != "Rint") {
fHists.clear();
}
return kSuccess;
}
ORDataProcessor::EReturnCode ORHistWriter::EndProcessing()
{
if(string(gDirectory->GetName()) != "Rint") return kSuccess;
map<int, TH1*>::iterator iHist;
for (iHist = fHists.begin(); iHist != fHists.end(); iHist++) {
delete iHist->second;
}
return kSuccess;
}