AliRoot Core  v5-06-30 (35d6c57)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliTPCCalibViewer Class Reference

Class for viewing/visualizing TPC calibration data. More...

#include <AliTPCCalibViewer.h>

Inheritance diagram for AliTPCCalibViewer:

Public Member Functions

 AliTPCCalibViewer ()
 
 AliTPCCalibViewer (const AliTPCCalibViewer &c)
 
 AliTPCCalibViewer (TTree *const tree)
 
 AliTPCCalibViewer (const char *fileName, const char *treeName="calPads")
 
AliTPCCalibVieweroperator= (const AliTPCCalibViewer &param)
 
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)
 
const char * 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
 
AliTPCCalPadGetCalPadOld (const char *desiredData, const char *cuts="", const char *calPadName="NoName") const
 
AliTPCCalPadGetCalPad (const char *desiredData, const char *cuts="", const char *calPadName="NoName") const
 
AliTPCCalROCGetCalROC (const char *desiredData, UInt_t sector, const char *cuts="") const
 
TObjArrayGetArrayOfCalPads ()
 
TObjArrayGetListOfVariables (Bool_t printList=kFALSE)
 
TObjArrayGetListOfNormalizationVariables (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...
 
TObjArrayfListOfObjectsToBeDeleted
 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...
 

Detailed Description

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.

Constructor & Destructor Documentation

AliTPCCalibViewer::AliTPCCalibViewer ( )

Definition at line 75 of file AliTPCCalibViewer.cxx.

AliTPCCalibViewer::AliTPCCalibViewer ( const AliTPCCalibViewer c)

dummy AliTPCCalibViewer copy constructor not yet working!!!

Definition at line 91 of file AliTPCCalibViewer.cxx.

AliTPCCalibViewer::AliTPCCalibViewer ( TTree *const  tree)

Constructor that initializes the calibration viewer

Definition at line 112 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 131 of file AliTPCCalibViewer.cxx.

AliTPCCalibViewer::~AliTPCCalibViewer ( )
virtual

AliTPCCalibViewer destructor all objects will be deleted, the file will be closed, the pictures will disappear

Definition at line 170 of file AliTPCCalibViewer.cxx.

Member Function Documentation

const char * 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":

  • 'Save' correct items:
    • form <replaceString>, take <variable>'s first char, add <removeString>, add rest of <variable>, e.g. "C!#EQmean" (<removeString> = "!#")
    • For each statistical information in "listOfNormalizationVariables":
      • ReplaceAll <variable><statistical_Information> with <replaceString><statistical_Information>
    • ReplaceAll <variable><abbreviation> with <replaceString><abbreviation>, e.g. "CEQmean~" -> "C!#EQmean~"
    • ReplaceAll <variable><appendStr> with <replaceString><appendStr>, e.g. "CEQmean.fElements" -> "C!#EQmean.fElements"
  • Do actual replacing:
    • ReplaceAll <variable> with <variable><fAbbreviation>, e.g. "CEQmean" -> "CEQmean~"
  • Undo saving:
    • For each statistical information in "listOfNormalizationVariables":
      • ReplaceAll <replaceString><statistical_Information> with <variable><statistical_Information>
    • ReplaceAll <replaceString><abbreviation> with <variable><abbreviation>, e.g. "C!#EQmean~" -> "CEQmean~"
    • ReplaceAll <replaceString><appendStr> with <variable><appendStr>, e.g. "C!#EQmean.fElements" -> "CEQmean.fElements"

Now all the missing "~" should be added.

Definition at line 207 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoFit(), AliTPCCalibViewerGUI::GetCutString(), and AliTPCCalibViewerGUI::GetDrawString().

TFriendElement* AliTPCCalibViewer::AddFriend ( const char *  treename,
const char *  filename 
)
inline

Definition at line 78 of file AliTPCCalibViewer.h.

Referenced by AddFiles(), and AddReferenceTree().

TFriendElement* AliTPCCalibViewer::AddFriend ( TTree *  tree,
const char *  alias,
Bool_t  warn = kFALSE 
)
inline

Definition at line 79 of file AliTPCCalibViewer.h.

TFriendElement* AliTPCCalibViewer::AddFriend ( const char *  treename,
TFile *  file 
)
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 1598 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoLoadTree(), and AliTPCCalibViewerGUItime::SetGuiTree().

void AliTPCCalibViewer::CreateObjectList ( const Char_t *  filename,
TObjArray calibObjects 
)
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 2143 of file AliTPCCalibViewer.cxx.

Referenced by MakeTree().

void AliTPCCalibViewer::Delete ( Option_t *  option = "")
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 193 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::~AliTPCCalibViewerGUI().

virtual void AliTPCCalibViewer::Draw ( Option_t *  opt = "")
inlinevirtual

Definition at line 43 of file AliTPCCalibViewer.h.

Referenced by Fit().

virtual Long64_t AliTPCCalibViewer::Draw ( const char *  varexp,
const TCut &  selection,
Option_t *  option = "",
Long64_t  nentries = 1000000000,
Long64_t  firstentry = 0 
)
inlinevirtual

Definition at line 44 of file AliTPCCalibViewer.h.

virtual Long64_t AliTPCCalibViewer::Draw ( const char *  varexp,
const char *  selection,
Option_t *  option = "",
Long64_t  nentries = 1000000000,
Long64_t  firstentry = 0 
)
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 527 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoDraw().

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 542 of file AliTPCCalibViewer.cxx.

void AliTPCCalibViewer::DrawLines ( TH1F *  cutHistoMean,
TVectorF  nsigma,
TLegend *  legend,
Int_t  color,
Bool_t  pm 
) const
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 1083 of file AliTPCCalibViewer.cxx.

Referenced by Integrate(), IntegrateOld(), SigmaCut(), and SigmaCutNew().

void AliTPCCalibViewer::DrawLines ( TGraph *  graph,
TVectorF  nsigma,
TLegend *  legend,
Int_t  color,
Bool_t  pm 
) const
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 1133 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 311 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoDraw(), 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 403 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 420 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoDraw(), 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 480 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 1631 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoFit(), and AliTPCCalibViewerGUI::GetDrawString().

void AliTPCCalibViewer::FormatHistoLabels ( TH1 *  histo) const

formats title and axis labels of histo removes '.fElements'

Definition at line 496 of file AliTPCCalibViewer.cxx.

Referenced by EasyDraw(), and EasyDraw1D().

TString& AliTPCCalibViewer::GetAbbreviation ( )
inline

Definition at line 38 of file AliTPCCalibViewer.h.

Referenced by AliTPCCalibViewerGUI::GetDrawString().

TString& AliTPCCalibViewer::GetAppendString ( )
inline

Definition at line 39 of file AliTPCCalibViewer.h.

TObjArray * AliTPCCalibViewer::GetArrayOfCalPads ( )

Returns a TObjArray with all AliTPCCalPads that are stored in the tree

  • this takes a while -

Definition at line 1609 of file AliTPCCalibViewer.cxx.

Int_t AliTPCCalibViewer::GetBin ( Float_t  value,
Int_t  nbins,
Double_t  binLow,
Double_t  binUp 
)
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 1188 of file AliTPCCalibViewer.cxx.

Referenced by 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

  • this takes a while -

Definition at line 1445 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoExport(), 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

  • this takes a while -

Definition at line 1420 of file AliTPCCalibViewer.cxx.

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 1477 of file AliTPCCalibViewer.cxx.

Referenced by 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 1566 of file AliTPCCalibViewer.cxx.

Referenced by AddAbbreviations(), 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 1494 of file AliTPCCalibViewer.cxx.

Referenced by AddAbbreviations(), GetArrayOfCalPads(), and AliTPCCalibViewerGUI::Initialize().

Double_t AliTPCCalibViewer::GetLTM ( Int_t  n,
const Double_t *const  array,
Double_t *const  sigma = 0,
Double_t  fraction = 0.9 
)
static

returns the LTM and sigma

Definition at line 1206 of file AliTPCCalibViewer.cxx.

Referenced by DrawHisto1D(), Integrate(), IntegrateOld(), and SigmaCut().

TTree* AliTPCCalibViewer::GetTree ( ) const
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 974 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoDraw(), 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 868 of file AliTPCCalibViewer.cxx.

TH1F * AliTPCCalibViewer::Integrate ( TH1F *const  histogram,
Float_t  mean = 0,
Float_t  sigma = 0,
Float_t  sigmaMax = 0,
Float_t  sigmaStep = -1 
)
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 \]

