]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
new design: digits are not data members of AliPHOS anymore (use IndexToObject)
authorschutz <schutz@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 4 Apr 2001 13:15:31 +0000 (13:15 +0000)
committerschutz <schutz@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 4 Apr 2001 13:15:31 +0000 (13:15 +0000)
PHOS/AliPHOSEmcRecPoint.cxx
PHOS/AliPHOSEmcRecPoint.h

index ac68121bdc5338720df6ed16622cd7835b5ef826..e08adf6005150fab298e48012877480fd225d2c0 100644 (file)
 
 // --- AliRoot header files ---
 
-#include "AliGenerator.h"
+ #include "AliGenerator.h"
 #include "AliPHOSGeometry.h"
 #include "AliPHOSEmcRecPoint.h"
 #include "AliRun.h"
-#include "AliPHOSIndexToObject.h"
 
 ClassImp(AliPHOSEmcRecPoint)
 
 //____________________________________________________________________________
-AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(Float_t W0, Float_t LocMaxCut)
-  : AliPHOSRecPoint()
+AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
 {
   // ctor
 
   fMulDigit   = 0 ;  
   fAmp   = 0. ;   
-  fEnergyList = new Float_t[fMaxDigit]; 
+  fCoreEnergy = 0 ; 
+  fEnergyList = 0 ;
   fGeom = AliPHOSGeometry::GetInstance() ;
-  fDelta     =  ((AliPHOSGeometry *) fGeom)->GetCrystalSize(0) ; 
-  fW0        = W0 ;          
-  fLocMaxCut = LocMaxCut ; 
   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
   
 }
@@ -73,6 +69,9 @@ void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
   // Adds a digit to the RecPoint
   //  and accumulates the total amplitude and the multiplicity 
   
+  if(fEnergyList == 0)
+    fEnergyList =  new Float_t[fMaxDigit]; 
+
   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
     fMaxDigit*=2 ; 
     Int_t * tempo = new ( Int_t[fMaxDigit] ) ; 
@@ -103,6 +102,8 @@ void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
   fEnergyList[fMulDigit]   = Energy ;
   fMulDigit++ ; 
   fAmp += Energy ; 
+
+  EvalPHOSMod(&digit) ;
 }
 
 //____________________________________________________________________________
