ROOT logo
/**************************************************************************
 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 *                                                                        *
 * Author: The ALICE Off-line Project.                                    *
 * Contributors are mentioned in the code where appropriate.              *
 *                                                                        *
 * Permission to use, copy, modify and distribute this software and its   *
 * documentation strictly for non-commercial purposes is hereby granted   *
 * without fee, provided that the above copyright notice appears in all   *
 * copies and that both the copyright notice and this permission notice   *
 * appear in the supporting documentation. The authors make no claims     *
 * about the suitability of this software for any purpose. It is          *
 * provided "as is" without express or implied warranty.                  *
 **************************************************************************/

/* $Id$ */

//--------------------------------------------------------------------
//          Options for the TPC Reconstruction in rec.C
//
//  4 options can be set to change the input for TPC reconstruction
//  which overwrites the usage of fUseHLTClusters of the AliTPCRecoParam
//
//  1) useRAW        - use RAW, if not present -> do nothing
//  2) useRAWorHLT   - use RAW, if not present -> use HLT clusters
//  3) useHLT        - use HLT clusters, if not present -> do nothing
//  4) useHLTorRAW   - use HLT clusters, if not present -> use RAW
//
//  -> The current default is useHLTorRAW
//--------------------------------------------------------------------

///////////////////////////////////////////////////////////////////////////////
//                                                                           //
// class for TPC reconstruction                                              //
//                                                                           //
///////////////////////////////////////////////////////////////////////////////

#include <TObject.h>
#include <TString.h>
#include <TObjString.h>
#include <TObjArray.h>
#include <TFile.h>

#include <AliLog.h>
#include <AliPID.h>
#include <AliESDpid.h>
#include <AliTPCPIDResponse.h>
#include "AliTPCReconstructor.h"
#include "AliRunLoader.h"
#include "AliRun.h"
#include "AliRawReader.h"
#include "AliTPCclusterer.h"
#include "AliTPCtracker.h"
#include "AliTPCParam.h"
#include "AliTPCParamSR.h"
#include "AliTPCcalibDB.h"
#include "AliTracker.h"
#include "AliMagF.h"

ClassImp(AliTPCReconstructor)


Int_t    AliTPCReconstructor::fgStreamLevel     = 0;        // stream (debug) level
AliTPCAltroEmulator *  AliTPCReconstructor::fAltroEmulator=0;    // ALTRO emulator

AliTPCReconstructor::AliTPCReconstructor():
AliReconstructor(),
fClusterer(NULL),
fArrSplines(NULL)
{
  //
  // default constructor
  //
  //
  //
  AliTPCcalibDB * calib = AliTPCcalibDB::Instance();
  const AliMagF * field = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
  calib->SetExBField(field);
  AliTPCParam* param = GetTPCParam();
  if (!param) {
    AliWarning("Loading default TPC parameters !");
    param = new AliTPCParamSR;
  }
  fClusterer = new AliTPCclusterer(param);
}

AliTPCReconstructor::AliTPCReconstructor(const AliTPCReconstructor& /*rec*/):
AliReconstructor(),
fClusterer(NULL),
fArrSplines(NULL)
{
  //
  // Dummy copu constructor
  //
}

AliTPCReconstructor& AliTPCReconstructor::operator=(const AliTPCReconstructor&){
  //
  // dummy operator
  //
  return *this;
}

//_____________________________________________________________________________
AliTPCReconstructor::~AliTPCReconstructor()
{
  if (fClusterer)   delete fClusterer;
  delete fArrSplines;
  delete AliTPCcalibDB::Instance();
}

//_____________________________________________________________________________
void AliTPCReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const {
  // single event local reconstruction
  // of TPC data
  fClusterer->SetInput(digitsTree);
  fClusterer->SetOutput(clustersTree);
  fClusterer->Digits2Clusters();
}

//_____________________________________________________________________________
void AliTPCReconstructor::Reconstruct(AliRawReader* rawReader, TTree* clustersTree) const {
  // single event local reconstruction
  // of TPC data starting from raw data

  fClusterer->SetOutput(clustersTree);
  fClusterer->Digits2Clusters(rawReader);
}