AliTPCCalibViewer_cxx_f68daaf.png
Picture from ROOT macro

Definition at line 1330 of file AliTPCCalibViewer.cxx.

TH1F * AliTPCCalibViewer::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

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 1351 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 889 of file AliTPCCalibViewer.cxx.

void AliTPCCalibViewer::MakeCalPadAliases ( TTree *  tree)
static

make standard aliases for standard calibration trees

Definition at line 2099 of file AliTPCCalibViewer.cxx.

Referenced by MakeTree().

void AliTPCCalibViewer::MakeTree ( const char *  fileName,
TObjArray array,
const char *  mapFileName = 0,
AliTPCCalPad *const  outlierPad = 0,
Float_t  ltmFraction = 0.9 
)
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 1855 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalPad::AddFriend(), AliTPCPreprocessorOnline::DumpToFile(), MakeTree(), and UnitTestAliTPCCalPadTree().

void AliTPCCalibViewer::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

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 2117 of file AliTPCCalibViewer.cxx.

void AliTPCCalibViewer::MakeTreeWithObjects ( const char *  fileName,
const TObjArray *const  array,
const char *  mapFileName = 0 
)
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 1707 of file AliTPCCalibViewer.cxx.

AliTPCCalibViewer & AliTPCCalibViewer::operator= ( const AliTPCCalibViewer param)

