]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - MFT/AliMFTClusterFinder.cxx
Reverting previous commit and going back to rev.64119Subversion Repositories:
[u/mrichter/AliRoot.git] / MFT / AliMFTClusterFinder.cxx
index 165ac720320b31f707d6311a4fc2b274dc057655..11ecb667283e9f7a23d05ff31a9c86ea6c1da4fc 100644 (file)
 #include "AliMFTSegmentation.h"
 #include "TTree.h"
 #include "TMath.h"
+#include "AliMFTConstants.h"
 #include "AliMFTClusterFinder.h"
 
+const Double_t AliMFTClusterFinder::fCutForAvailableDigits = AliMFTConstants::fCutForAvailableDigits;
+const Double_t AliMFTClusterFinder::fCutForAttachingDigits = AliMFTConstants::fCutForAttachingDigits;
+
+
 ClassImp(AliMFTClusterFinder)
 
 //====================================================================================================================================================
@@ -38,7 +43,8 @@ ClassImp(AliMFTClusterFinder)
 AliMFTClusterFinder::AliMFTClusterFinder() : 
   TObject(),
   fDigitsInCluster(0),
-  fCurrentDig(0),
+  fCurrentDigit(0),
+  fCurrentCluster(0),
   fSegmentation(0),
   fNPlanes(0)
 {
@@ -47,19 +53,40 @@ AliMFTClusterFinder::AliMFTClusterFinder() :
   
   for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) fClustersPerPlane[iPlane] = NULL;
   fDigitsInCluster = new TClonesArray("AliMFTDigit", fNMaxDigitsPerCluster);
+  fDigitsInCluster -> SetOwner(kTRUE);
 
 }
 
 //====================================================================================================================================================
 
 AliMFTClusterFinder::~AliMFTClusterFinder() {
-
+  
   AliDebug(1, "Deleting AliMFTClusterFinder...");
+  
+  for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) {
+    if (fClustersPerPlane[iPlane]) fClustersPerPlane[iPlane]->Delete(); delete fClustersPerPlane[iPlane]; fClustersPerPlane[iPlane] = 0x0;
+  }
+  if (fDigitsInCluster) fDigitsInCluster->Delete(); delete fDigitsInCluster; fDigitsInCluster = NULL;
+  delete fSegmentation; fSegmentation=NULL;
+  
+  AliDebug(1, "... done!");
+  
+}
 
-  for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) delete fClustersPerPlane[iPlane];
+//====================================================================================================================================================
 
+void AliMFTClusterFinder::Clear(const Option_t* /*opt*/) {
+  
+  AliDebug(1, "Clearing AliMFTClusterFinder...");
+  
+  for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) {
+    if(fClustersPerPlane[iPlane]) fClustersPerPlane[iPlane]->Clear("C"); 
+  }
+  if(fDigitsInCluster) fDigitsInCluster->Clear("C");   
+  fSegmentation->Clear("C"); 
+  
   AliDebug(1, "... done!");
-
+  
 }
 
 //====================================================================================================================================================
@@ -80,7 +107,7 @@ void AliMFTClusterFinder::StartEvent() {
   AliDebug(1, "Starting Event...");
   
   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
-    fClustersPerPlane[iPlane]->Clear();
+    fClustersPerPlane[iPlane]->Delete();
   }
 
   AliDebug(1, "... done!");
