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

///////////////////////////////////////////////////////////////////////////////
//
//     Data container for relative ITS-TPC alignment analysis
//     Holds an array of AliRelAlignerKalman objects
//     and takes care of merging when processing data in parallel
//
//     Origin: Mikolaj Krzewicki, Nikhef, Mikolaj.Krzewicki@cern.ch
//
//////////////////////////////////////////////////////////////////////////////

#include <TNamed.h>
#include <TGraphErrors.h>
#include <TAxis.h>
#include <TTree.h>
#include <TMath.h>
#include <TObjArray.h>
#include <TCollection.h>
#include <AliESDEvent.h>
#include <AliRelAlignerKalman.h>
#include "AliRelAlignerKalmanArray.h"
#include "TList.h"
#include "TBrowser.h"

ClassImp(AliRelAlignerKalmanArray)

//______________________________________________________________________________
AliRelAlignerKalmanArray::AliRelAlignerKalmanArray():
    TNamed("alignerArray","array of aligners"),
    fT0(0),
    fTimebinWidth(0),
    fSize(0),
    fOutRejSigmaOnMerge(10.),
    fOutRejSigmaOnSmooth(1.),
    fAlignerTemplate(),
    fPArray(NULL),
    fListOfGraphs(NULL)
{
  //ctor
}

//______________________________________________________________________________
AliRelAlignerKalmanArray::AliRelAlignerKalmanArray(Int_t t0, Int_t tend, Int_t slotwidth):
    TNamed("alignerArray","array of aligners"),
    fT0(t0),
    fTimebinWidth(slotwidth),
    fSize(0),
    fOutRejSigmaOnMerge(10.),
    fOutRejSigmaOnSmooth(1.),
    fAlignerTemplate(),
    fPArray(NULL),
    fListOfGraphs(new TList)
{
  //ctor
  if (slotwidth==0) fSize = 1;
  else fSize = (tend-t0)/slotwidth;
  fPArray = new AliRelAlignerKalman* [fSize];
  if (fPArray) for (Int_t i=0;i<fSize;i++){fPArray[i]=NULL;}//fill with zeros
  else fSize=0;
  fListOfGraphs->SetName("graphs");
  fListOfGraphs->SetOwner(kTRUE);
}

//______________________________________________________________________________
AliRelAlignerKalmanArray::AliRelAlignerKalmanArray( const AliRelAlignerKalmanArray& in):
    TNamed(in.GetName(), in.GetTitle()),
    fT0(in.fT0),
    fTimebinWidth(in.fTimebinWidth),
    fSize(in.fSize),
    fOutRejSigmaOnMerge(in.fOutRejSigmaOnMerge),
    fOutRejSigmaOnSmooth(in.fOutRejSigmaOnSmooth),
    fAlignerTemplate(in.fAlignerTemplate),
    fPArray(NULL),
    fListOfGraphs(new TList)
{
  //copy ctor
  fPArray = new AliRelAlignerKalman* [fSize];
  if (!fPArray) {fSize=0;return;} //if fail
  for (Int_t i=0;i<fSize;i++)
  {
    if (in.fPArray[i]) fPArray[i]=new AliRelAlignerKalman(*(in.fPArray[i]));
    else fPArray[i]=NULL;
  }
  fListOfGraphs->SetName("graphs");
  fListOfGraphs->SetOwner(kTRUE);
}

//______________________________________________________________________________
AliRelAlignerKalmanArray::~AliRelAlignerKalmanArray()
{
  //dtor
  ClearContents();
  delete [] fPArray;
  delete fListOfGraphs;
}

//______________________________________________________________________________
void AliRelAlignerKalmanArray::SetupArray(Int_t t0, Int_t tend, Int_t slotwidth)
{
  //setup array - clears old content
  ClearContents();
  fT0=t0;
  fTimebinWidth=slotwidth;
  Int_t tmpsize;
  if (slotwidth==0) tmpsize = 1;
  else tmpsize = (tend-t0)/slotwidth;
  if (tmpsize!=fSize)
  {
    delete [] fPArray;
    fPArray=new AliRelAlignerKalman* [tmpsize];
    if (fPArray) fSize=tmpsize;
    else fSize=0;
  }
  for (Int_t i=0;i<fSize;i++){fPArray[i]=NULL;}//fill with zeros
}

//______________________________________________________________________________
AliRelAlignerKalman* AliRelAlignerKalmanArray::GetAlignerTemplate()
{
  //get aligner template
  return &fAlignerTemplate;
}

//______________________________________________________________________________
Long64_t AliRelAlignerKalmanArray::Merge( TCollection* list )
{
  //Merge all the arrays
  //tlihe merge is vertical, meaning matching entries in tree are merged

  AliRelAlignerKalmanArray *arrayFromList;
  if (!list) return 0;
  TIter next(list);
  while ( (arrayFromList = dynamic_cast<AliRelAlignerKalmanArray*>(next())) )
  {
    if (fT0 != arrayFromList->fT0) continue;
    if (fTimebinWidth != arrayFromList->fTimebinWidth) continue;
    if (fSize != arrayFromList->fSize) continue;

    for (Int_t i=0; i<GetSize(); i++)
    {
      AliRelAlignerKalman* a1 = fPArray[i];
      AliRelAlignerKalman* a2 = arrayFromList->fPArray[i];
      if (a1 && a2)
      {
        a1->SetRejectOutliers(kFALSE);
        a1->SetOutRejSigma(fOutRejSigmaOnMerge); a1->Merge(a2);
      }
      else
        if (!a1 && a2) fPArray[i] = new AliRelAlignerKalman(*a2);
    }
  }
  fListOfGraphs->Delete();
  return 0;
}

