1 #ifndef __JGIZMOTOOLKIT__
2 #define __JGIZMOTOOLKIT__
18 #include "TGraphErrors.h"
20 #include "TGraph2DErrors.h"
21 #include "TMultiGraph.h"
29 #include "TIterator.h"
31 #include "TMethodCall.h"
46 namespace JPP {
using namespace JGIZMO; }
72 static const char*
const Add() {
return "Add"; }
73 static const char*
const add() {
return "add"; }
74 static const char*
const Subtract() {
return "Subtract"; }
75 static const char*
const subtract() {
return "subtract"; }
76 static const char*
const Multiply() {
return "Multiply"; }
77 static const char*
const multiply() {
return "multiply"; }
78 static const char*
const Divide() {
return "Divide"; }
79 static const char*
const divide() {
return "divide"; }
80 static const char*
const efficiency() {
return "efficiency"; }
81 static const char*
const stdev() {
return "stdev"; }
82 static const char*
const sqrt() {
return "sqrt"; }
83 static const char*
const Replace() {
return "Replace"; }
90 static const char*
const TIMESTAMP =
"#splitline{}{#splitline{%d:%m:%y}{ %H:%M}}%F1970-01-01 00:00:00";
103 inline TFile*
getFile(
const std::string& file_name,
const std::string& option =
"exist")
107 gErrorIgnoreLevel = kError;
113 if (i == zmap.end() || i->second == NULL || !i->second->IsOpen()) {
115 TFile* file = TFile::Open(file_name.c_str(), option.c_str());
117 zmap[file_name] = file;
140 if (in == NULL || !in->IsOpen()) {
145 return in->GetDirectory(
id.getDirectory());
163 const TRegexp regexp(
id.getObjectName());
165 TIter iter(dir->GetListOfKeys());
167 for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
169 const TString tag(key->GetName());
173 if (tag.Index(regexp) != -1) {
174 return key->ReadObj();
192 const TH1* h1 =
dynamic_cast<const TH1*
>(object);
196 if (h1->GetSumw2N()) {
197 for (Int_t i = 1; i <= h1->GetNbinsX(); ++i) {
198 if (h1->GetBinError(i) != 0.0) {
206 const TGraphErrors*
g1 =
dynamic_cast<const TGraphErrors*
>(object);
210 for (Int_t i = 0; i != g1->GetN(); ++i) {
211 if (g1->GetEY()[i] != 0.0) {
216 return g1->GetN() > 1;
220 const TGraph*
g1 =
dynamic_cast<const TGraph*
>(object);
223 return g1->GetN() > 1;
247 TString buffer(text);
249 if (
object != NULL) {
251 TClass* p = TClass::GetClass(object->ClassName());
257 TIterator* iter = p->GetListOfAllPublicMethods()->MakeIterator();
258 TMethod* method = NULL;
260 for (TMethod* p; (p = (TMethod*) iter->Next()) != NULL; ) {
261 if (buffer.Index(p->GetName()) != -1) {
262 if (method == NULL || strlen(p->GetName()) > strlen(method->GetName())) {
268 if (method == NULL) {
272 for (Ssiz_t index; (index = buffer.Index(method->GetName())) != -1; ) {
274 const TRegexp fp(
" *([^)]*)");
277 Ssiz_t pos = buffer.Index(fp, &len, index);
281 if (pos == -1 || pos != index + (Ssiz_t) strlen(method->GetName())) {
283 TMethodCall(p, method->GetName(), NULL).Execute(
object, value);
285 len = strlen(method->GetName());
289 TMethodCall(p, method->GetName(), NULL).Execute(
object, TString(buffer(pos + 1, len - 2)), value);
291 len += strlen(method->GetName());
294 buffer.Replace(index, len, TString::Format(
"%20.10e", value));
300 return TFormula(
"/tmp", buffer.Data()).Eval(0.0);
320 return getResult(TString(text.c_str()),
object);
334 const char* regexp(
"p[0-9]* *=");
336 TString buffer(text.c_str());
338 buffer = buffer(TRegexp(regexp));
339 buffer = buffer(1, buffer.Length() - 2);
341 if (!buffer.IsDigit()) {
345 return buffer.Atoi();
365 const char* regexp(
"=.*");
367 TString buffer(text.c_str());
369 buffer = buffer(TRegexp(regexp));
370 buffer = buffer(1, buffer.Length() - 1);
372 return getResult(std::string(buffer),
object);
395 const char* regexp(
"=.*");
397 TString buffer(text.c_str());
399 buffer = buffer(TRegexp(regexp));
400 buffer = buffer(1, buffer.Length() - 1);
403 istringstream
is((std::string(buffer)));
407 for (
int i = index; is >> value && i > 0; --i) {}
412 THROW(
JParseError,
"Text des not contain a number at given position " << buffer <<
' ' << index);
425 const Int_t
first = axis->GetFirst();
426 const Int_t last = axis->GetLast();
428 const Double_t
xmin = axis->GetBinLowEdge(first);
429 const Double_t
xmax = axis->GetBinLowEdge(last) + axis->GetBinWidth(last);
431 const Int_t
N = axis->GetNbins();
432 Double_t buffer[N+1];
434 buffer[0] = TMath::Power(10.0, axis->GetBinLowEdge(1));
436 for (Int_t i = 1; i <=
N; ++i) {
437 buffer[i] = TMath::Power(10.0, axis->GetBinLowEdge(i) + axis->GetBinWidth(i));
440 axis->Set(N, buffer);
442 if (axis->TestBit(TAxis::kAxisRange)) {
443 axis->SetRangeUser(TMath::Power(10.0, xmin), TMath::Power(10.0, xmax));
458 const TRegexp regexp[] = {
459 TString(
"^") + TString(parameter) + TString(
"[^a-zA-Z_]"),
460 TString(
"[^a-zA-Z_]") + TString(parameter) + TString(
"[^a-zA-Z_]"),
461 TString(
"[^a-zA-Z_]") + TString(parameter) + TString(
"$")
464 const TString replacement = TString(
"log10(") + TString(parameter) + TString(
")");
466 TString buffer(formula);
468 if (buffer.Length() == 1 && buffer[0] == parameter) {
470 buffer = replacement;
474 for (Ssiz_t pos = 0, i; pos < buffer.Length(); pos += replacement.Length()) {
475 if ((i = buffer.Index(regexp[0], pos)) != -1) { buffer.Replace((pos = i + 0), 1, replacement); }
476 else if ((i = buffer.Index(regexp[1], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
477 else if ((i = buffer.Index(regexp[2], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
492 inline void copy(
const TF1& from, TF1& to)
494 static_cast<TAttLine&
> (to) = static_cast<const TAttLine&> (from);
495 static_cast<TAttFill&
> (to) = static_cast<const TAttFill&> (from);
496 static_cast<TAttMarker&
>(to) = static_cast<const TAttMarker&>(from);
498 to.SetParameters(from.GetParameters());
500 to.SetNpx(from.GetNpx());
510 inline void copy(
const TF2& from, TF2& to)
512 copy(static_cast<const TF1&>(from), static_cast<TF1&>(to));
514 to.SetNpy(from.GetNpy());
549 fn.SetRange(f1->GetXmin(),
554 f1->SetRange(TMath::Power(10.0,f1->GetXmin()),
555 TMath::Power(10.0,f1->GetXmax()));
569 TString buffer = f1->GetExpFormula();
571 buffer =
"pow(10.0, " + buffer +
")";
573 TF1 fn(f1->GetName(), buffer);
595 fn.SetRange(f2->GetXmin(),
602 f2->SetRange(TMath::Power(10.0,f2->GetXmin()),
604 TMath::Power(10.0,f2->GetXmax()),
623 fn.SetRange(f2->GetXmin(),
630 f2->SetRange(f2->GetXmin(),
631 TMath::Power(10.0,f2->GetYmin()),
633 TMath::Power(10.0,f2->GetYmax()));
647 setLogarithmicX<TF1>(h1->GetListOfFunctions());
665 setLogarithmicX<TF2>(h2->GetListOfFunctions());
681 setLogarithmicY<TF2>(h2->GetListOfFunctions());
697 setLogarithmicX<TF1>(g1->GetListOfFunctions());
699 for (Int_t i = 0; i != g1->GetN(); ++i) {
700 g1->GetX()[i] =
pow(10.0, g1->GetX()[i]);
715 setLogarithmicY<TF1>(g1->GetListOfFunctions());
717 for (Int_t i = 0; i != g1->GetN(); ++i) {
718 g1->GetY()[i] =
pow(10.0, g1->GetY()[i]);
733 setLogarithmicX<TF1>(g1->GetListOfFunctions());
735 for (Int_t i = 0; i != g1->GetN(); ++i) {
736 g1->GetEX()[i] =
pow(10.0, g1->GetX()[i] + g1->GetEX()[i]) -
pow(10.0, g1->GetX()[i]);
737 g1->GetX() [i] =
pow(10.0, g1->GetX()[i]);
752 setLogarithmicY<TF1>(g1->GetListOfFunctions());
754 for (Int_t i = 0; i != g1->GetN(); ++i) {
755 g1->GetEY()[i] =
pow(10.0, g1->GetY()[i] + g1->GetEY()[i]) -
pow(10.0, g1->GetY()[i]);
756 g1->GetY() [i] =
pow(10.0, g1->GetY()[i]);
771 setLogarithmicX<TF2>(g2->GetListOfFunctions());
773 for (Int_t i = 0; i != g2->GetN(); ++i) {
774 g2->GetX()[i] =
pow(10.0, g2->GetX()[i]);
789 setLogarithmicY<TF2>(g2->GetListOfFunctions());
791 for (Int_t i = 0; i != g2->GetN(); ++i) {
792 g2->GetY()[i] =
pow(10.0, g2->GetY()[i]);
807 setLogarithmicX<TF2>(g2->GetListOfFunctions());
809 for (Int_t i = 0; i != g2->GetN(); ++i) {
810 g2->GetEX()[i] =
pow(10.0, g2->GetX()[i] + g2->GetEX()[i]) -
pow(10.0, g2->GetX()[i]);
811 g2->GetX() [i] =
pow(10.0, g2->GetX()[i]);
826 setLogarithmicY<TF2>(g2->GetListOfFunctions());
828 for (Int_t i = 0; i != g2->GetN(); ++i) {
829 g2->GetEY()[i] =
pow(10.0, g2->GetY()[i] + g2->GetEY()[i]) -
pow(10.0, g2->GetY()[i]);
830 g2->GetY() [i] =
pow(10.0, g2->GetY()[i]);
844 setLogarithmicX<TGraph>(gn->GetListOfGraphs());
857 setLogarithmicY<TGraph>(gn->GetListOfGraphs());
870 line->SetX1(
pow(10.0, line->GetX1()));
871 line->SetX2(
pow(10.0, line->GetX2()));
884 line->SetY1(
pow(10.0, line->GetY1()));
885 line->SetY2(
pow(10.0, line->GetY2()));
920 for (TIter i(list);
T* p =
dynamic_cast<T*
>(i.Next()); ) {
934 for (TIter i(list);
T* p =
dynamic_cast<T*
>(i.Next()); ) {
952 inline void convertToPDF(TH1& h1,
const std::string& option =
"NW",
const double factor = 1.0)
956 const bool normalise = (option.find(
'N') != string::npos || option.find(
'n') != string::npos);
957 const bool bin_width = (option.find(
'W') != string::npos || option.find(
'w') != string::npos);
958 const bool use_error = (option.find(
'E') != string::npos || option.find(
'e') != string::npos);
966 for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
967 W += h1.GetBinContent(i);
973 for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
975 const Double_t
w = W * (bin_width ? h1.GetXaxis()->GetBinWidth(i) : 1.0);
977 h1.SetBinContent(i, h1.GetBinContent(i) * factor /
w);
980 h1.SetBinError(i, h1.GetBinError(i) * factor /
w);
1001 inline void convertToPDF(TH2& h2,
const std::string& option =
"NXYW",
const double factor = 1.0)
1003 using namespace std;
1005 const bool normalise = (option.find(
'N') != string::npos || option.find(
'n') != string::npos);
1006 const bool X = (option.find(
'X') != string::npos || option.find(
'x') != string::npos);
1007 const bool Y = (option.find(
'Y') != string::npos || option.find(
'y') != string::npos);
1008 const bool bin_width = (option.find(
'W') != string::npos || option.find(
'w') != string::npos);
1009 const bool use_error = (option.find(
'E') != string::npos || option.find(
'e') != string::npos);
1019 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1020 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1021 W += h2.GetBinContent(ix,iy);
1028 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1029 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1031 const Double_t
w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) * h2.GetYaxis()->GetBinWidth(iy) : 1.0);
1033 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor /
w);
1036 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor /
w);
1044 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1050 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1051 W += h2.GetBinContent(ix,iy);
1057 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1059 const Double_t
w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) : 1.0);
1061 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor /
w);
1064 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor /
w);
1072 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
1078 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1079 W += h2.GetBinContent(ix,iy);
1085 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
1087 const Double_t
w = W * (bin_width ? h2.GetYaxis()->GetBinWidth(iy) : 1.0);
1089 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) /
w);
1092 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) /
w);
1116 inline void convertToPDF(TH3& h3,
const std::string& option =
"NXYW",
const double factor = 1.0)
1118 using namespace std;
1120 const bool normalise = (option.find(
'N') != string::npos || option.find(
'n') != string::npos);
1121 const bool X = (option.find(
'X') != string::npos || option.find(
'x') != string::npos);
1122 const bool Y = (option.find(
'Y') != string::npos || option.find(
'y') != string::npos);
1123 const bool Z = (option.find(
'Z') != string::npos || option.find(
'z') != string::npos);
1124 const bool bin_width = (option.find(
'W') != string::npos || option.find(
'w') != string::npos);
1125 const bool use_error = (option.find(
'E') != string::npos || option.find(
'e') != string::npos);
1135 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1136 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1137 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1138 W += h3.GetBinContent(ix,iy,iz);
1146 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1147 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1148 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1150 const Double_t
w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetYaxis()->GetBinWidth(iy) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1152 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor /
w);
1155 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor /
w);
1162 }
else if (
X &&
Z) {
1164 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1170 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1171 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1172 W += h3.GetBinContent(ix,iy,iz);
1179 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1180 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1182 const Double_t
w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1184 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor /
w);
1187 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor /
w);
1194 }
else if (
Y &&
Z) {
1196 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1202 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1203 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1204 W += h3.GetBinContent(ix,iy,iz);
1211 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1212 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1214 const Double_t
w = W * (bin_width ? h3.GetYaxis()->GetBinWidth(iy) * h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1216 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor /
w);
1219 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor /
w);
1226 }
else if (
X &&
Y) {
1228 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1234 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1235 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1236 W += h3.GetBinContent(ix,iy,iz);
1243 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1244 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1246 const Double_t
w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) * h3.GetYaxis()->GetBinWidth(iy) : 1.0);
1248 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor /
w);
1251 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor /
w);
1260 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1261 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1267 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1268 W += h3.GetBinContent(ix,iy,iz);
1274 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1276 const Double_t
w = W * (bin_width ? h3.GetXaxis()->GetBinWidth(ix) : 1.0);
1278 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) * factor /
w);
1281 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) * factor /
w);
1290 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1291 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1297 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1298 W += h3.GetBinContent(ix,iy,iz);
1304 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1306 const Double_t
w = W * (bin_width ? h3.GetYaxis()->GetBinWidth(iy) : 1.0);
1308 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) /
w);
1311 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) /
w);
1320 for (Int_t ix = 1; ix <= h3.GetXaxis()->GetNbins(); ++ix) {
1321 for (Int_t iy = 1; iy <= h3.GetYaxis()->GetNbins(); ++iy) {
1327 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1328 W += h3.GetBinContent(ix,iy,iz);
1334 for (Int_t iz = 1; iz <= h3.GetZaxis()->GetNbins(); ++iz) {
1336 const Double_t
w = W * (bin_width ? h3.GetZaxis()->GetBinWidth(iz) : 1.0);
1338 h3.SetBinContent(ix, iy, iz, h3.GetBinContent(ix,iy,iz) /
w);
1341 h3.SetBinError(ix, iy, iz, h3.GetBinError(ix,iy,iz) /
w);
1358 using namespace std;
1360 Double_t ymin = numeric_limits<Double_t>::max();
1361 Double_t ymax = numeric_limits<Double_t>::lowest();
1363 for (Int_t i = 0; i != g1.GetN(); ++i) {
1365 const Double_t y = g1.GetY()[i];
1367 if (y > ymax) { ymax = y; }
1368 if (y < ymin) { ymin = y; }
1371 g1.SetMinimum(ymin);
1372 g1.SetMaximum(ymax);
1383 using namespace std;
1385 Double_t zmin = numeric_limits<Double_t>::max();
1386 Double_t zmax = numeric_limits<Double_t>::lowest();
1388 for (Int_t i = 0; i != g2.GetN(); ++i) {
1390 const Double_t z = g2.GetZ()[i];
1392 if (z > zmax) { zmax = z; }
1393 if (z < zmin) { zmin = z; }
1396 g2.SetMinimum(zmin);
1397 g2.SetMaximum(zmax);
1417 double dx = (xmax -
xmin) * 0.1;
1419 if (xmin > dx || xmin < 0.0)
1449 using namespace JPP;
1451 if (axis->GetNbins() == (int) memo.size()) {
1453 for (
int i = 0; i != axis->GetNbins(); ++i) {
1455 const JPMTPhysicalAddress& address = memo[i];
1457 axis->SetBinLabel(i + 1, address.toString().c_str());
1480 using namespace JPP;
1484 if (axis ==
"x" || axis ==
"X")
1485 pax = h1.GetXaxis();
1486 else if (axis ==
"y" || axis ==
"Y")
1487 pax = h1.GetYaxis();
1488 else if (axis ==
"z" || axis ==
"Z")
1489 pax = h1.GetZaxis();
1493 if (pax->GetNbins() == (int) memo.size()) {
1495 for (
int i = 0; i != pax->GetNbins(); ++i) {
1499 pax->SetBinLabel(i + 1, address.toString().c_str());
1502 h1.LabelsOption(
"a", axis.c_str());
1519 return (key != NULL && TClass::GetClass(key->GetClassName())->IsTObject());
static const char *const sqrt()
Set contents to signed difference between squares.
TObject * getObject(const JRootObjectID &id)
Get first TObject with given identifier.
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
int getParameter(const std::string &text)
Get parameter number from text string.
static const char *const Divide()
ROOT TH1::Divide.
static const char *const Subtract()
ROOT TH1::Subtract.
Exception for a functional operation.
static const char *const SAME_AS_INPUT()
Set name of output histogram to name of input histogram.
#define THROW(JException_t, A)
Marco for throwing exception with std::ostream compatible message.
static const char *const Replace()
Set contents to associated function.
then JShowerPostfit f $INPUT_FILE o $OUTPUT_FILE N
void setLimits(TGraph &g1)
Set limits of TGraph.
static const char *const subtract()
Subtract contents with lookup bin in second histogram.
Auxiliary class to handle file name, ROOT directory and object name.
then fatal Wrong number of arguments fi set_variable STRING $argv[1] set_variable DETECTORXY_TXT $WORKDIR $DETECTORXY_TXT tail read X Y CHI2 RMS printf optimum n $X $Y $CHI2 $RMS awk v Y
static const char *const stdev()
Set contents to standard deviation.
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
void setLogarithmic(TAxis *axis)
Make histogram axis logarithmic (e.g. after using log10()).
const JPolynome f1(1.0, 2.0, 3.0)
Function.
bool isTAttLine(const TObject *object)
Get drawing option of TH1.
static const char *const efficiency()
Divide contents and multiply errors with inefficiency.
static const char *const Add()
ROOT TH1::Add.
Lookup table for PMT addresses in optical module.
static const char *const multiply()
Multiply contents with lookup bin in second histogram.
TFile * getFile(const std::string &file_name, const std::string &option="exist")
Get TFile pointer corresponding to give file name.
static const char *const add()
Add contents with lookup bin in second histogram.
Auxiliary data structure for JOpera1D.cc and JOpera2D.cc applications.
void setLogarithmicX(TList *list)
Make x-axis of objects in list logarithmic (e.g. after using log10()).
Double_t getResult(const TString &text, TObject *object=NULL)
Get result of given textual formula.
do set_variable OUTPUT_DIRECTORY $WORKDIR T
void convertToPDF(TH1 &h1, const std::string &option="NW", const double factor=1.0)
Convert 1D histogram to PDF.
void setAxisLabels(TAxis *axis, const JModuleAddressMap &memo)
Set axis with PMT address labels.
T pow(const T &x, const double y)
Power .
static const char *const Multiply()
ROOT TH1::Multiply.
static const char *const SAME_AS_OPERATION()
Set name of output histogram to name of operation.
void setLogarithmicY(TList *list)
Make y-axis of objects in list logarithmic (e.g. after using log10()).
then cat $TRIPOD_INITIAL<< EOF1 256877.5 4743716.7-2438.42 256815.5 4743395.0-2435.53 257096.2 4743636.0-2439.5EOFfiif[[!-f $DETECTOR]];then JEditDetector-a $DETECTOR_INITIAL-s"-1 addz -6.9"-o $DETECTOR--!eval`JPrintDetector-a $DETECTOR-O SUMMARY`for STRING in ${STRINGS[*]};do set_variable MODULE`getModule-a $DETECTOR-L"$STRING 0"`JEditDetector-a $DETECTOR-M"$MODULE setz -2.9"-o $DETECTOR--!donefiif[[!-f $TRIPOD]];then cp-p $TRIPOD_INITIAL $TRIPODfiJAcoustics.sh $DETECTOR_IDcat > acoustics_trigger_parameters txt<< EOFQ=0.0;TMax_s=0.020;quantile=0.9;numberOfHits=90;EOFJAcousticsEventBuilder.sh $DETECTOR $RUNS[*]INPUT_FILES=(`ls KM3NeT_ ${(l:8::0::0:) DETECTOR_ID}_0 *${^RUNS}_event.root`) cd $WORKDIRif[!$HOMEDIR-ef $WORKDIR];then cp-p $HOMEDIR/$DETECTOR $WORKDIR cp-p $HOMEDIR/$TRIPOD $WORKDIR cp-p $HOMEDIR/${^INPUT_FILES}$WORKDIR cp-p $HOMEDIR/{acoustics_fit_parameters, acoustics_trigger_parameters, disable, hydrophone, mechanics, sound_velocity, tripod, waveform}.txt $WORKDIRfisource $JPP_DIR/examples/JAcoustics/acoustics-fit-toolkit.shJConvertDetectorFormat-a $DETECTOR-o $HOMEDIR/detector_backup.datxJDetachPMTs-a $DETECTOR-o $DETECTORtimer_startinitialise stage_1B 0.002 0.1 0 > &stage log
void setRange(double &xmin, double &xmax, const bool logx)
Set axis range.
const JPMTAddressTranslator & getAddressTranslator(const int tdc) const
Get PMT address translator.
static const char *const divide()
Divide contents with lookup bin in second histogram.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Exception for parsing value.
bool isTObject(const TKey *key)
Check if given key corresponds to a TObject.
no fit printf nominal n $STRING awk v X
TDirectory * getDirectory(const JRootObjectID &id)
Get TDirectory pointer.
std::string getFilename(const std::string &file_name)
Get file name part, i.e. part after last JEEP::PATHNAME_SEPARATOR if any.
do set_variable MODULE getModule a $WORKDIR detector_a datx L $STRING JEditDetector a $WORKDIR detector_a datx M $MODULE setz o $WORKDIR detector_a datx JEditDetector a $WORKDIR detector_b datx M $MODULE setz o $WORKDIR detector_b datx done echo Output stored at $WORKDIR detector_a datx and $WORKDIR tripod_a txt JDrawDetector2D a $WORKDIR detector_a datx a $WORKDIR detector_b datx L BL o detector $FORMAT $BATCH JDrawDetector2D T $WORKDIR tripod_a txt T $WORKDIR tripod_b txt L BL o tripod $FORMAT $BATCH JCompareDetector a $WORKDIR detector_a datx b $WORKDIR detector_b datx o $WORKDIR abc root &dev null for KEY in X Y Z
Exception for accessing a value in a collection that is outside of its range.
TString getLogarithmic(const TString &formula, const char parameter)
Make given parameter in formula logarithmic (e.g. after using log10()).
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
static const char *const TIMESTAMP
Time stamp of earliest UTC time.
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
Double_t g1(const Double_t x)
Function.