@@ -133,6 +134,9 @@ Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
 {
   // Compares two RecPoints according to their position in the PHOS modules
 
+  Float_t delta = 1 ; //Width of "Sorting row". If you changibg this 
+                      //value (what is senseless) change as vell delta in
+                      //AliPHOSTrackSegmentMakerv* and other RecPoints...
   Int_t rv ; 
 
   AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ; 
@@ -142,16 +146,16 @@ Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
   Int_t phosmod2 = clu->GetPHOSMod() ;
 
   TVector3 locpos1; 
-  this->GetLocalPosition(locpos1) ;
+  GetLocalPosition(locpos1) ;
   TVector3 locpos2;  
   clu->GetLocalPosition(locpos2) ;  
 
   if(phosmod1 == phosmod2 ) {
-    Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/fDelta)-(Int_t)TMath::Ceil(locpos2.X()/fDelta) ;
+    Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
     if (rowdif> 0) 
-      rv = -1 ;
+      rv = 1 ;
      else if(rowdif < 0) 
-       rv = 1 ;
+       rv = -1 ;
     else if(locpos1.Z()>locpos2.Z()) 
       rv = -1 ;
     else 
@@ -171,118 +175,119 @@ Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
 //______________________________________________________________________________
 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
 {
-  // Execute action corresponding to one event
-  //  This member function is called when a AliPHOSRecPoint is clicked with the locator
-  //
-  //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
-  //  and switched off when the mouse button is released.
-  //
+//   Commented by Dmitri Peressounko: there is no possibility to ensure, 
+//   that AliPHOSIndexToObject keeps the correct information.
+
+//   // Execute action corresponding to one event
+//   //  This member function is called when a AliPHOSRecPoint is clicked with the locator
+//   //
+//   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
+//   //  and switched off when the mouse button is released.
+//   //
 
-  //   static Int_t pxold, pyold;
+//   //   static Int_t pxold, pyold;
 
-  AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
+//   AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
   
-  static TGraph *  digitgraph = 0 ;
+//   static TGraph *  digitgraph = 0 ;
   
-  if (!gPad->IsEditable()) return;
+//   if (!gPad->IsEditable()) return;
   
-  TH2F * histo = 0 ;
-  TCanvas * histocanvas ; 
+//   TH2F * histo = 0 ;
+//   TCanvas * histocanvas ; 
   
-  switch (event) {
+//   switch (event) {
     
-  case kButton1Down: {
-    AliPHOSDigit * digit ;
-    AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
-    Int_t iDigit;
-    Int_t relid[4] ;
+//   case kButton1Down: {
+//     AliPHOSDigit * digit ;
+//     AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
+//     Int_t iDigit;
+//     Int_t relid[4] ;
     
-    const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ; 
-    Float_t * xi = new Float_t[kMulDigit] ; 
-    Float_t * zi = new Float_t[kMulDigit] ; 
+//     const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ; 
+//     Float_t * xi = new Float_t[kMulDigit] ; 
+//     Float_t * zi = new Float_t[kMulDigit] ; 
     
-    // create the histogram for the single cluster 
-    // 1. gets histogram boundaries
-    Float_t ximax = -999. ; 
-    Float_t zimax = -999. ; 
-    Float_t ximin = 999. ; 
-    Float_t zimin = 999. ;
+//     // create the histogram for the single cluster 
+//     // 1. gets histogram boundaries
+//     Float_t ximax = -999. ; 
+//     Float_t zimax = -999. ; 
+//     Float_t ximin = 999. ; 
+//     Float_t zimin = 999. ;
     
-    for(iDigit=0; iDigit<kMulDigit; iDigit++) {
-      digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
-      phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
-      phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
-      if ( xi[iDigit] > ximax )
-       ximax = xi[iDigit] ; 
-      if ( xi[iDigit] < ximin )
-       ximin = xi[iDigit] ; 
-      if ( zi[iDigit] > zimax )
-       zimax = zi[iDigit] ; 
-      if ( zi[iDigit] < zimin )
-       zimin = zi[iDigit] ;     
-    }
-    ximax += phosgeom->GetCrystalSize(0) / 2. ;
-    zimax += phosgeom->GetCrystalSize(2) / 2. ;
-    ximin -= phosgeom->GetCrystalSize(0) / 2. ;
-    zimin -= phosgeom->GetCrystalSize(2) / 2. ;
-    Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
-    Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
+//     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
+//       digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
+//       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
+//       phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
+//       if ( xi[iDigit] > ximax )
+//     ximax = xi[iDigit] ; 
+//       if ( xi[iDigit] < ximin )
+//     ximin = xi[iDigit] ; 
+//       if ( zi[iDigit] > zimax )
+//     zimax = zi[iDigit] ; 
+//       if ( zi[iDigit] < zimin )
+//     zimin = zi[iDigit] ;     
+//     }
+//     ximax += phosgeom->GetCrystalSize(0) / 2. ;
+//     zimax += phosgeom->GetCrystalSize(2) / 2. ;
+//     ximin -= phosgeom->GetCrystalSize(0) / 2. ;
+//     zimin -= phosgeom->GetCrystalSize(2) / 2. ;
+//     Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
+//     Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
     
-    // 2. gets the histogram title
+//     // 2. gets the histogram title
     
-    Text_t title[100] ; 
-    sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
+//     Text_t title[100] ; 
+//     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
     
-    if (!histo) {
-      delete histo ; 
-      histo = 0 ; 
-    }
-    histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
+//     if (!histo) {
+//       delete histo ; 
+//       histo = 0 ; 
+//     }
+//     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
     
-    Float_t x, z ; 
-    for(iDigit=0; iDigit<kMulDigit; iDigit++) {
-      digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
-      phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
-      phosgeom->RelPosInModule(relid, x, z);
-      histo->Fill(x, z, fEnergyList[iDigit] ) ;
-    }
+//     Float_t x, z ; 
+//     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
+//       digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
+//       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
+//       phosgeom->RelPosInModule(relid, x, z);
+//       histo->Fill(x, z, fEnergyList[iDigit] ) ;
+//     }
     
-    if (!digitgraph) {
-      digitgraph = new TGraph(kMulDigit,xi,zi);
-      digitgraph-> SetMarkerStyle(5) ; 
-      digitgraph-> SetMarkerSize(1.) ;
-      digitgraph-> SetMarkerColor(1) ;
-      digitgraph-> Paint("P") ;
-    }
+//     if (!digitgraph) {
+//       digitgraph = new TGraph(kMulDigit,xi,zi);
+//       digitgraph-> SetMarkerStyle(5) ; 
+//       digitgraph-> SetMarkerSize(1.) ;
+//       digitgraph-> SetMarkerColor(1) ;
+//       digitgraph-> Paint("P") ;
+//     }
     
-    Print() ;
-    histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ; 
-    histocanvas->Draw() ; 
-    histo->Draw("lego1") ; 
+//     Print() ;
+//     histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ; 
+//     histocanvas->Draw() ; 
+//     histo->Draw("lego1") ; 
     
-    delete[] xi ; 
-    delete[] zi ; 
+//     delete[] xi ; 
+//     delete[] zi ; 
     
-    break;
-  }
+//     break;
+//   }
   
-  case kButton1Up: 
-    if (digitgraph) {
-      delete digitgraph  ;
-      digitgraph = 0 ;
-    }
-    break;
+//   case kButton1Up: 
+//     if (digitgraph) {
+//       delete digitgraph  ;
+//       digitgraph = 0 ;
+//     }
+//     break;
   
-   }
+//    }
 }
 
 //____________________________________________________________________________
-Float_t  AliPHOSEmcRecPoint::GetDispersion() const
+void  AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
 {
   // Calculates the dispersion of the shower at the origine of the RecPoint
 
-  AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
-
   Float_t d    = 0 ;
   Float_t wtot = 0 ;
 
@@ -296,32 +301,30 @@ Float_t  AliPHOSEmcRecPoint::GetDispersion() const
   
   Int_t iDigit;
   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
-    digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
+    digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
     Int_t relid[4] ;
     Float_t xi ;
     Float_t zi ;
     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
     phosgeom->RelPosInModule(relid, xi, zi);
-    Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
+    Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
     d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ; 
     wtot+=w ;
   }
 
+
   d /= wtot ;
 
-  return TMath::Sqrt(d) ;
+  fDispersion = TMath::Sqrt(d) ;
 }
 //______________________________________________________________________________