assignment operator - dummy not yet working!!!

Definition at line 153 of file AliTPCCalibViewer.cxx.

void AliTPCCalibViewer::SetAbbreviation ( const Char_t *  abr)
inline

Definition at line 40 of file AliTPCCalibViewer.h.

void AliTPCCalibViewer::SetAppendString ( const Char_t *  str)
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 664 of file AliTPCCalibViewer.cxx.

Referenced by AliTPCCalibViewerGUI::DoDraw(), 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 690 of file AliTPCCalibViewer.cxx.

TH1F * AliTPCCalibViewer::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

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 1248 of file AliTPCCalibViewer.cxx.

TH1F * AliTPCCalibViewer::SigmaCut ( TH1F *const  histogram,
Float_t  mean,
Float_t  sigma,
Float_t  sigmaMax,
Float_t  sigmaStep = -1,
Bool_t  pm = kFALSE 
)
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 \]

AliTPCCalibViewer_cxx_0eb126d.png
Picture from ROOT macro

Definition at line 1224 of file AliTPCCalibViewer.cxx.

TH1F * AliTPCCalibViewer::SigmaCut ( Int_t  n,
const Double_t *  array,
Double_t  mean,
Double_t  sigma,
Int_t  nbins,
const Double_t *  xbins,
Double_t  sigmaMax 
)
static

SigmaCut for variable binsize NOT YET IMPLEMENTED !!!

Definition at line 1320 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 776 of file AliTPCCalibViewer.cxx.

Member Data Documentation

TString AliTPCCalibViewer::fAbbreviation
protected
TString AliTPCCalibViewer::fAppendString
protected
TFile* AliTPCCalibViewer::fFile
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().

TObjArray* AliTPCCalibViewer::fListOfObjectsToBeDeleted
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().

TTree* AliTPCCalibViewer::fTree
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().

Bool_t AliTPCCalibViewer::fTreeMustBeDeleted
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().


The documentation for this class was generated from the following files: