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

////////////////////////////////////////////////////////////////////////////
//                                                                        //
//  Creates and handles digits from TRD hits                              //
//                                                                        //
//  Authors: C. Blume (blume@ikf.uni-frankfurt.de)                        //
//           C. Lippmann                                                  //
//           B. Vulpescu                                                  //
//                                                                        //
//  The following effects are included:                                   //
//      - Diffusion                                                       //
//      - ExB effects                                                     //
//      - Gas gain including fluctuations                                 //
//      - Pad-response (simple Gaussian approximation)                    //
//      - Time-response                                                   //
//      - Electronics noise                                               //
//      - Electronics gain                                                //
//      - Digitization                                                    //
//      - Zero suppression                                                //
//                                                                        //
////////////////////////////////////////////////////////////////////////////

#include <TGeoManager.h>
#include <TList.h>
#include <TMath.h>
#include <TRandom.h>
#include <TTree.h>

#include "AliRun.h"
#include "AliMC.h"
#include "AliRunLoader.h"
#include "AliLoader.h"
#include "AliConfig.h"
#include "AliDigitizationInput.h"
#include "AliRunLoader.h"
#include "AliLoader.h"
#include "AliLog.h"

#include "AliTRD.h"
#include "AliTRDhit.h"
#include "AliTRDdigitizer.h"
#include "AliTRDarrayDictionary.h"
#include "AliTRDarrayADC.h"
#include "AliTRDarraySignal.h"
#include "AliTRDdigitsManager.h"
#include "AliTRDgeometry.h"
#include "AliTRDpadPlane.h"
#include "AliTRDcalibDB.h"
#include "AliTRDSimParam.h"
#include "AliTRDCommonParam.h"
#include "AliTRDfeeParam.h"
#include "AliTRDmcmSim.h"
#include "AliTRDdigitsParam.h"

#include "Cal/AliTRDCalROC.h"
#include "Cal/AliTRDCalDet.h"
#include "Cal/AliTRDCalOnlineGainTableROC.h"

ClassImp(AliTRDdigitizer)

//_____________________________________________________________________________
AliTRDdigitizer::AliTRDdigitizer()
  :AliDigitizer()
  ,fRunLoader(0)
  ,fDigitsManager(0)
  ,fSDigitsManager(0)
  ,fSDigitsManagerList(0)
  ,fTRD(0)
  ,fGeo(0)
  ,fMcmSim(new AliTRDmcmSim)
  ,fEvent(0)
  ,fMasks(0)
  ,fCompress(kTRUE)
  ,fSDigits(kFALSE)
  ,fMergeSignalOnly(kFALSE)
{
  //
  // AliTRDdigitizer default constructor
  //
  
}

//_____________________________________________________________________________
AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)              
  :AliDigitizer(name,title)
  ,fRunLoader(0)
  ,fDigitsManager(0)
  ,fSDigitsManager(0)
  ,fSDigitsManagerList(0)
  ,fTRD(0)
  ,fGeo(0)
  ,fMcmSim(new AliTRDmcmSim)
  ,fEvent(0)
  ,fMasks(0)
  ,fCompress(kTRUE)
  ,fSDigits(kFALSE)
  ,fMergeSignalOnly(kFALSE)
{
  //
  // AliTRDdigitizer constructor
  //

}

//_____________________________________________________________________________
AliTRDdigitizer::AliTRDdigitizer(AliDigitizationInput* digInput
                               , const Text_t *name, const Text_t *title)
  :AliDigitizer(digInput,name,title)
  ,fRunLoader(0)
  ,fDigitsManager(0)
  ,fSDigitsManager(0)
  ,fSDigitsManagerList(0)
  ,fTRD(0)
  ,fGeo(0)
  ,fMcmSim(new AliTRDmcmSim)
  ,fEvent(0)
  ,fMasks(0)
  ,fCompress(kTRUE)
  ,fSDigits(kFALSE)
  ,fMergeSignalOnly(kFALSE)
{
  //
  // AliTRDdigitizer constructor
  //

}

//_____________________________________________________________________________
AliTRDdigitizer::AliTRDdigitizer(AliDigitizationInput* digInput)
  :AliDigitizer(digInput,"AliTRDdigitizer","TRD digitizer")
  ,fRunLoader(0)
  ,fDigitsManager(0)
  ,fSDigitsManager(0)
  ,fSDigitsManagerList(0)
  ,fTRD(0)
  ,fGeo(0)
  ,fMcmSim(new AliTRDmcmSim)
  ,fEvent(0)
  ,fMasks(0)
  ,fCompress(kTRUE)
  ,fSDigits(kFALSE)
  ,fMergeSignalOnly(kFALSE)
{
  //
  // AliTRDdigitizer constructor
  //

}

//_____________________________________________________________________________
AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
  :AliDigitizer(d)
  ,fRunLoader(0)
  ,fDigitsManager(0)
  ,fSDigitsManager(0)
  ,fSDigitsManagerList(0)
  ,fTRD(0)
  ,fGeo(0)
  ,fMcmSim(new AliTRDmcmSim)
  ,fEvent(0)
  ,fMasks(0)
  ,fCompress(d.fCompress)
  ,fSDigits(d.fSDigits)
  ,fMergeSignalOnly(d.fMergeSignalOnly)
{
  //
  // AliTRDdigitizer copy constructor
  //

}

//_____________________________________________________________________________
AliTRDdigitizer::~AliTRDdigitizer()
{
  //
  // AliTRDdigitizer destructor
  //

  delete fDigitsManager;
  fDigitsManager      = 0;

  // s-digitsmanager will be deleted via list
  fSDigitsManager     = 0;
  if (fSDigitsManagerList) {
    fSDigitsManagerList->Delete();
    delete fSDigitsManagerList;
  }
  fSDigitsManagerList = 0;

  delete [] fMasks;
  fMasks       = 0;

  delete fMcmSim;
  fMcmSim = 0;

  delete fGeo;
  fGeo = 0;

}

//_____________________________________________________________________________
AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
{
  //
  // Assignment operator
  //

  if (this != &d) {
    ((AliTRDdigitizer &) d).Copy(*this);
  }

  return *this;

}

//_____________________________________________________________________________
void AliTRDdigitizer::Copy(TObject &d) const
{
  //
  // Copy function
  //

  ((AliTRDdigitizer &) d).fRunLoader          = 0;
  ((AliTRDdigitizer &) d).fDigitsManager      = 0;
  ((AliTRDdigitizer &) d).fSDigitsManager     = 0;
  ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
  ((AliTRDdigitizer &) d).fTRD                = 0;
  ((AliTRDdigitizer &) d).fGeo                = 0;
  ((AliTRDdigitizer &) d).fEvent              = 0;
  ((AliTRDdigitizer &) d).fMasks              = 0;
  ((AliTRDdigitizer &) d).fCompress           = fCompress;
  ((AliTRDdigitizer &) d).fSDigits            = fSDigits;
  ((AliTRDdigitizer &) d).fMergeSignalOnly    = fMergeSignalOnly;

}