-Float_t AliPHOSEmcRecPoint::CoreEnergy()
+void AliPHOSEmcRecPoint::EvalCoreEnergy(TClonesArray * digits)
 {
   //This function calculates energy in the core, 
   //i.e. within radius rad = 3cm. Beyond this radius
   //in accoradnce with shower profile energy deposition 
   // should be less than 2%
 
-  AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
-
-  Float_t eCore = 0 ;
   Float_t coreRadius = 3 ;
 
   TVector3 locpos;
@@ -334,7 +337,7 @@ Float_t AliPHOSEmcRecPoint::CoreEnergy()
   
   Int_t iDigit;
   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
-    digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
+    digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
     Int_t relid[4] ;
     Float_t xi ;
     Float_t zi ;
@@ -342,19 +345,16 @@ Float_t AliPHOSEmcRecPoint::CoreEnergy()
     phosgeom->RelPosInModule(relid, xi, zi);    
     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
     if(distance < coreRadius)
-      eCore += fEnergyList[iDigit] ;
+      fCoreEnergy += fEnergyList[iDigit] ;
   }
 
-return eCore ;
 }
 
 //____________________________________________________________________________
-void  AliPHOSEmcRecPoint::GetElipsAxis(Float_t * lambda)
+void  AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
 {
   // Calculates the axis of the shower ellipsoid
 
-  AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
-
   Double_t wtot = 0. ;
   Double_t x    = 0.;
   Double_t z    = 0.;
@@ -367,13 +367,13 @@ void  AliPHOSEmcRecPoint::GetElipsAxis(Float_t * lambda)
   Int_t iDigit;
 
   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
-    digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
+    digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
     Int_t relid[4] ;
     Float_t xi ;
     Float_t zi ;
     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
     phosgeom->RelPosInModule(relid, xi, zi);
-    Double_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
+    Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
     dxx  += w * xi * xi ;
     x    += w * xi ;
     dzz  += w * zi * zi ;
@@ -403,30 +403,31 @@ void  AliPHOSEmcRecPoint::GetElipsAxis(Float_t * lambda)
 //   dxz = dxz/(CosX*CosZ) ;
 
 
-  lambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
-  if(lambda[0] > 0)
-    lambda[0] = TMath::Sqrt(lambda[0]) ;
+  fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
+  if(fLambda[0] > 0)
+    fLambda[0] = TMath::Sqrt(fLambda[0]) ;
 
-  lambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
-  if(lambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
-    lambda[1] = TMath::Sqrt(lambda[1]) ;
+  fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
+  if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
+    fLambda[1] = TMath::Sqrt(fLambda[1]) ;
   else
-    lambda[1]= 0. ;
+    fLambda[1]= 0. ;
 }
 
 //____________________________________________________________________________
-void AliPHOSEmcRecPoint::EvalAll(void )
+void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
 {
-  AliPHOSRecPoint::EvalAll() ;
-  EvalLocalPosition() ;
+  AliPHOSRecPoint::EvalAll(logWeight,digits) ;
+  EvalLocalPosition(logWeight, digits) ;
+  EvalElipsAxis(logWeight, digits) ;
+  EvalDispersion(logWeight, digits) ;
+  EvalCoreEnergy(digits);
 }
 //____________________________________________________________________________
-void AliPHOSEmcRecPoint::EvalLocalPosition(void )
+void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
 {
   // Calculates the center of gravity in the local PHOS-module coordinates 
 
-  AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
-
   Float_t wtot = 0. ;
  
   Int_t relid[4] ;
@@ -441,13 +442,13 @@ void AliPHOSEmcRecPoint::EvalLocalPosition(void )
   Int_t iDigit;
 
   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
-    digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) );
+    digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
 
     Float_t xi ;
     Float_t zi ;
     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
     phosgeom->RelPosInModule(relid, xi, zi);