@@ -97,119 +124,71 @@ void AliMFTClusterFinder::DigitsToClusters(const TObjArray *pDigitList) {
   AliDebug(1, Form("nPlanes = %d", fNPlanes));
 
   StartEvent(); 
+  Bool_t isDigAvailableForNewCluster = kTRUE;
   
+  TClonesArray *myDigitList = 0;
+
   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
 
     AliDebug(1, Form("Plane %02d", iPlane));
 
-    TClonesArray *myDigitList = (TClonesArray*) pDigitList->At(iPlane);
-    //    myDigitList->Sort();
+    myDigitList = (TClonesArray*) pDigitList->At(iPlane);
 
     AliDebug(1, Form("myDigitList->GetEntries() = %d", myDigitList->GetEntries()));
 
-    while (myDigitList->GetEntries()) {
-
-      fDigitsInCluster->Clear();
-      
-      Bool_t clusterUpdated=kTRUE;
-
-      //---------------------------------------------------------------------------------------------------------
-
-      while (clusterUpdated) {    // repeat the loop on the digits until no new compatible digit is found
-
-       clusterUpdated = kFALSE;
-
-       for (Int_t iDig=0; iDig<myDigitList->GetEntries(); iDig++) {
-         fCurrentDig = (AliMFTDigit*) myDigitList->At(iDig);
-         if (fDigitsInCluster->GetEntries()<fNMaxDigitsPerCluster) {
-           if (fDigitsInCluster->GetEntries()==0) {
-             new ((*fDigitsInCluster)[fDigitsInCluster->GetEntries()]) AliMFTDigit(*fCurrentDig);
-             myDigitList->Remove(fCurrentDig);
-             clusterUpdated = kTRUE;
-           }
-           else if (IsCurrentDigitCompatible()) {
-             new ((*fDigitsInCluster)[fDigitsInCluster->GetEntries()]) AliMFTDigit(*fCurrentDig);
-             myDigitList->Remove(fCurrentDig);
-             clusterUpdated = kTRUE;
-           }
-         }
-       }
-       
-       if (clusterUpdated) myDigitList->Compress();
-
-      }
-
-      //---------------------------------------------------------------------------------------------------------
-
-      AliDebug(1, "Building new cluster");
-
-      BuildNewCluster(iPlane);  // here the new cluster is built
-
-    }
-
-  }
-
-}
-
-//====================================================================================================================================================
-
-Bool_t AliMFTClusterFinder::IsCurrentDigitCompatible() {
-
-  // where it is decided if the current digit (fCurrentDig) is compatible with the current digits array (fDigitsInCluster)
-
-  for (Int_t iDigit=0; iDigit<fDigitsInCluster->GetEntries(); iDigit++) {
-    AliMFTDigit *tmpDig = (AliMFTDigit*) fDigitsInCluster->At(iDigit);
-    Int_t distX = TMath::Abs(tmpDig->GetPixelX() - fCurrentDig->GetPixelX());
-    Int_t distY = TMath::Abs(tmpDig->GetPixelY() - fCurrentDig->GetPixelY());
-    if (distX<=1 &&  distY<=1) return kTRUE;
-  }
-
-  return kFALSE;
+    Int_t cycleOverDigits = 0;
+    Double_t myCutForAvailableDigits = fCutForAvailableDigits;
     
-}
+    while (myDigitList->GetEntries()) {
 
-//====================================================================================================================================================
+      for (Int_t iDig=0; iDig<myDigitList->GetEntries(); iDig++) {
 
-void AliMFTClusterFinder::BuildNewCluster(Int_t plane) {
+       AliDebug(1, Form("Check %d: Digit %5d of %5d", cycleOverDigits, iDig, myDigitList->GetEntries()));
 
-  // where a new cluster is built, starting from the array of compatible digits (fDigitsInCluster)
+       fCurrentDigit = (AliMFTDigit*) myDigitList->At(iDig);
+       isDigAvailableForNewCluster = kTRUE;
 
-  AliDebug(1, Form("Starting cluster building from %d digits", fDigitsInCluster->GetEntries()));
+       for (Int_t iCluster=0; iCluster<fClustersPerPlane[iPlane]->GetEntries(); iCluster++) {
+         fCurrentCluster = (AliMFTCluster*) fClustersPerPlane[iPlane]->At(iCluster);
+         AliDebug(2, Form("Distance between cluster and digit = %f",fCurrentCluster->GetDistanceFromPixel(fCurrentDigit))); 
+         if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < fCutForAttachingDigits) { 
+           fCurrentCluster->AddPixel(fCurrentDigit); 
+           myDigitList->Remove(fCurrentDigit);
+           myDigitList->Compress();
+           iDig--;
+           isDigAvailableForNewCluster = kFALSE;
+           break; 
+         }
+         if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < myCutForAvailableDigits) isDigAvailableForNewCluster=kFALSE;
+       }
 
-  AliMFTCluster *newCluster = new AliMFTCluster();
+       if (isDigAvailableForNewCluster) {
+         AliMFTCluster *newCluster = new AliMFTCluster();
+         newCluster->AddPixel(fCurrentDigit);
+         myDigitList->Remove(fCurrentDigit);
+         myDigitList->Compress();
+         iDig--;
+         new ((*fClustersPerPlane[iPlane])[fClustersPerPlane[iPlane]->GetEntries()]) AliMFTCluster(*newCluster);
+               delete newCluster;
+       }
 
-  Double_t xCenters[fNMaxDigitsPerCluster] = {0};
-  Double_t yCenters[fNMaxDigitsPerCluster] = {0};
-  Double_t nElectrons = 0.;
+      }   // end of cycle over the digits
 
