![]() |
AliRoot Core
3dc7879 (3dc7879)
|
Class for viewing/visualizing TPC calibration data. More...
#include <AliTPCCalibViewer.h>
Public Member Functions | |
AliTPCCalibViewer () | |
AliTPCCalibViewer (const AliTPCCalibViewer &c) | |
AliTPCCalibViewer (TTree *const tree) | |
AliTPCCalibViewer (const char *fileName, const char *treeName="calPads") | |
AliTPCCalibViewer & | operator= (const AliTPCCalibViewer ¶m) |
virtual | ~AliTPCCalibViewer () |
virtual void | Delete (Option_t *option="") |
TString & | GetAbbreviation () |
TString & | GetAppendString () |
void | SetAbbreviation (const Char_t *abr) |
void | SetAppendString (const Char_t *str) |
virtual void | Draw (Option_t *opt="") |
virtual Long64_t | Draw (const char *varexp, const TCut &selection, Option_t *option="", Long64_t nentries=1000000000, Long64_t firstentry=0) |
virtual Long64_t | Draw (const char *varexp, const char *selection, Option_t *option="", Long64_t nentries=1000000000, Long64_t firstentry=0) |
TString | AddAbbreviations (const Char_t *c, Bool_t printDrawCommand=kFALSE) |
Int_t | EasyDraw (const char *drawCommand, const char *sector, const char *cuts=0, const char *drawOptions=0, Bool_t writeDrawCommand=kFALSE) const |
Int_t | EasyDraw (const char *drawCommand, Int_t sector, const char *cuts=0, const char *drawOptions=0, Bool_t writeDrawCommand=kFALSE) const |
Int_t | EasyDraw1D (const char *drawCommand, const char *sector, const char *cuts=0, const char *drawOptions=0, Bool_t writeDrawCommand=kFALSE) const |
Int_t | EasyDraw1D (const char *drawCommand, Int_t sector, const char *cuts=0, const char *drawOptions=0, Bool_t writeDrawCommand=kFALSE) const |
void | FormatHistoLabels (TH1 *histo) const |
Int_t | DrawHisto1D (const char *drawCommand, Int_t sector, const char *cuts=0, const char *sigmas="2;4;6", Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE) const |
Int_t | DrawHisto1D (const char *drawCommand, const char *sector, const char *cuts=0, const char *sigmas="2;4;6", Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE) const |
Int_t | SigmaCut (const char *drawCommand, Int_t sector, const char *cuts=0, Float_t sigmaMax=5, Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE, Bool_t pm=kFALSE, const char *sigmas="", Float_t sigmaStep=-1) const |
Int_t | SigmaCutNew (const char *drawCommand, const char *sector, const char *cuts=0, Float_t sigmaMax=5, Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE, Bool_t pm=kFALSE, const char *sigmas="", Float_t sigmaStep=-1) const |
Int_t | SigmaCut (const char *drawCommand, const char *sector, const char *cuts=0, Float_t sigmaMax=5, Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE, Bool_t pm=kFALSE, const char *sigmas="", Float_t sigmaStep=-1) const |
Int_t | Integrate (const char *drawCommand, const char *sector, const char *cuts=0, Float_t sigmaMax=5, Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE, const char *sigmas="", Float_t sigmaStep=-1) const |
Int_t | Integrate (const char *drawCommand, Int_t sector, const char *cuts=0, Float_t sigmaMax=5, Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE, const char *sigmas="", Float_t sigmaStep=-1) const |
Int_t | IntegrateOld (const char *drawCommand, const char *sector, const char *cuts=0, Float_t sigmaMax=5, Bool_t plotMean=kTRUE, Bool_t plotMedian=kTRUE, Bool_t plotLTM=kTRUE, const char *sigmas="", Float_t sigmaStep=-1) const |
AliTPCCalPad * | GetCalPadOld (const char *desiredData, const char *cuts="", const char *calPadName="NoName") const |
AliTPCCalPad * | GetCalPad (const char *desiredData, const char *cuts="", const char *calPadName="NoName") const |
AliTPCCalROC * | GetCalROC (const char *desiredData, UInt_t sector, const char *cuts="") const |
TObjArray * | GetArrayOfCalPads () |
TObjArray * | GetListOfVariables (Bool_t printList=kFALSE) |
TObjArray * | GetListOfNormalizationVariables (Bool_t printList=kFALSE) const |
TFriendElement * | AddReferenceTree (const char *filename, const char *treename="calPads", const char *refname="R") |
TFriendElement * | AddFriend (const char *treename, const char *filename) |
TFriendElement * | AddFriend (TTree *tree, const char *alias, Bool_t warn=kFALSE) |
TFriendElement * | AddFriend (const char *treename, TFile *file) |
TTree * | GetTree () const |
TString * | Fit (const char *drawCommand, const char *formula, const char *cuts, Double_t &chi2, TVectorD &fitParam, TMatrixD &covMatrix) |
Static Public Member Functions | |
static void | MakeTreeWithObjects (const char *fileName, const TObjArray *const array, const char *mapFileName=0) |
static void | MakeTree (const char *fileName, TObjArray *array, const char *mapFileName=0, AliTPCCalPad *const outlierPad=0, Float_t ltmFraction=0.9) |
static void | MakeTree (const char *outPutFileName, const Char_t *inputFileName, AliTPCCalPad *outlierPad=0, Float_t ltmFraction=0.9, const char *mapFileName="$ALICE_ROOT/TPC/Calib/MapCalibrationObjects.root") |
static void | CreateObjectList (const Char_t *filename, TObjArray *calibObjects) |
static void | MakeCalPadAliases (TTree *tree) |
static Double_t | GetLTM (Int_t n, const Double_t *const array, Double_t *const sigma=0, Double_t fraction=0.9) |
static Int_t | GetBin (Float_t value, Int_t nbins, Double_t binLow, Double_t binUp) |
static TH1F * | SigmaCut (Int_t n, const Float_t *array, Float_t mean, Float_t sigma, Int_t nbins, Float_t binLow, Float_t binUp, Float_t sigmaMax, Float_t sigmaStep=-1, Bool_t pm=kFALSE) |
static TH1F * | SigmaCut (TH1F *const histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep=-1, Bool_t pm=kFALSE) |
static TH1F * | Integrate (TH1F *const histogram, Float_t mean=0, Float_t sigma=0, Float_t sigmaMax=0, Float_t sigmaStep=-1) |
static TH1F * | Integrate (Int_t n, const Float_t *const array, Int_t nbins, Float_t binLow, Float_t binUp, Float_t mean=0, Float_t sigma=0, Float_t sigmaMax=0, Float_t sigmaStep=-1) |
static TH1F * | SigmaCut (Int_t n, const Double_t *array, Double_t mean, Double_t sigma, Int_t nbins, const Double_t *xbins, Double_t sigmaMax) |
Protected Member Functions | |
void | DrawLines (TH1F *cutHistoMean, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const |
void | DrawLines (TGraph *graph, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const |
Protected Attributes | |
TTree * | fTree |
tree containing visualization data (e.g. written by AliTPCCalPad::MakeTree(...) More... | |
TFile * | fFile |
file that contains a calPads tree (e.g. written by AliTPCCalPad::MakeTree(...) More... | |
TObjArray * | fListOfObjectsToBeDeleted |
Objects, that will be deleted when the destructor ist called. More... | |
Bool_t | fTreeMustBeDeleted |
decides weather the tree must be deleted in destructor or not More... | |
TString | fAbbreviation |
the abreviation for '.fElements' More... | |
TString | fAppendString |
'.fElements', stored in a TStrig More... | |
Class for viewing/visualizing TPC calibration data.
TPC calibration viewer/visualization class use Tree for visualization.
base on TTree functionality for visualization
Create a list of AliTPCCalPads, arrange them in an TObjArray. Pass this TObjArray to MakeTree and create the calibration Tree While craating this tree some statistical information are calculated Open the viewer with this Tree: AliTPCCalibViewer v("CalibTree.root") Have fun! EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
If you like to click, we recommand you the AliTPCCalibViewerGUI
THE DOCUMENTATION IS STILL NOT COMPLETED !!!!
ROOT includes
Definition at line 28 of file AliTPCCalibViewer.h.
AliTPCCalibViewer::AliTPCCalibViewer | ( | ) |
Definition at line 76 of file AliTPCCalibViewer.cxx.
AliTPCCalibViewer::AliTPCCalibViewer | ( | const AliTPCCalibViewer & | c | ) |
dummy AliTPCCalibViewer copy constructor not yet working!!!
Definition at line 92 of file AliTPCCalibViewer.cxx.
AliTPCCalibViewer::AliTPCCalibViewer | ( | TTree *const | tree | ) |
Constructor that initializes the calibration viewer
Definition at line 113 of file AliTPCCalibViewer.cxx.
AliTPCCalibViewer::AliTPCCalibViewer | ( | const char * | fileName, |
const char * | treeName = "calPads" |
||
) |
Constructor to initialize the calibration viewer the file 'fileName' contains the tree 'treeName'
Definition at line 132 of file AliTPCCalibViewer.cxx.
|
virtual |
AliTPCCalibViewer destructor all objects will be deleted, the file will be closed, the pictures will disappear
Definition at line 171 of file AliTPCCalibViewer.cxx.
TString AliTPCCalibViewer::AddAbbreviations | ( | const Char_t * | c, |
Bool_t | printDrawCommand = kFALSE |
||
) |
Replace all "<variable>" with "<variable><fAbbreviation>" (Adds forgotten "~") but take care on the statistical information, like "CEQmean_Mean" and also take care on correct given variables, like "CEQmean~"
For each variable out of "listOfVariables":
Now all the missing "~" should be added.
Definition at line 208 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoFit(), Draw(), AliTPCCalibViewerGUI::GetCutString(), and AliTPCCalibViewerGUI::GetDrawString().
|
inline |
Definition at line 78 of file AliTPCCalibViewer.h.
Referenced by AddFiles(), and AddReferenceTree().
|
inline |
Definition at line 79 of file AliTPCCalibViewer.h.
|
inline |
Definition at line 80 of file AliTPCCalibViewer.h.
TFriendElement * AliTPCCalibViewer::AddReferenceTree | ( | const char * | filename, |
const char * | treename = "calPads" , |
||
const char * | refname = "R" |
||
) |
add a reference tree to the current tree by default the treename is 'calPads' and the reference treename is 'R'
Definition at line 1633 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoLoadTree(), Draw(), and AliTPCCalibViewerGUItime::SetGuiTree().
|
static |
Function to create a TObjArray out of a given file the file must be in the following format: (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
type path dependingOnType
type == "CE": dependingOnType = CETmean CEQmean CETrms
type == "Pulser": dependingOnType = PulserTmean PulsterQmean PulserTrms
type == "Pedestals": dependingOnType = Pedestals Noise
type == "CalPad": dependingOnType = NameInFile NameToWriteToFile
Definition at line 2178 of file AliTPCCalibViewer.cxx.
Referenced by Draw(), and MakeTree().
|
virtual |
Should be called from AliTPCCalibViewerGUI class only. If you use Delete() do not call the destructor. All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
Definition at line 194 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::~AliTPCCalibViewerGUI().
|
inlinevirtual |
Definition at line 43 of file AliTPCCalibViewer.h.
Referenced by Fit().
|
inlinevirtual |
Definition at line 44 of file AliTPCCalibViewer.h.
|
inlinevirtual |
Definition at line 45 of file AliTPCCalibViewer.h.
Int_t AliTPCCalibViewer::DrawHisto1D | ( | const char * | drawCommand, |
Int_t | sector, | ||
const char * | cuts = 0 , |
||
const char * | sigmas = "2;4;6" , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE |
||
) | const |
Easy drawing of data, in principle the same as EasyDraw1D Difference: A line for the mean / median / LTM is drawn in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';' example: sigmas = "2; 4; 6;" at \(2 \sigma\), \(4 \sigma\) and \(6 \sigma\) a line is drawn. "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
Definition at line 562 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoDraw(), and Draw().
Int_t AliTPCCalibViewer::DrawHisto1D | ( | const char * | drawCommand, |
const char * | sector, | ||
const char * | cuts = 0 , |
||
const char * | sigmas = "2;4;6" , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE |
||
) | const |
Easy drawing of data, in principle the same as EasyDraw1D Difference: A line for the mean / median / LTM is drawn in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';' example: sigmas = "2; 4; 6;" at \(2 \sigma\), \(4 \sigma\) and \(6 \sigma\) a line is drawn. "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
Definition at line 577 of file AliTPCCalibViewer.cxx.
|
protected |
Private function for SigmaCut(...) and Integrate(...) Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
Definition at line 1118 of file AliTPCCalibViewer.cxx.
Referenced by Integrate(), IntegrateOld(), SigmaCut(), and SigmaCutNew().
|
protected |
Private function for SigmaCut(...) and Integrate(...) Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
Definition at line 1168 of file AliTPCCalibViewer.cxx.
Int_t AliTPCCalibViewer::EasyDraw | ( | const char * | drawCommand, |
const char * | sector, | ||
const char * | cuts = 0 , |
||
const char * | drawOptions = 0 , |
||
Bool_t | writeDrawCommand = kFALSE |
||
) | const |
easy drawing of data, use '~' for abbreviation of '.fElements' example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0") sector: sector-number - only the specified sector will be drwawn 'A'/'C' or 'a'/'c' - side A/C will be drawn 'ALL' - whole TPC will be drawn, projected on one side cuts: specifies cuts drawOptions: draw options like 'same' writeDrawCommand: write the command, that is passed to TTree::Draw
Definition at line 312 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoDraw(), Draw(), and EasyDraw().
Int_t AliTPCCalibViewer::EasyDraw | ( | const char * | drawCommand, |
Int_t | sector, | ||
const char * | cuts = 0 , |
||
const char * | drawOptions = 0 , |
||
Bool_t | writeDrawCommand = kFALSE |
||
) | const |
easy drawing of data, use '~' for abbreviation of '.fElements' example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0") sector: sector-number - only the specified sector will be drwawn cuts: specifies cuts drawOptions: draw options like 'same' writeDrawCommand: write the command, that is passed to TTree::Draw
Definition at line 438 of file AliTPCCalibViewer.cxx.
Int_t AliTPCCalibViewer::EasyDraw1D | ( | const char * | drawCommand, |
const char * | sector, | ||
const char * | cuts = 0 , |
||
const char * | drawOptions = 0 , |
||
Bool_t | writeDrawCommand = kFALSE |
||
) | const |
easy drawing of data, use '~' for abbreviation of '.fElements' example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0") sector: sector-number - the specified sector will be drwawn 'A'/'C' or 'a'/'c' - side A/C will be drawn 'ALL' - whole TPC will be drawn, projected on one side cuts: specifies cuts drawOptions: draw options like 'same' writeDrawCommand: write the command, that is passed to TTree::Draw
Definition at line 455 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoDraw(), Draw(), DrawHisto1D(), EasyDraw1D(), GetCalPadOld(), GetCalROC(), Integrate(), IntegrateOld(), SigmaCut(), and SigmaCutNew().
Int_t AliTPCCalibViewer::EasyDraw1D | ( | const char * | drawCommand, |
Int_t | sector, | ||
const char * | cuts = 0 , |
||
const char * | drawOptions = 0 , |
||
Bool_t | writeDrawCommand = kFALSE |
||
) | const |
easy drawing of data, use '~' for abbreviation of '.fElements' example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0") sector: sector-number - the specified sector will be drwawn cuts: specifies cuts drawOptions: draw options like 'same' writeDrawCommand: write the command, that is passed to TTree::Draw
Definition at line 515 of file AliTPCCalibViewer.cxx.
TString * AliTPCCalibViewer::Fit | ( | const char * | drawCommand, |
const char * | formula, | ||
const char * | cuts, | ||
Double_t & | chi2, | ||
TVectorD & | fitParam, | ||
TMatrixD & | covMatrix | ||
) |
fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts returns chi2, fitParam and covMatrix returns TString with fitted formula
Definition at line 1666 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoFit(), AliTPCCalibViewerGUI::GetDrawString(), and GetTree().
void AliTPCCalibViewer::FormatHistoLabels | ( | TH1 * | histo | ) | const |
formats title and axis labels of histo removes '.fElements'
Definition at line 531 of file AliTPCCalibViewer.cxx.
Referenced by Draw(), EasyDraw(), and EasyDraw1D().
|
inline |
Definition at line 38 of file AliTPCCalibViewer.h.
Referenced by AliTPCCalibViewerGUI::GetDrawString().
|
inline |
Definition at line 39 of file AliTPCCalibViewer.h.
TObjArray * AliTPCCalibViewer::GetArrayOfCalPads | ( | ) |
Returns a TObjArray with all AliTPCCalPads that are stored in the tree
Definition at line 1644 of file AliTPCCalibViewer.cxx.
Referenced by Draw().
|
static |
Returns the 'bin' for 'value' The interval between 'binLow' and 'binUp' is divided into 'nbins' equidistant bins avoid index out of bounds error: 'if (bin < binLow) bin = binLow' and vice versa
\[ GetBin(value) = \frac{nbins - 1}{binUp - binLow} \upoint (value - binLow) +1 \]
Definition at line 1223 of file AliTPCCalibViewer.cxx.
Referenced by GetTree(), Integrate(), and SigmaCut().
AliTPCCalPad * AliTPCCalibViewer::GetCalPad | ( | const char * | desiredData, |
const char * | cuts = "" , |
||
const char * | calPadName = "NoName" |
||
) | const |
creates a AliTPCCalPad out of the 'desiredData' the functionality of EasyDraw1D is used calPadName specifies the name of the created AliTPCCalPad
Definition at line 1480 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoExport(), Draw(), GetArrayOfCalPads(), MakeFit(), and RebuildData().
AliTPCCalPad * AliTPCCalibViewer::GetCalPadOld | ( | const char * | desiredData, |
const char * | cuts = "" , |
||
const char * | calPadName = "NoName" |
||
) | const |
creates a AliTPCCalPad out of the 'desiredData' the functionality of EasyDraw1D is used calPadName specifies the name of the created AliTPCCalPad
Definition at line 1455 of file AliTPCCalibViewer.cxx.
Referenced by Draw().
AliTPCCalROC * AliTPCCalibViewer::GetCalROC | ( | const char * | desiredData, |
UInt_t | sector, | ||
const char * | cuts = "" |
||
) | const |
creates a AliTPCCalROC out of the desiredData the functionality of EasyDraw1D is used sector specifies the sector of the created AliTPCCalROC
Definition at line 1512 of file AliTPCCalibViewer.cxx.
Referenced by Draw(), and MakeTreeWithObjects().
TObjArray * AliTPCCalibViewer::GetListOfNormalizationVariables | ( | Bool_t | printList = kFALSE | ) | const |
produces a list of available variables for normalization in the tree printList: print the list to the screen, after the scan is done
Definition at line 1601 of file AliTPCCalibViewer.cxx.
Referenced by AddAbbreviations(), Draw(), and AliTPCCalibViewerGUI::Initialize().
TObjArray * AliTPCCalibViewer::GetListOfVariables | ( | Bool_t | printList = kFALSE | ) |
scan the tree - produces a list of available variables in the tree printList: print the list to the screen, after the scan is done
Definition at line 1529 of file AliTPCCalibViewer.cxx.
Referenced by AddAbbreviations(), Draw(), GetArrayOfCalPads(), and AliTPCCalibViewerGUI::Initialize().
|
static |
returns the LTM and sigma
Definition at line 1241 of file AliTPCCalibViewer.cxx.
Referenced by DrawHisto1D(), GetTree(), Integrate(), IntegrateOld(), and SigmaCut().
|
inline |
Definition at line 81 of file AliTPCCalibViewer.h.
Referenced by LoadViewer(), MakeFit(), MakeTree(), RebuildData(), and startGUI().
Int_t AliTPCCalibViewer::Integrate | ( | const char * | drawCommand, |
const char * | sector, | ||
const char * | cuts = 0 , |
||
Float_t | sigmaMax = 5 , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE , |
||
const char * | sigmas = "" , |
||
Float_t | sigmaStep = -1 |
||
) | const |
Creates an integrated histogram \(S(t, \mu, \sigma)\), out of the input distribution distribution \(f(x, \mu, \sigma)\), given in "histogram" "mean" and "sigma" are \(\mu\) and \(\sigma\) of the distribution in "histogram", to be specified by the user sigmaMax: up to which sigma around the mean/median/LTM you want to integrate if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used The actual work is done on the array.
\[ f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx \]
Definition at line 1009 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoDraw(), Draw(), GetTree(), Integrate(), and IntegrateOld().
Int_t AliTPCCalibViewer::Integrate | ( | const char * | drawCommand, |
Int_t | sector, | ||
const char * | cuts = 0 , |
||
Float_t | sigmaMax = 5 , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE , |
||
const char * | sigmas = "" , |
||
Float_t | sigmaStep = -1 |
||
) | const |
Creates an integrated histogram \(S(t, \mu, \sigma)\), out of the input distribution distribution \(f(x, \mu, \sigma)\), given in "histogram" "mean" and "sigma" are \(\mu\) and \(\sigma\) of the distribution in "histogram", to be specified by the user sigmaMax: up to which sigma around the mean/median/LTM you want to integrate if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used The actual work is done on the array.
\[ f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx \]
Definition at line 903 of file AliTPCCalibViewer.cxx.
|
static |
Creates an integrated histogram \(S(t, \mu, \sigma)\), out of the input distribution distribution \(f(x, \mu, \sigma)\), given in "histogram" "mean" and "sigma" are \(\mu\) and \(\sigma\) of the distribution in "histogram", to be specified by the user sigmaMax: up to which sigma around the mean/median/LTM you want to integrate if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used The actual work is done on the array.
\[ f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx \]
Definition at line 1365 of file AliTPCCalibViewer.cxx.
|
static |
Creates an integrated histogram \(S(t, \mu, \sigma)\), out of the input distribution distribution \(f(x, \mu, \sigma)\), given in "histogram" "mean" and "sigma" are \(\mu\) and \(\sigma\) of the distribution in "histogram", to be specified by the user sigmaMax: up to which sigma around the mean/median/LTM you want to integrate if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used Here the actual work is done.
Definition at line 1386 of file AliTPCCalibViewer.cxx.
Int_t AliTPCCalibViewer::IntegrateOld | ( | const char * | drawCommand, |
const char * | sector, | ||
const char * | cuts = 0 , |
||
Float_t | sigmaMax = 5 , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE , |
||
const char * | sigmas = "" , |
||
Float_t | sigmaStep = -1 |
||
) | const |
Creates an integrated histogram \(S(t, \mu, \sigma)\), out of the input distribution distribution \(f(x, \mu, \sigma)\), given in "histogram" "mean" and "sigma" are \(\mu\) and \(\sigma\) of the distribution in "histogram", to be specified by the user sigmaMax: up to which sigma around the mean/median/LTM you want to integrate if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used The actual work is done on the array.
\[ f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx \]
Definition at line 924 of file AliTPCCalibViewer.cxx.
Referenced by Draw().
|
static |
make standard aliases for standard calibration trees
Definition at line 2134 of file AliTPCCalibViewer.cxx.
Referenced by GetTree(), and MakeTree().
|
static |
Write a tree with all available information if mapFileName is speciefied, the Map information are also written to the tree pads specified in outlierPad are not used for calculating statistics The following statistical information on the basis of a ROC are calculated: "_Median", "_Mean", "_LTM", "_RMS_LTM" "_Median_OutlierCutted", "_Mean_OutlierCutted", "_RMS_OutlierCutted", "_LTM_OutlierCutted", "_RMS_LTM_OutlierCutted" The following position variables are available: "row", "pad", "lx", "ly", "gx", "gy", "rpad", "channel"
The tree out of this function is the basis for the AliTPCCalibViewer and the AliTPCCalibViewerGUI.
Definition at line 1890 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalPad::AddFriend(), Draw(), AliTPCPreprocessorOnline::DumpToFile(), AliTPCCalPad::DumpUnitTestTrees(), MakeTree(), and UnitTestAliTPCCalPadTree().
|
static |
Function to create a calibration Tree with all available information. See also documentation to MakeTree the file "inputFileName" must be in the following format (see also CreateObjectList): (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
type path dependingOnType
type == "CE": dependingOnType = CETmean CEQmean CETrms
type == "Pulser": dependingOnType = PulserTmean PulsterQmean PulserTrms
type == "Pedestals": dependingOnType = Pedestals Noise
type == "CalPad": dependingOnType = NameInFile NameToWriteToFile
Definition at line 2152 of file AliTPCCalibViewer.cxx.
|
static |
Write tree with all available information im mapFileName is speciefied, the Map information are also written to the tree AliTPCCalPad-Objects are written directly to the tree, so that they can be accessd later on (does not work!!!)
Definition at line 1742 of file AliTPCCalibViewer.cxx.
Referenced by Draw().
AliTPCCalibViewer & AliTPCCalibViewer::operator= | ( | const AliTPCCalibViewer & | param | ) |
assignment operator - dummy not yet working!!!
Definition at line 154 of file AliTPCCalibViewer.cxx.
|
inline |
Definition at line 40 of file AliTPCCalibViewer.h.
|
inline |
Definition at line 41 of file AliTPCCalibViewer.h.
Int_t AliTPCCalibViewer::SigmaCut | ( | const char * | drawCommand, |
Int_t | sector, | ||
const char * | cuts = 0 , |
||
Float_t | sigmaMax = 5 , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE , |
||
Bool_t | pm = kFALSE , |
||
const char * | sigmas = "" , |
||
Float_t | sigmaStep = -1 |
||
) | const |
Creates a histogram \(S(t, \mu, \sigma)\), where you can see, how much of the data are inside sigma-intervals around the mean value The data of the distribution \(f(x, \mu, \sigma)\) are given in 'array', 'n' specifies the length of the array 'mean' and 'sigma' are \(\mu\) and \(\sigma\) of the distribution in 'array', to be specified by the user 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, \(t \sigma\)) sigmaStep: the binsize of the generated histogram
\[ f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = (\int_{\mu}^{\mu + t \sigma} f(x, \mu, \sigma) dx + \int_{\mu}^{\mu - t \sigma} f(x, \mu, \sigma) dx) / (\int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx) \]
Creates a histogram, where you can see, how much of the data are inside sigma-intervals around the mean/median/LTM with drawCommand, sector and cuts you specify your input data, see EasyDraw sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma) sigmaStep: the binsize of the generated histogram plotMean/plotMedian/plotLTM: specifies where to put the center
Definition at line 699 of file AliTPCCalibViewer.cxx.
Referenced by AliTPCCalibViewerGUI::DoDraw(), Draw(), GetTree(), and SigmaCut().
Int_t AliTPCCalibViewer::SigmaCut | ( | const char * | drawCommand, |
const char * | sector, | ||
const char * | cuts = 0 , |
||
Float_t | sigmaMax = 5 , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE , |
||
Bool_t | pm = kFALSE , |
||
const char * | sigmas = "" , |
||
Float_t | sigmaStep = -1 |
||
) | const |
Creates a histogram, where you can see, how much of the data are inside sigma-intervals around the mean/median/LTM with drawCommand, sector and cuts you specify your input data, see EasyDraw sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma) sigmaStep: the binsize of the generated histogram plotMean/plotMedian/plotLTM: specifies where to put the center
Definition at line 725 of file AliTPCCalibViewer.cxx.
|
static |
Creates a histogram \(S(t, \mu, \sigma)\), where you can see, how much of the data are inside sigma-intervals around the mean value The data of the distribution \(f(x, \mu, \sigma)\) are given in 'array', 'n' specifies the length of the array 'mean' and 'sigma' are \(\mu\) and \(\sigma\) of the distribution in 'array', to be specified by the user 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, \(t \sigma\)) sigmaStep: the binsize of the generated histogram Here the actual work is done.
Definition at line 1283 of file AliTPCCalibViewer.cxx.
|
static |
Creates a cumulative histogram \(S(t, \mu, \sigma)\), where you can see, how much of the data are inside sigma-intervals around the mean value The data of the distribution \(f(x, \mu, \sigma)\) are given in 'histogram' 'mean' and 'sigma' are \(\mu\) and \(\sigma\) of the distribution in 'histogram', to be specified by the user sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, \(t \sigma\)) sigmaStep: the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used pm: Decide weather \(t > 0\) (first case) or \(t\) arbitrary (secound case) The actual work is done on the array.
\[ f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = (\int_{\mu}^{\mu + t \sigma} f(x, \mu, \sigma) dx + \int_{\mu}^{\mu - t \sigma} f(x, \mu, \sigma) dx) / (\int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx), for t > 0 or f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{\mu}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx \]
Definition at line 1259 of file AliTPCCalibViewer.cxx.
|
static |
SigmaCut for variable binsize NOT YET IMPLEMENTED !!!
Definition at line 1355 of file AliTPCCalibViewer.cxx.
Int_t AliTPCCalibViewer::SigmaCutNew | ( | const char * | drawCommand, |
const char * | sector, | ||
const char * | cuts = 0 , |
||
Float_t | sigmaMax = 5 , |
||
Bool_t | plotMean = kTRUE , |
||
Bool_t | plotMedian = kTRUE , |
||
Bool_t | plotLTM = kTRUE , |
||
Bool_t | pm = kFALSE , |
||
const char * | sigmas = "" , |
||
Float_t | sigmaStep = -1 |
||
) | const |
Creates a histogram, where you can see, how much of the data are inside sigma-intervals around the mean/median/LTM with drawCommand, sector and cuts you specify your input data, see EasyDraw sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma) sigmaStep: the binsize of the generated histogram plotMean/plotMedian/plotLTM: specifies where to put the center
Definition at line 811 of file AliTPCCalibViewer.cxx.
Referenced by Draw().
|
protected |
the abreviation for '.fElements'
Definition at line 104 of file AliTPCCalibViewer.h.
Referenced by AddAbbreviations(), AliTPCCalibViewer(), EasyDraw(), EasyDraw1D(), Fit(), GetAbbreviation(), GetArrayOfCalPads(), GetCalPadOld(), GetCalROC(), operator=(), and SetAbbreviation().
|
protected |
'.fElements', stored in a TStrig
Definition at line 105 of file AliTPCCalibViewer.h.
Referenced by AddAbbreviations(), AliTPCCalibViewer(), EasyDraw(), EasyDraw1D(), Fit(), FormatHistoLabels(), GetAppendString(), GetListOfNormalizationVariables(), operator=(), and SetAppendString().
|
protected |
file that contains a calPads tree (e.g. written by AliTPCCalPad::MakeTree(...)
Definition at line 101 of file AliTPCCalibViewer.h.
Referenced by AliTPCCalibViewer(), Delete(), and ~AliTPCCalibViewer().
|
protected |
Objects, that will be deleted when the destructor ist called.
Definition at line 102 of file AliTPCCalibViewer.h.
Referenced by AddReferenceTree(), AliTPCCalibViewer(), Delete(), operator=(), and ~AliTPCCalibViewer().
|
protected |
tree containing visualization data (e.g. written by AliTPCCalPad::MakeTree(...)
Definition at line 100 of file AliTPCCalibViewer.h.
Referenced by AddFriend(), AliTPCCalibViewer(), Delete(), Draw(), DrawHisto1D(), EasyDraw(), EasyDraw1D(), Fit(), GetCalPad(), GetCalPadOld(), GetCalROC(), GetListOfVariables(), GetTree(), Integrate(), IntegrateOld(), operator=(), SigmaCut(), SigmaCutNew(), and ~AliTPCCalibViewer().
|
protected |
decides weather the tree must be deleted in destructor or not
Definition at line 103 of file AliTPCCalibViewer.h.
Referenced by AliTPCCalibViewer(), Delete(), operator=(), and ~AliTPCCalibViewer().