-    Float_t w = TMath::Max( 0., fW0 + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
+    Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
     x    += xi * w ;
     z    += zi * w ;
     wtot += w ;
@@ -521,13 +522,12 @@ Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
 }
 
 //____________________________________________________________________________
-Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy) const
+Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy,
+                                              Float_t locMaxCut,TClonesArray * digits) const
 { 
   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
   //  energy difference between two local maxima
 
-  AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
-
   AliPHOSDigit * digit ;
   AliPHOSDigit * digitN ;
   
@@ -535,28 +535,28 @@ Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEn
   Int_t iDigitN ;
   Int_t iDigit ;
 
-  for(iDigit = 0; iDigit < fMulDigit; iDigit++){
-    maxAt[iDigit] = (Int_t) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
-  }
+  for(iDigit = 0; iDigit < fMulDigit; iDigit++)
+    maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit])  ;
+
   
   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
     if(maxAt[iDigit] != -1) {
       digit = (AliPHOSDigit *) maxAt[iDigit] ;
           
       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {       
-       digitN = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigitN]) ) ; 
+       digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ; 
        
        if ( AreNeighbours(digit, digitN) ) {
          if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
            maxAt[iDigitN] = -1 ;
            // but may be digit too is not local max ?
-           if(fEnergyList[iDigit] < fEnergyList[iDigitN] + fLocMaxCut) 
+           if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
              maxAt[iDigit] = -1 ;
          }
          else {
            maxAt[iDigit] = -1 ;
            // but may be digitN too is not local max ?
-           if(fEnergyList[iDigit] > fEnergyList[iDigitN] - fLocMaxCut) 
+           if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
              maxAt[iDigitN] = -1 ; 
          } 
        } // if Areneighbours
@@ -576,39 +576,6 @@ Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEn
 }
 
 