//______________________________________________________________________________
Int_t AliRelAlignerKalmanArray::Timebin( UInt_t timestamp ) const
{
  //calculate binnumber for given timestamp
  if (fTimebinWidth==0) return 0;
  Int_t slot = (timestamp-fT0)/fTimebinWidth; //it's all integers!
  return slot;
}

//______________________________________________________________________________
AliRelAlignerKalman* AliRelAlignerKalmanArray::GetAligner(UInt_t ts)
{
  //get the aligner for specified timestamp
  Int_t tb = Timebin(ts);
  if (tb<0) return NULL;
  if (tb>=fSize) return NULL;
  if (!fPArray) return NULL;
  if (!fPArray[tb])
  {
    fPArray[tb] = new AliRelAlignerKalman( fAlignerTemplate );
    fPArray[tb]->SetTimeStamp(fT0+tb*fTimebinWidth);
  }
  return fPArray[tb];
}

//______________________________________________________________________________
AliRelAlignerKalman* AliRelAlignerKalmanArray::GetAligner(AliESDEvent* event)
{
  //get the aligner for this event and set relevant info
  AliRelAlignerKalman* a = GetAligner(event->GetTimeStamp());
  if (a) a->SetRunNumber(event->GetRunNumber());
  if (a) a->SetMagField(event->GetMagneticField());
  return a;
}

//______________________________________________________________________________
AliRelAlignerKalmanArray& AliRelAlignerKalmanArray::operator=(const AliRelAlignerKalmanArray& in)
{
  //assignment operator
  if(&in == this) return *this;
  TNamed::operator=(in);

  if (fSize!=in.fSize)
  {
    //if sizes different, delete curent and make a new one with new size
    ClearContents();
    delete [] fPArray;
    fPArray = new AliRelAlignerKalman* [in.fSize];
  }

  fOutRejSigmaOnMerge=in.fOutRejSigmaOnMerge;
  fOutRejSigmaOnSmooth=in.fOutRejSigmaOnSmooth;

  if (!fPArray) fSize=0;
  else fSize = in.fSize;
  fTimebinWidth = in.fTimebinWidth;
  fT0 = in.fT0;
  for (Int_t i=0;i<fSize;i++)
  {
    if (in.fPArray[i]) fPArray[i] = new AliRelAlignerKalman(*(in.fPArray[i]));
    else fPArray[i]=NULL;
  }
  return *this;
}

//______________________________________________________________________________
void AliRelAlignerKalmanArray::ClearContents()
{
  //clear contents of array
  for (Int_t i=0;i<fSize;i++)
  {
    if (fPArray[i]) delete fPArray[i];
  }
}

//______________________________________________________________________________
AliRelAlignerKalman* AliRelAlignerKalmanArray::At( Int_t i ) const
{
  //mimic TObjArray::At( Int_t i )
  if (i>=fSize) {printf("AliRelAlignerKalmanArray::At index %i out of bounds, size=%i\n",i,fSize); return NULL;}
  if (i<0) return NULL;
  return fPArray[i];
}

//______________________________________________________________________________
AliRelAlignerKalman* AliRelAlignerKalmanArray::operator[](Int_t i) const
{
  //mimic TObjArray::operator[](Int_t)
  if (i>=fSize) {printf("AliRelAlignerKalmanArray::operator[] index %i out of bounds, size=%i\n",i,fSize); return NULL;}
  if (i<0) return NULL;
  return fPArray[i];
}

//______________________________________________________________________________
AliRelAlignerKalman*& AliRelAlignerKalmanArray::operator[](Int_t i)
{
  //mimic TObjArray::operator[](Int_t) can be used as lvalue
  if (i>=fSize||i<0) {Error("operator[]", "index %i out of bounds, size=%i\n",i,fSize); return fPArray[0];}
  return fPArray[i];
}

//______________________________________________________________________________
AliRelAlignerKalman* AliRelAlignerKalmanArray::Last() const
{
  //return aligner in last filled slot
  if (fSize==0) return NULL;
  return fPArray[fSize-1];
}

//______________________________________________________________________________
void AliRelAlignerKalmanArray::Print(Option_t* option) const
{
  // print some info
  TString optionStr(option);
  printf( "%s: t0: %i, tend: %i, width: %i, size: %i, entries: %i\n",
          GetName(),
          fT0, (fT0+(fSize)*fTimebinWidth), fTimebinWidth,
          fSize, GetEntries() );
  if (optionStr.Contains("a"))
    for (Int_t i=0; i<fSize; i++)
    {
      AliRelAlignerKalman* al = fPArray[i];
      if (!al) continue;
      al->Print();
    }
}

