]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - TOF/AliTOF.cxx
Functions renamed to get a prefix PHOS
[u/mrichter/AliRoot.git] / TOF / AliTOF.cxx
index af9b7f8e4526d28dfeb3ccb16741ed5fd0f4a393..93407fdb9214e3642d28105dfef7fedd1c935969 100644 (file)
 
 /*
 $Log$
+Revision 1.38  2002/06/13 08:43:46  vicinanz
+Merging added and test macro
+
+Revision 1.37  2002/05/03 07:34:19  vicinanz
+Updated SDigitizer; Added AliTOFanalyzeSDigits.C macro
+
+Revision 1.36  2002/04/19 14:40:51  vicinanz
+Updated SDigitizer
+
+Revision 1.35  2002/03/21 13:52:53  vicinanz
+Minor changes to AliTOF constructor
+
+Revision 1.34  2002/02/20 13:41:38  hristov
+Default arguments set only in the header file
+
+Revision 1.33  2002/02/19 10:39:38  vicinanz
+t0 classes added and material update (steel added)
+
+Revision 1.31  2001/11/22 11:22:51  hristov
+Updated version of TOF digitization, N^2 problem solved (J.Chudoba)
+
+Revision 1.30  2001/10/21 18:30:39  hristov
+Several pointers were set to zero in the default constructors to avoid memory management problems
+
+Revision 1.29  2001/10/17 14:19:24  hristov
+delete replaced by delete []
+
+Revision 1.28  2001/10/05 12:02:01  vicinanz
+Minor improvements on Merger and SDigitizer
+
+Revision 1.27  2001/10/02 13:03:13  vicinanz
+Minor improvements on the code
+
+Revision 1.26  2001/09/27 10:39:20  vicinanz
+SDigitizer and Merger added
+
+Revision 1.25  2001/09/07 08:37:40  hristov
+Pointers initialised to 0 in the default constructor
+
 Revision 1.24  2001/09/05 16:31:00  hristov
 The deletion of TOF folders temporarily commented out
 
@@ -96,10 +135,14 @@ Introduction of the Copyright and cvs Log
 
 #include <iostream.h>
 #include <strstream.h>
+#include <fstream.h>
+#include <stdlib.h>
 
 #include "AliTOF.h"
 #include "AliTOFhit.h"
+#include "AliTOFhitT0.h"
 #include "AliTOFdigit.h"
+#include "AliTOFSDigit.h"
 #include "AliTOFRawSector.h"
 #include "AliTOFRoc.h"
 #include "AliTOFRawDigit.h"
@@ -128,18 +171,23 @@ AliTOF::AliTOF()
   //
   // Default constructor
   //
-  fFGeom = 0;
-  fDTask = 0;
-  fReTask = 0;
+  fFGeom = 0x0;
+  fDTask = 0x0;
+  fReTask = 0x0;
   fIshunt   = 0;
   fSDigits       = 0 ;
+  fNSDigits = 0;
   fDigits        = 0 ;
+  fReconParticles = 0x0;
   fName="TOF";
+  fMerger = 0x0;
+/* fp
   CreateTOFFolders();
+*/
 }
  
 //_____________________________________________________________________________
-AliTOF::AliTOF(const char *name, const char *title)
+AliTOF::AliTOF(const char *name, const char *title, Option_t *option)
        : AliDetector(name,title)
 {
   //
@@ -149,15 +197,24 @@ AliTOF::AliTOF(const char *name, const char *title)
   //
 
   // Initialization of hits, sdigits and digits array
-  //
-  fHits   = new TClonesArray("AliTOFhit",  405);
+  // added option for time zero analysis
+  if (strstr(option,"tzero")){
+    fHits   = new TClonesArray("AliTOFhitT0",  1000);
+    cout << "tzero option requires AliTOFv4T0 as TOF version (check Your Config.C)" << endl;
+  }else{
+    fHits   = new TClonesArray("AliTOFhit",  1000);
+  }
   gAlice->AddHitList(fHits);
   fIshunt  = 0;
+  fSDigits       = new TClonesArray("AliTOFSDigit",  1000);
+  fDigits        = new TClonesArray("AliTOFdigit",  1000);
+  fNSDigits = 0;
 
-  fSDigits       = new TClonesArray("AliTOFdigit",  405);
-
-  fDigits        = new TClonesArray("AliTOFdigit",  405);
-
+  fFGeom = 0x0;
+  fDTask = 0x0;
+  fReTask = 0x0;
+  fReconParticles = 0x0;
+  fMerger = 0x0;
 
   //
   // Digitization parameters
@@ -175,7 +232,7 @@ AliTOF::AliTOF(const char *name, const char *title)
   fZlenC   = 177.5;//cm length of module C
   fZlenB   = 141.0;//cm length of module B
   fZlenA   = 106.0;//cm length of module A
-  fZtof    = 370.5;//cm total semi-length of TOF detector
+  fZtof    = 371.5;//cm total semi-length of TOF detector
 
 // Strip Parameters
   fStripLn = 122.0;//cm  Strip Length
@@ -208,9 +265,10 @@ AliTOF::AliTOF(const char *name, const char *title)
                       // (TARODA)
   fNTdc       = 32;   // number of Tdc (Time to Digital Converter)
   fNPadXRoc   = (Int_t)fPadXSector/fNRoc; // number of pads for each ROC
-
+  /* fp 25 Sept 2001
   // Create TOF Folder Structure
-  CreateTOFFolders(); 
+  CreateTOFFolders();
+  */ 
 }
 
 //_____________________________________________________________________________