-  for (Int_t iDigit=0; iDigit<fDigitsInCluster->GetEntries(); iDigit++) {
-    AliMFTDigit *tmpDig = (AliMFTDigit*) fDigitsInCluster->At(iDigit);
-    xCenters[iDigit] = tmpDig->GetPixelCenterX();
-    yCenters[iDigit] = tmpDig->GetPixelCenterY();
-    nElectrons      += tmpDig->GetNElectrons();
-    for (Int_t iTrack=0; iTrack<tmpDig->GetNMCTracks(); iTrack++) newCluster->AddMCLabel(tmpDig->GetMCLabel(iTrack));
-  }
+      if (cycleOverDigits) myCutForAvailableDigits -= 0.5;
+      cycleOverDigits++;
 
-  newCluster -> SetX(TMath::Mean(fDigitsInCluster->GetEntries(), xCenters));
-  newCluster -> SetY(TMath::Mean(fDigitsInCluster->GetEntries(), yCenters));
-  newCluster -> SetZ(((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelCenterZ());
+    }   // no more digits to check in current plane!
 
-  Double_t minErrX = ((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelWidthX() / TMath::Sqrt(12.);
-  Double_t minErrY = ((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelWidthY() / TMath::Sqrt(12.);
-  Double_t minErrZ = ((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelWidthZ() / TMath::Sqrt(12.);
-  newCluster -> SetErrX( TMath::Max(TMath::RMS(fDigitsInCluster->GetEntries(), xCenters), minErrX) );
-  newCluster -> SetErrY( TMath::Max(TMath::RMS(fDigitsInCluster->GetEntries(), yCenters), minErrY) );
-  newCluster -> SetErrZ( minErrZ );
-    
-  newCluster -> SetNElectrons(nElectrons);
-  newCluster -> SetPlane(plane);
+    for (Int_t iCluster=0; iCluster<fClustersPerPlane[iPlane]->GetEntries(); iCluster++) {
+      fCurrentCluster = (AliMFTCluster*) fClustersPerPlane[iPlane]->At(iCluster);
+      fCurrentCluster -> TerminateCluster();
+    }
 
-  newCluster -> SetSize(fDigitsInCluster->GetEntries());
+    AliDebug(1, Form("Found %d clusters in plane %02d", fClustersPerPlane[iPlane]->GetEntries(), iPlane));
 
-  AliDebug(1, Form("Adding cluster #%02d to tree (%f, %f, %f)", 
-                  fClustersPerPlane[plane]->GetEntries(), newCluster->GetX(), newCluster->GetY(), newCluster->GetZ()));
+    myDigitList -> Delete();
 
-  new ((*fClustersPerPlane[plane])[fClustersPerPlane[plane]->GetEntries()]) AliMFTCluster(*newCluster);
+  }  // end of cycle over the planes
 
 }
 
@@ -226,10 +205,9 @@ void AliMFTClusterFinder::MakeClusterBranch(TTree *treeCluster) {
   if (treeCluster) {
     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
       AliDebug(1, Form("Setting Branch Plane_%02d for treeCluster",iPlane));
-      TBranch *branch = treeCluster->GetBranch(Form("Plane_%02d",iPlane));
-      if (branch) continue;
+      if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) continue;
       AliDebug(1, Form("Branch Plane_%02d does not exist, creating!",iPlane));
-      branch = treeCluster->Branch(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane])); 
+      treeCluster->Branch(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane])); 
     }
   }
 
@@ -266,6 +244,8 @@ void AliMFTClusterFinder::CreateClusters() {
   for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
     AliDebug(1, Form("plane %02d", iPlane));
     fClustersPerPlane[iPlane] = new TClonesArray("AliMFTCluster");
+    fClustersPerPlane[iPlane] -> SetOwner(kTRUE);
+
   }
 
 }