//______________________________________________________________________________
Int_t AliRelAlignerKalmanArray::GetEntries() const
{
  //get number of filled slots
  if (!fPArray) return 0;
  Int_t entries=0;
  for (Int_t i=0; i<fSize; i++)
  {
    if (fPArray[i]) entries++;
  }
  return entries;
}

//______________________________________________________________________________
void AliRelAlignerKalmanArray::FillTree( TTree* tree ) const
{
  AliRelAlignerKalman* al = NULL;
  tree->Branch("aligner","AliRelAlignerKalman",&al);
  //fill a tree with filled slots
  for (Int_t i=0; i<fSize; i++)
  {
    al = fPArray[i];
    if (al) tree->Fill();
  }
}

//______________________________________________________________________________
TGraphErrors* AliRelAlignerKalmanArray::MakeGraph(Int_t iparam) const
{
  //return a graph for the iparam-th parameter
  if (iparam>8 || iparam<0)
  {
    printf("wrong parameter number. must be from 0-8");
    return NULL;
  }

  Int_t n=GetEntries();
  TVectorD vx(n);
  TVectorD vy(n);
  TVectorD vex(n);
  TVectorD vey(n);
  Int_t entry=0;
  for (Int_t i=0; i<fSize; i++)
  {
    AliRelAlignerKalman* al = fPArray[i];
    if (!al) continue;
    vx(entry) = al->GetTimeStamp();
    vy(entry) = al->GetStateArr()[iparam];
    TMatrixDSym* cm = al->GetStateCov();
    vey(entry) = TMath::Sqrt((*cm)(iparam,iparam));
    entry++;
  }

  TGraphErrors* graph = new TGraphErrors(vx,vy,vex,vey);
  TString graphtitle;
  TString graphtitley;
  switch (iparam)
  {
  case 0:
    graphtitle="rotation \\psi";
    graphtitley="\\psi [rad]";
    break;
  case 1:
    graphtitle="rotation \\theta";
    graphtitley="\\theta [rad]";
    break;
  case 2:
    graphtitle="rotation \\phi";
    graphtitley="\\phi [rad]";
    break;
  case 3:
    graphtitle="shift x";
    graphtitley="x [cm]";
    break;
  case 4:
    graphtitle="shift y";
    graphtitley="y [cm]";
    break;
  case 5:
    graphtitle="shift z";
    graphtitley="z [cm]";
    break;
  case 6:
    graphtitle="TPC vd correction";
    graphtitley="correction factor []";
    break;
  case 7:
    graphtitle="TPC t0 correction";
    graphtitley="t0 correction [\\micros]";
    break;
  case 8:
    graphtitle="TPC dv/dy";
    graphtitley="dv/dy [cm/\\micros/m]";
    break;
  }
  graph->SetName(graphtitle);
  graph->SetTitle(graphtitle);
  TAxis* xas = graph->GetXaxis();
  TAxis* yas = graph->GetYaxis();
  xas->SetTitle("time");
  xas->SetTimeDisplay(1);
  yas->SetTitle(graphtitley);
  return graph;
}

//______________________________________________________________________________
AliRelAlignerKalmanArray* AliRelAlignerKalmanArray::MakeSmoothArray() const
{
  //return a smoothed version of the data
  AliRelAlignerKalmanArray* outputarr = new AliRelAlignerKalmanArray(fT0,
                                        fT0+fSize*fTimebinWidth,fTimebinWidth);

  AliRelAlignerKalman* tmpaligner = outputarr->GetAlignerTemplate();
  tmpaligner->SetOutRejSigma(fOutRejSigmaOnSmooth);
  //copy the first filled slot
  Int_t n=0;
  while (!fPArray[n]) {n++;}
  if (n==fSize) return NULL;
  *tmpaligner   = *fPArray[n];
  if (fPArray[n]->GetNUpdates()>10)
  (*outputarr)[n] = new AliRelAlignerKalman(*(fPArray[n]));
  //for the rest of slots use merge
  for (Int_t i=n+1;i<fSize;i++)
  {
    if (!fPArray[i]) continue;
    PropagateToTime(tmpaligner, fPArray[i]->GetTimeStamp());
    if (tmpaligner->Merge(fPArray[i]))
    (*outputarr)[i] = new AliRelAlignerKalman(*tmpaligner);
    else
    (*outputarr)[i] = new AliRelAlignerKalman(*(fPArray[i]));
  }
  
  tmpaligner->Reset(); //clean up

  return outputarr;
}

//______________________________________________________________________________
void AliRelAlignerKalmanArray::PropagateToTime(AliRelAlignerKalman* al, UInt_t timestamp ) const
{
  //propagate an aligner in time
  TMatrixDSym* cov = al->GetStateCov();
  TMatrixDSym corr(TMatrixDSym::kZero,*cov);
  UInt_t dt = (timestamp>al->GetTimeStamp())?
    timestamp-al->GetTimeStamp():al->GetTimeStamp()-timestamp;
  //the propagation matrix
  //to be added to the covariance matrix of kalman filter
  corr(6,6) = dt*1e-8;  //vdrift
  (*cov) += corr;
}

