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.                  *
 **************************************************************************/

////////////////////////////////////////////////////////////////////////////////
//                                                                            //
// AliTPCComposedCorrection class                                             //
//                                                                            //
// This class is creating a correction that is composed out of smaller        //
// corrections.                                                               //
// There are two ways the sub-corrections can be combined into this one:      //
// 1. kParallel: All corrections are applied at the given position x and      //
//    the dx terms are summed up (this commutes).                             //
// 2. kQueue: The corrections are called in order. The first one at the       //
//    given position x resulting in dx1, the second one is called at          //
//    the corrected position (x+dx1) resulting in dx2, the third one          //
//    is then called at position (x+dx1+dx2) and so forth. dx=dx1+dx2+...     //
//    is returned.                                                            //
// 3. kQueueResidual: like kQueue with the exception that in case of          //
//    a positive weight the 'Distortion' is called and in case of a negative  //
//    weight the 'Correction' is called, where the absolute of the weight     //
//    will be applied to the correction
// For the inverse of the correction this is taken into account by reversing  //
// the order the corrections are applied in the kQueue case (no issue for     //
// kParallel).                                                                //
//                                                                            //
// date: 27/04/2010                                                           //
// Authors: Magnus Mager, Stefan Rossegger, Jim Thomas                       //
//                                                                            //
// Example usage:                                                             //
//                                                                            //
//  AliMagF mag("mag","mag");                                                 //
//  AliTPCExBBShape exb;             // B field shape distortions             //
//  exb.SetBField(&mag);                                                      //
//                                                                            //
//  AliTPCExBTwist twist;            // ExB Twist distortions                 //
//  twist.SetXTwist(0.001);                                                   //
//                                                                            //
//   TObjArray cs;  cs.Add(&exb); cs.Add(&twist);                             //
//                                                                            //
//  AliTPCComposedCorrection cc;                                              //
//  cc.SetCorrections(&cs);                                                   //
//  cc.SetOmegaTauT1T2(wt,T1,T2);                                             //
//  cc.Print("DA");                                                           //               
//  cc.CreateHistoDRPhiinZR(0,100,100)->Draw("surf2");                        //
////////////////////////////////////////////////////////////////////////////////


#include <TCollection.h>
#include <TTimeStamp.h>
#include <TIterator.h>
#include <TMath.h>
#include "AliLog.h"

#include "AliTPCComposedCorrection.h"


AliTPCComposedCorrection::AliTPCComposedCorrection() 
  : AliTPCCorrection("composed_correction",
		     "composition of corrections"),
    fCorrections(0),
    fMode(kParallel),
    fWeights(0)  // weights of corrections
{
  //
  // default constructor
  //
}

AliTPCComposedCorrection::AliTPCComposedCorrection(TCollection *corrections,
						   AliTPCComposedCorrection::CompositionType mode)
  : AliTPCCorrection("composed_correction",
		     "composition of corrections"),
    fCorrections(corrections),
    fMode(mode),
    fWeights(0) //weights of correction
{
  //
  // Constructor that defines the set of corrections, this one is composed of.
  //
}

AliTPCComposedCorrection::~AliTPCComposedCorrection() {
  // 
  // destructor
  //
  if (!fCorrections) {
    AliInfo("No Correction-models were set: can not delete them");
  } else {
    TIterator *i=fCorrections->MakeIterator();
    AliTPCCorrection *c;
    while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
      delete c;
    }
    delete i;
  }
  if (fWeights) delete fWeights;
}


