]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - JETAN/AliUA1JetFinderV1.cxx
New class for dE/dx analysis (comparison with Bethe-Bloch, check of response function...
[u/mrichter/AliRoot.git] / JETAN / AliUA1JetFinderV1.cxx
index acd180ca3cc55fc6f52d1f1c52a33432d2bc6f5e..6dafbdd15e411be2ab191561122ae69b547363ff 100644 (file)
 #include "AliJetReader.h"
 #include "AliJet.h"
 #include "AliAODJet.h"
+#include "AliLog.h"
 
 
 ClassImp(AliUA1JetFinderV1)
 
-////////////////////////////////////////////////////////////////////////
-
-AliUA1JetFinderV1::AliUA1JetFinderV1()
+/////////////////////////////////////////////////////////////////////
 
+AliUA1JetFinderV1::AliUA1JetFinderV1() :
+    AliJetFinder(),
+    fLego(0)
 {
   // Constructor
-  fHeader = 0x0;
-  fLego   = 0x0;
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -55,6 +55,8 @@ AliUA1JetFinderV1::~AliUA1JetFinderV1()
 
 {
   // destructor
+  delete fLego;
+  fLego = 0;
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -97,16 +99,17 @@ void AliUA1JetFinderV1::FindJets()
     etaT[i] = lv->Eta();
     phiT[i] = ((lv->Phi() < 0) ? (lv->Phi()) + 2 * TMath::Pi() : lv->Phi());
     if (fReader->GetCutFlag(i) != 1) continue;
-    fLego->Fill(etaT[i], phiT[i], ptT[i]);
+    fLego ->Fill(etaT[i], phiT[i], ptT[i]);
     hPtTotal->Fill(ptT[i]);
     etbgTotal+= ptT[i];
   }
+
   fJets->SetNinput(nIn);
 
   // calculate total energy and fluctuation in map
   Double_t meanpt = hPtTotal->GetMean();
-  Double_t ptRMS = hPtTotal->GetRMS();
-  Double_t npart = hPtTotal->GetEntries();
+  Double_t ptRMS  = hPtTotal->GetRMS();
+  Double_t npart  = hPtTotal->GetEntries();
   Double_t dEtTotal = (TMath::Sqrt(npart))*TMath::Sqrt(meanpt * meanpt + ptRMS*ptRMS);
 
   // arrays to hold jets
@@ -159,11 +162,32 @@ void AliUA1JetFinderV1::FindJets()
      etbgTotal = etbgTotalN; // update with new background estimation
   } //end while
 
+  // add tracks to the jet if it wasn't yet done                                                         
+  if (header->GetBackgMode() == 0){
+    Float_t rc= header->GetRadius();
+    for(Int_t jpart = 0; jpart < nIn; jpart++){ // loop for all particles in array                     
+      for(Int_t ijet=0; ijet<nj; ijet++){
+        Float_t deta = etaT[jpart] - etaJet[ijet];
+        Float_t dphi = phiT[jpart] - phiJet[ijet];
+       if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
+       if (dphi >  TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
+        Float_t dr = TMath::Sqrt(deta * deta + dphi * dphi);
+       if(dr <= rc){ // particles inside this cone                                                    
+          injet[jpart] = ijet;
+          break;
+        }
+      }// end jets loop                                                                                
+    } //end particle loop                                                                              
+  }
+  
   // add jets to list
   Int_t* idxjets = new Int_t[nj];
   Int_t nselectj = 0;
-  printf("Found %d jets \n", nj);
+//  printf("Found %d jets \n", nj);
   
+  TRefArray *refs = 0;
+  Bool_t fromAod = !strcmp(fReader->ClassName(),"AliJetAODReader");
+  if (fromAod) refs = fReader->GetReferences();
   for(Int_t kj=0; kj<nj; kj++){
      if ((etaJet[kj] > (header->GetJetEtaMax())) ||
           (etaJet[kj] < (header->GetJetEtaMin())) ||
@@ -175,13 +199,21 @@ void AliUA1JetFinderV1::FindJets()
       en = TMath::Sqrt(px * px + py * py + pz * pz);
       fJets->AddJet(px, py, pz, en);
       AliAODJet jet(px, py, pz, en);
-      jet.Print("");
+
+      if (fromAod){
+        for(Int_t jpart = 0; jpart < nIn; jpart++) // loop for all particles in array
+          if (injet[jpart] == kj && fReader->GetCutFlag(jpart) == 1)
+                   jet.AddTrack(refs->At(jpart));  // check if the particle belongs to the jet and add the ref
+      }
+      
+      //jet.Print("");
       
       AddJet(jet);
       
       idxjets[nselectj] = kj;
       nselectj++;
-  }
+  } //end particle loop
+
   //add signal percentage and total signal  in AliJets for analysis tool
   Float_t* percentage  = new Float_t[nselectj];
   Int_t* ncells      = new Int_t[nselectj];
@@ -215,22 +247,22 @@ void AliUA1JetFinderV1::FindJets()
 
 
   //delete
-  delete ptT;
-  delete etaT;
-  delete phiT;
-  delete injet;
+  delete [] ptT;
+  delete [] etaT;
+  delete [] phiT;
+  delete [] injet;
   delete hPtTotal;
-  delete etaJet;
-  delete phiJet;
-  delete etJet;
-  delete etsigJet;
-  delete etallJet;
-  delete ncellsJet;
-  delete multJet;
-  delete idxjets;
-  delete percentage;
-  delete ncells;
-  delete mult;
+  delete [] etaJet;
+  delete [] phiJet;
+  delete [] etJet;
+  delete [] etsigJet;
+  delete [] etallJet;
+  delete [] ncellsJet;
+  delete [] multJet;
+  delete [] idxjets;
+  delete [] percentage;
+  delete [] ncells;
+  delete [] mult;
 
 
 }
@@ -243,19 +275,26 @@ void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t&
 {
 
    //dump lego
-  // check enough space! *to be done*
   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
-  Float_t etCell[60000];   //! Cell Energy
-  Float_t etaCell[60000];  //! Cell eta
-  Float_t phiCell[60000];  //! Cell phi
-  Int_t   flagCell[60000]; //! Cell flag
+  const Int_t nBinsMax = 120000; // we use a fixed array not to fragment memory
+  
+  const Int_t nBinEta = header->GetLegoNbinEta();
+  const Int_t nBinPhi = header->GetLegoNbinPhi();
+  if((nBinPhi*nBinEta)>nBinsMax){
+    AliError("Too many bins of the ETA-PHI histogram");
+  }
+
+  Float_t etCell[nBinsMax];   //! Cell Energy
+  Float_t etaCell[nBinsMax];  //! Cell eta
+  Float_t phiCell[nBinsMax];  //! Cell phi
+  Short_t   flagCell[nBinsMax]; //! Cell flag
 
   Int_t nCell = 0;
   TAxis* xaxis = fLego->GetXaxis();
   TAxis* yaxis = fLego->GetYaxis();
   Float_t e = 0.0;
-  for (Int_t i = 1; i <= header->GetLegoNbinEta(); i++) {
-      for (Int_t j = 1; j <= header->GetLegoNbinPhi(); j++) {
+  for (Int_t i = 1; i <= nBinEta; i++) {
+      for (Int_t j = 1; j <= nBinPhi; j++) {
               e = fLego->GetBinContent(i,j);
               if (e < 0.0) continue; // don't include this cells
               Float_t eta  = xaxis->GetBinCenter(i);
@@ -263,7 +302,7 @@ void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t&
               etCell[nCell]  = e;
               etaCell[nCell] = eta;
               phiCell[nCell] = phi;
-          flagCell[nCell] = 0; //default
+              flagCell[nCell] = 0; //default
               nCell++;
       }
   }
@@ -289,22 +328,23 @@ void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t&
   Int_t * index  = new Int_t[nCell];
   TMath::Sort(nCell, etCell, index);
   // variable used in centroide loop
-  Float_t eta = 0.0;
-  Float_t phi = 0.0;
-  Float_t eta0 = 0.0;
-  Float_t phi0 = 0.0;
-  Float_t etab = 0.0;
-  Float_t phib = 0.0;
-  Float_t etas = 0.0;
-  Float_t phis = 0.0;
-  Float_t ets = 0.0;
-  Float_t deta = 0.0;
-  Float_t dphi = 0.0;
-  Float_t dr = 0.0;
-  Float_t etsb = 0.0;
+  Float_t eta   = 0.0;
+  Float_t phi   = 0.0;
+  Float_t eta0  = 0.0;
+  Float_t phi0  = 0.0;
+  Float_t etab  = 0.0;
+  Float_t phib  = 0.0;
+  Float_t etas  = 0.0;
+  Float_t phis  = 0.0;
+  Float_t ets   = 0.0;
+  Float_t deta  = 0.0;
+  Float_t dphi  = 0.0;
+  Float_t dr    = 0.0;
+  Float_t etsb  = 0.0;
   Float_t etasb = 0.0;
   Float_t phisb = 0.0;
-
+  Float_t dphib = 0.0;
+  
 
   for(Int_t icell = 0; icell < nCell; icell++){
         Int_t jcell = index[icell];
@@ -324,21 +364,21 @@ void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t&
         phisb = 0.0;
         for(Int_t kcell =0; kcell < nCell; kcell++){
             Int_t lcell = index[kcell];
-            if(lcell == jcell) continue; // cell itself
-            if(flagCell[lcell] != 0) continue; // cell used before
-            if(etCell[lcell] > etCell[jcell]) continue;
+            if(lcell == jcell)                 continue; // cell itself
+            if(flagCell[lcell] != 0)           continue; // cell used before
+            if(etCell[lcell] > etCell[jcell])  continue; // can this happen
             //calculate dr
             deta = etaCell[lcell] - eta;
-                dphi = phiCell[lcell] - phi;
-                if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
-                if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
-                dr = TMath::Sqrt(deta * deta + dphi * dphi);
+           dphi = TMath::Abs(phiCell[lcell] - phi);
+           if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
+           dr = TMath::Sqrt(deta * deta + dphi * dphi);
             if(dr <= rc){
                // calculate offset from initiate cell
                deta = etaCell[lcell] - eta0;
                dphi = phiCell[lcell] - phi0;
-               if (dphi < -TMath::Pi()) dphi= -dphi - 2.0 * TMath::Pi();
-                   if (dphi > TMath::Pi()) dphi = 2.0 * TMath::Pi() - dphi;
+               if (dphi < - TMath::Pi()) dphi=  dphi + 2.0 * TMath::Pi();
+              if (dphi >   TMath::Pi()) dphi = dphi - 2.0 * TMath::Pi();
+              
                etas = etas + etCell[lcell]*deta;
                phis = phis + etCell[lcell]*dphi;
                ets = ets + etCell[lcell];
@@ -346,24 +386,26 @@ void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t&
                eta = eta0 + etas/ets;
                phi = phi0 + phis/ets;
                // if cone does not move much, just go to next step
-               dr = TMath::Sqrt((eta-etab)*(eta-etab) + (phi-phib)*(phi-phib));
+              dphib = TMath::Abs(phi - phib);
+              if (dphib > TMath::Pi()) dphib = 2. * TMath::Pi() - dphib;
+              dr = TMath::Sqrt((eta-etab)*(eta-etab) + dphib * dphib);
                if(dr <= minmove) break;
                // cone should not move more than max_mov
                dr = TMath::Sqrt((etas/ets)*(etas/ets) + (phis/ets)*(phis/ets));
                if(dr > maxmove){
-                  eta = etab;
-                  phi = phib;
-                  ets = etsb;
-                  etas = etasb;
-                  phis = phisb;
-               }else{ // store this loop information
-                 etab=eta;
-                 phib=phi;
-                 etsb = ets;
-                 etasb = etas;
-                 phisb = phis;
+                  eta = etab;
+                  phi = phib;
+                  ets = etsb;
+                  etas = etasb;
+                  phis = phisb;
+              } else { // store this loop information
+                  etab  = eta;
+                  phib  = phi;
+                  etsb  = ets;
+                  etasb = etas;
+                  phisb = phis;
                }
-            }
+            } // inside cone
         }//end of cells loop looking centroide
 
         //avoid cones overloap (to be implemented in the future)
@@ -428,8 +470,8 @@ void AliUA1JetFinderV1::RunAlgoritm(Float_t etbgTotal, Double_t dEtTotal, Int_t&
 
 
   //delete
-  delete index;
-  delete idx;
+  delete[] index;
+  delete[] idx;
 
 }
 ////////////////////////////////////////////////////////////////////////
@@ -820,11 +862,12 @@ void AliUA1JetFinderV1::Init()
 {
   // initializes some variables
   AliUA1JetHeaderV1* header = (AliUA1JetHeaderV1*) fHeader;
-   // book lego
   fLego = new
     TH2F("legoH","eta-phi",
         header->GetLegoNbinEta(), header->GetLegoEtaMin(),
         header->GetLegoEtaMax(),  header->GetLegoNbinPhi(),
         header->GetLegoPhiMin(),  header->GetLegoPhiMax());
+  // Do not store in current dir
+  fLego->SetDirectory(0);
 
 }