//_____________________________________________________________________________
AliTracker* AliTPCReconstructor::CreateTracker() const
{
// create a TPC tracker

  AliTPCParam* param = GetTPCParam();
  if (!param) {
    AliWarning("Loading default TPC parameters !");
    param = new AliTPCParamSR;
  }
  param->ReadGeoMatrices();
  
  AliTPCtracker* tracker = new AliTPCtracker(param);

  ParseOptions(tracker);

  return tracker;
}

//_____________________________________________________________________________
void AliTPCReconstructor::FillESD(TTree */*digitsTree*/, TTree */*clustersTree*/,
				  AliESDEvent* /*esd*/) const
{
// make PID
/*  Now done in AliESDpid
  Double_t parTPC[] = {50., 0.07, 5.};  // MIP nnormalized to channel 50 -MI
  AliTPCpidESD tpcPID(parTPC);
  tpcPID.MakePID(esd);
*/
}


//_____________________________________________________________________________
AliTPCParam* AliTPCReconstructor::GetTPCParam() const
{
// get the TPC parameters

  AliTPCParam* param = AliTPCcalibDB::Instance()->GetParameters();

  return param;
}

//_____________________________________________________________________________
void AliTPCReconstructor::SetSplinesFromOADB(const char* tmplt, AliESDpid *esdPID)
{
  //
  //  load splines from the OADB using 'template'
  //

  // only load splines if not already set
  if (!fArrSplines) {
    fArrSplines=new TObjArray(Int_t(AliPID::kSPECIES));
    fArrSplines->SetOwner();
    TString stemplate(tmplt);

    TString fileNamePIDresponse("$ALICE_ROOT/OADB/COMMON/PID/data/TPCPIDResponse.root");
    TFile f(fileNamePIDresponse.Data());

    TObjArray *arrPidResponseMaster=0x0;

    if (f.IsOpen() && !f.IsZombie()){
      arrPidResponseMaster=dynamic_cast<TObjArray*>(f.Get("TPCPIDResponse"));
    }
    f.Close();

    if (!arrPidResponseMaster){
      AliError("PID response array not found, cannot assign proper splines");
      return;
    }

    for (Int_t ispec=0; ispec<AliPID::kSPECIES; ++ispec)
    {
      Int_t ispec2=ispec;
      if (ispec==Int_t(AliPID::kMuon)) ispec2=Int_t(AliPID::kPion);

      TString particle=AliPID::ParticleName(ispec2);
      particle.ToUpper();

      TString splineName;
      splineName.Form(stemplate.Data(),particle.Data());
      TObject *spline=arrPidResponseMaster->FindObject(splineName.Data());
      if (!spline) {
        AliError(Form("No spline found for '%s'", splineName.Data()));
        continue;
      };
      AliInfo(Form("Adding Response function %d:%s",ispec,splineName.Data()));
      fArrSplines->AddAt(spline->Clone(), ispec);
    }    
    arrPidResponseMaster->Delete();
    delete arrPidResponseMaster;
    if (fArrSplines->GetEntries()!=Int_t(AliPID::kSPECIES)) {
      AliError("Splines not found for all species, cannot use proper PID");
      delete fArrSplines;
      fArrSplines=NULL;
      return;
    }
  }

  for (Int_t ispec=0; ispec<AliPID::kSPECIES; ++ispec)
  {
    esdPID->GetTPCResponse().SetResponseFunction( (AliPID::EParticleType)ispec, fArrSplines->UncheckedAt(ispec) );
  }

  esdPID->GetTPCResponse().SetUseDatabase(kTRUE);
}