@@ -252,7 +310,7 @@ void AliTOF::CreateTOFFolders()
   fReTask = new TTask(GetName(), tempo);
   aliceRe->Add(fReTask) ;
 
-  delete tempo ;
+  delete [] tempo ;
  
   // creates the TOF geometry  folder
   geomF->AddFolder("TOF", "Geometry for TOF") ;
@@ -261,13 +319,39 @@ void AliTOF::CreateTOFFolders()
 //_____________________________________________________________________________
 AliTOF::~AliTOF()
 {
-  // remove the alice folder 
+  // dtor:
+  // it remove also the alice folder 
   // and task that TOF creates instead of AliRun
   /* PH Temporarily commented because of problems
   TFolder * alice = (TFolder*)gROOT->GetListOfBrowsables()->FindObject("FPAlice") ;
   delete alice;
   alice = 0;
   */
+  if (fHits)
+    {
+      fHits->Delete ();
+      delete fHits;
+      fHits = 0;
+    }
+  if (fDigits)
+    {
+      fDigits->Delete ();
+      delete fDigits;
+      fDigits = 0;
+    }
+  if (fSDigits)
+    {
+      fSDigits->Delete ();
+      delete fSDigits;
+      fSDigits = 0;
+    }
+  if (fReconParticles)
+    {
+      fReconParticles->Delete ();
+      delete fReconParticles;
+      fReconParticles = 0;
+    }
+
 }
 
 //_____________________________________________________________________________
@@ -280,7 +364,18 @@ void AliTOF::AddHit(Int_t track, Int_t *vol, Float_t *hits)
   TClonesArray &lhits = *fHits;
   new(lhits[fNhits++]) AliTOFhit(fIshunt, track, vol, hits);
 }
+
+//_____________________________________________________________________________
+void AliTOF::AddT0Hit(Int_t track, Int_t *vol, Float_t *hits)
+{
+  //
+  // Add a TOF hit
+  // new with placement used
+  //
+  TClonesArray &lhits = *fHits;
+  new(lhits[fNhits++]) AliTOFhitT0(fIshunt, track, vol, hits);
+}
+
 //_____________________________________________________________________________
 void AliTOF::AddDigit(Int_t *tracks, Int_t *vol, Float_t *digits)
 {
@@ -292,6 +387,55 @@ void AliTOF::AddDigit(Int_t *tracks, Int_t *vol, Float_t *digits)
   new (ldigits[fNdigits++]) AliTOFdigit(tracks, vol, digits);
 }
 
+//___________________________________________
+void AliTOF::AddSDigit(Int_t tracknum, Int_t *vol, Float_t *digits)
+{
+     
+//
+// Add a TOF sdigit
+//
+        
+  TClonesArray &lSDigits = *fSDigits;   
+  new(lSDigits[fNSDigits++]) AliTOFSDigit(tracknum, vol, digits);
+}
+
+//_____________________________________________________________________________
+void AliTOF::SetTreeAddress ()
+{
+  // Set branch address for the Hits and Digits Tree.
+  char branchname[30];
+  AliDetector::SetTreeAddress ();
+
+  TBranch *branch;
+  TTree *treeD = gAlice->TreeD ();
+
+
+  if (treeD)
+    {
+      if (fDigits)
+       {
+         branch = treeD->GetBranch (branchname);
+         if (branch)
+           branch->SetAddress (&fDigits);
+       }
+
+    }
+//  if (fSDigits)
+    //  fSDigits->Clear ();
+
+  if (gAlice->TreeS () && fSDigits)
+    {
+      branch = gAlice->TreeS ()->GetBranch ("TOF");
+      if (branch)
+       branch->SetAddress (&fSDigits);
+    }
+
+  if (gAlice->TreeR() && fReconParticles) 
+    {
+      branch = gAlice->TreeR()->GetBranch("TOF"); 
+      if (branch) branch->SetAddress(&fReconParticles) ;
+    }   
+}
 
 //_____________________________________________________________________________
 void AliTOF::CreateGeometry()