Bool_t AliTPCComposedCorrection::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
  //
  // Add correction  - better name needed (left/right) - for the moment I assumme they commute
  // Why not to just use array of corrections - CPU consideration
  // Assumptions:
  //  - origin of distortion/correction are additive
  //  - corrections/distortion are small   and they commute
  //  - only correction ot the same type supported 
  const Int_t knCorr=100;
  if (corr==NULL) {
    AliError("Zerro pointer - correction");
    return kFALSE;
  }
  if (!fCorrections) fCorrections= new TObjArray(knCorr);
  // 1.) Case of  Composed correction
  AliTPCComposedCorrection * corrC = dynamic_cast<AliTPCComposedCorrection *>(corr);
  if (corrC){
    Int_t ncorrs= corrC->fCorrections->GetEntries();
    Bool_t isOK=kTRUE;
    for (Int_t icorr=0; icorr<ncorrs; icorr++){
      isOK&=AddCorrectionCompact(corrC->GetSubCorrection(icorr),weight*(*((corrC->fWeights)))[icorr]);
    }
    return isOK;
  }
  // 2.) Find subcorrection of the same type
  AliTPCCorrection * toAdd=0;
  Int_t ncorr=fCorrections->GetEntries();
  for (Int_t icorr=0; icorr<ncorr; icorr++){
    if (GetSubCorrection(icorr)==NULL)  continue;
    if (GetSubCorrection(icorr)->IsA()==corr->IsA())  toAdd=GetSubCorrection(icorr);
  }
  // 3.) create of givent type  if does not exist 
  if (toAdd==NULL){
    toAdd= (AliTPCCorrection*)((corr->IsA())->New());
    fCorrections->Add(toAdd);
  }
  // 4.) add to object of given type 
  return toAdd->AddCorrectionCompact(corr, weight);
}


AliTPCCorrection * AliTPCComposedCorrection::GetSubCorrection(Int_t ipos){
  //
  //
  //
  TObjArray *arr = (TObjArray*)fCorrections;
  return (AliTPCCorrection *)arr->At(ipos);
}

AliTPCCorrection * AliTPCComposedCorrection::GetSubCorrection(const char *cname){
  //
  //
  //
  TCollection *arr = fCorrections;
  return (AliTPCCorrection *)arr->FindObject(cname);
}



void AliTPCComposedCorrection::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
  //
  // This applies all correction and the specified manner (see general
  // class description for details).
  //

  if (!fCorrections) {
    AliInfo("No Corrections-models were set: can not calculate distortions");
    return;
  }
  TIterator *i=fCorrections->MakeIterator();
  AliTPCCorrection *c;
  Int_t weightIndex=0;
  switch (fMode) {
  case kParallel:
    Float_t dxi[3];
    for (int j=0;j<3;++j) dx[j]=0.;
    while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
      c->GetCorrection(x,roc,dxi);
      Double_t w=1;
      if (fWeights) w=(*fWeights)[weightIndex++];
      for (Int_t j=0;j<3;++j) dx[j]+=w*dxi[j];
    }
    break;
  case kQueue:
    Float_t xi[3];
    for (Int_t j=0;j<3;++j) xi[j]=x[j];
    while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
      c->GetCorrection(xi,roc,dx);
      Double_t w=1;
      if (fWeights) w=(*fWeights)[weightIndex++];
      for (Int_t j=0;j<3;++j) xi[j]+=w*dx[j];
    }
    for (Int_t j=0;j<3;++j) dx[j]=xi[j]-x[j];
    break;
  case kQueueResidual:
    //TODO: for the moment assume inverse of distortion
    //      check if this is what is desired
    GetDistortion(x,roc,dx);
    for (Int_t j=0;j<3;++j) dx[j]*=-1.;
    break;
  }
  delete i;
}