-// //____________________________________________________________________________
-// AliPHOSEmcRecPoint& AliPHOSEmcRecPoint::operator = (AliPHOSEmcRecPoint Clu) 
-// {
-//   int * dl = Clu.GetDigitsList() ; 
-//  if(fDigitsList) 
-//     delete fDigitsList  ;
-
-//   AliPHOSDigit * digit ;
-//   Int_t iDigit;
-
-//   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
-//     digit = (AliPHOSDigit *) dl[iDigit];
-//     AddDigit(*digit) ;
-//   }
-
-//   fAmp       = Clu.GetEnergy() ;
-//   fGeom      = Clu.GetGeom() ;
-//   TVector3 locpos;
-//   Clu.GetLocalPosition(locpos) ;
-//   fLocPos    = locpos;
-//   fMulDigit       = Clu.GetMultiplicity() ;
-//   fMaxDigit  = Clu.GetMaximumMultiplicity() ;
-//   fPHOSMod   = Clu.GetPHOSMod() ;
-//   fW0        = Clu.GetLogWeightCut() ; 
-//   fDelta     = Clu.GetDelta() ;
-//   fLocMaxCut = Clu.GetLocMaxCut() ;
-  
-//   delete dl ; 
-//   return *this ;
-// }
 
 //____________________________________________________________________________
 void AliPHOSEmcRecPoint::Print(Option_t * option) 
@@ -617,51 +584,21 @@ void AliPHOSEmcRecPoint::Print(Option_t * option)
   
   cout << "AliPHOSEmcRecPoint: " << endl ;
 
-  AliPHOSDigit * digit ; 
   Int_t iDigit;
-  AliPHOSGeometry * phosgeom =  (AliPHOSGeometry *) fGeom ;
-
-  Float_t xi ;
-  Float_t zi ;
-  Int_t relid[4] ; 
+  cout << " digits # = " ;
+  for(iDigit=0; iDigit<fMulDigit; iDigit++)
+    cout << fDigitsList[iDigit] << "  " ;  
+  cout << endl ;
+  
+  cout << " Energies = " ;
+  for(iDigit=0; iDigit<fMulDigit; iDigit++) 
+    cout  << fEnergyList[iDigit] << "  ";
+  cout << endl ;
   
-  AliPHOSIndexToObject * please =  AliPHOSIndexToObject::GetInstance() ; 
-
-  for(iDigit=0; iDigit<fMulDigit; iDigit++) {
-    digit = please->GimeDigit( fDigitsList[iDigit] ) ; 
-    phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
-    phosgeom->RelPosInModule(relid, xi, zi);
-    cout << " Id = " << digit->GetId() ;  
-    cout << "   module  = " << relid[0] ;  
-    cout << "   x  = " << xi ;  
-    cout << "   z  = " << zi ;  
-    cout << "   Energy = " << fEnergyList[iDigit] << endl ;
-  }
   cout << "       Multiplicity    = " << fMulDigit  << endl ;
   cout << "       Cluster Energy  = " << fAmp << endl ;
   cout << "       Stored at position " << GetIndexInList() << endl ; 
  
 }
-//______________________________________________________________________________
-// void AliPHOSEmcRecPoint::Streamer(TBuffer &R__b)
-// {
-//    // Stream an object of class AliPHOSEmcRecPoint.
-
-//    if (R__b.IsReading()) {
-//       Version_t R__v = R__b.ReadVersion(); if (R__v) { }
-//       AliPHOSRecPoint::Streamer(R__b);
-//       R__b >> fDelta;
-//       fEnergyList = new Float_t[fMulDigit] ; 
-//       R__b.ReadFastArray(fEnergyList, fMulDigit);
-//       R__b >> fLocMaxCut;
-//       R__b >> fW0;
-//    } else {
-//       R__b.WriteVersion(AliPHOSEmcRecPoint::IsA());
-//       AliPHOSRecPoint::Streamer(R__b);
-//       R__b << fDelta;
-//       R__b.WriteFastArray(fEnergyList, fMulDigit);
-//       R__b << fLocMaxCut;
-//       R__b << fW0;
-//    }
-// }
  
+  
index aae870eb72b292c0d04563c956186bb696cedb75..5b4dc1873f45f083ba130ce6084562a7e60fa65a 100644 (file)
@@ -13,7 +13,7 @@
 
 // --- ROOT system ---
 
-#include "TObject.h"
+//#include "TObject.h"
 #include "TArrayI.h"
  
 // --- Standard library ---