//_____________________________________________________________________________
void AliTPCReconstructor::GetPidSettings(AliESDpid *esdPID)
{
  //
  // Get TPC pid splines. They should be written to the OCDB during the CPass
  // the splines themselves are owned by the OCDB object
  //

  // parse options
  TString allopt(GetOption());
  TObjArray *optArray=allopt.Tokenize(";");

  // defines whether the pid was set via a specific option in the rec.C
  Bool_t pidSetInOptions = kFALSE;
  
  for (Int_t iopt=0; iopt<optArray->GetEntriesFast(); ++iopt){
    if (!optArray->At(iopt)) continue;
    TString option(static_cast<TObjString*>(optArray->At(iopt))->GetString().Strip(TString::kBoth,' '));

    if (!option.BeginsWith("PID.")) continue;

    // remove 'PID.' identifyer
    option.Remove(0,4);

    // parse PID type
    if (option.BeginsWith("Static=")){
      option.Remove(0,option.First('=')+1);
      if (option.Contains("LHC13b2_fix_PID")) {
        esdPID->GetTPCResponse().SetBetheBlochParameters(0.0320981, 19.9768, 2.52666e-16, 2.72123, 6.08092);
        esdPID->GetTPCResponse().SetMip(53.4968);
        pidSetInOptions=kTRUE;
      }
      
    } else if (option.BeginsWith("OADB=")) {
      option.Remove(0,option.First('=')+1);
      AliInfo(Form("Setting splines From OADB using template: '%s'",option.Data()));
      SetSplinesFromOADB(option, esdPID);
      pidSetInOptions=kTRUE;
    } else if (option.BeginsWith("OCDB=")){
      option.Remove(0,option.First('=')+1);
      // not yet implemented
    }
    
  }

  delete optArray;

  //
  // Initialisation of BB parameters from the OCDB.
  // They are stored in the AliTPCParam
  //
  if (!pidSetInOptions) {
    AliTPCParam* param = AliTPCcalibDB::Instance()->GetParameters();
    if (param) {
      TVectorD *paramBB=param->GetBetheBlochParameters();
      if (paramBB){
        esdPID->GetTPCResponse().SetBetheBlochParameters((*paramBB)(0),(*paramBB)(1),(*paramBB)(2),(*paramBB)(3),(*paramBB)(4));
        AliInfo(Form("Setting BB parameters from OCDB (AliTPCParam): %.2g, %.2g, %.2g, %.2g, %.2g",
                     (*paramBB)(0),(*paramBB)(1),(*paramBB)(2),(*paramBB)(3),(*paramBB)(4)));
      } else {
        AliError("Couldn't get BB parameters from OCDB, the old default values will be used instead");
      }
    } else {
      AliError("Couldn't get TPC parameters");
    }
  }
  
/*
  AliTPCcalibDB * calib = AliTPCcalibDB::Instance();
  
  //Get pid splines array
  TObjArray *arrSplines=calib->GetPidResponse();
  if (!arrSplines) return;
  AliTPCPIDResponse &tpcPID=esdPID->GetTPCResponse();
  tpcPID.SetUseDatabase(kTRUE);

  // check if parametrisations are already set.
  // since this is uniq for one run, we don't have to reload them
  if (tpcPID.GetResponseFunction(AliPID::kPion)) return;

  // get the default object
  TObject *defaultPID=arrSplines->At(AliPID::kUnknown);
  
  // loop over all particle species and set the response functions
  for (Int_t ispec=0; ispec<AliPID::kUnknown; ++ispec){
    TObject *pidSpline=arrSplines->At(ispec);
    if (!pidSpline) pidSpline=defaultPID;
    tpcPID.SetResponseFunction((AliPID::EParticleType)ispec,pidSpline);
  }
 */ 
}

