]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/Cal/AliTRDCalPIDLQ.cxx
o small fix
[u/mrichter/AliRoot.git] / TRD / Cal / AliTRDCalPIDLQ.cxx
CommitLineData
44dbae42 1/**************************************************************************
e7d6a389 2* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3* *
4* Author: The ALICE Off-line Project. *
5* Contributors are mentioned in the code where appropriate. *
6* *
7* Permission to use, copy, modify and distribute this software and its *
8* documentation strictly for non-commercial purposes is hereby granted *
9* without fee, provided that the above copyright notice appears in all *
10* copies and that both the copyright notice and this permission notice *
11* appear in the supporting documentation. The authors make no claims *
12* about the suitability of this software for any purpose. It is *
13* provided "as is" without express or implied warranty. *
14**************************************************************************/
44dbae42 15
16/* $Id$ */
17
18//////////////////////////////////////////////////////////////////////////
19// //
20// Container for the distributions of dE/dx and the time bin of the //
21// max. cluster for electrons and pions //
22// //
23// Authors: //
24// Prashant Shukla <shukla@pi0.physi.uni-heidelberg.de> //
25// Alex Bercuci <a.bercuci@gsi.de> //
26// //
27//////////////////////////////////////////////////////////////////////////
28
44dbae42 29#include <TFile.h>
30#include <TROOT.h>
9ded305e 31#include <TSystem.h>
44dbae42 32
33#include "AliLog.h"
34#include "AliPID.h"
35
9ded305e 36#include "TKDPDF.h"
11d80e40 37
44dbae42 38#include "AliTRDCalPIDLQ.h"
39#include "AliTRDcalibDB.h"
40
41ClassImp(AliTRDCalPIDLQ)
42
2e32a5ae 43Float_t AliTRDCalPIDLQ::fgTrackSegLength[kNLength] = { 3.7, 3.9, 4.2, 5.0 };
44dbae42 44
45//_________________________________________________________________________
46AliTRDCalPIDLQ::AliTRDCalPIDLQ()
47 :AliTRDCalPID("pid", "LQ PID references for TRD")
48{
49 //
50 // The Default constructor
51 //
52
1fcc78b7 53 //Init();
44dbae42 54
55}
56
57//_________________________________________________________________________
58AliTRDCalPIDLQ::AliTRDCalPIDLQ(const Text_t *name, const Text_t *title)
59 :AliTRDCalPID(name,title)
60{
61 //
62 // The main constructor
63 //
64
65 Init();
66
67}
68
44dbae42 69//_________________________________________________________________________
70Bool_t AliTRDCalPIDLQ::LoadReferences(Char_t *refFile)
71{
72 //
73 // Read the TRD dEdx histograms.
74 //
75
9ded305e 76 if(gSystem->AccessPathName(refFile, kReadPermission)){
77 AliError(Form("File %s.root not readable", refFile));
78 return kFALSE;
e7d6a389 79 }
9ded305e 80 if(!TFile::Open(refFile)){
81 AliError(Form("File %s corrupted", refFile));
44dbae42 82 return kFALSE;
83 }
9ded305e 84 TObjArray *pdf(NULL);
11d80e40 85 if (!( pdf = dynamic_cast<TObjArray*>(gFile->Get("PDF_LQ")))) {
9ded305e 86 AliError("PID data not available");
87 return kFALSE;
44dbae42 88 }
9ded305e 89 fModel=(TObjArray*)pdf->Clone("PDF");
90 gFile->Close();
44dbae42 91 return kTRUE;
44dbae42 92}
93
94//_________________________________________________________________________
11d80e40 95TObject* AliTRDCalPIDLQ::GetModel(Int_t ip, Int_t is, Int_t slices) const
44dbae42 96{
97 //
98 // Returns one selected dEdx histogram
99 //
100
11d80e40 101 if (is < 0 || is >= AliPID::kSPECIES) return NULL;
9ded305e 102 if(ip<0 || ip>= kNMom ) return NULL;
44dbae42 103
11d80e40 104 AliDebug(2, Form("Retrive dEdx distribution for %s @ p=%5.2f [GeV/c].", AliPID::ParticleShortName(is), fgTrackMomentum[ip]));
105 return fModel->At(GetModelID(ip, is, slices));
106}
107
108//_________________________________________________________________________
109Int_t AliTRDCalPIDLQ::GetNrefs()
110{
111// Returns the number of PDF distribution
112 return AliTRDCalPID::kNMom*AliPID::kSPECIES*2;
44dbae42 113}
114
115//_________________________________________________________________________
1fcc78b7 116Double_t AliTRDCalPIDLQ::GetProbability(Int_t spec, Float_t mom
117 , const Float_t * const dedx
11d80e40 118 , Float_t length
119 , Int_t slices) const
44dbae42 120{
11d80e40 121//
122// Core function of AliTRDCalPID class for calculating the
123// likelihood for species "spec" (see AliTRDtrack::kNspecie) of a
124// given momentum "mom" and a given dE/dx (slice "dedx") yield per
125// layer
126//
44dbae42 127
e7d6a389 128 if (spec < 0 || spec >= AliPID::kSPECIES) return 0.;
11d80e40 129
130 Bool_t k2D(slices==2);
131 Double_t x[]={0., 0.};
132 if(!CookdEdx(dedx, x, k2D)) return 0.;
e7d6a389 133
e7d6a389 134 // find the interval in momentum and track segment length which applies for this data
9ded305e 135/* Int_t ilength = 1;
2e32a5ae 136 while(ilength<kNLength-1 && length>fgTrackSegLength[ilength]){
e7d6a389 137 ilength++;
9ded305e 138 }*/
e7d6a389 139 Int_t imom = 1;
2e32a5ae 140 while(imom<kNMom-1 && mom>fgTrackMomentum[imom]) imom++;
44dbae42 141
9ded305e 142 Double_t p[2], e[2], r;
143 TKDPDF *pdf(NULL);
144
11d80e40 145 AliDebug(2, Form("Looking %cD for %s@p=%6.4f[GeV/c] log(dEdx)={%7.2f %7.2f}[a.u.] l=%4.2f[cm].", k2D?'2':'1', AliPID::ParticleShortName(spec), mom, x[0], x[1], length));
146 if(!(pdf = dynamic_cast<TKDPDF*>(fModel->At(GetModelID(imom-1, spec, slices))))) {
147 AliError(Form("%cD Ref data @ idx[%d] not found in DB.", k2D?'2':'1', GetModelID(imom-1, spec, slices)));
e7d6a389 148 return 0.;
149 }
a4f67853 150 if(!pdf->GetSize()) pdf->Bootstrap();
9ded305e 151 pdf->Eval(x, r, e[0], kFALSE);
152 p[0]=TMath::Abs(r); // conversion from interpolation to PDF
153 AliDebug(2, Form("LQ=%6.3f+-%5.2f%% @ %4.1f[GeV/c]", p[0], 1.E2*e[0]/p[0], fgTrackMomentum[imom-1]));
154
11d80e40 155 if(!(pdf = dynamic_cast<TKDPDF*>(fModel->At(GetModelID(imom, spec, slices))))){
156 AliError(Form("%cD Ref data @ idx[%d] not found in DB.", k2D?'2':'1', GetModelID(imom, spec, slices)));
9ded305e 157 return p[0];
e7d6a389 158 }
a4f67853 159 if(!pdf->GetSize()) pdf->Bootstrap();
9ded305e 160 pdf->Eval(x, r, e[1], kFALSE);
161 p[1]=TMath::Abs(r); // conversion from interpolation to PDF
162 AliDebug(2, Form("LQ=%6.3f+-%5.2f%% @ %4.1f[GeV/c]", p[1], 1.E2*e[1]/p[1], fgTrackMomentum[imom]));
e7d6a389 163
9ded305e 164 // return interpolation over momentum binning
165 if(mom < fgTrackMomentum[0]) return p[0];
166 else if(mom > fgTrackMomentum[kNMom-1]) return p[1];
167 else{
168 Double_t lmom[2] = {fgTrackMomentum[imom-1], fgTrackMomentum[imom]};
169 return p[0] + (p[1] - p[0])*(mom - lmom[0])/(lmom[1] - lmom[0]);
170 }
44dbae42 171}
172
173//_________________________________________________________________________
174void AliTRDCalPIDLQ::Init()
175{
9ded305e 176//
177// PID LQ list initialization
178//
11d80e40 179 fModel = new TObjArray(AliPID::kSPECIES * kNMom * 2);
180 fModel->SetOwner();
44dbae42 181}