@@ -27,11 +27,7 @@ class AliPHOSEmcRecPoint : public AliPHOSRecPoint  {
 
 public:
 
-  AliPHOSEmcRecPoint(){
-   // default ctor
-    fEnergyList = 0;  
-  } ;                    
-  AliPHOSEmcRecPoint(Float_t W0, Float_t LocMaxCut) ;
+  AliPHOSEmcRecPoint() ;
   AliPHOSEmcRecPoint(const AliPHOSEmcRecPoint & rp) {
     // cpy ctor requested by Coding Convention 
     // but not yet needed
@@ -40,26 +36,29 @@ public:
  
   virtual ~AliPHOSEmcRecPoint() ;  
 
-  virtual void  AddDigit(AliPHOSDigit & digit, Float_t Energy) ;  // add a digit to the digits list  
+  virtual void  AddDigit(AliPHOSDigit & digit, Float_t Energy) ;          // add a digit to the digits list  
   Int_t       Compare(const TObject * obj) const;                         // method for sorting  
-  Float_t     CoreEnergy() ;
-  void        EvalAll() ;
-  void        EvalLocalPosition() ;                                // computes the position in the PHOS module 
-  Float_t     GetDelta () const {     return fDelta ; }    
-  Float_t     GetDispersion() const ;                              // computes the dispersion of the shower
-  void        GetElipsAxis(Float_t * lambda) ;                     // computes the axis of shower ellipsoide
-  Float_t *   GetEnergiesList() const {    return fEnergyList ;}   // gets the list of energies making this recpoint
-  virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py) ; 
-  Float_t     GetLocMaxCut ()const  {    return fLocMaxCut ; }     // gets the cut of the local maximum search  
-  Float_t     GetLogWeightCut ()const { return fW0 ; }             // gets the logarythmic weight for the 
-                                                                   // center of gravity calculation
+
+  virtual void  EvalAll(Float_t logWeight,TClonesArray * digits) ;
+          void  EvalCoreEnergy(TClonesArray * digits) ;             
+  virtual void  EvalLocalPosition(Float_t logWeight,TClonesArray * digits) ;// computes the position in the PHOS module 
+  virtual void  EvalDispersion(Float_t logWeight,TClonesArray * digits) ;   // computes the dispersion of the shower
+  virtual void  EvalElipsAxis(Float_t logWeight, TClonesArray * digits );   // computes the axis of shower ellipsoide
+  virtual void  ExecuteEvent(Int_t event, Int_t px, Int_t py) ; 
+
+  Float_t         GetCoreEnergy()const {return fCoreEnergy ;}
+  virtual Float_t GetDispersion()const {return fDispersion ;}
+  virtual void    GetElipsAxis(Float_t * lambda)const { lambda[0] = fLambda[0] ;
+                                                        lambda[1] = fLambda[1] ; }
+  Float_t *   GetEnergiesList() const {return fEnergyList ;}       // gets the list of energies making this recpoint
   Float_t     GetMaximalEnergy(void) const ;                       // get the highest energy in the cluster
   Int_t       GetMaximumMultiplicity() const {return fMaxDigit ;}  // gets the maximum number of digits allowed
   Int_t       GetMultiplicity(void) const { return fMulDigit ; }   // gets the number of digits making this recpoint
   Int_t       GetMultiplicityAtLevel(const Float_t level) const ;  // computes multiplicity of digits with 
                                                                    // energy above relative level
-  Int_t       GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy) const ; // searches for the local maxima 
+  virtual Int_t GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy,
+                                    Float_t locMaxCut,TClonesArray * digits ) const ; 
+                                                                   // searches for the local maxima 
   Bool_t      IsEmc(void) const { return kTRUE ; }                 // true if the recpoint is in EMC
   Bool_t      IsSortable() const {return kTRUE ; }                 // says that emcrecpoints are sortable objects 
   void        Print(Option_t * opt = "void") ; 
@@ -71,14 +70,14 @@ public:
     return *this ; 
   }
 
- private:
+ protected:
 
-  Bool_t AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const ;
+  virtual Bool_t AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const ;
 
-  Float_t  fDelta ;          // parameter used to sort the clusters    
-  Float_t  *fEnergyList ;    //[fMulDigit] energy of digits
-  Float_t  fLocMaxCut ;      // minimum energy difference to distinguish two maxima 
-  Float_t  fW0 ;             // logarithmic weight factor for center of gravity calculation
+  Float_t fCoreEnergy ;
+  Float_t fLambda[2] ;        //
+  Float_t fDispersion ;
+  Float_t *fEnergyList ;    //[fMulDigit] energy of digits
   
   ClassDef(AliPHOSEmcRecPoint,1)  // EMC RecPoint (cluster)