]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/Cal/AliTRDCalPIDNN.cxx
Coding rules
[u/mrichter/AliRoot.git] / TRD / Cal / AliTRDCalPIDNN.cxx
1 /**************************************************************************
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  **************************************************************************/
15
16 /* $Id$ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 // Container of the distributions for the neural network                  //
21 //                                                                        //
22 // Author:                                                                //
23 // Alex Wilk <wilka@uni-muenster.de>                                      //
24 //                                                                        //
25 ////////////////////////////////////////////////////////////////////////////
26
27 #include <TFile.h>
28 #include <TROOT.h>
29 #include <TMultiLayerPerceptron.h>
30
31 #include "AliPID.h"
32 #include "AliLog.h"
33
34 #include "AliTRDgeometry.h"
35 #include "AliTRDCalPIDNN.h"
36 #include "AliTRDcalibDB.h"
37
38 ClassImp(AliTRDCalPIDNN)
39
40 //_________________________________________________________________________
41 AliTRDCalPIDNN::AliTRDCalPIDNN()
42   :AliTRDCalPID("pid", "NN PID references for TRD")
43 {
44   //
45   //  The Default constructor
46   //
47
48   Init();
49
50 }
51
52 //_________________________________________________________________________
53 AliTRDCalPIDNN::AliTRDCalPIDNN(const Text_t *name, const Text_t *title) 
54   :AliTRDCalPID(name,title)
55 {
56   //
57   //  The main constructor
58   //
59   
60   Init();
61
62 }
63
64 // //_____________________________________________________________________________
65 // AliTRDCalPIDNN::AliTRDCalPIDNN(const AliTRDCalPIDNN &c) 
66 //   :TNamed(c)
67 //   ,fMeanChargeRatio(c.fMeanChargeRatio)
68 //   ,fModel(0x0)
69 // {
70 //   //
71 //   // Copy constructor
72 //   //
73 // 
74 //   if (this != &c) ((AliTRDCalPIDNN &) c).Copy(*this);
75 //   
76 // }
77
78 //_________________________________________________________________________
79 AliTRDCalPIDNN::~AliTRDCalPIDNN()
80 {
81   //
82   // Destructor
83   //
84
85   //if (fModel) delete fModel;
86   
87 }
88
89 //_________________________________________________________________________
90 Bool_t AliTRDCalPIDNN::LoadReferences(Char_t *refFile)
91 {
92   //
93   // Read the TRD Neural Networks
94   //
95
96   // Read NN Root file  
97   TFile *nnFile = TFile::Open(refFile,"READ");
98   if (!nnFile || !nnFile->IsOpen()) {
99     AliError(Form("Opening TRD histgram file %s failed",refFile));
100     return kFALSE;
101   }
102   gROOT->cd();
103
104   // Read Networks
105   for (Int_t imom = 0; imom < kNMom; imom++) {
106     for (Int_t iplane = 0; iplane < AliTRDgeometry::kNlayer; iplane++) {
107       TMultiLayerPerceptron *nn = (TMultiLayerPerceptron *)
108          nnFile->Get(Form("NN_Mom%d_Plane%d",imom,iplane));
109       fModel->AddAt(nn,GetModelID(imom,0,iplane));
110     }
111   }
112
113   nnFile->Close();
114   delete nnFile;
115
116   return kTRUE;
117
118 }
119
120 //_________________________________________________________________________
121 TObject *AliTRDCalPIDNN::GetModel(Int_t ip, Int_t, Int_t iplane) const
122 {
123   //
124   // Returns one selected TMultiLayerPerceptron. iType not used.
125   //
126
127   if (ip<0 || ip>= kNMom) return 0x0;
128   
129   AliInfo(Form("Retrive MultiLayerPerceptron for %5.2f GeV/c for plane %d" 
130          ,fgTrackMomentum[ip]
131          ,iplane));
132   
133   return fModel->At(GetModelID(ip, 0, iplane));
134
135 }
136
137 //_________________________________________________________________________
138 Double_t AliTRDCalPIDNN::GetProbability(Int_t spec, Float_t mom
139                                       , const Float_t * const dedx
140                                       , Float_t, Int_t iplane) const
141 {
142   //
143   // Core function of AliTRDCalPID class for calculating the
144   // likelihood for species "spec" (see AliTRDtrack::kNspecie) of a
145   // given momentum "mom", a given dE/dx (slice "dedx") yield per
146   // layer in a given layer (iplane)
147   //
148
149   if (spec < 0 || spec >= AliPID::kSPECIES) return 0.;
150
151   // find the interval in momentum and track segment length which applies for this data
152   
153   Int_t imom = 1;
154   while (imom<AliTRDCalPID::kNMom-1 && mom>fgTrackMomentum[imom]) imom++;
155   Double_t lNN1, lNN2;
156   Double_t mom1 = fgTrackMomentum[imom-1], mom2 = fgTrackMomentum[imom];
157
158   TMultiLayerPerceptron *nn = 0x0;
159   if(!(nn = (TMultiLayerPerceptron *) fModel->At(GetModelID(imom-1, spec, iplane/*, ilength*/)))){
160     AliInfo(Form("Looking for mom(%f) plane(%d)", mom-1, iplane));
161     AliError(Form("NN id %d not found in DB.", GetModelID(imom-1, spec, iplane)));
162     return 0.;
163   }
164
165   Double_t ddedx[AliTRDCalPID::kNSlicesNN];
166
167   for (int inode=0; inode<AliTRDCalPID::kNSlicesNN; inode++) {
168     ddedx[inode] = (((Double_t) dedx[inode]/kMLPscale)*3)          // Bug fix! Needs new reference data or different calculation of dedx!!!!
169                  / (AliTRDcalibDB::Instance()->GetNumberOfTimeBins()/AliTRDCalPID::kNSlicesNN);
170   }
171
172   lNN1 = nn->Evaluate(spec, ddedx);
173   
174   if(!(nn = (TMultiLayerPerceptron*)fModel->At(GetModelID(imom, spec, iplane/*, ilength*/)))){
175     AliInfo(Form("Looking for mom(%f) plane(%d)", mom, iplane));
176     AliError(Form("NN id %d not found in DB.", GetModelID(imom, spec, iplane)));
177     return lNN1;
178   }
179   lNN2 = nn->Evaluate(spec, ddedx);
180   
181   // return interpolation over momentum binning
182   if      (mom < fgTrackMomentum[0]) {
183     return lNN1;
184   }
185   else if (mom > fgTrackMomentum[AliTRDCalPID::kNMom-1]) {
186     return lNN2;
187   }
188   else {
189     return lNN1 + (lNN2 - lNN1)*(mom - mom1)/(mom2 - mom1);
190   }
191   
192 }
193
194 //_________________________________________________________________________
195 void AliTRDCalPIDNN::Init()
196 {
197   //
198   // Initialization
199   //
200
201   fModel = new TObjArray(AliTRDgeometry::kNlayer * AliTRDCalPID::kNMom);
202   fModel->SetOwner();
203   
204 }
205
206 //_________________________________________________________________________
207 Int_t AliTRDCalPIDNN::GetModelID(Int_t mom, Int_t /*ii*/, Int_t plane) const
208 {
209   
210   // returns the ID of the NN distribution (66 MLPs, ID from 56 to 121)
211
212   return /*AliPID::kSPECIES * AliTRDCalPID::kNMom + */
213     plane * AliTRDCalPID::kNMom + mom;  
214   
215 }