@@ -361,20 +505,21 @@ void AliTOF::CreateMaterials()
 {
   //
   // Defines TOF materials for all versions
-  // Authors :   Maxim Martemianov, Boris Zagreev (ITEP) 
-  //            18/09/98 
-  // Revision: F. Pierella 5-3-2001
-  // Bologna University
+  // Revision: F. Pierella 18-VI-2002
   //
+
   Int_t   isxfld = gAlice->Field()->Integ();
   Float_t sxmgmx = gAlice->Field()->Max();
   //
-  //--- Quartz (SiO2) 
+  //--- Quartz (SiO2) to simulate float glass
+  //    density tuned to have correct float glass 
+  //    radiation length
   Float_t   aq[2] = { 28.0855,15.9994 };
   Float_t   zq[2] = { 14.,8. };
   Float_t   wq[2] = { 1.,2. };
-  Float_t   dq = 2.20;
+  Float_t   dq = 2.55; // std value: 2.2
   Int_t nq = -2;
+
   // --- Freon C2F4H2 (TOF-TDR pagg.)
   // Geant Manual CONS110-1, pag. 43 (Geant, Detector Description and Simulation Tool)
   Float_t afre[3]  = {12.011,18.998,1.007};
@@ -421,12 +566,14 @@ void AliTOF::CreateMaterials()
   Float_t wmatg10[4] = { .259,.288,.248,.205 };
   Float_t densg10  = 1.7;
   Int_t nlmatg10 = -4;
-  // --- DME 
-  Float_t adme[5] = { 12.,1.,16.,19.,79. };
-  Float_t zdme[5] = {  6.,1., 8., 9.,35. };
-  Float_t wmatdme[5] = { .4056,.0961,.2562,.1014,.1407 };
-  Float_t densdme  = .00205;
-  Int_t nlmatdme = 5;
+
+  // plexiglass CH2=C(CH3)CO2CH3
+  Float_t aplex[3] = { 12.,1.,16.};
+  Float_t zplex[3] = {  6.,1., 8.};
+  Float_t wmatplex[3] = {5.,8.,2.};
+  Float_t densplex  =1.16;
+  Int_t nplex = -3;
+
   // ---- ALUMINA (AL203) 
   Float_t aal[2] = { 27.,16.};
   Float_t zal[2] = { 13., 8.};
@@ -439,6 +586,12 @@ void AliTOF::CreateMaterials()
   Float_t wwa[2] = {  2.,  1. };
   Float_t dwa    = 1.0;
   Int_t nwa = -2;
+
+// stainless steel
+  Float_t asteel[4] = { 55.847,51.9961,58.6934,28.0855 };
+  Float_t zsteel[4] = { 26.,24.,28.,14. };
+  Float_t wsteel[4] = { .715,.18,.1,.005 };
+
   //
   //AliMaterial(0, "Vacuum$", 1e-16, 1e-16, 1e-16, 1e16, 1e16);
   AliMaterial( 1, "Air$",14.61,7.3,0.001205,30423.24,67500.);
@@ -446,15 +599,16 @@ void AliTOF::CreateMaterials()
   AliMaterial( 3, "C  $",  12.01,  6.0, 2.265,18.8, 74.4);
   AliMixture ( 4, "Polyethilene$", ape, zpe, dpe, npe, wpe);
   AliMixture ( 5, "G10$", ag10, zg10, densg10, nlmatg10, wmatg10);
-  AliMixture ( 6, "DME ", adme, zdme, densdme, nlmatdme, wmatdme);
+  AliMixture ( 6, "PLE$", aplex, zplex, densplex, nplex, wmatplex);
   AliMixture ( 7, "CO2$", ac, zc, dc, nc, wc);
   AliMixture ( 8, "ALUMINA$", aal, zal, densal, nlmatal, wmatal);
   AliMaterial( 9, "Al $", 26.98, 13., 2.7, 8.9, 37.2);
   AliMaterial(10, "C-TRD$", 12.01, 6., 2.265*18.8/69.282*15./100, 18.8, 74.4); // for 15%
   AliMixture (11, "Mylar$",  amy, zmy, dmy, nmy, wmy);
   AliMixture (12, "Freon$",  afre, zfre, densfre, nfre, wfre);
-  AliMixture (13, "Quartz$", aq, zq, dq, nq, wq);
+  AliMixture (13, "Glass$", aq, zq, dq, nq, wq);
   AliMixture (14, "Water$",  awa, zwa, dwa, nwa, wwa);
+  AliMixture (15, "STAINLESS STEEL$", asteel, zsteel, 7.88, 4, wsteel);
 
   Float_t epsil, stmin, deemax, stemax;
  
@@ -475,7 +629,7 @@ void AliTOF::CreateMaterials()
   AliMedium( 3, "C  $"  ,  3, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
   AliMedium( 4, "Pol$"  ,  4, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
   AliMedium( 5, "G10$"  ,  5, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
-  AliMedium( 6, "DME$"  ,  6, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
+  AliMedium( 6, "PLE$"  ,  6, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
   AliMedium( 7, "CO2$"  ,  7, 0, isxfld, sxmgmx, 10., -.01, -.1, .01, -.01);
   AliMedium( 8,"ALUMINA$", 8, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
   AliMedium( 9,"Al Frame$",9, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
@@ -486,6 +640,7 @@ void AliTOF::CreateMaterials()
   AliMedium(14, "Fre-S$", 12, 1, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
   AliMedium(15, "Glass$", 13, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
   AliMedium(16, "Water$", 14, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
+  AliMedium(17, "STEEL$", 15, 0, isxfld, sxmgmx, 10., stemax, deemax, epsil, stmin);
 }
 
 //_____________________________________________________________________________
@@ -496,7 +651,22 @@ Int_t AliTOF::DistancetoPrimitive(Int_t , Int_t ) const
   //
   return 9999;
 }
+
+//___________________________________________
+void AliTOF::ResetHits ()
+{
+  // Reset number of clusters and the cluster array for this detector
+  AliDetector::ResetHits ();
+}
+
+//____________________________________________
+void AliTOF::ResetDigits ()
+{
+  //
+  // Reset number of digits and the digits array for this detector
+  AliDetector::ResetDigits ();
+  //
+} 
 //_____________________________________________________________________________
 void AliTOF::Init()
 {
@@ -514,9 +684,9 @@ void AliTOF::MakeBranch(Option_t* option, const char *file)
  //
  // Initializes the Branches of the TOF inside the 
  // trees written for each event. 
- //  AliDetector::MakeBranch initializes just the 
+ // AliDetector::MakeBranch initializes just the 
  // Branch inside TreeH. Here we add the branches in 
- // TreeD and TreeS.
+ // TreeD, TreeS and TreeR.
  //
   AliDetector::MakeBranch(option,file);
 
@@ -526,6 +696,7 @@ void AliTOF::MakeBranch(Option_t* option, const char *file)
   
   const char *oD = strstr(option,"D");
   const char *oS = strstr(option,"S");
+  const char *oR = strstr(option,"R");
 
   if (oD)
   //
@@ -549,6 +720,17 @@ void AliTOF::MakeBranch(Option_t* option, const char *file)
                               branchname, &fSDigits,buffersize, file) ;
   }
 
+  if (oR)
+  //
+  // one branch for TOF reconstructed particles
+  //
+
+
+  if (fReconParticles && gAlice->TreeR() && oR){
+             MakeBranchInTree(gAlice->TreeR(),
+                              branchname, &fReconParticles,buffersize, file) ;
+  }
+
 }
 
 //____________________________________________________________________________
@@ -576,132 +758,63 @@ void AliTOF::FinishEvent()
 void AliTOF::SDigits2Digits()
 {
 //
-// Generate digits
+// Generate digits performing merging
 //
+  /*
     int nparticles = gAlice->GetNtrack();
     cout << "Particles       :" <<nparticles<<endl;
     if (nparticles > 0 ) {
+      
       AliTOF::Hits2Digits();
+      
     }
+  */
+  cout<<"AliTOF::SDigits2Digits"<<endl; 
+    if (fMerger) {
+      fMerger->Init();
+      cout<<"AliTOF::SDigits2Digits Init"<<endl; 
+      fMerger->Digitise();
+      cout<<"AliTOF::SDigits2Digits Digitise() "<<endl; 
+     }
 }
 
-//____________________________________________________________________________
-void AliTOF::Hits2Digits()
+//---------------------------------------------------------------------
+void   AliTOF::SetMerger(AliTOFMerger* merger)
 {
-//
-// Starting from the Hits Tree (TreeH), this
-// function writes the TOF Digits Branch in the Tree (TreeD) storing 
-// the digits informations.
-// It has to be called just at the end of an event or 
-// at the end of a whole run.
-// It could  also be called by AliTOF::Finish Event()
-// Just for MC events. 
-//
-// Called by the ROOT script Hits2Digits.C
-//
-// Simulation of detector response.
-
-   Int_t ver = this->IsVersion();
-   if(ver==0) return; // no digits for AliTOFv0
-
-  Int_t nhits = 0;       // total number of hits for the current track
-  Int_t evNumber = 0;    // evnumber
-  Int_t    tracks[3];    // track info
-  Int_t    vol[5];       // dummy location for digit
-  Float_t  digit[2];     // TOF digit variables
-  
-
-  // Get pointers to Alice detectors and Hits containers
-  AliDetector* TOF  = gAlice->GetDetector("TOF");
-
-
-  TTree*  tD = gAlice->TreeD();
-
-  TTree* tH = gAlice->TreeH(); // pointer to the hits tree
-  Stat_t ntracks = tH->GetEntries();
+// Set pointer to merger
+    fMerger = merger;
+}
 
-  cout << "Total number of processed tracks in event " << gAlice->GetEvNumber() <<
-  " :" << ntracks << endl;
+//---------------------------------------------------------------------
+AliTOFMerger*  AliTOF::Merger()
+{
+// Return pointer to merger
+    return fMerger;
+}
 
-  // do nothing if no tracked particles
 
-  if( ntracks > 0){
+//---------------------------------------------------------------------
 
-  // ptr to the current TOF hit
-  AliTOFhit* tofHit;
+void AliTOF::Hits2SDigits()
+{
+//
+// Use the TOF SDigitizer to make TOF SDigits
+//
+//
+  //#ifdef DEBUG
+  cout<<"ALiTOF::Hits2SDigits> start...\n";
+  //#endif
   
-  // Start loop on tracks in the hits containers
-  // check for the total number of processed hits
-  Int_t totnhits   =0;
-  Int_t totndigits =0;
+  //char * fileSDigits = 0 ;
+  char * fileHeader = 0;
+  AliTOFSDigitizer * sd = new AliTOFSDigitizer(fileHeader) ;
 
-  if(TOF) {
+  sd->Exec("") ;
+  sd->Print("");
 
-  for (Int_t track=0; track<ntracks;track++) {
+  delete sd ;
   
-      // loop on all hits for the current track
-
-      for(tofHit=(AliTOFhit*)TOF->FirstHit(track); tofHit; tofHit=(AliTOFhit*)TOF->NextHit()) {
-        ++nhits;  
-        ++totnhits;
-
-        vol[0] = tofHit->GetSector();
-        vol[1] = tofHit->GetPlate();
-        vol[2] = tofHit->GetPadx();
-        vol[3] = tofHit->GetPadz();
-        vol[4] = tofHit->GetStrip();
-
-        // 95% of efficiency to be inserted here
-        // edge effect to be inserted here
-        // cross talk  to be inserted here
-
-        Float_t idealtime = tofHit->GetTof(); // unit s
-        idealtime *= 1.E+12;  // conversion from s to ps
-                              // fTimeRes is given usually in ps
-        Float_t tdctime   = gRandom->Gaus(idealtime, fTimeRes);        
-        digit[0] = tdctime;
-
-        // typical Landau Distribution to be inserted here
-        // instead of Gaussian Distribution
-        Float_t idealcharge = tofHit->GetEdep();
-        Float_t adccharge = gRandom->Gaus(idealcharge, fChrgRes);
-        digit[1] = adccharge;
-        Int_t tracknum = tofHit -> GetTrack();
-        tracks[0] = tracknum;
-       tracks[1] = 0;
-       tracks[2] = 0;
-
-        Bool_t overlap = CheckOverlap(vol, digit, tracknum);
-        if(!overlap) 
-        AddDigit(tracks, vol, digit);
-        ++totndigits;
-        } // end loop on hits for the current track
-
-     } // end loop on ntracks
-
-  // some statistics concerning digitization
-  cout << "Total number of processed TOF hits: " << totnhits   << endl;
-  cout << "Total number of created TOF digits: " << totndigits << endl;
-
-  } // close if TOF switched ON
-
-} // close if( ntracks > 0)
-
-
-// fill and write the branch
-
-   evNumber = gAlice->GetEvNumber();
-   char hname[30];
-   sprintf(hname,"TreeD%d",evNumber);
-
-   tD->Fill();
-
-   tD->Write(hname,TObject::kOverwrite);
-
-   // reset tree
-   gAlice->TreeD()->Reset();
 }
-
 //___________________________________________________________________________
 Bool_t AliTOF::CheckOverlap(Int_t* vol, Float_t* digit,Int_t Track)
 {
@@ -709,30 +822,51 @@ Bool_t AliTOF::CheckOverlap(Int_t* vol, Float_t* digit,Int_t Track)
 // Checks if 2 or more hits belong to the same pad.
 // In this case the data assigned to the digit object
 // are the ones of the first hit in order of Time.
-//
-// Called only by Hits2Digits.
+// 2 hits from the same track on the same pad are collected.
+// Called only by Hits2SDigits.
+// This procedure has to be optimized in the next TOF release.
 //
 
-        Bool_t overlap = 0;
+        Bool_t overlap = kFALSE;
         Int_t  vol2[5];
 
-        for (Int_t ndig=0; ndig<fNdigits; ndig++){
-          AliTOFdigit* currentDigit = (AliTOFdigit*)(fDigits->UncheckedAt(ndig));
+        for (Int_t ndig=0; ndig<fSDigits->GetEntries(); ndig++){
+          AliTOFdigit* currentDigit = (AliTOFdigit*)(fSDigits->UncheckedAt(ndig));
            currentDigit->GetLocation(vol2);
-           Bool_t idem=1;
+           Bool_t idem= kTRUE;
+          // check on digit volume
            for (Int_t i=0;i<=4;i++){
-               if (vol[i]!=vol2[i]) idem=0;}
-           if (idem){
+              if (!idem) break;
+               if (vol[i]!=vol2[i]) idem=kFALSE;}
+
+           if (idem){  // same pad fired
              Float_t tdc2 = digit[0];
               Float_t tdc1 = currentDigit->GetTdc();
-              if (tdc1>tdc2){
-                  currentDigit->SetTdc(tdc2); 
-                  currentDigit->SetAdc(digit[1]);
-              }
-              currentDigit->AddTrack(Track);
-              overlap = 1;
-           }
-        }
+
+             // we separate two digits on the same pad if
+             // they are separated in time by at least 25 ns
+             // remember that tdc time is given in ps
+
+              if (TMath::Abs(tdc1-tdc2)<25000){
+                 // in case of overlap we take the earliest
+                 if (tdc1>tdc2){
+                   currentDigit->SetTdc(tdc2); 
+                   currentDigit->SetAdc(digit[1]);
+                 }
+                 else {
+                  currentDigit->SetTdc(tdc1);
+                  currentDigit->SetAdc(digit[1]);
+                 }
+                  currentDigit->AddTrack(Track); // add track number in the track array
+                  overlap = kTRUE;
+                 return overlap;
+              } else 
+               overlap= kFALSE;
+
+           } // close if (idem) -> two digits on the same TOF pad
+
+        } // end loop on existing sdigits
+
         return overlap;
 }
 
@@ -882,3 +1016,20 @@ void AliTOF::Raw2Digits(Int_t evNumber)
   tD->Write(0,TObject::kOverwrite);
 } 
 
+////////////////////////////////////////////////////////////////////////
+void AliTOF::RecreateSDigitsArray() {
+//
+// delete TClonesArray fSDigits and create it again
+//  needed for backward compatability with PPR test production
+//
+  delete fSDigits;
+  fSDigits       = new TClonesArray("AliTOFSDigit",  1000);
+}
+////////////////////////////////////////////////////////////////////////
+void AliTOF::CreateSDigitsArray() {
+//
+// create TClonesArray fSDigits
+//  needed for backward compatability with PPR test production
+//
+  fSDigits       = new TClonesArray("AliTOFSDigit",  1000);
+}