37{
40
41 struct {
42 double XMIN = -0.2*PI;
43 double XMAX = +0.2*PI;
44 double FACTOR = 1.0;
45 double YMIN = 10.0;
46 } parameters;
47
50 string detectorFile;
52 string option;
54
55 try {
56
58
63
64 JParser<> zap(
"Program to fit histogram in output of JMuonCompass.cc.");
65
67 zap[
'f'] =
make_field(inputFile,
"input files (output from JMuonCompass).");
69 zap[
'a'] =
make_field(detectorFile,
"detector file.");
70 zap[
'A'] =
make_field(overwriteDetector,
"overwrite detector file provided through '-a' with corrected orientations;"\
71 "\nUse option -AA to also correct compass calibrations.");
72 zap[
'O'] =
make_field(option,
"ROOT fit option, see TH1::Fit.") =
"";
74
75 zap(argc, argv);
76 }
77 catch(const exception &error) {
78 FATAL(error.what() << endl);
79 }
80
81
82 gErrorIgnoreLevel = kFatal;
83
84
86
87 try {
89 }
92 }
93
95
96 if (option.find('S') == string::npos) { option += 'S'; }
98
99 TH2D* h2 = NULL;
100
101 for (vector<string>::const_iterator i = inputFile.begin(); i != inputFile.end(); ++i) {
102
103 NOTICE(
"Processing " << *i << endl);
104
105 TFile in(i->c_str(), "exist");
106
107 if (!in.IsOpen()) {
108 FATAL(
"File " << *i <<
" not opened." << endl);
109 }
110
111 TH2D* p = dynamic_cast<TH2D*>(in.Get(h2_t));
112
113 if (p == NULL) {
114 FATAL(
"File " << *i <<
" has no histogram <" << h2_t <<
">." << endl);
115 }
116
117 if (h2 == NULL)
118 h2 = (TH2D*) p->Clone();
119 else
120 h2->Add(p);
121
122 h2->SetDirectory(0);
123
124 in.Close();
125 }
126
127 if (h2 == NULL) {
128 FATAL(
"No histogram <" << h2_t <<
">." << endl);
129 }
130
131
132
133 struct result_type {
134
135 result_type() :
136 value(0.0),
137 error(0.0)
138 {}
139
140 result_type(double value,
141 double error) :
142 value(value),
143 error(error)
144 {}
145
146 double value;
147 double error;
148 };
149
151
152 map_type zmap;
153
154
155
156
157 const TAxis* x_axis = h2->GetXaxis();
158 const TAxis* y_axis = h2->GetYaxis();
159
160 TH1D h0("h0", NULL, x_axis->GetNbins(), -0.5, x_axis->GetNbins() - 0.5);
161
164
165 TH2D H2d("detector", NULL,
166 strings.size(), -0.5, strings.size() - 0.5,
168
169 for (Int_t ix = 1; ix <= H2d.GetXaxis()->GetNbins(); ++ix) {
170 H2d.GetXaxis()->SetBinLabel(ix,
MAKE_CSTRING(strings.at(ix-1)));
171 }
172 for (Int_t iy = 1; iy <= H2d.GetYaxis()->GetNbins(); ++iy) {
174 }
175
176 TH2D* H2c = (TH2D*) H2d.Clone("chi2");
177 TH2D* H2q = (TH2D*) H2d.Clone("status");
178
179
181
182 for (Int_t ix = 1; ix <= x_axis->GetNbins(); ++ix) {
183
184 const int index = ix - 1;
186
188 continue;
189 }
190
191 DEBUG(
"Bin " << setw(4) << ix <<
" -> " << setw(8) << module.
getID() << endl);
192
193 TH1D h1(
MAKE_CSTRING(module.
getID() <<
".1D"), NULL, y_axis->GetNbins(), y_axis->GetXmin(), y_axis->GetXmax());
194
195
196
197 Double_t ymin = numeric_limits<double>::max();
198 Double_t ymax = numeric_limits<double>::lowest();
199 Double_t x0 = 0.0;
200
201 for (Int_t iy = 1; iy <= y_axis->GetNbins(); ++iy) {
202
203 h1.SetBinContent(iy, h2->GetBinContent(ix,iy));
204 h1.SetBinError (iy, h2->GetEntries() * 1.0e-7 * parameters.FACTOR);
205
206 const Double_t
x = h1.GetBinCenter (iy);
207 const Double_t
y = h1.GetBinContent(iy);
208
209 if (y < ymin) {
211 }
212
213 if (y > ymax) { ymax =
y; }
214 if (y < ymin) { ymin =
y; }
215 }
216
217 if (ymax - ymin >= parameters.YMIN) {
218
219 Double_t
xmin = x0 + parameters.XMIN;
220 Double_t
xmax = x0 + parameters.XMAX;
221
222 Double_t Rx = 0.0;
223
224 if (xmin < y_axis->GetXmin()) { Rx = +PI; }
225 if (xmax > y_axis->GetXmax()) { Rx = -PI; }
226
227
228
229 if (Rx != 0.0) {
230
231 NOTICE(
"Module " << setw(8) << module.
getID() <<
" offset " << Rx/PI <<
" [pi]" << endl);
232
233 x0 += Rx;
236
237 h1.Reset();
238
239 const Double_t dx = (y_axis->GetXmax() - y_axis->GetXmin()) / y_axis->GetNbins();
240
241 for (Int_t iy = 1; iy <= y_axis->GetNbins(); ++iy) {
242
243 Int_t i = iy + (Int_t) (Rx / dx);
244
245 if (i < 1) { i += (Int_t) (2.0*PI / dx); }
246 if (i > y_axis->GetNbins()) { i -= (Int_t) (2.0*PI / dx); }
247
248 h1.SetBinContent(i, h2->GetBinContent(ix,iy));
249 h1.SetBinError (i, h2->GetEntries() * 1.0e-7 * parameters.FACTOR);
250 }
251 }
252
253 if (xmin < y_axis->GetXmin()) {
xmin = y_axis->GetXmin();
xmax =
xmin + 0.8*PI; }
254 if (xmax > y_axis->GetXmax()) {
xmax = y_axis->GetXmax();
xmin =
xmax - 0.8*PI; }
255
256 TF1*
f1 =
new TF1(
"f1",
"[0] - [1]*cos(x - [2])");
257
258 f1->SetParameter(0, 0.5 * (ymin + ymax));
259 f1->SetParameter(1, 0.5 * (ymax - ymin));
260 f1->SetParameter(2, x0);
261
262 f1->SetParError(0, 0.01*(ymin+ymax));
263 f1->SetParError(1, 0.01*(ymin+ymax));
264 f1->SetParError(2, 0.05);
265
266 DEBUG(
"Start values:" << endl);
267
268 for (
int i = 0; i !=
f1->GetNpar(); ++i) {
269 DEBUG(left << setw(12) <<
f1->GetParName (i) <<
' ' <<
271 }
272
273
274
275 TFitResultPtr
result = h1.Fit(f1, option.c_str(),
"same", xmin, xmax);
276
277 const bool status =
result->IsValid();
278
279
280
281 if (Rx != 0.0) {
282 f1->SetParameter(2,
f1->GetParameter(2) - Rx);
283 }
284
286 cout << "Slice: "
287 << setw(4) << ix << ' '
288 << setw(16) << h1.GetName() << ' '
289 <<
FIXED(7,3) <<
f1->GetParameter(2) <<
" +/- "
290 <<
FIXED(7,3) <<
f1->GetParError(2) <<
" [rad] "
294 << (status ? "" : "failed") << endl;
295 }
296
297 if (status) {
298 zmap[index] = result_type(
f1->GetParameter(2),
f1->GetParError(2));
299 }
300
301 struct {
302 const int ix;
303 const int iy;
304 } H2 = {
305 strings.getIndex(module.
getString()) + 1,
306 module .getFloor()
307 };
308
309 if (status) {
310 H2d .SetBinContent(H2.ix, H2.iy,
f1->GetParameter(2));
311 }
313 H2c->SetBinContent(H2.ix, H2.iy,
result->Chi2() /
result->Ndf());
314 }
315 {
316 H2q->SetBinContent(H2.ix, H2.iy, status ? +1.0 : -1.0);
317 }
318
319 h1.Write();
320
322 }
323 }
324
325 if (zmap.empty()) {
326 NOTICE(
"No calibration results." << endl);
327 }
328
329 for (map_type::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
330 h0.SetBinContent(i->first, i->second.value);
331 h0.SetBinError (i->first, i->second.error);
332 }
333
334 if (overwriteDetector) {
335
336 NOTICE(
"Store calibration data on file " << detectorFile << endl);
337
339
340 for (map_type::const_iterator i = zmap.begin(); i != zmap.end(); ++i) {
341
343 const JVector3D center =
module.getPosition();
344
345 module.sub(center);
346
347 module.rotate(JRotation3Z(i->second.value));
348
349 module.add(center);
350
351 if (overwriteDetector > 1) {
352 module.setQuaternion(module.getQuaternion() * JQuaternion3Z(i->second.value));
353 }
354 }
355
357 }
358
359 h0 .Write();
360 h2 ->Write();
361 H2d .Write();
362 H2c->Write();
363 H2q->Write();
364
365 out.Close();
366}
#define DEBUG(A)
Message macros.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
#define MAKE_CSTRING(A)
Make C-string.
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
int getFloor() const
Get floor number.
int getString() const
Get string number.
Router for direct addressing of module data in detector data structure.
Data structure for a composite optical module.
Utility class to parse parameter values.
Data structure for vector in three dimensions.
int getID() const
Get identifier.
Auxiliary class to handle multiple boolean-like I/O.
Utility class to parse command line options.
const JPolynome f1(1.0, 2.0, 3.0)
Function.
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
void store(const std::string &file_name, const JDetector &detector)
Store detector to output file.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Auxiliary data structure for floating point format specification.
Router for mapping of string identifier to index.
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Auxiliary data structure for floating point format specification.