//_____________________________________________________________________________
void AliTPCReconstructor::ParseOptions( AliTPCtracker* tracker ) const
{
// parse options from rec.C and set in clusterer and tracker
  
  TString option = GetOption();
  
  Int_t useHLTClusters = 3;

  if (option.Contains("use")) {
    
    AliInfo(Form("Overide TPC RecoParam with option %s",option.Data()));
    
    if (option.Contains("useRAW")) {
      useHLTClusters = 1;
      if (option.Contains("useRAWorHLT"))
	useHLTClusters = 2;
    }
    else if (option.Contains("useHLT")) {
      useHLTClusters = 3;
      if (option.Contains("useHLTorRAW"))
	useHLTClusters = 4;
    }
  }
  else {
    const AliTPCRecoParam* param = GetRecoParam();
    useHLTClusters = param->GetUseHLTClusters();
  }

  AliInfo(Form("Usage of HLT clusters in TPC reconstruction : %d", useHLTClusters));

  fClusterer->SetUseHLTClusters(useHLTClusters);
  tracker->SetUseHLTClusters(useHLTClusters);

  return;
}
 AliTPCReconstructor.cxx:1
 AliTPCReconstructor.cxx:2
 AliTPCReconstructor.cxx:3
 AliTPCReconstructor.cxx:4
 AliTPCReconstructor.cxx:5
 AliTPCReconstructor.cxx:6
 AliTPCReconstructor.cxx:7
 AliTPCReconstructor.cxx:8
 AliTPCReconstructor.cxx:9
 AliTPCReconstructor.cxx:10
 AliTPCReconstructor.cxx:11
 AliTPCReconstructor.cxx:12
 AliTPCReconstructor.cxx:13
 AliTPCReconstructor.cxx:14
 AliTPCReconstructor.cxx:15
 AliTPCReconstructor.cxx:16
 AliTPCReconstructor.cxx:17
 AliTPCReconstructor.cxx:18
 AliTPCReconstructor.cxx:19
 AliTPCReconstructor.cxx:20
 AliTPCReconstructor.cxx:21
 AliTPCReconstructor.cxx:22
 AliTPCReconstructor.cxx:23
 AliTPCReconstructor.cxx:24
 AliTPCReconstructor.cxx:25
 AliTPCReconstructor.cxx:26
 AliTPCReconstructor.cxx:27
 AliTPCReconstructor.cxx:28
 AliTPCReconstructor.cxx:29
 AliTPCReconstructor.cxx:30
 AliTPCReconstructor.cxx:31
 AliTPCReconstructor.cxx:32
 AliTPCReconstructor.cxx:33
 AliTPCReconstructor.cxx:34
 AliTPCReconstructor.cxx:35
 AliTPCReconstructor.cxx:36
 AliTPCReconstructor.cxx:37
 AliTPCReconstructor.cxx:38
 AliTPCReconstructor.cxx:39
 AliTPCReconstructor.cxx:40
 AliTPCReconstructor.cxx:41
 AliTPCReconstructor.cxx:42
 AliTPCReconstructor.cxx:43
 AliTPCReconstructor.cxx:44
 AliTPCReconstructor.cxx:45
 AliTPCReconstructor.cxx:46
 AliTPCReconstructor.cxx:47
 AliTPCReconstructor.cxx:48
 AliTPCReconstructor.cxx:49
 AliTPCReconstructor.cxx:50
 AliTPCReconstructor.cxx:51
 AliTPCReconstructor.cxx:52
 AliTPCReconstructor.cxx:53
 AliTPCReconstructor.cxx:54
 AliTPCReconstructor.cxx:55
 AliTPCReconstructor.cxx:56
 AliTPCReconstructor.cxx:57
 AliTPCReconstructor.cxx:58
 AliTPCReconstructor.cxx:59
 AliTPCReconstructor.cxx:60
 AliTPCReconstructor.cxx:61
 AliTPCReconstructor.cxx:62
 AliTPCReconstructor.cxx:63
 AliTPCReconstructor.cxx:64
 AliTPCReconstructor.cxx:65
 AliTPCReconstructor.cxx:66
 AliTPCReconstructor.cxx:67
 AliTPCReconstructor.cxx:68
 AliTPCReconstructor.cxx:69
 AliTPCReconstructor.cxx:70
 AliTPCReconstructor.cxx:71
 AliTPCReconstructor.cxx:72
 AliTPCReconstructor.cxx:73
 AliTPCReconstructor.cxx:74
 AliTPCReconstructor.cxx:75
 AliTPCReconstructor.cxx:76
 AliTPCReconstructor.cxx:77
 AliTPCReconstructor.cxx:78
 AliTPCReconstructor.cxx:79
 AliTPCReconstructor.cxx:80
 AliTPCReconstructor.cxx:81
 AliTPCReconstructor.cxx:82
 AliTPCReconstructor.cxx:83
 AliTPCReconstructor.cxx:84
 AliTPCReconstructor.cxx:85
 AliTPCReconstructor.cxx:86
 AliTPCReconstructor.cxx:87
 AliTPCReconstructor.cxx:88
 AliTPCReconstructor.cxx:89
 AliTPCReconstructor.cxx:90
 AliTPCReconstructor.cxx:91
 AliTPCReconstructor.cxx:92
 AliTPCReconstructor.cxx:93
 AliTPCReconstructor.cxx:94
 AliTPCReconstructor.cxx:95
 AliTPCReconstructor.cxx:96
 AliTPCReconstructor.cxx:97
 AliTPCReconstructor.cxx:98
 AliTPCReconstructor.cxx:99
 AliTPCReconstructor.cxx:100
 AliTPCReconstructor.cxx:101
 AliTPCReconstructor.cxx:102
 AliTPCReconstructor.cxx:103
 AliTPCReconstructor.cxx:104
 AliTPCReconstructor.cxx:105
 AliTPCReconstructor.cxx:106
 AliTPCReconstructor.cxx:107
 AliTPCReconstructor.cxx:108
 AliTPCReconstructor.cxx:109
 AliTPCReconstructor.cxx:110
 AliTPCReconstructor.cxx:111
 AliTPCReconstructor.cxx:112
 AliTPCReconstructor.cxx:113
 AliTPCReconstructor.cxx:114
 AliTPCReconstructor.cxx:115
 AliTPCReconstructor.cxx:116
 AliTPCReconstructor.cxx:117
 AliTPCReconstructor.cxx:118
 AliTPCReconstructor.cxx:119
 AliTPCReconstructor.cxx:120
 AliTPCReconstructor.cxx:121
 AliTPCReconstructor.cxx:122
 AliTPCReconstructor.cxx:123
 AliTPCReconstructor.cxx:124
 AliTPCReconstructor.cxx:125
 AliTPCReconstructor.cxx:126
 AliTPCReconstructor.cxx:127
 AliTPCReconstructor.cxx:128
 AliTPCReconstructor.cxx:129
 AliTPCReconstructor.cxx:130
 AliTPCReconstructor.cxx:131
 AliTPCReconstructor.cxx:132
 AliTPCReconstructor.cxx:133
 AliTPCReconstructor.cxx:134
 AliTPCReconstructor.cxx:135
 AliTPCReconstructor.cxx:136
 AliTPCReconstructor.cxx:137
 AliTPCReconstructor.cxx:138
 AliTPCReconstructor.cxx:139
 AliTPCReconstructor.cxx:140
 AliTPCReconstructor.cxx:141
 AliTPCReconstructor.cxx:142
 AliTPCReconstructor.cxx:143
 AliTPCReconstructor.cxx:144
 AliTPCReconstructor.cxx:145
 AliTPCReconstructor.cxx:146
 AliTPCReconstructor.cxx:147
 AliTPCReconstructor.cxx:148
 AliTPCReconstructor.cxx:149
 AliTPCReconstructor.cxx:150
 AliTPCReconstructor.cxx:151
 AliTPCReconstructor.cxx:152
 AliTPCReconstructor.cxx:153
 AliTPCReconstructor.cxx:154
 AliTPCReconstructor.cxx:155
 AliTPCReconstructor.cxx:156
 AliTPCReconstructor.cxx:157
 AliTPCReconstructor.cxx:158
 AliTPCReconstructor.cxx:159
 AliTPCReconstructor.cxx:160
 AliTPCReconstructor.cxx:161
 AliTPCReconstructor.cxx:162
 AliTPCReconstructor.cxx:163
 AliTPCReconstructor.cxx:164
 AliTPCReconstructor.cxx:165
 AliTPCReconstructor.cxx:166
 AliTPCReconstructor.cxx:167
 AliTPCReconstructor.cxx:168
 AliTPCReconstructor.cxx:169
 AliTPCReconstructor.cxx:170
 AliTPCReconstructor.cxx:171
 AliTPCReconstructor.cxx:172
 AliTPCReconstructor.cxx:173
 AliTPCReconstructor.cxx:174
 AliTPCReconstructor.cxx:175
 AliTPCReconstructor.cxx:176
 AliTPCReconstructor.cxx:177
 AliTPCReconstructor.cxx:178
 AliTPCReconstructor.cxx:179
 AliTPCReconstructor.cxx:180
 AliTPCReconstructor.cxx:181
 AliTPCReconstructor.cxx:182
 AliTPCReconstructor.cxx:183
 AliTPCReconstructor.cxx:184
 AliTPCReconstructor.cxx:185
 AliTPCReconstructor.cxx:186
 AliTPCReconstructor.cxx:187
 AliTPCReconstructor.cxx:188
 AliTPCReconstructor.cxx:189
 AliTPCReconstructor.cxx:190
 AliTPCReconstructor.cxx:191
 AliTPCReconstructor.cxx:192
 AliTPCReconstructor.cxx:193
 AliTPCReconstructor.cxx:194
 AliTPCReconstructor.cxx:195
 AliTPCReconstructor.cxx:196
 AliTPCReconstructor.cxx:197
 AliTPCReconstructor.cxx:198
 AliTPCReconstructor.cxx:199
 AliTPCReconstructor.cxx:200
 AliTPCReconstructor.cxx:201
 AliTPCReconstructor.cxx:202
 AliTPCReconstructor.cxx:203
 AliTPCReconstructor.cxx:204
 AliTPCReconstructor.cxx:205
 AliTPCReconstructor.cxx:206
 AliTPCReconstructor.cxx:207
 AliTPCReconstructor.cxx:208
 AliTPCReconstructor.cxx:209
 AliTPCReconstructor.cxx:210
 AliTPCReconstructor.cxx:211
 AliTPCReconstructor.cxx:212
 AliTPCReconstructor.cxx:213
 AliTPCReconstructor.cxx:214
 AliTPCReconstructor.cxx:215
 AliTPCReconstructor.cxx:216
 AliTPCReconstructor.cxx:217
 AliTPCReconstructor.cxx:218
 AliTPCReconstructor.cxx:219
 AliTPCReconstructor.cxx:220
 AliTPCReconstructor.cxx:221
 AliTPCReconstructor.cxx:222
 AliTPCReconstructor.cxx:223
 AliTPCReconstructor.cxx:224
 AliTPCReconstructor.cxx:225
 AliTPCReconstructor.cxx:226
 AliTPCReconstructor.cxx:227
 AliTPCReconstructor.cxx:228
 AliTPCReconstructor.cxx:229
 AliTPCReconstructor.cxx:230
 AliTPCReconstructor.cxx:231
 AliTPCReconstructor.cxx:232
 AliTPCReconstructor.cxx:233
 AliTPCReconstructor.cxx:234
 AliTPCReconstructor.cxx:235
 AliTPCReconstructor.cxx:236
 AliTPCReconstructor.cxx:237
 AliTPCReconstructor.cxx:238
 AliTPCReconstructor.cxx:239
 AliTPCReconstructor.cxx:240
 AliTPCReconstructor.cxx:241
 AliTPCReconstructor.cxx:242
 AliTPCReconstructor.cxx:243
 AliTPCReconstructor.cxx:244
 AliTPCReconstructor.cxx:245
 AliTPCReconstructor.cxx:246
 AliTPCReconstructor.cxx:247
 AliTPCReconstructor.cxx:248
 AliTPCReconstructor.cxx:249
 AliTPCReconstructor.cxx:250
 AliTPCReconstructor.cxx:251
 AliTPCReconstructor.cxx:252
 AliTPCReconstructor.cxx:253
 AliTPCReconstructor.cxx:254
 AliTPCReconstructor.cxx:255
 AliTPCReconstructor.cxx:256
 AliTPCReconstructor.cxx:257
 AliTPCReconstructor.cxx:258
 AliTPCReconstructor.cxx:259
 AliTPCReconstructor.cxx:260
 AliTPCReconstructor.cxx:261
 AliTPCReconstructor.cxx:262
 AliTPCReconstructor.cxx:263
 AliTPCReconstructor.cxx:264
 AliTPCReconstructor.cxx:265
 AliTPCReconstructor.cxx:266
 AliTPCReconstructor.cxx:267
 AliTPCReconstructor.cxx:268
 AliTPCReconstructor.cxx:269
 AliTPCReconstructor.cxx:270
 AliTPCReconstructor.cxx:271
 AliTPCReconstructor.cxx:272
 AliTPCReconstructor.cxx:273
 AliTPCReconstructor.cxx:274
 AliTPCReconstructor.cxx:275
 AliTPCReconstructor.cxx:276
 AliTPCReconstructor.cxx:277
 AliTPCReconstructor.cxx:278
 AliTPCReconstructor.cxx:279
 AliTPCReconstructor.cxx:280
 AliTPCReconstructor.cxx:281
 AliTPCReconstructor.cxx:282
 AliTPCReconstructor.cxx:283
 AliTPCReconstructor.cxx:284
 AliTPCReconstructor.cxx:285
 AliTPCReconstructor.cxx:286
 AliTPCReconstructor.cxx:287
 AliTPCReconstructor.cxx:288
 AliTPCReconstructor.cxx:289
 AliTPCReconstructor.cxx:290
 AliTPCReconstructor.cxx:291
 AliTPCReconstructor.cxx:292
 AliTPCReconstructor.cxx:293
 AliTPCReconstructor.cxx:294
 AliTPCReconstructor.cxx:295
 AliTPCReconstructor.cxx:296
 AliTPCReconstructor.cxx:297
 AliTPCReconstructor.cxx:298
 AliTPCReconstructor.cxx:299
 AliTPCReconstructor.cxx:300
 AliTPCReconstructor.cxx:301
 AliTPCReconstructor.cxx:302
 AliTPCReconstructor.cxx:303
 AliTPCReconstructor.cxx:304
 AliTPCReconstructor.cxx:305
 AliTPCReconstructor.cxx:306
 AliTPCReconstructor.cxx:307
 AliTPCReconstructor.cxx:308
 AliTPCReconstructor.cxx:309
 AliTPCReconstructor.cxx:310
 AliTPCReconstructor.cxx:311
 AliTPCReconstructor.cxx:312
 AliTPCReconstructor.cxx:313
 AliTPCReconstructor.cxx:314
 AliTPCReconstructor.cxx:315
 AliTPCReconstructor.cxx:316
 AliTPCReconstructor.cxx:317
 AliTPCReconstructor.cxx:318
 AliTPCReconstructor.cxx:319
 AliTPCReconstructor.cxx:320
 AliTPCReconstructor.cxx:321
 AliTPCReconstructor.cxx:322
 AliTPCReconstructor.cxx:323
 AliTPCReconstructor.cxx:324
 AliTPCReconstructor.cxx:325
 AliTPCReconstructor.cxx:326
 AliTPCReconstructor.cxx:327
 AliTPCReconstructor.cxx:328
 AliTPCReconstructor.cxx:329
 AliTPCReconstructor.cxx:330
 AliTPCReconstructor.cxx:331
 AliTPCReconstructor.cxx:332
 AliTPCReconstructor.cxx:333
 AliTPCReconstructor.cxx:334
 AliTPCReconstructor.cxx:335
 AliTPCReconstructor.cxx:336
 AliTPCReconstructor.cxx:337
 AliTPCReconstructor.cxx:338
 AliTPCReconstructor.cxx:339
 AliTPCReconstructor.cxx:340
 AliTPCReconstructor.cxx:341
 AliTPCReconstructor.cxx:342
 AliTPCReconstructor.cxx:343
 AliTPCReconstructor.cxx:344
 AliTPCReconstructor.cxx:345
 AliTPCReconstructor.cxx:346
 AliTPCReconstructor.cxx:347
 AliTPCReconstructor.cxx:348
 AliTPCReconstructor.cxx:349
 AliTPCReconstructor.cxx:350
 AliTPCReconstructor.cxx:351
 AliTPCReconstructor.cxx:352
 AliTPCReconstructor.cxx:353
 AliTPCReconstructor.cxx:354
 AliTPCReconstructor.cxx:355
 AliTPCReconstructor.cxx:356
 AliTPCReconstructor.cxx:357
 AliTPCReconstructor.cxx:358
 AliTPCReconstructor.cxx:359
 AliTPCReconstructor.cxx:360
 AliTPCReconstructor.cxx:361
 AliTPCReconstructor.cxx:362
 AliTPCReconstructor.cxx:363