1 #ifndef __JGIZMOTOOLKIT__
2 #define __JGIZMOTOOLKIT__
17 #include "TGraphErrors.h"
19 #include "TMultiGraph.h"
25 #include "TIterator.h"
27 #include "TMethodCall.h"
42 namespace JPP {
using namespace JGIZMO; }
67 static const char*
const Add() {
return "Add"; }
68 static const char*
const add() {
return "add"; }
69 static const char*
const Subtract() {
return "Subtract"; }
70 static const char*
const subtract() {
return "subtract"; }
71 static const char*
const Multiply() {
return "Multiply"; }
72 static const char*
const multiply() {
return "multiply"; }
73 static const char*
const Divide() {
return "Divide"; }
74 static const char*
const divide() {
return "divide"; }
75 static const char*
const efficiency() {
return "efficiency"; }
76 static const char*
const stdev() {
return "stdev"; }
77 static const char*
const sqrt() {
return "sqrt"; }
78 static const char*
const Replace() {
return "Replace"; }
85 static const char*
const TIMESTAMP =
"#splitline{}{#splitline{%d:%m:%y}{ %H:%M}}%F1970-01-01 00:00:00";
98 inline TFile*
getFile(
const std::string& file_name,
const std::string& option =
"exist")
102 gErrorIgnoreLevel = kError;
108 if (i == zmap.end() || i->second == NULL || !i->second->IsOpen()) {
110 TFile*
file = TFile::Open(file_name.c_str(), option.c_str());
112 zmap[file_name] =
file;
135 if (in == NULL || !in->IsOpen()) {
140 return in->GetDirectory(
id.getDirectory());
158 const TRegexp regexp(
id.getObjectName());
160 TIter iter(dir->GetListOfKeys());
162 for (TKey* key; (key = (TKey*) iter.Next()) != NULL; ) {
164 const TString tag(key->GetName());
168 if (tag.Index(regexp) != -1) {
169 return key->ReadObj();
187 const TH1*
h1 =
dynamic_cast<const TH1*
>(object);
191 if (h1->GetSumw2N()) {
192 for (Int_t i = 1; i <= h1->GetNbinsX(); ++i) {
193 if (h1->GetBinError(i) != 0.0) {
201 const TGraphErrors*
g1 =
dynamic_cast<const TGraphErrors*
>(object);
205 for (Int_t i = 0; i != g1->GetN(); ++i) {
206 if (g1->GetEY()[i] != 0.0) {
211 return g1->GetN() > 1;
215 const TGraph*
g1 =
dynamic_cast<const TGraph*
>(object);
218 return g1->GetN() > 1;
242 TString buffer(text);
244 if (
object != NULL) {
246 TClass* p = TClass::GetClass(object->ClassName());
252 TIterator* iter = p->GetListOfAllPublicMethods()->MakeIterator();
253 TMethod* method = NULL;
255 for (TMethod* p; (p = (TMethod*) iter->Next()) != NULL; ) {
256 if (buffer.Index(p->GetName()) != -1) {
257 if (method == NULL || strlen(p->GetName()) > strlen(method->GetName())) {
263 if (method == NULL) {
267 for (Ssiz_t index; (index = buffer.Index(method->GetName())) != -1; ) {
269 const TRegexp fp(
" *([^)]*)");
272 Ssiz_t pos = buffer.Index(fp, &len, index);
276 if (pos == -1 || pos != index + (Ssiz_t) strlen(method->GetName())) {
278 TMethodCall(p, method->GetName(), NULL).Execute(
object, value);
280 len = strlen(method->GetName());
284 TMethodCall(p, method->GetName(), NULL).Execute(
object, TString(buffer(pos + 1, len - 2)), value);
286 len += strlen(method->GetName());
289 buffer.Replace(index, len, TString::Format(
"%15.5e", value));
295 return TFormula(
"/tmp", buffer.Data()).Eval(0.0);
315 return getResult(TString(text.c_str()),
object);
329 const char* regexp(
"p[0-9]* *=");
331 TString buffer(text.c_str());
333 buffer = buffer(TRegexp(regexp));
334 buffer = buffer(1, buffer.Length() - 2);
336 if (!buffer.IsDigit()) {
340 return buffer.Atoi();
360 const char* regexp(
"=.*");
362 TString buffer(text.c_str());
364 buffer = buffer(TRegexp(regexp));
365 buffer = buffer(1, buffer.Length() - 1);
367 return getResult(std::string(buffer),
object);
390 const char* regexp(
"=.*");
392 TString buffer(text.c_str());
394 buffer = buffer(TRegexp(regexp));
395 buffer = buffer(1, buffer.Length() - 1);
398 istringstream
is((std::string(buffer)));
402 for (
int i = index; is >> value && i > 0; --i) {}
407 THROW(
JParseError,
"Text des not contain a number at given position " << buffer <<
' ' << index);
420 const Int_t
first = axis->GetFirst();
421 const Int_t last = axis->GetLast();
423 const Double_t xmin = axis->GetBinLowEdge(first);
424 const Double_t xmax = axis->GetBinLowEdge(last) + axis->GetBinWidth(last);
426 const Int_t
N = axis->GetNbins();
427 Double_t buffer[N+1];
429 buffer[0] = TMath::Power(10.0, axis->GetBinLowEdge(1));
431 for (Int_t i = 1; i <=
N; ++i) {
432 buffer[i] = TMath::Power(10.0, axis->GetBinLowEdge(i) + axis->GetBinWidth(i));
435 axis->Set(N, buffer);
437 if (axis->TestBit(TAxis::kAxisRange)) {
438 axis->SetRangeUser(TMath::Power(10.0, xmin), TMath::Power(10.0, xmax));
453 const TRegexp regexp[] = {
454 TString(
"^") + TString(parameter) + TString(
"[^a-zA-Z_]"),
455 TString(
"[^a-zA-Z_]") + TString(parameter) + TString(
"[^a-zA-Z_]"),
456 TString(
"[^a-zA-Z_]") + TString(parameter) + TString(
"$")
459 const TString replacement = TString(
"log10(") + TString(parameter) + TString(
")");
461 TString buffer(formula);
463 if (buffer.Length() == 1 && buffer[0] == parameter) {
465 buffer = replacement;
469 for (Ssiz_t pos = 0, i; pos < buffer.Length(); pos += replacement.Length()) {
470 if ((i = buffer.Index(regexp[0], pos)) != -1) { buffer.Replace((pos = i + 0), 1, replacement); }
471 else if ((i = buffer.Index(regexp[1], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
472 else if ((i = buffer.Index(regexp[2], pos)) != -1) { buffer.Replace((pos = i + 1), 1, replacement); }
486 inline void copy(
const TF1& from, TF1& to)
488 static_cast<TAttLine&
> (to) = static_cast<const TAttLine&> (from);
489 static_cast<TAttFill&
> (to) = static_cast<const TAttFill&> (from);
490 static_cast<TAttMarker&
>(to) = static_cast<const TAttMarker&>(from);
492 to.SetParameters(from.GetParameters());
494 to.SetNpx(from.GetNpx());
504 inline void copy(
const TF2& from, TF2& to)
506 copy(static_cast<const TF1&>(from), static_cast<TF1&>(to));
508 to.SetNpy(from.GetNpy());
525 fn.SetRange(f1->GetXmin(),
530 f1->SetRange(TMath::Power(10.0,f1->GetXmin()),
531 TMath::Power(10.0,f1->GetXmax()));
549 fn.SetRange(f2->GetXmin(),
556 f2->SetRange(TMath::Power(10.0,f2->GetXmin()),
558 TMath::Power(10.0,f2->GetXmax()),
577 fn.SetRange(f2->GetXmin(),
584 f2->SetRange(f2->GetXmin(),
585 TMath::Power(10.0,f2->GetYmin()),
587 TMath::Power(10.0,f2->GetYmax()));
601 for (TIter iter(h1->GetListOfFunctions()); TF1* f1 =
dynamic_cast<TF1*
>(iter.Next()); ) {
621 for (TIter iter(h2->GetListOfFunctions()); TF2* f2 =
dynamic_cast<TF2*
>(iter.Next()); ) {
639 for (TIter iter(h2->GetListOfFunctions()); TF2* f2 =
dynamic_cast<TF2*
>(iter.Next()); ) {
657 for (TIter iter(g1->GetListOfFunctions()); TF1* f1 =
dynamic_cast<TF1*
>(iter.Next()); ) {
661 for (Int_t i = 0; i != g1->GetN(); ++i) {
662 g1->GetX()[i] =
pow(10.0, g1->GetX()[i]);
677 for (TIter iter(g2->GetListOfFunctions()); TF2* f2 =
dynamic_cast<TF2*
>(iter.Next()); ) {
681 for (Int_t i = 0; i != g2->GetN(); ++i) {
682 g2->GetX()[i] =
pow(10.0, g2->GetX()[i]);
697 for (TIter iter(g2->GetListOfFunctions()); TF2* f2 =
dynamic_cast<TF2*
>(iter.Next()); ) {
701 for (Int_t i = 0; i != g2->GetN(); ++i) {
702 g2->GetY()[i] =
pow(10.0, g2->GetY()[i]);
717 for (TIter i1(m1->GetListOfGraphs()); TGraph*
g1 =
dynamic_cast<TGraph*
>(i1()); ) {
736 inline void convertToPDF(TH1&
h1,
const std::string& option =
"NW",
const double factor = 1.0)
740 const bool normalise = (option.find(
'N') != string::npos || option.find(
'n') != string::npos);
741 const bool bin_width = (option.find(
'W') != string::npos || option.find(
'w') != string::npos);
742 const bool use_error = (option.find(
'E') != string::npos || option.find(
'e') != string::npos);
750 for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
751 W += h1.GetBinContent(i);
757 for (Int_t i = 1; i <= h1.GetXaxis()->GetNbins(); ++i) {
759 const Double_t
w = W * (bin_width ? h1.GetXaxis()->GetBinWidth(i) : 1.0);
761 h1.SetBinContent(i, h1.GetBinContent(i) * factor /
w);
764 h1.SetBinError(i, h1.GetBinError(i) * factor /
w);
785 inline void convertToPDF(TH2& h2,
const std::string& option =
"NXYW",
const double factor = 1.0)
789 const bool normalise = (option.find(
'N') != string::npos || option.find(
'n') != string::npos);
790 const bool X = (option.find(
'X') != string::npos || option.find(
'x') != string::npos);
791 const bool Y = (option.find(
'Y') != string::npos || option.find(
'y') != string::npos);
792 const bool bin_width = (option.find(
'W') != string::npos || option.find(
'w') != string::npos);
793 const bool use_error = (option.find(
'E') != string::npos || option.find(
'e') != string::npos);
803 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
804 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
805 W += h2.GetBinContent(ix,iy);
812 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
813 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
815 const Double_t
w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) * h2.GetYaxis()->GetBinWidth(iy) : 1.0);
817 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor /
w);
820 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor /
w);
828 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
834 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
835 W += h2.GetBinContent(ix,iy);
841 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
843 const Double_t
w = W * (bin_width ? h2.GetXaxis()->GetBinWidth(ix) : 1.0);
845 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) * factor /
w);
848 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) * factor /
w);
856 for (Int_t ix = 1; ix <= h2.GetXaxis()->GetNbins(); ++ix) {
862 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
863 W += h2.GetBinContent(ix,iy);
869 for (Int_t iy = 1; iy <= h2.GetYaxis()->GetNbins(); ++iy) {
871 const Double_t
w = W * (bin_width ? h2.GetYaxis()->GetBinWidth(iy) : 1.0);
873 h2.SetBinContent(ix, iy, h2.GetBinContent(ix,iy) /
w);
876 h2.SetBinError(ix, iy, h2.GetBinError(ix,iy) /
w);
894 Double_t ymin = numeric_limits<Double_t>::max();
895 Double_t ymax = numeric_limits<Double_t>::lowest();
897 for (Int_t i = 0; i != g1.GetN(); ++i) {
899 const Double_t y = g1.GetY()[i];
901 if (y > ymax) { ymax = y; }
902 if (y < ymin) { ymin = y; }
919 Double_t zmin = numeric_limits<Double_t>::max();
920 Double_t zmax = numeric_limits<Double_t>::lowest();
922 for (Int_t i = 0; i != g2.GetN(); ++i) {
924 const Double_t z = g2.GetZ()[i];
926 if (z > zmax) { zmax = z; }
927 if (z < zmin) { zmin = z; }
951 double dx = (xmax - xmin) * 0.1;
953 if (xmin > dx || xmin < 0.0)
985 if (axis->GetNbins() == (int) memo.size()) {
987 for (
int i = 0; i != axis->GetNbins(); ++i) {
989 const JPMTPhysicalAddress& address = memo[i];
991 axis->SetBinLabel(i + 1, address.toString().c_str());
1014 using namespace JPP;
1018 if (axis ==
"x" || axis ==
"X")
1019 pax = h1.GetXaxis();
1020 else if (axis ==
"y" || axis ==
"Y")
1021 pax = h1.GetYaxis();
1022 else if (axis ==
"z" || axis ==
"Z")
1023 pax = h1.GetZaxis();
1027 if (pax->GetNbins() == (int) memo.size()) {
1029 for (
int i = 0; i != pax->GetNbins(); ++i) {
1033 pax->SetBinLabel(i + 1, address.toString().c_str());
1036 h1.LabelsOption(
"a", axis.c_str());
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.
static const char *const SAME_AS_INPUT()
Set name of output histogram to name of input histogram.
void setLogarithmicX(TF1 *f1)
Make parameter x of function logarithmic (e.g. after filling with log10()).
#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.
then for HISTOGRAM in h0 h1
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 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"` typeset -Z 4 STRING JOpera1D -f hydrophone.root
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 filling with log10()).
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.
Double_t getResult(const TString &text, TObject *object=NULL)
Get result of given textual formula.
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 .
then break fi done getCenter read X Y Z let X
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 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.
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.
void setLogarithmicY(TF2 *f2)
Make parameter y of function logarithmic (e.g. after filling with log10()).
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 filling with 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 source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
static const char *const TIMESTAMP
Time stamp of earliest UTC time.
then usage $script< string identifier >< detectorfile > event file(toashort file)+" "\nNote that the event files and toashort files should be one-to-one related." fi if (( $
Double_t g1(const Double_t x)
Function.