//_____________________________________________________________________________
void AliTRDdigitizer::Digitize(const Option_t* option)
{
  //
  // Executes the merging
  //

  Int_t iInput;

  AliTRDdigitsManager *sdigitsManager;

  TString optionString = option;
  if (optionString.Contains("deb")) {
    AliLog::SetClassDebugLevel("AliTRDdigitizer",1);
    AliInfo("Called with debug option");
  }

  // The AliRoot file is already connected by the manager
  AliRunLoader *inrl = 0x0;
  
  if (gAlice) {
    AliDebug(1,"AliRun object found on file.");
  }
  else {
    inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(0));
    inrl->LoadgAlice();
    gAlice = inrl->GetAliRun();
    if (!gAlice) {
      AliError("Could not find AliRun object.");
      return;
    }
  }
                                                                           
  Int_t nInput = fDigInput->GetNinputs();
  fMasks       = new Int_t[nInput];
  for (iInput = 0; iInput < nInput; iInput++) {
    fMasks[iInput] = fDigInput->GetMask(iInput);
  }

  //
  // Initialization
  //

  AliRunLoader *orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());

  if (InitDetector()) {

    AliLoader *ogime = orl->GetLoader("TRDLoader");

    TTree *tree = 0;
    if (fSDigits) { 
      // If we produce SDigits
      tree = ogime->TreeS();
      if (!tree) {
	ogime->MakeTree("S");
	tree = ogime->TreeS();
      }
    }
    else {
      // If we produce Digits
      tree = ogime->TreeD();
      if (!tree) {
	ogime->MakeTree("D");
        tree = ogime->TreeD();
      }
    }

    MakeBranch(tree);

  }
 
  for (iInput = 0; iInput < nInput; iInput++) {

    AliDebug(1,Form("Add input stream %d",iInput));

    // Check if the input tree exists
    inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
    AliLoader *gime = inrl->GetLoader("TRDLoader");

    TTree *treees = gime->TreeS();
    if (treees == 0x0) {
      if (gime->LoadSDigits()) {
        AliError(Form("Error Occured while loading S. Digits for input %d.",iInput));
        return;
      }
      treees = gime->TreeS();
    }
    
    if (treees == 0x0) {
      AliError(Form("Input stream %d does not exist",iInput));
      return;
    } 

    // Read the s-digits via digits manager
    sdigitsManager = new AliTRDdigitsManager();
    sdigitsManager->SetSDigits(kTRUE);
    
    AliRunLoader *rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
    AliLoader *gimme = rl->GetLoader("TRDLoader");
    if (!gimme->TreeS()) 
      {
	gimme->LoadSDigits();
      }

    sdigitsManager->ReadDigits(gimme->TreeS());
   
    // Add the s-digits to the input list 
    AddSDigitsManager(sdigitsManager);

  }

  // Convert the s-digits to normal digits
  AliDebug(1,"Do the conversion");
  SDigits2Digits();

  // Store the digits
  AliDebug(1,"Write the digits");
  fDigitsManager->WriteDigits();

  // Write parameters
  orl->CdGAFile();

  // Clean up
  DeleteSDigitsManager();

  AliDebug(1,"Done");

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
{
  //
  // Opens a ROOT-file with TRD-hits and reads in the hit-tree
  //
  // Connect the AliRoot file containing Geometry, Kine, and Hits
  //  

  TString evfoldname = AliConfig::GetDefaultEventFolderName();

  fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
  if (!fRunLoader) {
    fRunLoader = AliRunLoader::Open(file,evfoldname,"UPDATE");
  }  
  if (!fRunLoader) {
    AliError(Form("Can not open session for file %s.",file));
    return kFALSE;
  }
   
  if (!fRunLoader->GetAliRun()) {
    fRunLoader->LoadgAlice();
  }
  gAlice = fRunLoader->GetAliRun();
  
  if (gAlice) {
    AliDebug(1,"AliRun object found on file.");
  }
  else {
    AliError("Could not find AliRun object.");
    return kFALSE;
  }

  fEvent = nEvent;

  AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
  if (!loader) {
    AliError("Can not get TRD loader from Run Loader");
    return kFALSE;
  }
  
  if (InitDetector()) {
    TTree *tree = 0;
    if (fSDigits) { 
      // If we produce SDigits
      tree = loader->TreeS();
      if (!tree) {
        loader->MakeTree("S");
        tree = loader->TreeS();
      }
    }
    else {
      // If we produce Digits
      tree = loader->TreeD();
      if (!tree) {
        loader->MakeTree("D");
        tree = loader->TreeD();
      }
    }
    return MakeBranch(tree);
  }
  else {
    return kFALSE;
  }

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::Open(AliRunLoader * const runLoader, Int_t nEvent)
{
  //
  // Opens a ROOT-file with TRD-hits and reads in the hit-tree
  //
  // Connect the AliRoot file containing Geometry, Kine, and Hits
  //  

  fRunLoader = runLoader;
  if (!fRunLoader) {
    AliError("RunLoader does not exist");
    return kFALSE;
  }
   
  if (!fRunLoader->GetAliRun()) {
    fRunLoader->LoadgAlice();
  }
  gAlice = fRunLoader->GetAliRun();
  
  if (gAlice) {
    AliDebug(1,"AliRun object found on file.");
  }
  else {
    AliError("Could not find AliRun object.");
    return kFALSE;
  }

  fEvent = nEvent;

  AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
  if (!loader) {
    AliError("Can not get TRD loader from Run Loader");
    return kFALSE;
  }
  
  if (InitDetector()) {
    TTree *tree = 0;
    if (fSDigits) { 
      // If we produce SDigits
      tree = loader->TreeS();
      if (!tree) {
        loader->MakeTree("S");
        tree = loader->TreeS();
      }
    }
    else {
      // If we produce Digits
      tree = loader->TreeD();
      if (!tree) {
        loader->MakeTree("D");
        tree = loader->TreeD();
      }
    }
    return MakeBranch(tree);
  }
  else {
    return kFALSE;
  }

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::InitDetector()
{
  //
  // Sets the pointer to the TRD detector and the geometry
  //

  // Get the pointer to the detector class and check for version 1
  fTRD = (AliTRD *) gAlice->GetDetector("TRD");
  if (!fTRD) {
    AliFatal("No TRD module found");
    exit(1);
  }
  if (fTRD->IsVersion() != 1) {
    AliFatal("TRD must be version 1 (slow simulator)");
    exit(1);
  }

  // Get the geometry
  fGeo = new AliTRDgeometry();

  // Create a digits manager
  if (fDigitsManager) {
    delete fDigitsManager;
  }
  fDigitsManager = new AliTRDdigitsManager();
  fDigitsManager->SetSDigits(fSDigits);
  fDigitsManager->CreateArrays();
  fDigitsManager->SetEvent(fEvent);

  // The list for the input s-digits manager to be merged
  if (fSDigitsManagerList) {
    fSDigitsManagerList->Delete();
  } 
  else {
    fSDigitsManagerList = new TList();
  }

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::MakeBranch(TTree *tree) const
{
  // 
  // Create the branches for the digits array
  //

  return fDigitsManager->MakeBranch(tree);

}

//_____________________________________________________________________________
void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
{
  //
  // Add a digits manager for s-digits to the input list.
  //

  fSDigitsManagerList->Add(man);

}

//_____________________________________________________________________________
void AliTRDdigitizer::DeleteSDigitsManager()
{
  //
  // Removes digits manager from the input list.
  //

  fSDigitsManagerList->Delete();

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::MakeDigits()
{
  //
  // Creates digits.
  //

  AliDebug(1,"Start creating digits");

  if (!fGeo) {
    AliError("No geometry defined");
    return kFALSE;
  }

  AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
  if (!calibration) {
    AliFatal("Could not get calibration object");
    return kFALSE;
  }

  const Int_t kNdet  = AliTRDgeometry::Ndet();

  Float_t **hits = new Float_t*[kNdet];
  Int_t    *nhit = new Int_t[kNdet];
  memset(nhit,0,kNdet*sizeof(Int_t));

  AliTRDarraySignal *signals = 0x0;

  // Check the number of time bins from simParam against OCDB,
  // if OCDB value is not supposed to be used.
  // As default, the value from OCDB is taken
  if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
    if (calibration->GetNumberOfTimeBinsDCS() != AliTRDSimParam::Instance()->GetNTimeBins()) {
      AliWarning(Form("Number of time bins is different to OCDB value [SIM=%d, OCDB=%d]"
                     ,AliTRDSimParam::Instance()->GetNTimeBins()
                     ,calibration->GetNumberOfTimeBinsDCS()));
    }
    // Save the values for the raw data headers
    fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
  }
  else {
    // Get the OCDB values
    Int_t nTB = calibration->GetNumberOfTimeBinsDCS();
    if (nTB < 0) { // Currently -1 gets returned for "undefined" and "mixed",
                   // one might go back to -1 undefined and -2 mixed?
      AliError("No useful DCS information available for this run! Using standard values.");
      // // We fall back to the standard OCDB object, 
      // // cache the current run number..
      // Long64_t run = calibration->GetRun();
      // calibration->SetRun(0);
      // nTB = calibration->GetNumberOfTimeBinsDCS();
      // // ..to set it again
      // calibration->SetRun(run);
      // // If there's no standard OCDB object, we can still fail
      // if (nTB < 0) {
      // 	AliFatal("No standard object found in the OCDB!");
      // }
      nTB = AliTRDSimParam::Instance()->GetNTimeBins();
    }
    // Save the values for the raw data headers
    fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(nTB);
  }

  // Save the values for the raw data headers
  fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
 
  // Sort all hits according to detector number
  if (!SortHits(hits,nhit)) {
    AliError("Sorting hits failed");
    delete [] hits;
    delete [] nhit;
    return kFALSE;
  }

  // Loop through all detectors
  for (Int_t det = 0; det < kNdet; det++) {

    // Detectors that are switched off, not installed, etc.
    if ((!calibration->IsChamberNoData(det))    &&
        ( fGeo->ChamberInGeometry(det))         &&
        (nhit[det] > 0)) {

      signals = new AliTRDarraySignal();

      // Convert the hits of the current detector to detector signals
      if (!ConvertHits(det,hits[det],nhit[det],signals)) {
	AliError(Form("Conversion of hits failed for detector=%d",det));
        delete [] hits;
        delete [] nhit;
        delete signals;
        signals = 0x0;
        return kFALSE;
      }

      // Convert the detector signals to digits or s-digits
      if (!ConvertSignals(det,signals)) {
	AliError(Form("Conversion of signals failed for detector=%d",det));
        delete [] hits;
        delete [] nhit;
        delete signals;
        signals = 0x0;
	return kFALSE;
      }

      // Delete the signals array
      delete signals;
      signals = 0x0;

    } // if: detector status

    delete [] hits[det];

  } // for: detector

  if (!fSDigits) {
    if (AliDataLoader *trklLoader 
          = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
      if (trklLoader->Tree())
        trklLoader->WriteData("OVERWRITE");
    }
  }
    
  delete [] hits;
  delete [] nhit;

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::SortHits(Float_t **hits, Int_t *nhit)
{
  //
  // Read all the hits and sorts them according to detector number
  // in the output array <hits>.
  //

  AliDebug(1,"Start sorting hits");

  const Int_t kNdet = AliTRDgeometry::Ndet();
  // Size of the hit vector
  const Int_t kNhit = 6;

  Float_t *xyz      = 0;
  Int_t    nhitTrk  = 0;

  Int_t   *lhit     = new Int_t[kNdet];
  memset(lhit,0,kNdet*sizeof(Int_t));

  for (Int_t det = 0; det < kNdet; det++) {
    hits[det] = 0x0;
  }

  AliLoader *gimme = fRunLoader->GetLoader("TRDLoader");
  if (!gimme->TreeH()) {
    gimme->LoadHits();
  }
  TTree *hitTree = gimme->TreeH();
  if (hitTree == 0x0) {
    AliError("Can not get TreeH");
    delete [] lhit;
    return kFALSE;
  }
  fTRD->SetTreeAddress();

  // Get the number of entries in the hit tree
  // (Number of primary particles creating a hit somewhere)
  Int_t nTrk = (Int_t) hitTree->GetEntries();
  AliDebug(1,Form("Found %d tracks",nTrk));

  // Loop through all the tracks in the tree
  for (Int_t iTrk = 0; iTrk < nTrk; iTrk++) {

    gAlice->GetMCApp()->ResetHits();
    hitTree->GetEvent(iTrk);

    if (!fTRD->Hits()) {
      AliError(Form("No hits array for track = %d",iTrk));
      continue;
    }

    // Number of hits for this track
    nhitTrk = fTRD->Hits()->GetEntriesFast();

    Int_t hitCnt = 0;
    // Loop through the TRD hits
    AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
    while (hit) {

      hitCnt++;

      // Don't analyze test hits
      if (((Int_t) hit->GetCharge()) != 0) {

        Int_t   trk  = hit->Track();
        Int_t   det  = hit->GetDetector();
        Int_t   q    = hit->GetCharge();
        Float_t x    = hit->X();
        Float_t y    = hit->Y();
        Float_t z    = hit->Z();
        Float_t time = hit->GetTime();

        if (nhit[det] == lhit[det]) {
          // Inititialization of new detector
          xyz        = new Float_t[kNhit*(nhitTrk+lhit[det])];
          if (hits[det]) {
            memcpy(xyz,hits[det],sizeof(Float_t)*kNhit*lhit[det]);
            delete [] hits[det];
	  }
          lhit[det] += nhitTrk;
          hits[det]  = xyz;
	}
        else {
          xyz        = hits[det];
	}
        xyz[nhit[det]*kNhit+0] = x;
        xyz[nhit[det]*kNhit+1] = y;
        xyz[nhit[det]*kNhit+2] = z;
        xyz[nhit[det]*kNhit+3] = q;
        xyz[nhit[det]*kNhit+4] = trk;
        xyz[nhit[det]*kNhit+5] = time;
        nhit[det]++;

      } // if: charge != 0

      hit = (AliTRDhit *) fTRD->NextHit();   

    } // for: hits of one track

  } // for: tracks

  delete [] lhit;

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::ConvertHits(Int_t det
                                  , const Float_t * const hits
                                  , Int_t nhit
                                  , AliTRDarraySignal *signals)
{
  //
  // Converts the detectorwise sorted hits to detector signals
  //

  AliDebug(1,Form("Start converting hits for detector=%d (nhits=%d)",det,nhit));

  // Number of pads included in the pad response
  const Int_t kNpad      = 3;
  // Number of track dictionary arrays
  const Int_t kNdict     = AliTRDdigitsManager::kNDict;
  // Size of the hit vector
  const Int_t kNhit      = 6;

  // Width of the amplification region
  const Float_t kAmWidth = AliTRDgeometry::AmThick();
  // Width of the drift region
  const Float_t kDrWidth = AliTRDgeometry::DrThick();
  // Drift + amplification region 
  const Float_t kDrMin   =          - 0.5 * kAmWidth;
  const Float_t kDrMax   = kDrWidth + 0.5 * kAmWidth;
  
  Int_t    iPad          = 0;
  Int_t    dict          = 0;
  Int_t    timeBinTRFend = 1;

  Double_t pos[3];
  Double_t loc[3];
  Double_t padSignal[kNpad];
  Double_t signalOld[kNpad];

  AliTRDarrayDictionary *dictionary[kNdict];

  AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
  AliTRDCommonParam *commonParam = AliTRDCommonParam::Instance();
  AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();

  if (!commonParam) {
    AliFatal("Could not get common parameterss");
    return kFALSE;
  }
  if (!simParam) {
    AliFatal("Could not get simulation parameters");
    return kFALSE;
  }  
  if (!calibration) {
    AliFatal("Could not get calibration object");  
    return kFALSE;
  }

  // Get the detector wise calibration objects
  AliTRDCalROC       *calVdriftROC      = 0;
  Float_t             calVdriftDetValue = 0.0;
  const AliTRDCalDet *calVdriftDet      = calibration->GetVdriftDet();  
  AliTRDCalROC       *calT0ROC          = 0;
  Float_t             calT0DetValue     = 0.0;
  const AliTRDCalDet *calT0Det          = calibration->GetT0Det();  
  Double_t            calExBDetValue    = 0.0;
  const AliTRDCalDet *calExBDet         = calibration->GetExBDet();

  if (simParam->TRFOn()) {
    timeBinTRFend = ((Int_t) (simParam->GetTRFhi() 
                  * commonParam->GetSamplingFrequency())) - 1;
  }

  Int_t   nTimeTotal   = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
  Float_t samplingRate = commonParam->GetSamplingFrequency();
  Float_t elAttachProp = simParam->GetElAttachProp() / 100.0; 

  AliTRDpadPlane *padPlane = fGeo->GetPadPlane(det);
  Int_t   layer   = fGeo->GetLayer(det);         //update
  Float_t row0    = padPlane->GetRow0ROC();
  Int_t   nRowMax = padPlane->GetNrows();
  Int_t   nColMax = padPlane->GetNcols();

  // Create a new array for the signals
  signals->Allocate(nRowMax,nColMax,nTimeTotal);

  // Create a new array for the dictionary
  for (dict = 0; dict < kNdict; dict++) {       
    dictionary[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
    dictionary[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
  }      

  // Loop through the hits in this detector
  for (Int_t hit = 0; hit < nhit; hit++) {

    pos[0]          = hits[hit*kNhit+0];
    pos[1]          = hits[hit*kNhit+1];
    pos[2]          = hits[hit*kNhit+2];
    Float_t q       = hits[hit*kNhit+3];
    Float_t hittime = hits[hit*kNhit+5];
    Int_t   track   = ((Int_t) hits[hit*kNhit+4]);

    Int_t   inDrift = 1;

    // Find the current volume with the geo manager
    gGeoManager->SetCurrentPoint(pos);
    gGeoManager->FindNode();      
    if (strstr(gGeoManager->GetPath(),"/UK")) {
      inDrift = 0;
    }

    // Get the calibration objects
    calVdriftROC      = calibration->GetVdriftROC(det);
    calVdriftDetValue = calVdriftDet->GetValue(det);
    calT0ROC          = calibration->GetT0ROC(det);
    calT0DetValue     = calT0Det->GetValue(det);
    calExBDetValue    = calExBDet->GetValue(det);

    // Go to the local coordinate system:
    // loc[0] - col  direction in amplification or driftvolume
    // loc[1] - row  direction in amplification or driftvolume
    // loc[2] - time direction in amplification or driftvolume
    gGeoManager->MasterToLocal(pos,loc);
    if (inDrift) {
      // Relative to middle of amplification region
      loc[2] = loc[2] - kDrWidth/2.0 - kAmWidth/2.0;
    } 

    // The driftlength [cm] (w/o diffusion yet !).
    // It is negative if the hit is between pad plane and anode wires.
    Double_t driftlength = -1.0 * loc[2];

    // Stupid patch to take care of TR photons that are absorbed
    // outside the chamber volume. A real fix would actually need
    // a more clever implementation of the TR hit generation
    if (q < 0.0) {
      if ((loc[1] < padPlane->GetRowEndROC()) ||
          (loc[1] > padPlane->GetRow0ROC())) {
        continue;
      }
      if ((driftlength < kDrMin) ||
          (driftlength > kDrMax)) {
        continue;
      }
    }

    // Get row and col of unsmeared electron to retrieve drift velocity
    // The pad row (z-direction)
    Int_t    rowE         = padPlane->GetPadRowNumberROC(loc[1]);
    if (rowE < 0) {
      continue;
    }
    Double_t rowOffset    = padPlane->GetPadRowOffsetROC(rowE,loc[1]);
    // The pad column (rphi-direction)
    Double_t offsetTilt   = padPlane->GetTiltOffset(rowOffset);
    Int_t    colE         = padPlane->GetPadColNumber(loc[0]+offsetTilt);
    if (colE < 0) {
      continue;	  
    }
    Double_t colOffset    = 0.0;

    // Normalized drift length
    Float_t  driftvelocity  = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
    Double_t absdriftlength = TMath::Abs(driftlength);
    if (commonParam->ExBOn()) {
      absdriftlength /= TMath::Sqrt(1.0 / (1.0 + calExBDetValue*calExBDetValue));
    }

    // Loop over all electrons of this hit
    // TR photons produce hits with negative charge
    Int_t nEl = ((Int_t) TMath::Abs(q));
    for (Int_t iEl = 0; iEl < nEl; iEl++) {

      // Now the real local coordinate system of the ROC
      // column direction: locC
      // row direction:    locR 
      // time direction:   locT
      // locR and locC are identical to the coordinates of the corresponding
      // volumina of the drift or amplification region.
      // locT is defined relative to the wire plane (i.e. middle of amplification
      // region), meaning locT = 0, and is negative for hits coming from the
      // drift region. 
      Double_t locC = loc[0];
      Double_t locR = loc[1];
      Double_t locT = loc[2];

      // Electron attachment
      if (simParam->ElAttachOn()) {
        if (gRandom->Rndm() < (absdriftlength * elAttachProp)) {
          continue;
        }
      }
          
      // Apply the diffusion smearing
      if (simParam->DiffusionOn()) {
        if (!(Diffusion(driftvelocity,absdriftlength,calExBDetValue,locR,locC,locT))) {
          continue;
	}
      }

      // Apply E x B effects (depends on drift direction)
      if (commonParam->ExBOn()) {
        locC = locC + calExBDetValue * driftlength;
      }

      // The electron position after diffusion and ExB in pad coordinates.
      // The pad row (z-direction)
      rowE       = padPlane->GetPadRowNumberROC(locR);
      if (rowE < 0) continue;
      rowOffset  = padPlane->GetPadRowOffsetROC(rowE,locR);

      // The pad column (rphi-direction)
      offsetTilt = padPlane->GetTiltOffset(rowOffset);
      colE       = padPlane->GetPadColNumber(locC+offsetTilt);
      if (colE < 0) continue;         
      colOffset  = padPlane->GetPadColOffset(colE,locC+offsetTilt);

      // Also re-retrieve drift velocity because col and row may have changed
      driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
      Float_t t0    = calT0DetValue     + calT0ROC->GetValue(colE,rowE);

      // Convert the position to drift time [mus], using either constant drift velocity or
      // time structure of drift cells (non-isochronity, GARFIELD calculation).
      // Also add absolute time of hits to take pile-up events into account properly
      Double_t drifttime;
      if (simParam->TimeStructOn()) {
	// Get z-position with respect to anode wire
        Double_t zz  =  row0 - locR + padPlane->GetAnodeWireOffset();
	zz -= ((Int_t)(2 * zz)) / 2.0;
        if (zz > 0.25) {
          zz  = 0.5 - zz;
        }
        // Use drift time map (GARFIELD)
        drifttime = commonParam->TimeStruct(driftvelocity,0.5*kAmWidth-1.0*locT,zz)
                  + hittime;
      } 
      else {
	// Use constant drift velocity
        drifttime = TMath::Abs(locT) / driftvelocity
                  + hittime;
      }

      // Apply the gas gain including fluctuations
      Double_t ggRndm = 0.0;
      do {
        ggRndm = gRandom->Rndm();
      } while (ggRndm <= 0);
      Double_t signal = -(simParam->GetGasGain()) * TMath::Log(ggRndm);

      // Apply the pad response 
      if (simParam->PRFOn()) {
        // The distance of the electron to the center of the pad 
        // in units of pad width
        Double_t dist = (colOffset - 0.5*padPlane->GetColSize(colE))
                      / padPlane->GetColSize(colE);
        // This is a fixed parametrization, i.e. not dependent on
        // calibration values !
        if (!(calibration->PadResponse(signal,dist,layer,padSignal))) continue;
      }
      else {
        padSignal[0] = 0.0;
        padSignal[1] = signal;
        padSignal[2] = 0.0;
      }

      // The time bin (always positive), with t0 distortion
      Double_t timeBinIdeal = drifttime * samplingRate + t0;
      // Protection 
      if (TMath::Abs(timeBinIdeal) > 2*nTimeTotal) {
        timeBinIdeal = 2 * nTimeTotal;
      }
      Int_t    timeBinTruncated = ((Int_t) timeBinIdeal);
      // The distance of the position to the middle of the timebin
      Double_t timeOffset       = ((Float_t) timeBinTruncated 
                                + 0.5 - timeBinIdeal) / samplingRate;
          
      // Sample the time response inside the drift region
      // + additional time bins before and after.
      // The sampling is done always in the middle of the time bin
      for (Int_t iTimeBin = TMath::Max(timeBinTruncated,0)
          ;iTimeBin < TMath::Min(timeBinTruncated+timeBinTRFend,nTimeTotal)
	  ;iTimeBin++) {

        // Apply the time response
        Double_t timeResponse = 1.0;
        Double_t crossTalk    = 0.0;
        Double_t time         = (iTimeBin - timeBinTruncated) / samplingRate + timeOffset;

        if (simParam->TRFOn()) {
          timeResponse = simParam->TimeResponse(time);
        }
        if (simParam->CTOn()) {
          crossTalk    = simParam->CrossTalk(time);
        }

        signalOld[0] = 0.0;
        signalOld[1] = 0.0;
        signalOld[2] = 0.0;

        for (iPad = 0; iPad < kNpad; iPad++) {

          Int_t colPos = colE + iPad - 1;
          if (colPos <        0) continue;
          if (colPos >= nColMax) break;

          // Add the signals
          signalOld[iPad] = signals->GetData(rowE,colPos,iTimeBin);

          if (colPos != colE) {
	    // Cross talk added to non-central pads
            signalOld[iPad] += padSignal[iPad] 
	                     * (timeResponse + crossTalk);
          } 
          else {
	    // W/o cross talk at central pad
            signalOld[iPad] += padSignal[iPad] 
	                     * timeResponse;
          }

          signals->SetData(rowE,colPos,iTimeBin,signalOld[iPad]);

          // Store the track index in the dictionary
          // Note: We store index+1 in order to allow the array to be compressed
	  // Note2: Taking out the +1 in track
          if (signalOld[iPad] > 0.0) { 
            for (dict = 0; dict < kNdict; dict++) {
	      Int_t oldTrack = dictionary[dict]->GetData(rowE,colPos,iTimeBin); 
	      if (oldTrack == track) break;
	      if (oldTrack ==  -1 ) {
		dictionary[dict]->SetData(rowE,colPos,iTimeBin,track);
	        break;
	      }
            }
          }

	} // Loop: pads

      } // Loop: time bins

    } // Loop: electrons of a single hit

  } // Loop: hits

  AliDebug(2,Form("Finished analyzing %d hits",nhit));

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::ConvertSignals(Int_t det, AliTRDarraySignal *signals)
{
  //
  // Convert signals to digits
  //

  AliDebug(1,Form("Start converting the signals for detector %d",det));

  if (fSDigits) {
    // Convert the signal array to s-digits
    if (!Signal2SDigits(det,signals)) {
      return kFALSE;
    }
  }
  else {
    // Convert the signal array to digits
    if (!Signal2ADC(det,signals)) {
      return kFALSE;
    }
    // Run digital processing for digits
    RunDigitalProcessing(det);
  }

  // Compress the arrays
  CompressOutputArrays(det);

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::Signal2ADC(Int_t det, AliTRDarraySignal *signals)
{
  //
  // Converts the sampled electron signals to ADC values for a given chamber
  //

  AliDebug(1,Form("Start converting signals to ADC values for detector=%d",det));

  AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
  if (!calibration) {
    AliFatal("Could not get calibration object");
    return kFALSE;
  }

  AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
  if (!simParam) {
    AliFatal("Could not get simulation parameters");
    return kFALSE;
  }  

  // Converts number of electrons to fC
  const Double_t kEl2fC = 1.602e-19 * 1.0e15; 

  // Coupling factor
  Double_t coupling     = simParam->GetPadCoupling() 
                        * simParam->GetTimeCoupling();
  // Electronics conversion factor
  Double_t convert      = kEl2fC 
                        * simParam->GetChipGain();
  // ADC conversion factor
  Double_t adcConvert   = simParam->GetADCoutRange()
                        / simParam->GetADCinRange();
  // The electronics baseline in mV
  Double_t baseline     = simParam->GetADCbaseline() 
                        / adcConvert;
  // The electronics baseline in electrons
  Double_t baselineEl   = baseline
                        / convert;

  Int_t row  = 0;
  Int_t col  = 0;
  Int_t time = 0;

  Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
  Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
  Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
  if (fSDigitsManager->GetDigitsParam()->GetNTimeBins(det)) {
    nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
  }
  else {
    AliFatal("Could not get number of time bins");
    return kFALSE;
  }

  // The gain factor calibration objects
  const AliTRDCalDet *calGainFactorDet      = calibration->GetGainFactorDet();  
  AliTRDCalROC       *calGainFactorROC      = 0x0;
  Float_t             calGainFactorDetValue = 0.0;

  AliTRDarrayADC     *digits                = 0x0;

  if (!signals) {
    AliError(Form("Signals array for detector %d does not exist\n",det));
    return kFALSE;
  }
  if (signals->HasData()) {
    // Expand the container if neccessary
    signals->Expand();   
  }
  else {
    // Create missing containers
    signals->Allocate(nRowMax,nColMax,nTimeTotal);      
  }

  // Get the container for the digits of this detector
  if (fDigitsManager->HasSDigits()) {
    AliError("Digits manager has s-digits");
    return kFALSE;
  }

  digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
  // Allocate memory space for the digits buffer
  if (!digits->HasData()) {
    digits->Allocate(nRowMax,nColMax,nTimeTotal);
  }

  // Get the calibration objects
  calGainFactorROC      = calibration->GetGainFactorROC(det);
  calGainFactorDetValue = calGainFactorDet->GetValue(det);

  // Create the digits for this chamber
  for (row  = 0; row  <  nRowMax; row++ ) {
    for (col  = 0; col  <  nColMax; col++ ) {

      // halfchamber masking
      Int_t iMcm            = (Int_t)(col/18); // current group of 18 col pads
      Int_t halfchamberside = (iMcm>3 ? 1 : 0); // 0=Aside, 1=Bside
      // Halfchambers that are switched off, masked by calibration
      if (calibration->IsHalfChamberNoData(det, halfchamberside)) 
	continue;

      // Check whether pad is masked
      // Bridged pads are not considered yet!!!
      if (calibration->IsPadMasked(det,col,row) || 
          calibration->IsPadNotConnected(det,col,row)) {
        continue;
      }

      // The gain factors
      Float_t padgain    = calGainFactorDetValue 
                         * calGainFactorROC->GetValue(col,row);
      if (padgain <= 0) {
        AliError(Form("Not a valid gain %f, %d %d %d",padgain,det,col,row));
      }

      for (time = 0; time < nTimeTotal; time++) {

	// Get the signal amplitude
        Float_t signalAmp = signals->GetData(row,col,time);
        // Pad and time coupling
        signalAmp *= coupling;
	// Gain factors
	signalAmp *= padgain;

        // Add the noise, starting from minus ADC baseline in electrons
        signalAmp  = TMath::Max((Double_t) gRandom->Gaus(signalAmp,simParam->GetNoise())
                               ,-baselineEl);

        // Convert to mV
        signalAmp *= convert;
        // Add ADC baseline in mV
        signalAmp += baseline;

 	// Convert to ADC counts. Set the overflow-bit fADCoutRange if the
	// signal is larger than fADCinRange
        Short_t adc  = 0;
        if (signalAmp >= simParam->GetADCinRange()) {
          adc = ((Short_t) simParam->GetADCoutRange());
	}
        else {
          adc = TMath::Nint(signalAmp * adcConvert);
	}

        // Saving all digits
	digits->SetData(row,col,time,adc);

      } // for: time

    } // for: col
  } // for: row

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::Signal2SDigits(Int_t det, AliTRDarraySignal *signals)
{
  //
  // Converts the sampled electron signals to s-digits
  //

  AliDebug(1,Form("Start converting signals to s-digits for detector=%d",det));

  AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
  if (!calibration) {
    AliFatal("Could not get calibration object");
    return kFALSE;
  }

  Int_t row  = 0;
  Int_t col  = 0;
  Int_t time = 0;

  Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
  Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
  Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);

  // Get the container for the digits of this detector
  if (!fDigitsManager->HasSDigits()) {
    AliError("Digits manager has no s-digits");
    return kFALSE;
  }

  AliTRDarraySignal *digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
  // Allocate memory space for the digits buffer
  if (!digits->HasData()) {
    digits->Allocate(nRowMax,nColMax,nTimeTotal);
  }

  // Create the sdigits for this chamber
  for (row  = 0; row  <  nRowMax; row++ ) {
    for (col  = 0; col  <  nColMax; col++ ) {

      // halfchamber masking
      Int_t iMcm            = (Int_t)(col/18); // current group of 18 col pads
      Int_t halfchamberside = (iMcm>3 ? 1 : 0); // 0=Aside, 1=Bside
      // Halfchambers that are switched off, masked by calibration
      if (calibration->IsHalfChamberNoData(det, halfchamberside))
	continue;

      for (time = 0; time < nTimeTotal; time++) {
        digits->SetData(row,col,time,signals->GetData(row,col,time));
      } // for: time
    } // for: col
  } // for: row

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::Digits2SDigits(AliTRDdigitsManager * const manDig
                                     , AliTRDdigitsManager * const manSDig)
{
  //
  // Converts digits into s-digits. Needed for embedding into real data.
  //

  AliDebug(1,"Start converting digits to s-digits");

  if (!fGeo) {
    fGeo = new AliTRDgeometry();
  }

  AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
  if (!calibration) {
    AliFatal("Could not get calibration object");
    return kFALSE;
  }

  AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
  if (!simParam) {
    AliFatal("Could not get simulation parameters");
    return kFALSE;
  }  

  // Converts number of electrons to fC
  const Double_t kEl2fC = 1.602e-19 * 1.0e15; 

  // Coupling factor
  Double_t coupling     = simParam->GetPadCoupling() 
                        * simParam->GetTimeCoupling();
  // Electronics conversion factor
  Double_t convert      = kEl2fC 
                        * simParam->GetChipGain();
  // ADC conversion factor
  Double_t adcConvert   = simParam->GetADCoutRange()
                        / simParam->GetADCinRange();
  // The electronics baseline in mV
  Double_t baseline     = simParam->GetADCbaseline() 
                        / adcConvert;
  // The electronics baseline in electrons
  //Double_t baselineEl   = baseline
  //                      / convert;

  // The gainfactor calibration objects
  // Not used since these digits are supposed to be from real raw data
  //const AliTRDCalDet *calGainFactorDet      = calibration->GetGainFactorDet();  
  //AliTRDCalROC       *calGainFactorROC      = 0;
  //Float_t             calGainFactorDetValue = 0.0;

  Int_t row  = 0;
  Int_t col  = 0;
  Int_t time = 0;

  for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {

    Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
    Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
    Int_t nTimeTotal = manDig->GetDigitsParam()->GetNTimeBins(det);

    // Get the calibration objects
    //calGainFactorROC      = calibration->GetGainFactorROC(det);
    //calGainFactorDetValue = calGainFactorDet->GetValue(det);

    // Get the digits
    AliTRDarrayADC *digits = (AliTRDarrayADC *) manDig->GetDigits(det);

    if (!manSDig->HasSDigits()) {
      AliError("SDigits manager has no s-digits");
      return kFALSE;
    }
    // Get the s-digits
    AliTRDarraySignal     *sdigits = (AliTRDarraySignal *)     manSDig->GetSDigits(det);
    AliTRDarrayDictionary *tracks0 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,0);
    AliTRDarrayDictionary *tracks1 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,1);
    AliTRDarrayDictionary *tracks2 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,2);
    // Allocate memory space for the digits buffer
    sdigits->Allocate(nRowMax,nColMax,nTimeTotal);
    tracks0->Allocate(nRowMax,nColMax,nTimeTotal);
    tracks1->Allocate(nRowMax,nColMax,nTimeTotal);
    tracks2->Allocate(nRowMax,nColMax,nTimeTotal);

    // Keep the digits param
    manSDig->GetDigitsParam()->SetNTimeBinsAll(manDig->GetDigitsParam()->GetNTimeBins(0));
    manSDig->GetDigitsParam()->SetADCbaselineAll(manDig->GetDigitsParam()->GetADCbaseline(0));

    if (digits->HasData()) {

      digits->Expand();

      // Create the sdigits for this chamber
      for (row  = 0; row  <  nRowMax; row++ ) {
        for (col  = 0; col  <  nColMax; col++ ) {

          // The gain factors
          //Float_t padgain = calGainFactorDetValue 
          //                * calGainFactorROC->GetValue(col,row);

          for (time = 0; time < nTimeTotal; time++) {

            Short_t  adcVal = digits->GetData(row,col,time);
            Double_t signal = (Double_t) adcVal;
            // ADC -> signal in mV
            signal /= adcConvert;
            // Subtract baseline in mV
            signal -= baseline;
            // Signal in mV -> signal in #electrons
            signal /= convert;
            // Gain factor
            //signal /= padgain; // Not needed for real data
            // Pad and time coupling
            signal /= coupling;

            sdigits->SetData(row,col,time,signal);
            tracks0->SetData(row,col,time,0);
            tracks1->SetData(row,col,time,0);
            tracks2->SetData(row,col,time,0);

          } // for: time

        } // for: col
      } // for: row
  
    } // if: has data

    sdigits->Compress(0);
    tracks0->Compress();
    tracks1->Compress();
    tracks2->Compress();

    // No compress just remove
    manDig->RemoveDigits(det);
    manDig->RemoveDictionaries(det);      

  } // for: det

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::SDigits2Digits()
{
  //
  // Merges the input s-digits and converts them to normal digits
  //

  if (!MergeSDigits()) {
    return kFALSE;
  }

  return ConvertSDigits();

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::MergeSDigits()
{
  //
  // Merges the input s-digits:
  //   - The amplitude of the different inputs are summed up.
  //   - Of the track IDs from the input dictionaries only one is
  //     kept for each input. This works for maximal 3 different merged inputs.
  //

  // Number of track dictionary arrays
  const Int_t kNDict = AliTRDdigitsManager::kNDict;

  AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
  if (!simParam) {
    AliFatal("Could not get simulation parameters");
    return kFALSE;
  }
  
  AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
  if (!calibration) {
    AliFatal("Could not get calibration object");
    return kFALSE;
  }
  
  Int_t iDict = 0;
  Int_t jDict = 0;

  AliTRDarraySignal     *digitsA;
  AliTRDarraySignal     *digitsB;
  AliTRDarrayDictionary *dictionaryA[kNDict];
  AliTRDarrayDictionary *dictionaryB[kNDict];

  AliTRDdigitsManager   *mergeSDigitsManager = 0x0;
  // Get the first s-digits
  fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
  if (!fSDigitsManager) { 
    AliError("No SDigits manager");
    return kFALSE;
  }

  // Loop through the other sets of s-digits
  mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(fSDigitsManager);

  if (mergeSDigitsManager) {
    AliDebug(1,Form("Merge %d input files.",fSDigitsManagerList->GetSize()));
  }
  else {
    AliDebug(1,"Only one input file.");
  }
  
  Int_t iMerge = 0;

  while (mergeSDigitsManager) {

    iMerge++;

    // Loop through the detectors
    for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {

      Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet);
      if (mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet) != nTimeTotal) {
        AliError(Form("Mismatch in the number of time bins [%d,%d] in detector %d"
                     ,nTimeTotal
  		     ,mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet)
		     ,iDet));
        return kFALSE;
      }

      Int_t nRowMax = fGeo->GetPadPlane(iDet)->GetNrows();
      Int_t nColMax = fGeo->GetPadPlane(iDet)->GetNcols();
	  
      // Loop through the pixels of one detector and add the signals
      digitsA = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(iDet);    
      digitsB = (AliTRDarraySignal *) mergeSDigitsManager->GetSDigits(iDet); 
      digitsA->Expand();  
      if (!digitsA->HasData()) continue;
      digitsB->Expand();    
      if (!digitsB->HasData()) continue;
	  
      for (iDict = 0; iDict < kNDict; iDict++) {
	dictionaryA[iDict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(iDet,iDict);
	dictionaryB[iDict] = (AliTRDarrayDictionary *) mergeSDigitsManager->GetDictionary(iDet,iDict);
	dictionaryA[iDict]->Expand();  
        dictionaryB[iDict]->Expand();
      }

      // Merge only detectors that contain a signal
      Bool_t doMerge = kTRUE;
      if (fMergeSignalOnly) {
        if (digitsA->GetOverThreshold(0) == 0) {                             
	  doMerge = kFALSE;
	}
      }
	  
      if (doMerge) {
	      
	AliDebug(1,Form("Merge detector %d of input no.%d",iDet,iMerge+1));
	      
	for (Int_t iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
	  for (Int_t iCol  = 0; iCol  <  nColMax;   iCol++ ) {
	    for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
         	
	      // Add the amplitudes of the summable digits 
	      Float_t ampA = digitsA->GetData(iRow,iCol,iTime);
	      Float_t ampB = digitsB->GetData(iRow,iCol,iTime);
	      ampA += ampB;
	      digitsA->SetData(iRow,iCol,iTime,ampA);

	      // Add the mask to the track id if defined.
	      for (iDict = 0; iDict < kNDict; iDict++) {
		Int_t trackB = dictionaryB[iDict]->GetData(iRow,iCol,iTime);
	        if ((fMasks) && (trackB > 0))  {
		  for (jDict = 0; jDict < kNDict; jDict++) { 
		    Int_t trackA = dictionaryA[iDict]->GetData(iRow,iCol,iTime); 
		    if (trackA == 0) {
		      trackA = trackB + fMasks[iMerge];
		      dictionaryA[iDict]->SetData(iRow,iCol,iTime,trackA);  
		    } // if:  track A == 0
		  } // for: jDict
		} // if:  fMasks and trackB > 0
	      } // for: iDict

	    } // for: iTime
	  } // for: iCol
        } // for: iRow

      } // if:  doMerge

      mergeSDigitsManager->RemoveDigits(iDet);
      mergeSDigitsManager->RemoveDictionaries(iDet);
  
      if (fCompress) {
        digitsA->Compress(0); 
        for (iDict = 0; iDict < kNDict; iDict++) {                                     
	  dictionaryA[iDict]->Compress();
        }
      }
      
    } // for: detectors    
      
    // The next set of s-digits
    mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(mergeSDigitsManager);
    
  } // while: mergeDigitsManagers
  
  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::ConvertSDigits()
{
  //
  // Converts s-digits to normal digits
  //

  AliTRDarraySignal *digitsIn = 0x0;

  if (!fSDigitsManager->HasSDigits()) {
    AliError("No s-digits in digits manager");
    return kFALSE;
  }

  // Loop through the detectors
  for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {

    // Get the merged s-digits (signals)
    digitsIn = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(det);
    if (!digitsIn->HasData()) {
      AliDebug(2,Form("No digits for det=%d",det));
      continue;
    }
    
    // Convert the merged sdigits to digits
    if (!Signal2ADC(det,digitsIn)) {
      continue;
    }

    // Copy the dictionary information to the output array
    if (!CopyDictionary(det)) {
      continue;
    }

    // Delete 
    fSDigitsManager->RemoveDigits(det);
    fSDigitsManager->RemoveDictionaries(det);

    // Run digital processing
    RunDigitalProcessing(det);

    // Compress the arrays
    CompressOutputArrays(det);

  } // for: detector numbers

  if (AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
    if (trklLoader->Tree())
      trklLoader->WriteData("OVERWRITE");
  }

  // Save the values for the raw data headers
  if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
    fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
  }
  else {
    fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
  }
  fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());

  return kTRUE;

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::CopyDictionary(Int_t det)
{
  //
  // Copies the dictionary information from the s-digits arrays
  // to the output arrays
  //

  AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
  if (!calibration) {
    AliFatal("Could not get calibration object");
    return kFALSE;
  }

  AliDebug(1,Form("Start copying dictionaries for detector=%d",det));

  const Int_t kNDict = AliTRDdigitsManager::kNDict;
  AliTRDarrayDictionary *dictionaryIn[kNDict];
  AliTRDarrayDictionary *dictionaryOut[kNDict];

  Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
  Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
  Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);

  Int_t row  = 0;
  Int_t col  = 0;
  Int_t time = 0;
  Int_t dict = 0;

  for (dict = 0; dict < kNDict; dict++) {

    dictionaryIn[dict]  = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(det,dict);
    dictionaryIn[dict]->Expand();
    dictionaryOut[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
    dictionaryOut[dict]->Allocate(nRowMax,nColMax,nTimeTotal);

    for (row = 0; row < nRowMax; row++) {
      for (col = 0; col < nColMax; col++) {
        for (time = 0; time < nTimeTotal; time++) {
	  Int_t track = dictionaryIn[dict]->GetData(row,col,time);
	  dictionaryOut[dict]->SetData(row,col,time,track);
	} // for: time
      } // for: col
    } // for: row
    
  } // for: dictionaries
  
  return kTRUE;

}

//_____________________________________________________________________________
void AliTRDdigitizer::CompressOutputArrays(Int_t det)
{
  //
  // Compress the output arrays
  //

  const Int_t kNDict = AliTRDdigitsManager::kNDict;
  AliTRDarrayDictionary *dictionary = 0x0;

  if (fCompress) {

    if (!fSDigits) {
      AliTRDarrayADC *digits = 0x0;  
      digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
      digits->Compress();
    }

    if (fSDigits) {
      AliTRDarraySignal *digits = 0x0; 
      digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
      digits->Compress(0);
    }

    for (Int_t dict = 0; dict < kNDict; dict++) {
      dictionary = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
      dictionary->Compress();
    }

  }

}

//_____________________________________________________________________________
Bool_t AliTRDdigitizer::WriteDigits() const
{
  //
  // Writes out the TRD-digits and the dictionaries
  //

  // Write parameters
  fRunLoader->CdGAFile();

  // Store the digits and the dictionary in the tree
  return fDigitsManager->WriteDigits();

}

//_____________________________________________________________________________
void AliTRDdigitizer::InitOutput(Int_t iEvent)
{
  //
  // Initializes the output branches
  //

  fEvent = iEvent;
   
  if (!fRunLoader) {
    AliError("Run Loader is NULL");
    return;  
  }

  AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
  if (!loader) {
    AliError("Can not get TRD loader from Run Loader");
    return;
  }

  TTree *tree = 0;
  
  if (fSDigits) { 
    // If we produce SDigits
    tree = loader->TreeS();
    if (!tree) {
      loader->MakeTree("S");
      tree = loader->TreeS();
    }
  }
  else {
    // If we produce Digits
    tree = loader->TreeD();
    if (!tree) {
      loader->MakeTree("D");
      tree = loader->TreeD();
    }
  }
  fDigitsManager->SetEvent(iEvent);
  fDigitsManager->MakeBranch(tree);

}
  
//_____________________________________________________________________________
Int_t AliTRDdigitizer::Diffusion(Float_t vdrift, Double_t absdriftlength
                               , Double_t exbvalue
                               , Double_t &lRow, Double_t &lCol, Double_t &lTime)
{
  //
  // Applies the diffusion smearing to the position of a single electron.
  // Depends on absolute drift length.
  //
  
  Float_t diffL = 0.0;
  Float_t diffT = 0.0;

  if (AliTRDCommonParam::Instance()->GetDiffCoeff(diffL,diffT,vdrift)) {

    Float_t driftSqrt = TMath::Sqrt(absdriftlength);
    Float_t sigmaT    = driftSqrt * diffT;
    Float_t sigmaL    = driftSqrt * diffL;
    lRow  = gRandom->Gaus(lRow ,sigmaT);
    if (AliTRDCommonParam::Instance()->ExBOn()) {
      lCol  = gRandom->Gaus(lCol ,sigmaT * 1.0 / (1.0 + exbvalue*exbvalue));
      lTime = gRandom->Gaus(lTime,sigmaL * 1.0 / (1.0 + exbvalue*exbvalue));
    }
    else {
      lCol  = gRandom->Gaus(lCol ,sigmaT);
      lTime = gRandom->Gaus(lTime,sigmaL);
    }

    return 1;

  }
  else {

    return 0;

  }

}
  
//_____________________________________________________________________________
void AliTRDdigitizer::RunDigitalProcessing(Int_t det)
{
  //
  // Run the digital processing in the TRAP
  //

  AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();

  AliTRDarrayADC *digits = fDigitsManager->GetDigits(det);
  if (!digits)
    return;

  //Call the methods in the mcm class using the temporary array as input  
  // process the data in the same order as in hardware
  for (Int_t side = 0; side <= 1; side++) {
    for(Int_t rob = side; rob < digits->GetNrow() / 2; rob += 2) {
      for(Int_t mcm = 0; mcm < 16; mcm++) {
	fMcmSim->Init(det, rob, mcm);
	fMcmSim->SetDataByPad(digits, fDigitsManager);
	fMcmSim->Filter();
	if (feeParam->GetTracklet()) {
	  fMcmSim->Tracklet();
	  fMcmSim->StoreTracklets();
	}
	fMcmSim->ZSMapping();
	fMcmSim->WriteData(digits);
      }
    }
  }
}

 AliTRDdigitizer.cxx:1
 AliTRDdigitizer.cxx:2
 AliTRDdigitizer.cxx:3
 AliTRDdigitizer.cxx:4
 AliTRDdigitizer.cxx:5
 AliTRDdigitizer.cxx:6
 AliTRDdigitizer.cxx:7
 AliTRDdigitizer.cxx:8
 AliTRDdigitizer.cxx:9
 AliTRDdigitizer.cxx:10
 AliTRDdigitizer.cxx:11
 AliTRDdigitizer.cxx:12
 AliTRDdigitizer.cxx:13
 AliTRDdigitizer.cxx:14
 AliTRDdigitizer.cxx:15
 AliTRDdigitizer.cxx:16
 AliTRDdigitizer.cxx:17
 AliTRDdigitizer.cxx:18
 AliTRDdigitizer.cxx:19
 AliTRDdigitizer.cxx:20
 AliTRDdigitizer.cxx:21
 AliTRDdigitizer.cxx:22
 AliTRDdigitizer.cxx:23
 AliTRDdigitizer.cxx:24
 AliTRDdigitizer.cxx:25
 AliTRDdigitizer.cxx:26
 AliTRDdigitizer.cxx:27
 AliTRDdigitizer.cxx:28
 AliTRDdigitizer.cxx:29
 AliTRDdigitizer.cxx:30
 AliTRDdigitizer.cxx:31
 AliTRDdigitizer.cxx:32
 AliTRDdigitizer.cxx:33
 AliTRDdigitizer.cxx:34
 AliTRDdigitizer.cxx:35
 AliTRDdigitizer.cxx:36
 AliTRDdigitizer.cxx:37
 AliTRDdigitizer.cxx:38
 AliTRDdigitizer.cxx:39
 AliTRDdigitizer.cxx:40
 AliTRDdigitizer.cxx:41
 AliTRDdigitizer.cxx:42
 AliTRDdigitizer.cxx:43
 AliTRDdigitizer.cxx:44
 AliTRDdigitizer.cxx:45
 AliTRDdigitizer.cxx:46
 AliTRDdigitizer.cxx:47
 AliTRDdigitizer.cxx:48
 AliTRDdigitizer.cxx:49
 AliTRDdigitizer.cxx:50
 AliTRDdigitizer.cxx:51
 AliTRDdigitizer.cxx:52
 AliTRDdigitizer.cxx:53
 AliTRDdigitizer.cxx:54
 AliTRDdigitizer.cxx:55
 AliTRDdigitizer.cxx:56
 AliTRDdigitizer.cxx:57
 AliTRDdigitizer.cxx:58
 AliTRDdigitizer.cxx:59
 AliTRDdigitizer.cxx:60
 AliTRDdigitizer.cxx:61
 AliTRDdigitizer.cxx:62
 AliTRDdigitizer.cxx:63
 AliTRDdigitizer.cxx:64
 AliTRDdigitizer.cxx:65
 AliTRDdigitizer.cxx:66
 AliTRDdigitizer.cxx:67
 AliTRDdigitizer.cxx:68
 AliTRDdigitizer.cxx:69
 AliTRDdigitizer.cxx:70
 AliTRDdigitizer.cxx:71
 AliTRDdigitizer.cxx:72
 AliTRDdigitizer.cxx:73
 AliTRDdigitizer.cxx:74
 AliTRDdigitizer.cxx:75
 AliTRDdigitizer.cxx:76
 AliTRDdigitizer.cxx:77
 AliTRDdigitizer.cxx:78
 AliTRDdigitizer.cxx:79
 AliTRDdigitizer.cxx:80
 AliTRDdigitizer.cxx:81
 AliTRDdigitizer.cxx:82
 AliTRDdigitizer.cxx:83
 AliTRDdigitizer.cxx:84
 AliTRDdigitizer.cxx:85
 AliTRDdigitizer.cxx:86
 AliTRDdigitizer.cxx:87
 AliTRDdigitizer.cxx:88
 AliTRDdigitizer.cxx:89
 AliTRDdigitizer.cxx:90
 AliTRDdigitizer.cxx:91
 AliTRDdigitizer.cxx:92
 AliTRDdigitizer.cxx:93
 AliTRDdigitizer.cxx:94
 AliTRDdigitizer.cxx:95
 AliTRDdigitizer.cxx:96
 AliTRDdigitizer.cxx:97
 AliTRDdigitizer.cxx:98
 AliTRDdigitizer.cxx:99
 AliTRDdigitizer.cxx:100
 AliTRDdigitizer.cxx:101
 AliTRDdigitizer.cxx:102
 AliTRDdigitizer.cxx:103
 AliTRDdigitizer.cxx:104
 AliTRDdigitizer.cxx:105
 AliTRDdigitizer.cxx:106
 AliTRDdigitizer.cxx:107
 AliTRDdigitizer.cxx:108
 AliTRDdigitizer.cxx:109
 AliTRDdigitizer.cxx:110
 AliTRDdigitizer.cxx:111
 AliTRDdigitizer.cxx:112
 AliTRDdigitizer.cxx:113
 AliTRDdigitizer.cxx:114
 AliTRDdigitizer.cxx:115
 AliTRDdigitizer.cxx:116
 AliTRDdigitizer.cxx:117
 AliTRDdigitizer.cxx:118
 AliTRDdigitizer.cxx:119
 AliTRDdigitizer.cxx:120
 AliTRDdigitizer.cxx:121
 AliTRDdigitizer.cxx:122
 AliTRDdigitizer.cxx:123
 AliTRDdigitizer.cxx:124
 AliTRDdigitizer.cxx:125
 AliTRDdigitizer.cxx:126
 AliTRDdigitizer.cxx:127
 AliTRDdigitizer.cxx:128
 AliTRDdigitizer.cxx:129
 AliTRDdigitizer.cxx:130
 AliTRDdigitizer.cxx:131
 AliTRDdigitizer.cxx:132
 AliTRDdigitizer.cxx:133
 AliTRDdigitizer.cxx:134
 AliTRDdigitizer.cxx:135
 AliTRDdigitizer.cxx:136
 AliTRDdigitizer.cxx:137
 AliTRDdigitizer.cxx:138
 AliTRDdigitizer.cxx:139
 AliTRDdigitizer.cxx:140
 AliTRDdigitizer.cxx:141
 AliTRDdigitizer.cxx:142
 AliTRDdigitizer.cxx:143
 AliTRDdigitizer.cxx:144
 AliTRDdigitizer.cxx:145
 AliTRDdigitizer.cxx:146
 AliTRDdigitizer.cxx:147
 AliTRDdigitizer.cxx:148
 AliTRDdigitizer.cxx:149
 AliTRDdigitizer.cxx:150
 AliTRDdigitizer.cxx:151
 AliTRDdigitizer.cxx:152
 AliTRDdigitizer.cxx:153
 AliTRDdigitizer.cxx:154
 AliTRDdigitizer.cxx:155
 AliTRDdigitizer.cxx:156
 AliTRDdigitizer.cxx:157
 AliTRDdigitizer.cxx:158
 AliTRDdigitizer.cxx:159
 AliTRDdigitizer.cxx:160
 AliTRDdigitizer.cxx:161
 AliTRDdigitizer.cxx:162
 AliTRDdigitizer.cxx:163
 AliTRDdigitizer.cxx:164
 AliTRDdigitizer.cxx:165
 AliTRDdigitizer.cxx:166
 AliTRDdigitizer.cxx:167
 AliTRDdigitizer.cxx:168
 AliTRDdigitizer.cxx:169
 AliTRDdigitizer.cxx:170
 AliTRDdigitizer.cxx:171
 AliTRDdigitizer.cxx:172
 AliTRDdigitizer.cxx:173
 AliTRDdigitizer.cxx:174
 AliTRDdigitizer.cxx:175
 AliTRDdigitizer.cxx:176
 AliTRDdigitizer.cxx:177
 AliTRDdigitizer.cxx:178
 AliTRDdigitizer.cxx:179
 AliTRDdigitizer.cxx:180
 AliTRDdigitizer.cxx:181
 AliTRDdigitizer.cxx:182
 AliTRDdigitizer.cxx:183
 AliTRDdigitizer.cxx:184
 AliTRDdigitizer.cxx:185
 AliTRDdigitizer.cxx:186
 AliTRDdigitizer.cxx:187
 AliTRDdigitizer.cxx:188
 AliTRDdigitizer.cxx:189
 AliTRDdigitizer.cxx:190
 AliTRDdigitizer.cxx:191
 AliTRDdigitizer.cxx:192
 AliTRDdigitizer.cxx:193
 AliTRDdigitizer.cxx:194
 AliTRDdigitizer.cxx:195
 AliTRDdigitizer.cxx:196
 AliTRDdigitizer.cxx:197
 AliTRDdigitizer.cxx:198
 AliTRDdigitizer.cxx:199
 AliTRDdigitizer.cxx:200
 AliTRDdigitizer.cxx:201
 AliTRDdigitizer.cxx:202
 AliTRDdigitizer.cxx:203
 AliTRDdigitizer.cxx:204
 AliTRDdigitizer.cxx:205
 AliTRDdigitizer.cxx:206
 AliTRDdigitizer.cxx:207
 AliTRDdigitizer.cxx:208
 AliTRDdigitizer.cxx:209
 AliTRDdigitizer.cxx:210
 AliTRDdigitizer.cxx:211
 AliTRDdigitizer.cxx:212
 AliTRDdigitizer.cxx:213
 AliTRDdigitizer.cxx:214
 AliTRDdigitizer.cxx:215
 AliTRDdigitizer.cxx:216
 AliTRDdigitizer.cxx:217
 AliTRDdigitizer.cxx:218
 AliTRDdigitizer.cxx:219
 AliTRDdigitizer.cxx:220
 AliTRDdigitizer.cxx:221
 AliTRDdigitizer.cxx:222
 AliTRDdigitizer.cxx:223
 AliTRDdigitizer.cxx:224
 AliTRDdigitizer.cxx:225
 AliTRDdigitizer.cxx:226
 AliTRDdigitizer.cxx:227
 AliTRDdigitizer.cxx:228
 AliTRDdigitizer.cxx:229
 AliTRDdigitizer.cxx:230
 AliTRDdigitizer.cxx:231
 AliTRDdigitizer.cxx:232
 AliTRDdigitizer.cxx:233
 AliTRDdigitizer.cxx:234
 AliTRDdigitizer.cxx:235
 AliTRDdigitizer.cxx:236
 AliTRDdigitizer.cxx:237
 AliTRDdigitizer.cxx:238
 AliTRDdigitizer.cxx:239
 AliTRDdigitizer.cxx:240
 AliTRDdigitizer.cxx:241
 AliTRDdigitizer.cxx:242
 AliTRDdigitizer.cxx:243
 AliTRDdigitizer.cxx:244
 AliTRDdigitizer.cxx:245
 AliTRDdigitizer.cxx:246
 AliTRDdigitizer.cxx:247
 AliTRDdigitizer.cxx:248
 AliTRDdigitizer.cxx:249
 AliTRDdigitizer.cxx:250
 AliTRDdigitizer.cxx:251
 AliTRDdigitizer.cxx:252
 AliTRDdigitizer.cxx:253
 AliTRDdigitizer.cxx:254
 AliTRDdigitizer.cxx:255
 AliTRDdigitizer.cxx:256
 AliTRDdigitizer.cxx:257
 AliTRDdigitizer.cxx:258
 AliTRDdigitizer.cxx:259
 AliTRDdigitizer.cxx:260
 AliTRDdigitizer.cxx:261
 AliTRDdigitizer.cxx:262
 AliTRDdigitizer.cxx:263
 AliTRDdigitizer.cxx:264
 AliTRDdigitizer.cxx:265
 AliTRDdigitizer.cxx:266
 AliTRDdigitizer.cxx:267
 AliTRDdigitizer.cxx:268
 AliTRDdigitizer.cxx:269
 AliTRDdigitizer.cxx:270
 AliTRDdigitizer.cxx:271
 AliTRDdigitizer.cxx:272
 AliTRDdigitizer.cxx:273
 AliTRDdigitizer.cxx:274
 AliTRDdigitizer.cxx:275
 AliTRDdigitizer.cxx:276
 AliTRDdigitizer.cxx:277
 AliTRDdigitizer.cxx:278
 AliTRDdigitizer.cxx:279
 AliTRDdigitizer.cxx:280
 AliTRDdigitizer.cxx:281
 AliTRDdigitizer.cxx:282
 AliTRDdigitizer.cxx:283
 AliTRDdigitizer.cxx:284
 AliTRDdigitizer.cxx:285
 AliTRDdigitizer.cxx:286
 AliTRDdigitizer.cxx:287
 AliTRDdigitizer.cxx:288
 AliTRDdigitizer.cxx:289
 AliTRDdigitizer.cxx:290
 AliTRDdigitizer.cxx:291
 AliTRDdigitizer.cxx:292
 AliTRDdigitizer.cxx:293
 AliTRDdigitizer.cxx:294
 AliTRDdigitizer.cxx:295
 AliTRDdigitizer.cxx:296
 AliTRDdigitizer.cxx:297
 AliTRDdigitizer.cxx:298
 AliTRDdigitizer.cxx:299
 AliTRDdigitizer.cxx:300
 AliTRDdigitizer.cxx:301
 AliTRDdigitizer.cxx:302
 AliTRDdigitizer.cxx:303
 AliTRDdigitizer.cxx:304
 AliTRDdigitizer.cxx:305
 AliTRDdigitizer.cxx:306
 AliTRDdigitizer.cxx:307
 AliTRDdigitizer.cxx:308
 AliTRDdigitizer.cxx:309
 AliTRDdigitizer.cxx:310
 AliTRDdigitizer.cxx:311
 AliTRDdigitizer.cxx:312
 AliTRDdigitizer.cxx:313
 AliTRDdigitizer.cxx:314
 AliTRDdigitizer.cxx:315
 AliTRDdigitizer.cxx:316
 AliTRDdigitizer.cxx:317
 AliTRDdigitizer.cxx:318
 AliTRDdigitizer.cxx:319
 AliTRDdigitizer.cxx:320
 AliTRDdigitizer.cxx:321
 AliTRDdigitizer.cxx:322
 AliTRDdigitizer.cxx:323
 AliTRDdigitizer.cxx:324
 AliTRDdigitizer.cxx:325
 AliTRDdigitizer.cxx:326
 AliTRDdigitizer.cxx:327
 AliTRDdigitizer.cxx:328
 AliTRDdigitizer.cxx:329
 AliTRDdigitizer.cxx:330
 AliTRDdigitizer.cxx:331
 AliTRDdigitizer.cxx:332
 AliTRDdigitizer.cxx:333
 AliTRDdigitizer.cxx:334
 AliTRDdigitizer.cxx:335
 AliTRDdigitizer.cxx:336
 AliTRDdigitizer.cxx:337
 AliTRDdigitizer.cxx:338
 AliTRDdigitizer.cxx:339
 AliTRDdigitizer.cxx:340
 AliTRDdigitizer.cxx:341
 AliTRDdigitizer.cxx:342
 AliTRDdigitizer.cxx:343
 AliTRDdigitizer.cxx:344
 AliTRDdigitizer.cxx:345
 AliTRDdigitizer.cxx:346
 AliTRDdigitizer.cxx:347
 AliTRDdigitizer.cxx:348
 AliTRDdigitizer.cxx:349
 AliTRDdigitizer.cxx:350
 AliTRDdigitizer.cxx:351
 AliTRDdigitizer.cxx:352
 AliTRDdigitizer.cxx:353
 AliTRDdigitizer.cxx:354
 AliTRDdigitizer.cxx:355
 AliTRDdigitizer.cxx:356
 AliTRDdigitizer.cxx:357
 AliTRDdigitizer.cxx:358
 AliTRDdigitizer.cxx:359
 AliTRDdigitizer.cxx:360
 AliTRDdigitizer.cxx:361
 AliTRDdigitizer.cxx:362
 AliTRDdigitizer.cxx:363
 AliTRDdigitizer.cxx:364
 AliTRDdigitizer.cxx:365
 AliTRDdigitizer.cxx:366
 AliTRDdigitizer.cxx:367
 AliTRDdigitizer.cxx:368
 AliTRDdigitizer.cxx:369
 AliTRDdigitizer.cxx:370
 AliTRDdigitizer.cxx:371
 AliTRDdigitizer.cxx:372
 AliTRDdigitizer.cxx:373
 AliTRDdigitizer.cxx:374
 AliTRDdigitizer.cxx:375
 AliTRDdigitizer.cxx:376
 AliTRDdigitizer.cxx:377
 AliTRDdigitizer.cxx:378
 AliTRDdigitizer.cxx:379
 AliTRDdigitizer.cxx:380
 AliTRDdigitizer.cxx:381
 AliTRDdigitizer.cxx:382
 AliTRDdigitizer.cxx:383
 AliTRDdigitizer.cxx:384
 AliTRDdigitizer.cxx:385
 AliTRDdigitizer.cxx:386
 AliTRDdigitizer.cxx:387
 AliTRDdigitizer.cxx:388
 AliTRDdigitizer.cxx:389
 AliTRDdigitizer.cxx:390
 AliTRDdigitizer.cxx:391
 AliTRDdigitizer.cxx:392
 AliTRDdigitizer.cxx:393
 AliTRDdigitizer.cxx:394
 AliTRDdigitizer.cxx:395
 AliTRDdigitizer.cxx:396
 AliTRDdigitizer.cxx:397
 AliTRDdigitizer.cxx:398
 AliTRDdigitizer.cxx:399
 AliTRDdigitizer.cxx:400
 AliTRDdigitizer.cxx:401
 AliTRDdigitizer.cxx:402
 AliTRDdigitizer.cxx:403
 AliTRDdigitizer.cxx:404
 AliTRDdigitizer.cxx:405
 AliTRDdigitizer.cxx:406
 AliTRDdigitizer.cxx:407
 AliTRDdigitizer.cxx:408
 AliTRDdigitizer.cxx:409
 AliTRDdigitizer.cxx:410
 AliTRDdigitizer.cxx:411
 AliTRDdigitizer.cxx:412
 AliTRDdigitizer.cxx:413
 AliTRDdigitizer.cxx:414
 AliTRDdigitizer.cxx:415
 AliTRDdigitizer.cxx:416
 AliTRDdigitizer.cxx:417
 AliTRDdigitizer.cxx:418
 AliTRDdigitizer.cxx:419
 AliTRDdigitizer.cxx:420
 AliTRDdigitizer.cxx:421
 AliTRDdigitizer.cxx:422
 AliTRDdigitizer.cxx:423
 AliTRDdigitizer.cxx:424
 AliTRDdigitizer.cxx:425
 AliTRDdigitizer.cxx:426
 AliTRDdigitizer.cxx:427
 AliTRDdigitizer.cxx:428
 AliTRDdigitizer.cxx:429
 AliTRDdigitizer.cxx:430
 AliTRDdigitizer.cxx:431
 AliTRDdigitizer.cxx:432
 AliTRDdigitizer.cxx:433
 AliTRDdigitizer.cxx:434
 AliTRDdigitizer.cxx:435
 AliTRDdigitizer.cxx:436
 AliTRDdigitizer.cxx:437
 AliTRDdigitizer.cxx:438
 AliTRDdigitizer.cxx:439
 AliTRDdigitizer.cxx:440
 AliTRDdigitizer.cxx:441
 AliTRDdigitizer.cxx:442
 AliTRDdigitizer.cxx:443
 AliTRDdigitizer.cxx:444
 AliTRDdigitizer.cxx:445
 AliTRDdigitizer.cxx:446
 AliTRDdigitizer.cxx:447
 AliTRDdigitizer.cxx:448
 AliTRDdigitizer.cxx:449
 AliTRDdigitizer.cxx:450
 AliTRDdigitizer.cxx:451
 AliTRDdigitizer.cxx:452
 AliTRDdigitizer.cxx:453
 AliTRDdigitizer.cxx:454
 AliTRDdigitizer.cxx:455
 AliTRDdigitizer.cxx:456
 AliTRDdigitizer.cxx:457
 AliTRDdigitizer.cxx:458
 AliTRDdigitizer.cxx:459
 AliTRDdigitizer.cxx:460
 AliTRDdigitizer.cxx:461
 AliTRDdigitizer.cxx:462
 AliTRDdigitizer.cxx:463
 AliTRDdigitizer.cxx:464
 AliTRDdigitizer.cxx:465
 AliTRDdigitizer.cxx:466
 AliTRDdigitizer.cxx:467
 AliTRDdigitizer.cxx:468
 AliTRDdigitizer.cxx:469
 AliTRDdigitizer.cxx:470
 AliTRDdigitizer.cxx:471
 AliTRDdigitizer.cxx:472
 AliTRDdigitizer.cxx:473
 AliTRDdigitizer.cxx:474
 AliTRDdigitizer.cxx:475
 AliTRDdigitizer.cxx:476
 AliTRDdigitizer.cxx:477
 AliTRDdigitizer.cxx:478
 AliTRDdigitizer.cxx:479
 AliTRDdigitizer.cxx:480
 AliTRDdigitizer.cxx:481
 AliTRDdigitizer.cxx:482
 AliTRDdigitizer.cxx:483
 AliTRDdigitizer.cxx:484
 AliTRDdigitizer.cxx:485
 AliTRDdigitizer.cxx:486
 AliTRDdigitizer.cxx:487
 AliTRDdigitizer.cxx:488
 AliTRDdigitizer.cxx:489
 AliTRDdigitizer.cxx:490
 AliTRDdigitizer.cxx:491
 AliTRDdigitizer.cxx:492
 AliTRDdigitizer.cxx:493
 AliTRDdigitizer.cxx:494
 AliTRDdigitizer.cxx:495
 AliTRDdigitizer.cxx:496
 AliTRDdigitizer.cxx:497
 AliTRDdigitizer.cxx:498
 AliTRDdigitizer.cxx:499
 AliTRDdigitizer.cxx:500
 AliTRDdigitizer.cxx:501
 AliTRDdigitizer.cxx:502
 AliTRDdigitizer.cxx:503
 AliTRDdigitizer.cxx:504
 AliTRDdigitizer.cxx:505
 AliTRDdigitizer.cxx:506
 AliTRDdigitizer.cxx:507
 AliTRDdigitizer.cxx:508
 AliTRDdigitizer.cxx:509
 AliTRDdigitizer.cxx:510
 AliTRDdigitizer.cxx:511
 AliTRDdigitizer.cxx:512
 AliTRDdigitizer.cxx:513
 AliTRDdigitizer.cxx:514
 AliTRDdigitizer.cxx:515
 AliTRDdigitizer.cxx:516
 AliTRDdigitizer.cxx:517
 AliTRDdigitizer.cxx:518
 AliTRDdigitizer.cxx:519
 AliTRDdigitizer.cxx:520
 AliTRDdigitizer.cxx:521
 AliTRDdigitizer.cxx:522
 AliTRDdigitizer.cxx:523
 AliTRDdigitizer.cxx:524
 AliTRDdigitizer.cxx:525
 AliTRDdigitizer.cxx:526
 AliTRDdigitizer.cxx:527
 AliTRDdigitizer.cxx:528
 AliTRDdigitizer.cxx:529
 AliTRDdigitizer.cxx:530
 AliTRDdigitizer.cxx:531
 AliTRDdigitizer.cxx:532
 AliTRDdigitizer.cxx:533
 AliTRDdigitizer.cxx:534
 AliTRDdigitizer.cxx:535
 AliTRDdigitizer.cxx:536
 AliTRDdigitizer.cxx:537
 AliTRDdigitizer.cxx:538
 AliTRDdigitizer.cxx:539
 AliTRDdigitizer.cxx:540
 AliTRDdigitizer.cxx:541
 AliTRDdigitizer.cxx:542
 AliTRDdigitizer.cxx:543
 AliTRDdigitizer.cxx:544
 AliTRDdigitizer.cxx:545
 AliTRDdigitizer.cxx:546
 AliTRDdigitizer.cxx:547
 AliTRDdigitizer.cxx:548
 AliTRDdigitizer.cxx:549
 AliTRDdigitizer.cxx:550
 AliTRDdigitizer.cxx:551
 AliTRDdigitizer.cxx:552
 AliTRDdigitizer.cxx:553
 AliTRDdigitizer.cxx:554
 AliTRDdigitizer.cxx:555
 AliTRDdigitizer.cxx:556
 AliTRDdigitizer.cxx:557
 AliTRDdigitizer.cxx:558
 AliTRDdigitizer.cxx:559
 AliTRDdigitizer.cxx:560
 AliTRDdigitizer.cxx:561
 AliTRDdigitizer.cxx:562
 AliTRDdigitizer.cxx:563
 AliTRDdigitizer.cxx:564
 AliTRDdigitizer.cxx:565
 AliTRDdigitizer.cxx:566
 AliTRDdigitizer.cxx:567
 AliTRDdigitizer.cxx:568
 AliTRDdigitizer.cxx:569
 AliTRDdigitizer.cxx:570
 AliTRDdigitizer.cxx:571
 AliTRDdigitizer.cxx:572
 AliTRDdigitizer.cxx:573
 AliTRDdigitizer.cxx:574
 AliTRDdigitizer.cxx:575
 AliTRDdigitizer.cxx:576
 AliTRDdigitizer.cxx:577
 AliTRDdigitizer.cxx:578
 AliTRDdigitizer.cxx:579
 AliTRDdigitizer.cxx:580
 AliTRDdigitizer.cxx:581
 AliTRDdigitizer.cxx:582
 AliTRDdigitizer.cxx:583
 AliTRDdigitizer.cxx:584
 AliTRDdigitizer.cxx:585
 AliTRDdigitizer.cxx:586
 AliTRDdigitizer.cxx:587
 AliTRDdigitizer.cxx:588
 AliTRDdigitizer.cxx:589
 AliTRDdigitizer.cxx:590
 AliTRDdigitizer.cxx:591
 AliTRDdigitizer.cxx:592
 AliTRDdigitizer.cxx:593
 AliTRDdigitizer.cxx:594
 AliTRDdigitizer.cxx:595
 AliTRDdigitizer.cxx:596
 AliTRDdigitizer.cxx:597
 AliTRDdigitizer.cxx:598
 AliTRDdigitizer.cxx:599
 AliTRDdigitizer.cxx:600
 AliTRDdigitizer.cxx:601
 AliTRDdigitizer.cxx:602
 AliTRDdigitizer.cxx:603
 AliTRDdigitizer.cxx:604
 AliTRDdigitizer.cxx:605
 AliTRDdigitizer.cxx:606
 AliTRDdigitizer.cxx:607
 AliTRDdigitizer.cxx:608
 AliTRDdigitizer.cxx:609
 AliTRDdigitizer.cxx:610
 AliTRDdigitizer.cxx:611
 AliTRDdigitizer.cxx:612
 AliTRDdigitizer.cxx:613
 AliTRDdigitizer.cxx:614
 AliTRDdigitizer.cxx:615
 AliTRDdigitizer.cxx:616
 AliTRDdigitizer.cxx:617
 AliTRDdigitizer.cxx:618
 AliTRDdigitizer.cxx:619
 AliTRDdigitizer.cxx:620
 AliTRDdigitizer.cxx:621
 AliTRDdigitizer.cxx:622
 AliTRDdigitizer.cxx:623
 AliTRDdigitizer.cxx:624
 AliTRDdigitizer.cxx:625
 AliTRDdigitizer.cxx:626
 AliTRDdigitizer.cxx:627
 AliTRDdigitizer.cxx:628
 AliTRDdigitizer.cxx:629
 AliTRDdigitizer.cxx:630
 AliTRDdigitizer.cxx:631
 AliTRDdigitizer.cxx:632
 AliTRDdigitizer.cxx:633
 AliTRDdigitizer.cxx:634
 AliTRDdigitizer.cxx:635
 AliTRDdigitizer.cxx:636
 AliTRDdigitizer.cxx:637
 AliTRDdigitizer.cxx:638
 AliTRDdigitizer.cxx:639
 AliTRDdigitizer.cxx:640
 AliTRDdigitizer.cxx:641
 AliTRDdigitizer.cxx:642
 AliTRDdigitizer.cxx:643
 AliTRDdigitizer.cxx:644
 AliTRDdigitizer.cxx:645
 AliTRDdigitizer.cxx:646
 AliTRDdigitizer.cxx:647
 AliTRDdigitizer.cxx:648
 AliTRDdigitizer.cxx:649
 AliTRDdigitizer.cxx:650
 AliTRDdigitizer.cxx:651
 AliTRDdigitizer.cxx:652
 AliTRDdigitizer.cxx:653
 AliTRDdigitizer.cxx:654
 AliTRDdigitizer.cxx:655
 AliTRDdigitizer.cxx:656
 AliTRDdigitizer.cxx:657
 AliTRDdigitizer.cxx:658
 AliTRDdigitizer.cxx:659
 AliTRDdigitizer.cxx:660
 AliTRDdigitizer.cxx:661
 AliTRDdigitizer.cxx:662
 AliTRDdigitizer.cxx:663
 AliTRDdigitizer.cxx:664
 AliTRDdigitizer.cxx:665
 AliTRDdigitizer.cxx:666
 AliTRDdigitizer.cxx:667
 AliTRDdigitizer.cxx:668
 AliTRDdigitizer.cxx:669
 AliTRDdigitizer.cxx:670
 AliTRDdigitizer.cxx:671
 AliTRDdigitizer.cxx:672
 AliTRDdigitizer.cxx:673
 AliTRDdigitizer.cxx:674
 AliTRDdigitizer.cxx:675
 AliTRDdigitizer.cxx:676
 AliTRDdigitizer.cxx:677
 AliTRDdigitizer.cxx:678
 AliTRDdigitizer.cxx:679
 AliTRDdigitizer.cxx:680
 AliTRDdigitizer.cxx:681
 AliTRDdigitizer.cxx:682
 AliTRDdigitizer.cxx:683
 AliTRDdigitizer.cxx:684
 AliTRDdigitizer.cxx:685
 AliTRDdigitizer.cxx:686
 AliTRDdigitizer.cxx:687
 AliTRDdigitizer.cxx:688
 AliTRDdigitizer.cxx:689
 AliTRDdigitizer.cxx:690
 AliTRDdigitizer.cxx:691
 AliTRDdigitizer.cxx:692
 AliTRDdigitizer.cxx:693
 AliTRDdigitizer.cxx:694
 AliTRDdigitizer.cxx:695
 AliTRDdigitizer.cxx:696
 AliTRDdigitizer.cxx:697
 AliTRDdigitizer.cxx:698
 AliTRDdigitizer.cxx:699
 AliTRDdigitizer.cxx:700
 AliTRDdigitizer.cxx:701
 AliTRDdigitizer.cxx:702
 AliTRDdigitizer.cxx:703
 AliTRDdigitizer.cxx:704
 AliTRDdigitizer.cxx:705
 AliTRDdigitizer.cxx:706
 AliTRDdigitizer.cxx:707
 AliTRDdigitizer.cxx:708
 AliTRDdigitizer.cxx:709
 AliTRDdigitizer.cxx:710
 AliTRDdigitizer.cxx:711
 AliTRDdigitizer.cxx:712
 AliTRDdigitizer.cxx:713
 AliTRDdigitizer.cxx:714
 AliTRDdigitizer.cxx:715
 AliTRDdigitizer.cxx:716
 AliTRDdigitizer.cxx:717
 AliTRDdigitizer.cxx:718
 AliTRDdigitizer.cxx:719
 AliTRDdigitizer.cxx:720
 AliTRDdigitizer.cxx:721
 AliTRDdigitizer.cxx:722
 AliTRDdigitizer.cxx:723
 AliTRDdigitizer.cxx:724
 AliTRDdigitizer.cxx:725
 AliTRDdigitizer.cxx:726
 AliTRDdigitizer.cxx:727
 AliTRDdigitizer.cxx:728
 AliTRDdigitizer.cxx:729
 AliTRDdigitizer.cxx:730
 AliTRDdigitizer.cxx:731
 AliTRDdigitizer.cxx:732
 AliTRDdigitizer.cxx:733
 AliTRDdigitizer.cxx:734
 AliTRDdigitizer.cxx:735
 AliTRDdigitizer.cxx:736
 AliTRDdigitizer.cxx:737
 AliTRDdigitizer.cxx:738
 AliTRDdigitizer.cxx:739
 AliTRDdigitizer.cxx:740
 AliTRDdigitizer.cxx:741
 AliTRDdigitizer.cxx:742
 AliTRDdigitizer.cxx:743
 AliTRDdigitizer.cxx:744
 AliTRDdigitizer.cxx:745
 AliTRDdigitizer.cxx:746
 AliTRDdigitizer.cxx:747
 AliTRDdigitizer.cxx:748
 AliTRDdigitizer.cxx:749
 AliTRDdigitizer.cxx:750
 AliTRDdigitizer.cxx:751
 AliTRDdigitizer.cxx:752
 AliTRDdigitizer.cxx:753
 AliTRDdigitizer.cxx:754
 AliTRDdigitizer.cxx:755
 AliTRDdigitizer.cxx:756
 AliTRDdigitizer.cxx:757
 AliTRDdigitizer.cxx:758
 AliTRDdigitizer.cxx:759
 AliTRDdigitizer.cxx:760
 AliTRDdigitizer.cxx:761
 AliTRDdigitizer.cxx:762
 AliTRDdigitizer.cxx:763
 AliTRDdigitizer.cxx:764
 AliTRDdigitizer.cxx:765
 AliTRDdigitizer.cxx:766
 AliTRDdigitizer.cxx:767
 AliTRDdigitizer.cxx:768
 AliTRDdigitizer.cxx:769
 AliTRDdigitizer.cxx:770
 AliTRDdigitizer.cxx:771
 AliTRDdigitizer.cxx:772
 AliTRDdigitizer.cxx:773
 AliTRDdigitizer.cxx:774
 AliTRDdigitizer.cxx:775
 AliTRDdigitizer.cxx:776
 AliTRDdigitizer.cxx:777
 AliTRDdigitizer.cxx:778
 AliTRDdigitizer.cxx:779
 AliTRDdigitizer.cxx:780
 AliTRDdigitizer.cxx:781
 AliTRDdigitizer.cxx:782
 AliTRDdigitizer.cxx:783
 AliTRDdigitizer.cxx:784
 AliTRDdigitizer.cxx:785
 AliTRDdigitizer.cxx:786
 AliTRDdigitizer.cxx:787
 AliTRDdigitizer.cxx:788
 AliTRDdigitizer.cxx:789
 AliTRDdigitizer.cxx:790
 AliTRDdigitizer.cxx:791
 AliTRDdigitizer.cxx:792
 AliTRDdigitizer.cxx:793
 AliTRDdigitizer.cxx:794
 AliTRDdigitizer.cxx:795
 AliTRDdigitizer.cxx:796
 AliTRDdigitizer.cxx:797
 AliTRDdigitizer.cxx:798
 AliTRDdigitizer.cxx:799
 AliTRDdigitizer.cxx:800
 AliTRDdigitizer.cxx:801
 AliTRDdigitizer.cxx:802
 AliTRDdigitizer.cxx:803
 AliTRDdigitizer.cxx:804
 AliTRDdigitizer.cxx:805
 AliTRDdigitizer.cxx:806
 AliTRDdigitizer.cxx:807
 AliTRDdigitizer.cxx:808
 AliTRDdigitizer.cxx:809
 AliTRDdigitizer.cxx:810
 AliTRDdigitizer.cxx:811
 AliTRDdigitizer.cxx:812
 AliTRDdigitizer.cxx:813
 AliTRDdigitizer.cxx:814
 AliTRDdigitizer.cxx:815
 AliTRDdigitizer.cxx:816
 AliTRDdigitizer.cxx:817
 AliTRDdigitizer.cxx:818
 AliTRDdigitizer.cxx:819
 AliTRDdigitizer.cxx:820
 AliTRDdigitizer.cxx:821
 AliTRDdigitizer.cxx:822
 AliTRDdigitizer.cxx:823
 AliTRDdigitizer.cxx:824
 AliTRDdigitizer.cxx:825
 AliTRDdigitizer.cxx:826
 AliTRDdigitizer.cxx:827
 AliTRDdigitizer.cxx:828
 AliTRDdigitizer.cxx:829
 AliTRDdigitizer.cxx:830
 AliTRDdigitizer.cxx:831
 AliTRDdigitizer.cxx:832
 AliTRDdigitizer.cxx:833
 AliTRDdigitizer.cxx:834
 AliTRDdigitizer.cxx:835
 AliTRDdigitizer.cxx:836
 AliTRDdigitizer.cxx:837
 AliTRDdigitizer.cxx:838
 AliTRDdigitizer.cxx:839
 AliTRDdigitizer.cxx:840
 AliTRDdigitizer.cxx:841
 AliTRDdigitizer.cxx:842
 AliTRDdigitizer.cxx:843
 AliTRDdigitizer.cxx:844
 AliTRDdigitizer.cxx:845
 AliTRDdigitizer.cxx:846
 AliTRDdigitizer.cxx:847
 AliTRDdigitizer.cxx:848
 AliTRDdigitizer.cxx:849
 AliTRDdigitizer.cxx:850
 AliTRDdigitizer.cxx:851
 AliTRDdigitizer.cxx:852
 AliTRDdigitizer.cxx:853
 AliTRDdigitizer.cxx:854
 AliTRDdigitizer.cxx:855
 AliTRDdigitizer.cxx:856
 AliTRDdigitizer.cxx:857
 AliTRDdigitizer.cxx:858
 AliTRDdigitizer.cxx:859
 AliTRDdigitizer.cxx:860
 AliTRDdigitizer.cxx:861
 AliTRDdigitizer.cxx:862
 AliTRDdigitizer.cxx:863
 AliTRDdigitizer.cxx:864
 AliTRDdigitizer.cxx:865
 AliTRDdigitizer.cxx:866
 AliTRDdigitizer.cxx:867
 AliTRDdigitizer.cxx:868
 AliTRDdigitizer.cxx:869
 AliTRDdigitizer.cxx:870
 AliTRDdigitizer.cxx:871
 AliTRDdigitizer.cxx:872
 AliTRDdigitizer.cxx:873
 AliTRDdigitizer.cxx:874
 AliTRDdigitizer.cxx:875
 AliTRDdigitizer.cxx:876
 AliTRDdigitizer.cxx:877
 AliTRDdigitizer.cxx:878
 AliTRDdigitizer.cxx:879
 AliTRDdigitizer.cxx:880
 AliTRDdigitizer.cxx:881
 AliTRDdigitizer.cxx:882
 AliTRDdigitizer.cxx:883
 AliTRDdigitizer.cxx:884
 AliTRDdigitizer.cxx:885
 AliTRDdigitizer.cxx:886
 AliTRDdigitizer.cxx:887
 AliTRDdigitizer.cxx:888
 AliTRDdigitizer.cxx:889
 AliTRDdigitizer.cxx:890
 AliTRDdigitizer.cxx:891
 AliTRDdigitizer.cxx:892
 AliTRDdigitizer.cxx:893
 AliTRDdigitizer.cxx:894
 AliTRDdigitizer.cxx:895
 AliTRDdigitizer.cxx:896
 AliTRDdigitizer.cxx:897
 AliTRDdigitizer.cxx:898
 AliTRDdigitizer.cxx:899
 AliTRDdigitizer.cxx:900
 AliTRDdigitizer.cxx:901
 AliTRDdigitizer.cxx:902
 AliTRDdigitizer.cxx:903
 AliTRDdigitizer.cxx:904
 AliTRDdigitizer.cxx:905
 AliTRDdigitizer.cxx:906
 AliTRDdigitizer.cxx:907
 AliTRDdigitizer.cxx:908
 AliTRDdigitizer.cxx:909
 AliTRDdigitizer.cxx:910
 AliTRDdigitizer.cxx:911
 AliTRDdigitizer.cxx:912
 AliTRDdigitizer.cxx:913
 AliTRDdigitizer.cxx:914
 AliTRDdigitizer.cxx:915
 AliTRDdigitizer.cxx:916
 AliTRDdigitizer.cxx:917
 AliTRDdigitizer.cxx:918
 AliTRDdigitizer.cxx:919
 AliTRDdigitizer.cxx:920
 AliTRDdigitizer.cxx:921
 AliTRDdigitizer.cxx:922
 AliTRDdigitizer.cxx:923
 AliTRDdigitizer.cxx:924
 AliTRDdigitizer.cxx:925
 AliTRDdigitizer.cxx:926
 AliTRDdigitizer.cxx:927
 AliTRDdigitizer.cxx:928
 AliTRDdigitizer.cxx:929
 AliTRDdigitizer.cxx:930
 AliTRDdigitizer.cxx:931
 AliTRDdigitizer.cxx:932
 AliTRDdigitizer.cxx:933
 AliTRDdigitizer.cxx:934
 AliTRDdigitizer.cxx:935
 AliTRDdigitizer.cxx:936
 AliTRDdigitizer.cxx:937
 AliTRDdigitizer.cxx:938
 AliTRDdigitizer.cxx:939
 AliTRDdigitizer.cxx:940
 AliTRDdigitizer.cxx:941
 AliTRDdigitizer.cxx:942
 AliTRDdigitizer.cxx:943
 AliTRDdigitizer.cxx:944
 AliTRDdigitizer.cxx:945
 AliTRDdigitizer.cxx:946
 AliTRDdigitizer.cxx:947
 AliTRDdigitizer.cxx:948
 AliTRDdigitizer.cxx:949
 AliTRDdigitizer.cxx:950
 AliTRDdigitizer.cxx:951
 AliTRDdigitizer.cxx:952
 AliTRDdigitizer.cxx:953
 AliTRDdigitizer.cxx:954
 AliTRDdigitizer.cxx:955
 AliTRDdigitizer.cxx:956
 AliTRDdigitizer.cxx:957
 AliTRDdigitizer.cxx:958
 AliTRDdigitizer.cxx:959
 AliTRDdigitizer.cxx:960
 AliTRDdigitizer.cxx:961
 AliTRDdigitizer.cxx:962
 AliTRDdigitizer.cxx:963
 AliTRDdigitizer.cxx:964
 AliTRDdigitizer.cxx:965
 AliTRDdigitizer.cxx:966
 AliTRDdigitizer.cxx:967
 AliTRDdigitizer.cxx:968
 AliTRDdigitizer.cxx:969
 AliTRDdigitizer.cxx:970
 AliTRDdigitizer.cxx:971
 AliTRDdigitizer.cxx:972
 AliTRDdigitizer.cxx:973
 AliTRDdigitizer.cxx:974
 AliTRDdigitizer.cxx:975
 AliTRDdigitizer.cxx:976
 AliTRDdigitizer.cxx:977
 AliTRDdigitizer.cxx:978
 AliTRDdigitizer.cxx:979
 AliTRDdigitizer.cxx:980
 AliTRDdigitizer.cxx:981
 AliTRDdigitizer.cxx:982
 AliTRDdigitizer.cxx:983
 AliTRDdigitizer.cxx:984
 AliTRDdigitizer.cxx:985
 AliTRDdigitizer.cxx:986
 AliTRDdigitizer.cxx:987
 AliTRDdigitizer.cxx:988
 AliTRDdigitizer.cxx:989
 AliTRDdigitizer.cxx:990
 AliTRDdigitizer.cxx:991
 AliTRDdigitizer.cxx:992
 AliTRDdigitizer.cxx:993
 AliTRDdigitizer.cxx:994
 AliTRDdigitizer.cxx:995
 AliTRDdigitizer.cxx:996
 AliTRDdigitizer.cxx:997
 AliTRDdigitizer.cxx:998
 AliTRDdigitizer.cxx:999
 AliTRDdigitizer.cxx:1000
 AliTRDdigitizer.cxx:1001
 AliTRDdigitizer.cxx:1002
 AliTRDdigitizer.cxx:1003
 AliTRDdigitizer.cxx:1004
 AliTRDdigitizer.cxx:1005
 AliTRDdigitizer.cxx:1006
 AliTRDdigitizer.cxx:1007
 AliTRDdigitizer.cxx:1008
 AliTRDdigitizer.cxx:1009
 AliTRDdigitizer.cxx:1010
 AliTRDdigitizer.cxx:1011
 AliTRDdigitizer.cxx:1012
 AliTRDdigitizer.cxx:1013
 AliTRDdigitizer.cxx:1014
 AliTRDdigitizer.cxx:1015
 AliTRDdigitizer.cxx:1016
 AliTRDdigitizer.cxx:1017
 AliTRDdigitizer.cxx:1018
 AliTRDdigitizer.cxx:1019
 AliTRDdigitizer.cxx:1020
 AliTRDdigitizer.cxx:1021
 AliTRDdigitizer.cxx:1022
 AliTRDdigitizer.cxx:1023
 AliTRDdigitizer.cxx:1024
 AliTRDdigitizer.cxx:1025
 AliTRDdigitizer.cxx:1026
 AliTRDdigitizer.cxx:1027
 AliTRDdigitizer.cxx:1028
 AliTRDdigitizer.cxx:1029
 AliTRDdigitizer.cxx:1030
 AliTRDdigitizer.cxx:1031
 AliTRDdigitizer.cxx:1032
 AliTRDdigitizer.cxx:1033
 AliTRDdigitizer.cxx:1034
 AliTRDdigitizer.cxx:1035
 AliTRDdigitizer.cxx:1036
 AliTRDdigitizer.cxx:1037
 AliTRDdigitizer.cxx:1038
 AliTRDdigitizer.cxx:1039
 AliTRDdigitizer.cxx:1040
 AliTRDdigitizer.cxx:1041
 AliTRDdigitizer.cxx:1042
 AliTRDdigitizer.cxx:1043
 AliTRDdigitizer.cxx:1044
 AliTRDdigitizer.cxx:1045
 AliTRDdigitizer.cxx:1046
 AliTRDdigitizer.cxx:1047
 AliTRDdigitizer.cxx:1048
 AliTRDdigitizer.cxx:1049
 AliTRDdigitizer.cxx:1050
 AliTRDdigitizer.cxx:1051
 AliTRDdigitizer.cxx:1052
 AliTRDdigitizer.cxx:1053
 AliTRDdigitizer.cxx:1054
 AliTRDdigitizer.cxx:1055
 AliTRDdigitizer.cxx:1056
 AliTRDdigitizer.cxx:1057
 AliTRDdigitizer.cxx:1058
 AliTRDdigitizer.cxx:1059
 AliTRDdigitizer.cxx:1060
 AliTRDdigitizer.cxx:1061
 AliTRDdigitizer.cxx:1062
 AliTRDdigitizer.cxx:1063
 AliTRDdigitizer.cxx:1064
 AliTRDdigitizer.cxx:1065
 AliTRDdigitizer.cxx:1066
 AliTRDdigitizer.cxx:1067
 AliTRDdigitizer.cxx:1068
 AliTRDdigitizer.cxx:1069
 AliTRDdigitizer.cxx:1070
 AliTRDdigitizer.cxx:1071
 AliTRDdigitizer.cxx:1072
 AliTRDdigitizer.cxx:1073
 AliTRDdigitizer.cxx:1074
 AliTRDdigitizer.cxx:1075
 AliTRDdigitizer.cxx:1076
 AliTRDdigitizer.cxx:1077
 AliTRDdigitizer.cxx:1078
 AliTRDdigitizer.cxx:1079
 AliTRDdigitizer.cxx:1080
 AliTRDdigitizer.cxx:1081
 AliTRDdigitizer.cxx:1082
 AliTRDdigitizer.cxx:1083
 AliTRDdigitizer.cxx:1084
 AliTRDdigitizer.cxx:1085
 AliTRDdigitizer.cxx:1086
 AliTRDdigitizer.cxx:1087
 AliTRDdigitizer.cxx:1088
 AliTRDdigitizer.cxx:1089
 AliTRDdigitizer.cxx:1090
 AliTRDdigitizer.cxx:1091
 AliTRDdigitizer.cxx:1092
 AliTRDdigitizer.cxx:1093
 AliTRDdigitizer.cxx:1094
 AliTRDdigitizer.cxx:1095
 AliTRDdigitizer.cxx:1096
 AliTRDdigitizer.cxx:1097
 AliTRDdigitizer.cxx:1098
 AliTRDdigitizer.cxx:1099
 AliTRDdigitizer.cxx:1100
 AliTRDdigitizer.cxx:1101
 AliTRDdigitizer.cxx:1102
 AliTRDdigitizer.cxx:1103
 AliTRDdigitizer.cxx:1104
 AliTRDdigitizer.cxx:1105
 AliTRDdigitizer.cxx:1106
 AliTRDdigitizer.cxx:1107
 AliTRDdigitizer.cxx:1108
 AliTRDdigitizer.cxx:1109
 AliTRDdigitizer.cxx:1110
 AliTRDdigitizer.cxx:1111
 AliTRDdigitizer.cxx:1112
 AliTRDdigitizer.cxx:1113
 AliTRDdigitizer.cxx:1114
 AliTRDdigitizer.cxx:1115
 AliTRDdigitizer.cxx:1116
 AliTRDdigitizer.cxx:1117
 AliTRDdigitizer.cxx:1118
 AliTRDdigitizer.cxx:1119
 AliTRDdigitizer.cxx:1120
 AliTRDdigitizer.cxx:1121
 AliTRDdigitizer.cxx:1122
 AliTRDdigitizer.cxx:1123
 AliTRDdigitizer.cxx:1124
 AliTRDdigitizer.cxx:1125
 AliTRDdigitizer.cxx:1126
 AliTRDdigitizer.cxx:1127
 AliTRDdigitizer.cxx:1128
 AliTRDdigitizer.cxx:1129
 AliTRDdigitizer.cxx:1130
 AliTRDdigitizer.cxx:1131
 AliTRDdigitizer.cxx:1132
 AliTRDdigitizer.cxx:1133
 AliTRDdigitizer.cxx:1134
 AliTRDdigitizer.cxx:1135
 AliTRDdigitizer.cxx:1136
 AliTRDdigitizer.cxx:1137
 AliTRDdigitizer.cxx:1138
 AliTRDdigitizer.cxx:1139
 AliTRDdigitizer.cxx:1140
 AliTRDdigitizer.cxx:1141
 AliTRDdigitizer.cxx:1142
 AliTRDdigitizer.cxx:1143
 AliTRDdigitizer.cxx:1144
 AliTRDdigitizer.cxx:1145
 AliTRDdigitizer.cxx:1146
 AliTRDdigitizer.cxx:1147
 AliTRDdigitizer.cxx:1148
 AliTRDdigitizer.cxx:1149
 AliTRDdigitizer.cxx:1150
 AliTRDdigitizer.cxx:1151
 AliTRDdigitizer.cxx:1152
 AliTRDdigitizer.cxx:1153
 AliTRDdigitizer.cxx:1154
 AliTRDdigitizer.cxx:1155
 AliTRDdigitizer.cxx:1156
 AliTRDdigitizer.cxx:1157
 AliTRDdigitizer.cxx:1158
 AliTRDdigitizer.cxx:1159
 AliTRDdigitizer.cxx:1160
 AliTRDdigitizer.cxx:1161
 AliTRDdigitizer.cxx:1162
 AliTRDdigitizer.cxx:1163
 AliTRDdigitizer.cxx:1164
 AliTRDdigitizer.cxx:1165
 AliTRDdigitizer.cxx:1166
 AliTRDdigitizer.cxx:1167
 AliTRDdigitizer.cxx:1168
 AliTRDdigitizer.cxx:1169
 AliTRDdigitizer.cxx:1170
 AliTRDdigitizer.cxx:1171
 AliTRDdigitizer.cxx:1172
 AliTRDdigitizer.cxx:1173
 AliTRDdigitizer.cxx:1174
 AliTRDdigitizer.cxx:1175
 AliTRDdigitizer.cxx:1176
 AliTRDdigitizer.cxx:1177
 AliTRDdigitizer.cxx:1178
 AliTRDdigitizer.cxx:1179
 AliTRDdigitizer.cxx:1180
 AliTRDdigitizer.cxx:1181
 AliTRDdigitizer.cxx:1182
 AliTRDdigitizer.cxx:1183
 AliTRDdigitizer.cxx:1184
 AliTRDdigitizer.cxx:1185
 AliTRDdigitizer.cxx:1186
 AliTRDdigitizer.cxx:1187
 AliTRDdigitizer.cxx:1188
 AliTRDdigitizer.cxx:1189
 AliTRDdigitizer.cxx:1190
 AliTRDdigitizer.cxx:1191
 AliTRDdigitizer.cxx:1192
 AliTRDdigitizer.cxx:1193
 AliTRDdigitizer.cxx:1194
 AliTRDdigitizer.cxx:1195
 AliTRDdigitizer.cxx:1196
 AliTRDdigitizer.cxx:1197
 AliTRDdigitizer.cxx:1198
 AliTRDdigitizer.cxx:1199
 AliTRDdigitizer.cxx:1200
 AliTRDdigitizer.cxx:1201
 AliTRDdigitizer.cxx:1202
 AliTRDdigitizer.cxx:1203
 AliTRDdigitizer.cxx:1204
 AliTRDdigitizer.cxx:1205
 AliTRDdigitizer.cxx:1206
 AliTRDdigitizer.cxx:1207
 AliTRDdigitizer.cxx:1208
 AliTRDdigitizer.cxx:1209
 AliTRDdigitizer.cxx:1210
 AliTRDdigitizer.cxx:1211
 AliTRDdigitizer.cxx:1212
 AliTRDdigitizer.cxx:1213
 AliTRDdigitizer.cxx:1214
 AliTRDdigitizer.cxx:1215
 AliTRDdigitizer.cxx:1216
 AliTRDdigitizer.cxx:1217
 AliTRDdigitizer.cxx:1218
 AliTRDdigitizer.cxx:1219
 AliTRDdigitizer.cxx:1220
 AliTRDdigitizer.cxx:1221
 AliTRDdigitizer.cxx:1222
 AliTRDdigitizer.cxx:1223
 AliTRDdigitizer.cxx:1224
 AliTRDdigitizer.cxx:1225
 AliTRDdigitizer.cxx:1226
 AliTRDdigitizer.cxx:1227
 AliTRDdigitizer.cxx:1228
 AliTRDdigitizer.cxx:1229
 AliTRDdigitizer.cxx:1230
 AliTRDdigitizer.cxx:1231
 AliTRDdigitizer.cxx:1232
 AliTRDdigitizer.cxx:1233
 AliTRDdigitizer.cxx:1234
 AliTRDdigitizer.cxx:1235
 AliTRDdigitizer.cxx:1236
 AliTRDdigitizer.cxx:1237
 AliTRDdigitizer.cxx:1238
 AliTRDdigitizer.cxx:1239
 AliTRDdigitizer.cxx:1240
 AliTRDdigitizer.cxx:1241
 AliTRDdigitizer.cxx:1242
 AliTRDdigitizer.cxx:1243
 AliTRDdigitizer.cxx:1244
 AliTRDdigitizer.cxx:1245
 AliTRDdigitizer.cxx:1246
 AliTRDdigitizer.cxx:1247
 AliTRDdigitizer.cxx:1248
 AliTRDdigitizer.cxx:1249
 AliTRDdigitizer.cxx:1250
 AliTRDdigitizer.cxx:1251
 AliTRDdigitizer.cxx:1252
 AliTRDdigitizer.cxx:1253
 AliTRDdigitizer.cxx:1254
 AliTRDdigitizer.cxx:1255
 AliTRDdigitizer.cxx:1256
 AliTRDdigitizer.cxx:1257
 AliTRDdigitizer.cxx:1258
 AliTRDdigitizer.cxx:1259
 AliTRDdigitizer.cxx:1260
 AliTRDdigitizer.cxx:1261
 AliTRDdigitizer.cxx:1262
 AliTRDdigitizer.cxx:1263
 AliTRDdigitizer.cxx:1264
 AliTRDdigitizer.cxx:1265
 AliTRDdigitizer.cxx:1266
 AliTRDdigitizer.cxx:1267
 AliTRDdigitizer.cxx:1268
 AliTRDdigitizer.cxx:1269
 AliTRDdigitizer.cxx:1270
 AliTRDdigitizer.cxx:1271
 AliTRDdigitizer.cxx:1272
 AliTRDdigitizer.cxx:1273
 AliTRDdigitizer.cxx:1274
 AliTRDdigitizer.cxx:1275
 AliTRDdigitizer.cxx:1276
 AliTRDdigitizer.cxx:1277
 AliTRDdigitizer.cxx:1278
 AliTRDdigitizer.cxx:1279
 AliTRDdigitizer.cxx:1280
 AliTRDdigitizer.cxx:1281
 AliTRDdigitizer.cxx:1282
 AliTRDdigitizer.cxx:1283
 AliTRDdigitizer.cxx:1284
 AliTRDdigitizer.cxx:1285
 AliTRDdigitizer.cxx:1286
 AliTRDdigitizer.cxx:1287
 AliTRDdigitizer.cxx:1288
 AliTRDdigitizer.cxx:1289
 AliTRDdigitizer.cxx:1290
 AliTRDdigitizer.cxx:1291
 AliTRDdigitizer.cxx:1292
 AliTRDdigitizer.cxx:1293
 AliTRDdigitizer.cxx:1294
 AliTRDdigitizer.cxx:1295
 AliTRDdigitizer.cxx:1296
 AliTRDdigitizer.cxx:1297
 AliTRDdigitizer.cxx:1298
 AliTRDdigitizer.cxx:1299
 AliTRDdigitizer.cxx:1300
 AliTRDdigitizer.cxx:1301
 AliTRDdigitizer.cxx:1302
 AliTRDdigitizer.cxx:1303
 AliTRDdigitizer.cxx:1304
 AliTRDdigitizer.cxx:1305
 AliTRDdigitizer.cxx:1306
 AliTRDdigitizer.cxx:1307
 AliTRDdigitizer.cxx:1308
 AliTRDdigitizer.cxx:1309
 AliTRDdigitizer.cxx:1310
 AliTRDdigitizer.cxx:1311
 AliTRDdigitizer.cxx:1312
 AliTRDdigitizer.cxx:1313
 AliTRDdigitizer.cxx:1314
 AliTRDdigitizer.cxx:1315
 AliTRDdigitizer.cxx:1316
 AliTRDdigitizer.cxx:1317
 AliTRDdigitizer.cxx:1318
 AliTRDdigitizer.cxx:1319
 AliTRDdigitizer.cxx:1320
 AliTRDdigitizer.cxx:1321
 AliTRDdigitizer.cxx:1322
 AliTRDdigitizer.cxx:1323
 AliTRDdigitizer.cxx:1324
 AliTRDdigitizer.cxx:1325
 AliTRDdigitizer.cxx:1326
 AliTRDdigitizer.cxx:1327
 AliTRDdigitizer.cxx:1328
 AliTRDdigitizer.cxx:1329
 AliTRDdigitizer.cxx:1330
 AliTRDdigitizer.cxx:1331
 AliTRDdigitizer.cxx:1332
 AliTRDdigitizer.cxx:1333
 AliTRDdigitizer.cxx:1334
 AliTRDdigitizer.cxx:1335
 AliTRDdigitizer.cxx:1336
 AliTRDdigitizer.cxx:1337
 AliTRDdigitizer.cxx:1338
 AliTRDdigitizer.cxx:1339
 AliTRDdigitizer.cxx:1340
 AliTRDdigitizer.cxx:1341
 AliTRDdigitizer.cxx:1342
 AliTRDdigitizer.cxx:1343
 AliTRDdigitizer.cxx:1344
 AliTRDdigitizer.cxx:1345
 AliTRDdigitizer.cxx:1346
 AliTRDdigitizer.cxx:1347
 AliTRDdigitizer.cxx:1348
 AliTRDdigitizer.cxx:1349
 AliTRDdigitizer.cxx:1350
 AliTRDdigitizer.cxx:1351
 AliTRDdigitizer.cxx:1352
 AliTRDdigitizer.cxx:1353
 AliTRDdigitizer.cxx:1354
 AliTRDdigitizer.cxx:1355
 AliTRDdigitizer.cxx:1356
 AliTRDdigitizer.cxx:1357
 AliTRDdigitizer.cxx:1358
 AliTRDdigitizer.cxx:1359
 AliTRDdigitizer.cxx:1360
 AliTRDdigitizer.cxx:1361
 AliTRDdigitizer.cxx:1362
 AliTRDdigitizer.cxx:1363
 AliTRDdigitizer.cxx:1364
 AliTRDdigitizer.cxx:1365
 AliTRDdigitizer.cxx:1366
 AliTRDdigitizer.cxx:1367
 AliTRDdigitizer.cxx:1368
 AliTRDdigitizer.cxx:1369
 AliTRDdigitizer.cxx:1370
 AliTRDdigitizer.cxx:1371
 AliTRDdigitizer.cxx:1372
 AliTRDdigitizer.cxx:1373
 AliTRDdigitizer.cxx:1374
 AliTRDdigitizer.cxx:1375
 AliTRDdigitizer.cxx:1376
 AliTRDdigitizer.cxx:1377
 AliTRDdigitizer.cxx:1378
 AliTRDdigitizer.cxx:1379
 AliTRDdigitizer.cxx:1380
 AliTRDdigitizer.cxx:1381
 AliTRDdigitizer.cxx:1382
 AliTRDdigitizer.cxx:1383
 AliTRDdigitizer.cxx:1384
 AliTRDdigitizer.cxx:1385
 AliTRDdigitizer.cxx:1386
 AliTRDdigitizer.cxx:1387
 AliTRDdigitizer.cxx:1388
 AliTRDdigitizer.cxx:1389
 AliTRDdigitizer.cxx:1390
 AliTRDdigitizer.cxx:1391
 AliTRDdigitizer.cxx:1392
 AliTRDdigitizer.cxx:1393
 AliTRDdigitizer.cxx:1394
 AliTRDdigitizer.cxx:1395
 AliTRDdigitizer.cxx:1396
 AliTRDdigitizer.cxx:1397
 AliTRDdigitizer.cxx:1398
 AliTRDdigitizer.cxx:1399
 AliTRDdigitizer.cxx:1400
 AliTRDdigitizer.cxx:1401
 AliTRDdigitizer.cxx:1402
 AliTRDdigitizer.cxx:1403
 AliTRDdigitizer.cxx:1404
 AliTRDdigitizer.cxx:1405
 AliTRDdigitizer.cxx:1406
 AliTRDdigitizer.cxx:1407
 AliTRDdigitizer.cxx:1408
 AliTRDdigitizer.cxx:1409
 AliTRDdigitizer.cxx:1410
 AliTRDdigitizer.cxx:1411
 AliTRDdigitizer.cxx:1412
 AliTRDdigitizer.cxx:1413
 AliTRDdigitizer.cxx:1414
 AliTRDdigitizer.cxx:1415
 AliTRDdigitizer.cxx:1416
 AliTRDdigitizer.cxx:1417
 AliTRDdigitizer.cxx:1418
 AliTRDdigitizer.cxx:1419
 AliTRDdigitizer.cxx:1420
 AliTRDdigitizer.cxx:1421
 AliTRDdigitizer.cxx:1422
 AliTRDdigitizer.cxx:1423
 AliTRDdigitizer.cxx:1424
 AliTRDdigitizer.cxx:1425
 AliTRDdigitizer.cxx:1426
 AliTRDdigitizer.cxx:1427
 AliTRDdigitizer.cxx:1428
 AliTRDdigitizer.cxx:1429
 AliTRDdigitizer.cxx:1430
 AliTRDdigitizer.cxx:1431
 AliTRDdigitizer.cxx:1432
 AliTRDdigitizer.cxx:1433
 AliTRDdigitizer.cxx:1434
 AliTRDdigitizer.cxx:1435
 AliTRDdigitizer.cxx:1436
 AliTRDdigitizer.cxx:1437
 AliTRDdigitizer.cxx:1438
 AliTRDdigitizer.cxx:1439
 AliTRDdigitizer.cxx:1440
 AliTRDdigitizer.cxx:1441
 AliTRDdigitizer.cxx:1442
 AliTRDdigitizer.cxx:1443
 AliTRDdigitizer.cxx:1444
 AliTRDdigitizer.cxx:1445
 AliTRDdigitizer.cxx:1446
 AliTRDdigitizer.cxx:1447
 AliTRDdigitizer.cxx:1448
 AliTRDdigitizer.cxx:1449
 AliTRDdigitizer.cxx:1450
 AliTRDdigitizer.cxx:1451
 AliTRDdigitizer.cxx:1452
 AliTRDdigitizer.cxx:1453
 AliTRDdigitizer.cxx:1454
 AliTRDdigitizer.cxx:1455
 AliTRDdigitizer.cxx:1456
 AliTRDdigitizer.cxx:1457
 AliTRDdigitizer.cxx:1458
 AliTRDdigitizer.cxx:1459
 AliTRDdigitizer.cxx:1460
 AliTRDdigitizer.cxx:1461
 AliTRDdigitizer.cxx:1462
 AliTRDdigitizer.cxx:1463
 AliTRDdigitizer.cxx:1464
 AliTRDdigitizer.cxx:1465
 AliTRDdigitizer.cxx:1466
 AliTRDdigitizer.cxx:1467
 AliTRDdigitizer.cxx:1468
 AliTRDdigitizer.cxx:1469
 AliTRDdigitizer.cxx:1470
 AliTRDdigitizer.cxx:1471
 AliTRDdigitizer.cxx:1472
 AliTRDdigitizer.cxx:1473
 AliTRDdigitizer.cxx:1474
 AliTRDdigitizer.cxx:1475
 AliTRDdigitizer.cxx:1476
 AliTRDdigitizer.cxx:1477
 AliTRDdigitizer.cxx:1478
 AliTRDdigitizer.cxx:1479
 AliTRDdigitizer.cxx:1480
 AliTRDdigitizer.cxx:1481
 AliTRDdigitizer.cxx:1482
 AliTRDdigitizer.cxx:1483
 AliTRDdigitizer.cxx:1484
 AliTRDdigitizer.cxx:1485
 AliTRDdigitizer.cxx:1486
 AliTRDdigitizer.cxx:1487
 AliTRDdigitizer.cxx:1488
 AliTRDdigitizer.cxx:1489
 AliTRDdigitizer.cxx:1490
 AliTRDdigitizer.cxx:1491
 AliTRDdigitizer.cxx:1492
 AliTRDdigitizer.cxx:1493
 AliTRDdigitizer.cxx:1494
 AliTRDdigitizer.cxx:1495
 AliTRDdigitizer.cxx:1496
 AliTRDdigitizer.cxx:1497
 AliTRDdigitizer.cxx:1498
 AliTRDdigitizer.cxx:1499
 AliTRDdigitizer.cxx:1500
 AliTRDdigitizer.cxx:1501
 AliTRDdigitizer.cxx:1502
 AliTRDdigitizer.cxx:1503
 AliTRDdigitizer.cxx:1504
 AliTRDdigitizer.cxx:1505
 AliTRDdigitizer.cxx:1506
 AliTRDdigitizer.cxx:1507
 AliTRDdigitizer.cxx:1508
 AliTRDdigitizer.cxx:1509
 AliTRDdigitizer.cxx:1510
 AliTRDdigitizer.cxx:1511
 AliTRDdigitizer.cxx:1512
 AliTRDdigitizer.cxx:1513
 AliTRDdigitizer.cxx:1514
 AliTRDdigitizer.cxx:1515
 AliTRDdigitizer.cxx:1516
 AliTRDdigitizer.cxx:1517
 AliTRDdigitizer.cxx:1518
 AliTRDdigitizer.cxx:1519
 AliTRDdigitizer.cxx:1520
 AliTRDdigitizer.cxx:1521
 AliTRDdigitizer.cxx:1522
 AliTRDdigitizer.cxx:1523
 AliTRDdigitizer.cxx:1524
 AliTRDdigitizer.cxx:1525
 AliTRDdigitizer.cxx:1526
 AliTRDdigitizer.cxx:1527
 AliTRDdigitizer.cxx:1528
 AliTRDdigitizer.cxx:1529
 AliTRDdigitizer.cxx:1530
 AliTRDdigitizer.cxx:1531
 AliTRDdigitizer.cxx:1532
 AliTRDdigitizer.cxx:1533
 AliTRDdigitizer.cxx:1534
 AliTRDdigitizer.cxx:1535
 AliTRDdigitizer.cxx:1536
 AliTRDdigitizer.cxx:1537
 AliTRDdigitizer.cxx:1538
 AliTRDdigitizer.cxx:1539
 AliTRDdigitizer.cxx:1540
 AliTRDdigitizer.cxx:1541
 AliTRDdigitizer.cxx:1542
 AliTRDdigitizer.cxx:1543
 AliTRDdigitizer.cxx:1544
 AliTRDdigitizer.cxx:1545
 AliTRDdigitizer.cxx:1546
 AliTRDdigitizer.cxx:1547
 AliTRDdigitizer.cxx:1548
 AliTRDdigitizer.cxx:1549
 AliTRDdigitizer.cxx:1550
 AliTRDdigitizer.cxx:1551
 AliTRDdigitizer.cxx:1552
 AliTRDdigitizer.cxx:1553
 AliTRDdigitizer.cxx:1554
 AliTRDdigitizer.cxx:1555
 AliTRDdigitizer.cxx:1556
 AliTRDdigitizer.cxx:1557
 AliTRDdigitizer.cxx:1558
 AliTRDdigitizer.cxx:1559
 AliTRDdigitizer.cxx:1560
 AliTRDdigitizer.cxx:1561
 AliTRDdigitizer.cxx:1562
 AliTRDdigitizer.cxx:1563
 AliTRDdigitizer.cxx:1564
 AliTRDdigitizer.cxx:1565
 AliTRDdigitizer.cxx:1566
 AliTRDdigitizer.cxx:1567
 AliTRDdigitizer.cxx:1568
 AliTRDdigitizer.cxx:1569
 AliTRDdigitizer.cxx:1570
 AliTRDdigitizer.cxx:1571
 AliTRDdigitizer.cxx:1572
 AliTRDdigitizer.cxx:1573
 AliTRDdigitizer.cxx:1574
 AliTRDdigitizer.cxx:1575
 AliTRDdigitizer.cxx:1576
 AliTRDdigitizer.cxx:1577
 AliTRDdigitizer.cxx:1578
 AliTRDdigitizer.cxx:1579
 AliTRDdigitizer.cxx:1580
 AliTRDdigitizer.cxx:1581
 AliTRDdigitizer.cxx:1582
 AliTRDdigitizer.cxx:1583
 AliTRDdigitizer.cxx:1584
 AliTRDdigitizer.cxx:1585
 AliTRDdigitizer.cxx:1586
 AliTRDdigitizer.cxx:1587
 AliTRDdigitizer.cxx:1588
 AliTRDdigitizer.cxx:1589
 AliTRDdigitizer.cxx:1590
 AliTRDdigitizer.cxx:1591
 AliTRDdigitizer.cxx:1592
 AliTRDdigitizer.cxx:1593
 AliTRDdigitizer.cxx:1594
 AliTRDdigitizer.cxx:1595
 AliTRDdigitizer.cxx:1596
 AliTRDdigitizer.cxx:1597
 AliTRDdigitizer.cxx:1598
 AliTRDdigitizer.cxx:1599
 AliTRDdigitizer.cxx:1600
 AliTRDdigitizer.cxx:1601
 AliTRDdigitizer.cxx:1602
 AliTRDdigitizer.cxx:1603
 AliTRDdigitizer.cxx:1604
 AliTRDdigitizer.cxx:1605
 AliTRDdigitizer.cxx:1606
 AliTRDdigitizer.cxx:1607
 AliTRDdigitizer.cxx:1608
 AliTRDdigitizer.cxx:1609
 AliTRDdigitizer.cxx:1610
 AliTRDdigitizer.cxx:1611
 AliTRDdigitizer.cxx:1612
 AliTRDdigitizer.cxx:1613
 AliTRDdigitizer.cxx:1614
 AliTRDdigitizer.cxx:1615
 AliTRDdigitizer.cxx:1616
 AliTRDdigitizer.cxx:1617
 AliTRDdigitizer.cxx:1618
 AliTRDdigitizer.cxx:1619
 AliTRDdigitizer.cxx:1620
 AliTRDdigitizer.cxx:1621
 AliTRDdigitizer.cxx:1622
 AliTRDdigitizer.cxx:1623
 AliTRDdigitizer.cxx:1624
 AliTRDdigitizer.cxx:1625
 AliTRDdigitizer.cxx:1626
 AliTRDdigitizer.cxx:1627
 AliTRDdigitizer.cxx:1628
 AliTRDdigitizer.cxx:1629
 AliTRDdigitizer.cxx:1630
 AliTRDdigitizer.cxx:1631
 AliTRDdigitizer.cxx:1632
 AliTRDdigitizer.cxx:1633
 AliTRDdigitizer.cxx:1634
 AliTRDdigitizer.cxx:1635
 AliTRDdigitizer.cxx:1636
 AliTRDdigitizer.cxx:1637
 AliTRDdigitizer.cxx:1638
 AliTRDdigitizer.cxx:1639
 AliTRDdigitizer.cxx:1640
 AliTRDdigitizer.cxx:1641
 AliTRDdigitizer.cxx:1642
 AliTRDdigitizer.cxx:1643
 AliTRDdigitizer.cxx:1644
 AliTRDdigitizer.cxx:1645
 AliTRDdigitizer.cxx:1646
 AliTRDdigitizer.cxx:1647
 AliTRDdigitizer.cxx:1648
 AliTRDdigitizer.cxx:1649
 AliTRDdigitizer.cxx:1650
 AliTRDdigitizer.cxx:1651
 AliTRDdigitizer.cxx:1652
 AliTRDdigitizer.cxx:1653
 AliTRDdigitizer.cxx:1654
 AliTRDdigitizer.cxx:1655
 AliTRDdigitizer.cxx:1656
 AliTRDdigitizer.cxx:1657
 AliTRDdigitizer.cxx:1658
 AliTRDdigitizer.cxx:1659
 AliTRDdigitizer.cxx:1660
 AliTRDdigitizer.cxx:1661
 AliTRDdigitizer.cxx:1662
 AliTRDdigitizer.cxx:1663
 AliTRDdigitizer.cxx:1664
 AliTRDdigitizer.cxx:1665
 AliTRDdigitizer.cxx:1666
 AliTRDdigitizer.cxx:1667
 AliTRDdigitizer.cxx:1668
 AliTRDdigitizer.cxx:1669
 AliTRDdigitizer.cxx:1670
 AliTRDdigitizer.cxx:1671
 AliTRDdigitizer.cxx:1672
 AliTRDdigitizer.cxx:1673
 AliTRDdigitizer.cxx:1674
 AliTRDdigitizer.cxx:1675
 AliTRDdigitizer.cxx:1676
 AliTRDdigitizer.cxx:1677
 AliTRDdigitizer.cxx:1678
 AliTRDdigitizer.cxx:1679
 AliTRDdigitizer.cxx:1680
 AliTRDdigitizer.cxx:1681
 AliTRDdigitizer.cxx:1682
 AliTRDdigitizer.cxx:1683
 AliTRDdigitizer.cxx:1684
 AliTRDdigitizer.cxx:1685
 AliTRDdigitizer.cxx:1686
 AliTRDdigitizer.cxx:1687
 AliTRDdigitizer.cxx:1688
 AliTRDdigitizer.cxx:1689
 AliTRDdigitizer.cxx:1690
 AliTRDdigitizer.cxx:1691
 AliTRDdigitizer.cxx:1692
 AliTRDdigitizer.cxx:1693
 AliTRDdigitizer.cxx:1694
 AliTRDdigitizer.cxx:1695
 AliTRDdigitizer.cxx:1696
 AliTRDdigitizer.cxx:1697
 AliTRDdigitizer.cxx:1698
 AliTRDdigitizer.cxx:1699
 AliTRDdigitizer.cxx:1700
 AliTRDdigitizer.cxx:1701
 AliTRDdigitizer.cxx:1702
 AliTRDdigitizer.cxx:1703
 AliTRDdigitizer.cxx:1704
 AliTRDdigitizer.cxx:1705
 AliTRDdigitizer.cxx:1706
 AliTRDdigitizer.cxx:1707
 AliTRDdigitizer.cxx:1708
 AliTRDdigitizer.cxx:1709
 AliTRDdigitizer.cxx:1710
 AliTRDdigitizer.cxx:1711
 AliTRDdigitizer.cxx:1712
 AliTRDdigitizer.cxx:1713
 AliTRDdigitizer.cxx:1714
 AliTRDdigitizer.cxx:1715
 AliTRDdigitizer.cxx:1716
 AliTRDdigitizer.cxx:1717
 AliTRDdigitizer.cxx:1718
 AliTRDdigitizer.cxx:1719
 AliTRDdigitizer.cxx:1720
 AliTRDdigitizer.cxx:1721
 AliTRDdigitizer.cxx:1722
 AliTRDdigitizer.cxx:1723
 AliTRDdigitizer.cxx:1724
 AliTRDdigitizer.cxx:1725
 AliTRDdigitizer.cxx:1726
 AliTRDdigitizer.cxx:1727
 AliTRDdigitizer.cxx:1728
 AliTRDdigitizer.cxx:1729
 AliTRDdigitizer.cxx:1730
 AliTRDdigitizer.cxx:1731
 AliTRDdigitizer.cxx:1732
 AliTRDdigitizer.cxx:1733
 AliTRDdigitizer.cxx:1734
 AliTRDdigitizer.cxx:1735
 AliTRDdigitizer.cxx:1736
 AliTRDdigitizer.cxx:1737
 AliTRDdigitizer.cxx:1738
 AliTRDdigitizer.cxx:1739
 AliTRDdigitizer.cxx:1740
 AliTRDdigitizer.cxx:1741
 AliTRDdigitizer.cxx:1742
 AliTRDdigitizer.cxx:1743
 AliTRDdigitizer.cxx:1744
 AliTRDdigitizer.cxx:1745
 AliTRDdigitizer.cxx:1746
 AliTRDdigitizer.cxx:1747
 AliTRDdigitizer.cxx:1748
 AliTRDdigitizer.cxx:1749
 AliTRDdigitizer.cxx:1750
 AliTRDdigitizer.cxx:1751
 AliTRDdigitizer.cxx:1752
 AliTRDdigitizer.cxx:1753
 AliTRDdigitizer.cxx:1754
 AliTRDdigitizer.cxx:1755
 AliTRDdigitizer.cxx:1756
 AliTRDdigitizer.cxx:1757
 AliTRDdigitizer.cxx:1758
 AliTRDdigitizer.cxx:1759
 AliTRDdigitizer.cxx:1760
 AliTRDdigitizer.cxx:1761
 AliTRDdigitizer.cxx:1762
 AliTRDdigitizer.cxx:1763
 AliTRDdigitizer.cxx:1764
 AliTRDdigitizer.cxx:1765
 AliTRDdigitizer.cxx:1766
 AliTRDdigitizer.cxx:1767
 AliTRDdigitizer.cxx:1768
 AliTRDdigitizer.cxx:1769
 AliTRDdigitizer.cxx:1770
 AliTRDdigitizer.cxx:1771
 AliTRDdigitizer.cxx:1772
 AliTRDdigitizer.cxx:1773
 AliTRDdigitizer.cxx:1774
 AliTRDdigitizer.cxx:1775
 AliTRDdigitizer.cxx:1776
 AliTRDdigitizer.cxx:1777
 AliTRDdigitizer.cxx:1778
 AliTRDdigitizer.cxx:1779
 AliTRDdigitizer.cxx:1780
 AliTRDdigitizer.cxx:1781
 AliTRDdigitizer.cxx:1782
 AliTRDdigitizer.cxx:1783
 AliTRDdigitizer.cxx:1784
 AliTRDdigitizer.cxx:1785
 AliTRDdigitizer.cxx:1786
 AliTRDdigitizer.cxx:1787
 AliTRDdigitizer.cxx:1788
 AliTRDdigitizer.cxx:1789
 AliTRDdigitizer.cxx:1790
 AliTRDdigitizer.cxx:1791
 AliTRDdigitizer.cxx:1792
 AliTRDdigitizer.cxx:1793
 AliTRDdigitizer.cxx:1794
 AliTRDdigitizer.cxx:1795
 AliTRDdigitizer.cxx:1796
 AliTRDdigitizer.cxx:1797
 AliTRDdigitizer.cxx:1798
 AliTRDdigitizer.cxx:1799
 AliTRDdigitizer.cxx:1800
 AliTRDdigitizer.cxx:1801
 AliTRDdigitizer.cxx:1802
 AliTRDdigitizer.cxx:1803
 AliTRDdigitizer.cxx:1804
 AliTRDdigitizer.cxx:1805
 AliTRDdigitizer.cxx:1806
 AliTRDdigitizer.cxx:1807
 AliTRDdigitizer.cxx:1808
 AliTRDdigitizer.cxx:1809
 AliTRDdigitizer.cxx:1810
 AliTRDdigitizer.cxx:1811
 AliTRDdigitizer.cxx:1812
 AliTRDdigitizer.cxx:1813
 AliTRDdigitizer.cxx:1814
 AliTRDdigitizer.cxx:1815
 AliTRDdigitizer.cxx:1816
 AliTRDdigitizer.cxx:1817
 AliTRDdigitizer.cxx:1818
 AliTRDdigitizer.cxx:1819
 AliTRDdigitizer.cxx:1820
 AliTRDdigitizer.cxx:1821
 AliTRDdigitizer.cxx:1822
 AliTRDdigitizer.cxx:1823
 AliTRDdigitizer.cxx:1824
 AliTRDdigitizer.cxx:1825
 AliTRDdigitizer.cxx:1826
 AliTRDdigitizer.cxx:1827
 AliTRDdigitizer.cxx:1828
 AliTRDdigitizer.cxx:1829
 AliTRDdigitizer.cxx:1830
 AliTRDdigitizer.cxx:1831
 AliTRDdigitizer.cxx:1832
 AliTRDdigitizer.cxx:1833
 AliTRDdigitizer.cxx:1834
 AliTRDdigitizer.cxx:1835
 AliTRDdigitizer.cxx:1836
 AliTRDdigitizer.cxx:1837
 AliTRDdigitizer.cxx:1838
 AliTRDdigitizer.cxx:1839
 AliTRDdigitizer.cxx:1840
 AliTRDdigitizer.cxx:1841
 AliTRDdigitizer.cxx:1842
 AliTRDdigitizer.cxx:1843
 AliTRDdigitizer.cxx:1844
 AliTRDdigitizer.cxx:1845
 AliTRDdigitizer.cxx:1846
 AliTRDdigitizer.cxx:1847
 AliTRDdigitizer.cxx:1848
 AliTRDdigitizer.cxx:1849
 AliTRDdigitizer.cxx:1850
 AliTRDdigitizer.cxx:1851
 AliTRDdigitizer.cxx:1852
 AliTRDdigitizer.cxx:1853
 AliTRDdigitizer.cxx:1854
 AliTRDdigitizer.cxx:1855
 AliTRDdigitizer.cxx:1856
 AliTRDdigitizer.cxx:1857
 AliTRDdigitizer.cxx:1858
 AliTRDdigitizer.cxx:1859
 AliTRDdigitizer.cxx:1860
 AliTRDdigitizer.cxx:1861
 AliTRDdigitizer.cxx:1862
 AliTRDdigitizer.cxx:1863
 AliTRDdigitizer.cxx:1864
 AliTRDdigitizer.cxx:1865
 AliTRDdigitizer.cxx:1866
 AliTRDdigitizer.cxx:1867
 AliTRDdigitizer.cxx:1868
 AliTRDdigitizer.cxx:1869
 AliTRDdigitizer.cxx:1870
 AliTRDdigitizer.cxx:1871
 AliTRDdigitizer.cxx:1872
 AliTRDdigitizer.cxx:1873
 AliTRDdigitizer.cxx:1874
 AliTRDdigitizer.cxx:1875
 AliTRDdigitizer.cxx:1876
 AliTRDdigitizer.cxx:1877
 AliTRDdigitizer.cxx:1878
 AliTRDdigitizer.cxx:1879
 AliTRDdigitizer.cxx:1880
 AliTRDdigitizer.cxx:1881
 AliTRDdigitizer.cxx:1882
 AliTRDdigitizer.cxx:1883
 AliTRDdigitizer.cxx:1884
 AliTRDdigitizer.cxx:1885
 AliTRDdigitizer.cxx:1886
 AliTRDdigitizer.cxx:1887
 AliTRDdigitizer.cxx:1888
 AliTRDdigitizer.cxx:1889
 AliTRDdigitizer.cxx:1890
 AliTRDdigitizer.cxx:1891
 AliTRDdigitizer.cxx:1892
 AliTRDdigitizer.cxx:1893
 AliTRDdigitizer.cxx:1894
 AliTRDdigitizer.cxx:1895
 AliTRDdigitizer.cxx:1896
 AliTRDdigitizer.cxx:1897
 AliTRDdigitizer.cxx:1898
 AliTRDdigitizer.cxx:1899
 AliTRDdigitizer.cxx:1900
 AliTRDdigitizer.cxx:1901
 AliTRDdigitizer.cxx:1902
 AliTRDdigitizer.cxx:1903
 AliTRDdigitizer.cxx:1904
 AliTRDdigitizer.cxx:1905
 AliTRDdigitizer.cxx:1906
 AliTRDdigitizer.cxx:1907
 AliTRDdigitizer.cxx:1908
 AliTRDdigitizer.cxx:1909
 AliTRDdigitizer.cxx:1910
 AliTRDdigitizer.cxx:1911
 AliTRDdigitizer.cxx:1912
 AliTRDdigitizer.cxx:1913
 AliTRDdigitizer.cxx:1914
 AliTRDdigitizer.cxx:1915
 AliTRDdigitizer.cxx:1916
 AliTRDdigitizer.cxx:1917
 AliTRDdigitizer.cxx:1918
 AliTRDdigitizer.cxx:1919
 AliTRDdigitizer.cxx:1920
 AliTRDdigitizer.cxx:1921
 AliTRDdigitizer.cxx:1922
 AliTRDdigitizer.cxx:1923
 AliTRDdigitizer.cxx:1924
 AliTRDdigitizer.cxx:1925
 AliTRDdigitizer.cxx:1926
 AliTRDdigitizer.cxx:1927
 AliTRDdigitizer.cxx:1928
 AliTRDdigitizer.cxx:1929
 AliTRDdigitizer.cxx:1930
 AliTRDdigitizer.cxx:1931
 AliTRDdigitizer.cxx:1932
 AliTRDdigitizer.cxx:1933
 AliTRDdigitizer.cxx:1934
 AliTRDdigitizer.cxx:1935
 AliTRDdigitizer.cxx:1936
 AliTRDdigitizer.cxx:1937
 AliTRDdigitizer.cxx:1938
 AliTRDdigitizer.cxx:1939
 AliTRDdigitizer.cxx:1940
 AliTRDdigitizer.cxx:1941
 AliTRDdigitizer.cxx:1942
 AliTRDdigitizer.cxx:1943
 AliTRDdigitizer.cxx:1944
 AliTRDdigitizer.cxx:1945
 AliTRDdigitizer.cxx:1946
 AliTRDdigitizer.cxx:1947
 AliTRDdigitizer.cxx:1948
 AliTRDdigitizer.cxx:1949
 AliTRDdigitizer.cxx:1950
 AliTRDdigitizer.cxx:1951
 AliTRDdigitizer.cxx:1952
 AliTRDdigitizer.cxx:1953
 AliTRDdigitizer.cxx:1954
 AliTRDdigitizer.cxx:1955
 AliTRDdigitizer.cxx:1956
 AliTRDdigitizer.cxx:1957
 AliTRDdigitizer.cxx:1958
 AliTRDdigitizer.cxx:1959
 AliTRDdigitizer.cxx:1960
 AliTRDdigitizer.cxx:1961
 AliTRDdigitizer.cxx:1962
 AliTRDdigitizer.cxx:1963
 AliTRDdigitizer.cxx:1964
 AliTRDdigitizer.cxx:1965
 AliTRDdigitizer.cxx:1966
 AliTRDdigitizer.cxx:1967
 AliTRDdigitizer.cxx:1968
 AliTRDdigitizer.cxx:1969
 AliTRDdigitizer.cxx:1970
 AliTRDdigitizer.cxx:1971
 AliTRDdigitizer.cxx:1972
 AliTRDdigitizer.cxx:1973
 AliTRDdigitizer.cxx:1974
 AliTRDdigitizer.cxx:1975
 AliTRDdigitizer.cxx:1976
 AliTRDdigitizer.cxx:1977
 AliTRDdigitizer.cxx:1978
 AliTRDdigitizer.cxx:1979
 AliTRDdigitizer.cxx:1980
 AliTRDdigitizer.cxx:1981
 AliTRDdigitizer.cxx:1982
 AliTRDdigitizer.cxx:1983
 AliTRDdigitizer.cxx:1984
 AliTRDdigitizer.cxx:1985
 AliTRDdigitizer.cxx:1986
 AliTRDdigitizer.cxx:1987
 AliTRDdigitizer.cxx:1988
 AliTRDdigitizer.cxx:1989
 AliTRDdigitizer.cxx:1990
 AliTRDdigitizer.cxx:1991
 AliTRDdigitizer.cxx:1992
 AliTRDdigitizer.cxx:1993
 AliTRDdigitizer.cxx:1994