//______________________________________________________________________________
void AliRelAlignerKalmanArray::Browse(TBrowser* b )
{
   if (!b) return;
   if (!fListOfGraphs) 
   {
     fListOfGraphs=new TList();
     fListOfGraphs->SetName("graphs");
     fListOfGraphs->SetOwner(kTRUE);
   }
   for (Int_t i=0;i<9;i++)
   {
     TGraphErrors* gr = dynamic_cast<TGraphErrors*>(fListOfGraphs->At(i));
     if (gr) continue;
     fListOfGraphs->AddAt(MakeGraph(i),i);
   }
   b->Add(fListOfGraphs);
}

 AliRelAlignerKalmanArray.cxx:1
 AliRelAlignerKalmanArray.cxx:2
 AliRelAlignerKalmanArray.cxx:3
 AliRelAlignerKalmanArray.cxx:4
 AliRelAlignerKalmanArray.cxx:5
 AliRelAlignerKalmanArray.cxx:6
 AliRelAlignerKalmanArray.cxx:7
 AliRelAlignerKalmanArray.cxx:8
 AliRelAlignerKalmanArray.cxx:9
 AliRelAlignerKalmanArray.cxx:10
 AliRelAlignerKalmanArray.cxx:11
 AliRelAlignerKalmanArray.cxx:12
 AliRelAlignerKalmanArray.cxx:13
 AliRelAlignerKalmanArray.cxx:14
 AliRelAlignerKalmanArray.cxx:15
 AliRelAlignerKalmanArray.cxx:16
 AliRelAlignerKalmanArray.cxx:17
 AliRelAlignerKalmanArray.cxx:18
 AliRelAlignerKalmanArray.cxx:19
 AliRelAlignerKalmanArray.cxx:20
 AliRelAlignerKalmanArray.cxx:21
 AliRelAlignerKalmanArray.cxx:22
 AliRelAlignerKalmanArray.cxx:23
 AliRelAlignerKalmanArray.cxx:24
 AliRelAlignerKalmanArray.cxx:25
 AliRelAlignerKalmanArray.cxx:26
 AliRelAlignerKalmanArray.cxx:27
 AliRelAlignerKalmanArray.cxx:28
 AliRelAlignerKalmanArray.cxx:29
 AliRelAlignerKalmanArray.cxx:30
 AliRelAlignerKalmanArray.cxx:31
 AliRelAlignerKalmanArray.cxx:32
 AliRelAlignerKalmanArray.cxx:33
 AliRelAlignerKalmanArray.cxx:34
 AliRelAlignerKalmanArray.cxx:35
 AliRelAlignerKalmanArray.cxx:36
 AliRelAlignerKalmanArray.cxx:37
 AliRelAlignerKalmanArray.cxx:38
 AliRelAlignerKalmanArray.cxx:39
 AliRelAlignerKalmanArray.cxx:40
 AliRelAlignerKalmanArray.cxx:41
 AliRelAlignerKalmanArray.cxx:42
 AliRelAlignerKalmanArray.cxx:43
 AliRelAlignerKalmanArray.cxx:44
 AliRelAlignerKalmanArray.cxx:45
 AliRelAlignerKalmanArray.cxx:46
 AliRelAlignerKalmanArray.cxx:47
 AliRelAlignerKalmanArray.cxx:48
 AliRelAlignerKalmanArray.cxx:49
 AliRelAlignerKalmanArray.cxx:50
 AliRelAlignerKalmanArray.cxx:51
 AliRelAlignerKalmanArray.cxx:52
 AliRelAlignerKalmanArray.cxx:53
 AliRelAlignerKalmanArray.cxx:54
 AliRelAlignerKalmanArray.cxx:55
 AliRelAlignerKalmanArray.cxx:56
 AliRelAlignerKalmanArray.cxx:57
 AliRelAlignerKalmanArray.cxx:58
 AliRelAlignerKalmanArray.cxx:59
 AliRelAlignerKalmanArray.cxx:60
 AliRelAlignerKalmanArray.cxx:61
 AliRelAlignerKalmanArray.cxx:62
 AliRelAlignerKalmanArray.cxx:63
 AliRelAlignerKalmanArray.cxx:64
 AliRelAlignerKalmanArray.cxx:65
 AliRelAlignerKalmanArray.cxx:66
 AliRelAlignerKalmanArray.cxx:67
 AliRelAlignerKalmanArray.cxx:68
 AliRelAlignerKalmanArray.cxx:69
 AliRelAlignerKalmanArray.cxx:70
 AliRelAlignerKalmanArray.cxx:71
 AliRelAlignerKalmanArray.cxx:72
 AliRelAlignerKalmanArray.cxx:73
 AliRelAlignerKalmanArray.cxx:74
 AliRelAlignerKalmanArray.cxx:75
 AliRelAlignerKalmanArray.cxx:76
 AliRelAlignerKalmanArray.cxx:77
 AliRelAlignerKalmanArray.cxx:78
 AliRelAlignerKalmanArray.cxx:79
 AliRelAlignerKalmanArray.cxx:80
 AliRelAlignerKalmanArray.cxx:81
 AliRelAlignerKalmanArray.cxx:82
 AliRelAlignerKalmanArray.cxx:83
 AliRelAlignerKalmanArray.cxx:84
 AliRelAlignerKalmanArray.cxx:85
 AliRelAlignerKalmanArray.cxx:86
 AliRelAlignerKalmanArray.cxx:87
 AliRelAlignerKalmanArray.cxx:88
 AliRelAlignerKalmanArray.cxx:89
 AliRelAlignerKalmanArray.cxx:90
 AliRelAlignerKalmanArray.cxx:91
 AliRelAlignerKalmanArray.cxx:92
 AliRelAlignerKalmanArray.cxx:93
 AliRelAlignerKalmanArray.cxx:94
 AliRelAlignerKalmanArray.cxx:95
 AliRelAlignerKalmanArray.cxx:96
 AliRelAlignerKalmanArray.cxx:97
 AliRelAlignerKalmanArray.cxx:98
 AliRelAlignerKalmanArray.cxx:99
 AliRelAlignerKalmanArray.cxx:100
 AliRelAlignerKalmanArray.cxx:101
 AliRelAlignerKalmanArray.cxx:102
 AliRelAlignerKalmanArray.cxx:103
 AliRelAlignerKalmanArray.cxx:104
 AliRelAlignerKalmanArray.cxx:105
 AliRelAlignerKalmanArray.cxx:106
 AliRelAlignerKalmanArray.cxx:107
 AliRelAlignerKalmanArray.cxx:108
 AliRelAlignerKalmanArray.cxx:109
 AliRelAlignerKalmanArray.cxx:110
 AliRelAlignerKalmanArray.cxx:111
 AliRelAlignerKalmanArray.cxx:112
 AliRelAlignerKalmanArray.cxx:113
 AliRelAlignerKalmanArray.cxx:114
 AliRelAlignerKalmanArray.cxx:115
 AliRelAlignerKalmanArray.cxx:116
 AliRelAlignerKalmanArray.cxx:117
 AliRelAlignerKalmanArray.cxx:118
 AliRelAlignerKalmanArray.cxx:119
 AliRelAlignerKalmanArray.cxx:120
 AliRelAlignerKalmanArray.cxx:121
 AliRelAlignerKalmanArray.cxx:122
 AliRelAlignerKalmanArray.cxx:123
 AliRelAlignerKalmanArray.cxx:124
 AliRelAlignerKalmanArray.cxx:125
 AliRelAlignerKalmanArray.cxx:126
 AliRelAlignerKalmanArray.cxx:127
 AliRelAlignerKalmanArray.cxx:128
 AliRelAlignerKalmanArray.cxx:129
 AliRelAlignerKalmanArray.cxx:130
 AliRelAlignerKalmanArray.cxx:131
 AliRelAlignerKalmanArray.cxx:132
 AliRelAlignerKalmanArray.cxx:133
 AliRelAlignerKalmanArray.cxx:134
 AliRelAlignerKalmanArray.cxx:135
 AliRelAlignerKalmanArray.cxx:136
 AliRelAlignerKalmanArray.cxx:137
 AliRelAlignerKalmanArray.cxx:138
 AliRelAlignerKalmanArray.cxx:139
 AliRelAlignerKalmanArray.cxx:140
 AliRelAlignerKalmanArray.cxx:141
 AliRelAlignerKalmanArray.cxx:142
 AliRelAlignerKalmanArray.cxx:143
 AliRelAlignerKalmanArray.cxx:144
 AliRelAlignerKalmanArray.cxx:145
 AliRelAlignerKalmanArray.cxx:146
 AliRelAlignerKalmanArray.cxx:147
 AliRelAlignerKalmanArray.cxx:148
 AliRelAlignerKalmanArray.cxx:149
 AliRelAlignerKalmanArray.cxx:150
 AliRelAlignerKalmanArray.cxx:151
 AliRelAlignerKalmanArray.cxx:152
 AliRelAlignerKalmanArray.cxx:153
 AliRelAlignerKalmanArray.cxx:154
 AliRelAlignerKalmanArray.cxx:155
 AliRelAlignerKalmanArray.cxx:156
 AliRelAlignerKalmanArray.cxx:157
 AliRelAlignerKalmanArray.cxx:158
 AliRelAlignerKalmanArray.cxx:159
 AliRelAlignerKalmanArray.cxx:160
 AliRelAlignerKalmanArray.cxx:161
 AliRelAlignerKalmanArray.cxx:162
 AliRelAlignerKalmanArray.cxx:163
 AliRelAlignerKalmanArray.cxx:164
 AliRelAlignerKalmanArray.cxx:165
 AliRelAlignerKalmanArray.cxx:166
 AliRelAlignerKalmanArray.cxx:167
 AliRelAlignerKalmanArray.cxx:168
 AliRelAlignerKalmanArray.cxx:169
 AliRelAlignerKalmanArray.cxx:170
 AliRelAlignerKalmanArray.cxx:171
 AliRelAlignerKalmanArray.cxx:172
 AliRelAlignerKalmanArray.cxx:173
 AliRelAlignerKalmanArray.cxx:174
 AliRelAlignerKalmanArray.cxx:175
 AliRelAlignerKalmanArray.cxx:176
 AliRelAlignerKalmanArray.cxx:177
 AliRelAlignerKalmanArray.cxx:178
 AliRelAlignerKalmanArray.cxx:179
 AliRelAlignerKalmanArray.cxx:180
 AliRelAlignerKalmanArray.cxx:181
 AliRelAlignerKalmanArray.cxx:182
 AliRelAlignerKalmanArray.cxx:183
 AliRelAlignerKalmanArray.cxx:184
 AliRelAlignerKalmanArray.cxx:185
 AliRelAlignerKalmanArray.cxx:186
 AliRelAlignerKalmanArray.cxx:187
 AliRelAlignerKalmanArray.cxx:188
 AliRelAlignerKalmanArray.cxx:189
 AliRelAlignerKalmanArray.cxx:190
 AliRelAlignerKalmanArray.cxx:191
 AliRelAlignerKalmanArray.cxx:192
 AliRelAlignerKalmanArray.cxx:193
 AliRelAlignerKalmanArray.cxx:194
 AliRelAlignerKalmanArray.cxx:195
 AliRelAlignerKalmanArray.cxx:196
 AliRelAlignerKalmanArray.cxx:197
 AliRelAlignerKalmanArray.cxx:198
 AliRelAlignerKalmanArray.cxx:199
 AliRelAlignerKalmanArray.cxx:200
 AliRelAlignerKalmanArray.cxx:201
 AliRelAlignerKalmanArray.cxx:202
 AliRelAlignerKalmanArray.cxx:203
 AliRelAlignerKalmanArray.cxx:204
 AliRelAlignerKalmanArray.cxx:205
 AliRelAlignerKalmanArray.cxx:206
 AliRelAlignerKalmanArray.cxx:207
 AliRelAlignerKalmanArray.cxx:208
 AliRelAlignerKalmanArray.cxx:209
 AliRelAlignerKalmanArray.cxx:210
 AliRelAlignerKalmanArray.cxx:211
 AliRelAlignerKalmanArray.cxx:212
 AliRelAlignerKalmanArray.cxx:213
 AliRelAlignerKalmanArray.cxx:214
 AliRelAlignerKalmanArray.cxx:215
 AliRelAlignerKalmanArray.cxx:216
 AliRelAlignerKalmanArray.cxx:217
 AliRelAlignerKalmanArray.cxx:218
 AliRelAlignerKalmanArray.cxx:219
 AliRelAlignerKalmanArray.cxx:220
 AliRelAlignerKalmanArray.cxx:221
 AliRelAlignerKalmanArray.cxx:222
 AliRelAlignerKalmanArray.cxx:223
 AliRelAlignerKalmanArray.cxx:224
 AliRelAlignerKalmanArray.cxx:225
 AliRelAlignerKalmanArray.cxx:226
 AliRelAlignerKalmanArray.cxx:227
 AliRelAlignerKalmanArray.cxx:228
 AliRelAlignerKalmanArray.cxx:229
 AliRelAlignerKalmanArray.cxx:230
 AliRelAlignerKalmanArray.cxx:231
 AliRelAlignerKalmanArray.cxx:232
 AliRelAlignerKalmanArray.cxx:233
 AliRelAlignerKalmanArray.cxx:234
 AliRelAlignerKalmanArray.cxx:235
 AliRelAlignerKalmanArray.cxx:236
 AliRelAlignerKalmanArray.cxx:237
 AliRelAlignerKalmanArray.cxx:238
 AliRelAlignerKalmanArray.cxx:239
 AliRelAlignerKalmanArray.cxx:240
 AliRelAlignerKalmanArray.cxx:241
 AliRelAlignerKalmanArray.cxx:242
 AliRelAlignerKalmanArray.cxx:243
 AliRelAlignerKalmanArray.cxx:244
 AliRelAlignerKalmanArray.cxx:245
 AliRelAlignerKalmanArray.cxx:246
 AliRelAlignerKalmanArray.cxx:247
 AliRelAlignerKalmanArray.cxx:248
 AliRelAlignerKalmanArray.cxx:249
 AliRelAlignerKalmanArray.cxx:250
 AliRelAlignerKalmanArray.cxx:251
 AliRelAlignerKalmanArray.cxx:252
 AliRelAlignerKalmanArray.cxx:253
 AliRelAlignerKalmanArray.cxx:254
 AliRelAlignerKalmanArray.cxx:255
 AliRelAlignerKalmanArray.cxx:256
 AliRelAlignerKalmanArray.cxx:257
 AliRelAlignerKalmanArray.cxx:258
 AliRelAlignerKalmanArray.cxx:259
 AliRelAlignerKalmanArray.cxx:260
 AliRelAlignerKalmanArray.cxx:261
 AliRelAlignerKalmanArray.cxx:262
 AliRelAlignerKalmanArray.cxx:263
 AliRelAlignerKalmanArray.cxx:264
 AliRelAlignerKalmanArray.cxx:265
 AliRelAlignerKalmanArray.cxx:266
 AliRelAlignerKalmanArray.cxx:267
 AliRelAlignerKalmanArray.cxx:268
 AliRelAlignerKalmanArray.cxx:269
 AliRelAlignerKalmanArray.cxx:270
 AliRelAlignerKalmanArray.cxx:271
 AliRelAlignerKalmanArray.cxx:272
 AliRelAlignerKalmanArray.cxx:273
 AliRelAlignerKalmanArray.cxx:274
 AliRelAlignerKalmanArray.cxx:275
 AliRelAlignerKalmanArray.cxx:276
 AliRelAlignerKalmanArray.cxx:277
 AliRelAlignerKalmanArray.cxx:278
 AliRelAlignerKalmanArray.cxx:279
 AliRelAlignerKalmanArray.cxx:280
 AliRelAlignerKalmanArray.cxx:281
 AliRelAlignerKalmanArray.cxx:282
 AliRelAlignerKalmanArray.cxx:283
 AliRelAlignerKalmanArray.cxx:284
 AliRelAlignerKalmanArray.cxx:285
 AliRelAlignerKalmanArray.cxx:286
 AliRelAlignerKalmanArray.cxx:287
 AliRelAlignerKalmanArray.cxx:288
 AliRelAlignerKalmanArray.cxx:289
 AliRelAlignerKalmanArray.cxx:290
 AliRelAlignerKalmanArray.cxx:291
 AliRelAlignerKalmanArray.cxx:292
 AliRelAlignerKalmanArray.cxx:293
 AliRelAlignerKalmanArray.cxx:294
 AliRelAlignerKalmanArray.cxx:295
 AliRelAlignerKalmanArray.cxx:296
 AliRelAlignerKalmanArray.cxx:297
 AliRelAlignerKalmanArray.cxx:298
 AliRelAlignerKalmanArray.cxx:299
 AliRelAlignerKalmanArray.cxx:300
 AliRelAlignerKalmanArray.cxx:301
 AliRelAlignerKalmanArray.cxx:302
 AliRelAlignerKalmanArray.cxx:303
 AliRelAlignerKalmanArray.cxx:304
 AliRelAlignerKalmanArray.cxx:305
 AliRelAlignerKalmanArray.cxx:306
 AliRelAlignerKalmanArray.cxx:307
 AliRelAlignerKalmanArray.cxx:308
 AliRelAlignerKalmanArray.cxx:309
 AliRelAlignerKalmanArray.cxx:310
 AliRelAlignerKalmanArray.cxx:311
 AliRelAlignerKalmanArray.cxx:312
 AliRelAlignerKalmanArray.cxx:313
 AliRelAlignerKalmanArray.cxx:314
 AliRelAlignerKalmanArray.cxx:315
 AliRelAlignerKalmanArray.cxx:316
 AliRelAlignerKalmanArray.cxx:317
 AliRelAlignerKalmanArray.cxx:318
 AliRelAlignerKalmanArray.cxx:319
 AliRelAlignerKalmanArray.cxx:320
 AliRelAlignerKalmanArray.cxx:321
 AliRelAlignerKalmanArray.cxx:322
 AliRelAlignerKalmanArray.cxx:323
 AliRelAlignerKalmanArray.cxx:324
 AliRelAlignerKalmanArray.cxx:325
 AliRelAlignerKalmanArray.cxx:326
 AliRelAlignerKalmanArray.cxx:327
 AliRelAlignerKalmanArray.cxx:328
 AliRelAlignerKalmanArray.cxx:329
 AliRelAlignerKalmanArray.cxx:330
 AliRelAlignerKalmanArray.cxx:331
 AliRelAlignerKalmanArray.cxx:332
 AliRelAlignerKalmanArray.cxx:333
 AliRelAlignerKalmanArray.cxx:334
 AliRelAlignerKalmanArray.cxx:335
 AliRelAlignerKalmanArray.cxx:336
 AliRelAlignerKalmanArray.cxx:337
 AliRelAlignerKalmanArray.cxx:338
 AliRelAlignerKalmanArray.cxx:339
 AliRelAlignerKalmanArray.cxx:340
 AliRelAlignerKalmanArray.cxx:341
 AliRelAlignerKalmanArray.cxx:342
 AliRelAlignerKalmanArray.cxx:343
 AliRelAlignerKalmanArray.cxx:344
 AliRelAlignerKalmanArray.cxx:345
 AliRelAlignerKalmanArray.cxx:346
 AliRelAlignerKalmanArray.cxx:347
 AliRelAlignerKalmanArray.cxx:348
 AliRelAlignerKalmanArray.cxx:349
 AliRelAlignerKalmanArray.cxx:350
 AliRelAlignerKalmanArray.cxx:351
 AliRelAlignerKalmanArray.cxx:352
 AliRelAlignerKalmanArray.cxx:353
 AliRelAlignerKalmanArray.cxx:354
 AliRelAlignerKalmanArray.cxx:355
 AliRelAlignerKalmanArray.cxx:356
 AliRelAlignerKalmanArray.cxx:357
 AliRelAlignerKalmanArray.cxx:358
 AliRelAlignerKalmanArray.cxx:359
 AliRelAlignerKalmanArray.cxx:360
 AliRelAlignerKalmanArray.cxx:361
 AliRelAlignerKalmanArray.cxx:362
 AliRelAlignerKalmanArray.cxx:363
 AliRelAlignerKalmanArray.cxx:364
 AliRelAlignerKalmanArray.cxx:365
 AliRelAlignerKalmanArray.cxx:366
 AliRelAlignerKalmanArray.cxx:367
 AliRelAlignerKalmanArray.cxx:368
 AliRelAlignerKalmanArray.cxx:369
 AliRelAlignerKalmanArray.cxx:370
 AliRelAlignerKalmanArray.cxx:371
 AliRelAlignerKalmanArray.cxx:372
 AliRelAlignerKalmanArray.cxx:373
 AliRelAlignerKalmanArray.cxx:374
 AliRelAlignerKalmanArray.cxx:375
 AliRelAlignerKalmanArray.cxx:376
 AliRelAlignerKalmanArray.cxx:377
 AliRelAlignerKalmanArray.cxx:378
 AliRelAlignerKalmanArray.cxx:379
 AliRelAlignerKalmanArray.cxx:380
 AliRelAlignerKalmanArray.cxx:381
 AliRelAlignerKalmanArray.cxx:382
 AliRelAlignerKalmanArray.cxx:383
 AliRelAlignerKalmanArray.cxx:384
 AliRelAlignerKalmanArray.cxx:385
 AliRelAlignerKalmanArray.cxx:386
 AliRelAlignerKalmanArray.cxx:387
 AliRelAlignerKalmanArray.cxx:388
 AliRelAlignerKalmanArray.cxx:389
 AliRelAlignerKalmanArray.cxx:390
 AliRelAlignerKalmanArray.cxx:391
 AliRelAlignerKalmanArray.cxx:392
 AliRelAlignerKalmanArray.cxx:393
 AliRelAlignerKalmanArray.cxx:394
 AliRelAlignerKalmanArray.cxx:395
 AliRelAlignerKalmanArray.cxx:396
 AliRelAlignerKalmanArray.cxx:397
 AliRelAlignerKalmanArray.cxx:398
 AliRelAlignerKalmanArray.cxx:399
 AliRelAlignerKalmanArray.cxx:400
 AliRelAlignerKalmanArray.cxx:401
 AliRelAlignerKalmanArray.cxx:402
 AliRelAlignerKalmanArray.cxx:403
 AliRelAlignerKalmanArray.cxx:404
 AliRelAlignerKalmanArray.cxx:405
 AliRelAlignerKalmanArray.cxx:406
 AliRelAlignerKalmanArray.cxx:407
 AliRelAlignerKalmanArray.cxx:408
 AliRelAlignerKalmanArray.cxx:409
 AliRelAlignerKalmanArray.cxx:410
 AliRelAlignerKalmanArray.cxx:411
 AliRelAlignerKalmanArray.cxx:412
 AliRelAlignerKalmanArray.cxx:413
 AliRelAlignerKalmanArray.cxx:414
 AliRelAlignerKalmanArray.cxx:415
 AliRelAlignerKalmanArray.cxx:416
 AliRelAlignerKalmanArray.cxx:417
 AliRelAlignerKalmanArray.cxx:418
 AliRelAlignerKalmanArray.cxx:419
 AliRelAlignerKalmanArray.cxx:420
 AliRelAlignerKalmanArray.cxx:421
 AliRelAlignerKalmanArray.cxx:422
 AliRelAlignerKalmanArray.cxx:423
 AliRelAlignerKalmanArray.cxx:424
 AliRelAlignerKalmanArray.cxx:425
 AliRelAlignerKalmanArray.cxx:426
 AliRelAlignerKalmanArray.cxx:427
 AliRelAlignerKalmanArray.cxx:428
 AliRelAlignerKalmanArray.cxx:429
 AliRelAlignerKalmanArray.cxx:430
 AliRelAlignerKalmanArray.cxx:431
 AliRelAlignerKalmanArray.cxx:432
 AliRelAlignerKalmanArray.cxx:433
 AliRelAlignerKalmanArray.cxx:434
 AliRelAlignerKalmanArray.cxx:435
 AliRelAlignerKalmanArray.cxx:436
 AliRelAlignerKalmanArray.cxx:437
 AliRelAlignerKalmanArray.cxx:438
 AliRelAlignerKalmanArray.cxx:439
 AliRelAlignerKalmanArray.cxx:440
 AliRelAlignerKalmanArray.cxx:441
 AliRelAlignerKalmanArray.cxx:442
 AliRelAlignerKalmanArray.cxx:443
 AliRelAlignerKalmanArray.cxx:444
 AliRelAlignerKalmanArray.cxx:445
 AliRelAlignerKalmanArray.cxx:446
 AliRelAlignerKalmanArray.cxx:447
 AliRelAlignerKalmanArray.cxx:448
 AliRelAlignerKalmanArray.cxx:449
 AliRelAlignerKalmanArray.cxx:450
 AliRelAlignerKalmanArray.cxx:451
 AliRelAlignerKalmanArray.cxx:452
 AliRelAlignerKalmanArray.cxx:453
 AliRelAlignerKalmanArray.cxx:454
 AliRelAlignerKalmanArray.cxx:455
 AliRelAlignerKalmanArray.cxx:456
 AliRelAlignerKalmanArray.cxx:457
 AliRelAlignerKalmanArray.cxx:458
 AliRelAlignerKalmanArray.cxx:459
 AliRelAlignerKalmanArray.cxx:460
 AliRelAlignerKalmanArray.cxx:461
 AliRelAlignerKalmanArray.cxx:462
 AliRelAlignerKalmanArray.cxx:463
 AliRelAlignerKalmanArray.cxx:464
 AliRelAlignerKalmanArray.cxx:465
 AliRelAlignerKalmanArray.cxx:466