void AliTPCComposedCorrection::GetDistortion(const Float_t x[],const Short_t roc,Float_t dx[]) {
  //
  // This applies all distortions and the specified manner (see general
  // class descxiption for details).
  //

  if (!fCorrections) {
    AliInfo("No Corrections-models were set: can not calculate distortions");
    return;
  }
  
  if (fMode==kQueueResidual && !fWeights) {
    AliInfo("kQueueResidual mode was selected but no weights were given. Switching to kQueue instead.");
    fMode=kQueue;
  }
  
  TIterator *i=fCorrections->MakeReverseIterator();
  AliTPCCorrection *c;
  Int_t weightIndex=0;
  switch (fMode) {
  case kParallel:
    Float_t dxi[3];
    for (int j=0;j<3;++j) dx[j]=0.;
    while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
      c->GetDistortion(x,roc,dxi);
      Double_t w=1;
      if (fWeights) w=(*fWeights)[weightIndex++];
      for (Int_t j=0;j<3;++j) dx[j]+=w*dxi[j];
    }
    break;
  case kQueue:
    Float_t xi[3];
    for (Int_t j=0;j<3;++j) xi[j]=x[j];
    while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
      c->GetDistortion(xi,roc,dx);
      Double_t w=1;
      if (fWeights) w=(*fWeights)[weightIndex++];
      for (Int_t j=0;j<3;++j) xi[j]+=w*dx[j];
    }
    for (Int_t j=0;j<3;++j) dx[j]=xi[j]-x[j];
    break;
  case kQueueResidual:
    Float_t xi2[3];
    for (Int_t j=0;j<3;++j) xi2[j]=x[j];
    while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
      Double_t w=(*fWeights)[weightIndex++];
      if (w>0) c->GetDistortion(xi2,roc,dx);
      else c->GetCorrection(xi2,roc,dx);
      for (Int_t j=0;j<3;++j) xi2[j]+=TMath::Abs(w)*dx[j];
    }
    for (Int_t j=0;j<3;++j) dx[j]=xi2[j]-x[j];
    break;
  }
  delete i;
}


void AliTPCComposedCorrection::Print(Option_t* option) const {
  //
  // Print function to check which correction classes are used 
  // option=="d" prints details regarding the setted magnitude 
  // option=="a" prints the C0 and C1 coefficents for calibration purposes
  //

  printf("Composed TPC spacepoint correction \"%s\" -- composed of:\n",GetTitle());
  TString opt = option; opt.ToLower();
  Int_t in=1;
  if (!fCorrections) {
    printf("   - composed correction is empty!\n");
    return;
  }
  TIterator *i=fCorrections->MakeIterator();
  AliTPCCorrection *c;
  while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
    if (opt.Contains("d")) {
      printf("\n");
      printf("%d. %s\t%s\n",in,c->GetTitle(), c->GetName());
      c->Print(option);
    } else {
      printf("%d. %s\t%s\n",in,c->GetTitle(), c->GetName());
    }
    ++in;
  }
  if (in==1) printf("  Info: The correction compound is empty: No corrections set\n");
  delete i;
}


void AliTPCComposedCorrection::Init() {
  //
  // Initialization funtion (not used at the moment)
  //
  if (!fCorrections) {
    AliInfo("No Correction-models were set");
    return;
  }
  TIterator *i=fCorrections->MakeIterator();
  AliTPCCorrection *c;
  while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) 
    c->Init();
  delete i;
  
}

void AliTPCComposedCorrection::Update(const TTimeStamp &timeStamp) {
  //
  // Update function 
  //
  if (!fCorrections) {
    AliInfo("No Correction-models were set");
    return;
  }

  TIterator *i=fCorrections->MakeIterator();
  AliTPCCorrection *c;
  while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) 
    c->Update(timeStamp);
  delete i;
 
}



void AliTPCComposedCorrection::SetOmegaTauT1T2(Float_t omegaTau,Float_t t1,Float_t t2) {
  //
  // Gives the possibility to set the OmegaTau plus Tensor corrections T1 and T2 (effective omega Tau)
  // to each subcorrection (since they might become event specific due to changing drift velocity)
  //
  // The omegaTau comes idealy from the Database, since it is a function of drift velocity, B and E field 
  // e.g. omegaTau  = -10.0 * Bz * vdrift / Ez ; // with Bz in kG and Ez in V/cm
  //      omegaTau  = -0.325 for Bz=5kG, Ez=400V/cm and vdrift = 2.6cm/muSec
  // The T1 and T2 tensors were measured in a dedicated calibration run
  //
  // Note: overwrites previously set values!
  // 

  if (!fCorrections) {
    AliInfo("No Correction-models were set");
    return;
  }

  TIterator *i=fCorrections->MakeIterator();
  AliTPCCorrection *c;
  while (0!=(c=dynamic_cast<AliTPCCorrection*>(i->Next()))) {
    c->SetOmegaTauT1T2(omegaTau,t1,t2);
  }
  delete i;
}

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