]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
single cell cluster is fixed
authorbnandi <bnandi@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 24 Jul 2006 11:15:28 +0000 (11:15 +0000)
committerbnandi <bnandi@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 24 Jul 2006 11:15:28 +0000 (11:15 +0000)
PMD/AliPMDClusteringV1.cxx
PMD/AliPMDSDigitsRead.C [new file with mode: 0644]

index ad48f8c24abac445cc010a760c4523552caa583e..2a15ac3a70db232a74f6777611f0deee72641388 100644 (file)
@@ -85,7 +85,7 @@ void AliPMDClusteringV1::DoClust(Int_t idet, Int_t ismn, Double_t celladc[48][96
 
   int i, i1, i2, j, nmx1, incr, id, jd;
   Int_t   celldataX[15], celldataY[15];
-  Float_t clusdata[7];
+  Float_t clusdata[6];
 
   Double_t  cutoff, ave;
 
@@ -209,7 +209,7 @@ void AliPMDClusteringV1::DoClust(Int_t idet, Int_t ismn, Double_t celladc[48][96
              celldataY[ihit] = fClTr[ihit][i1]%10000;
            }
        }
-
+      //printf("%d %f %f\n",idet,cluXC,cluYC );
       pmdcl = new AliPMDcluster(idet, ismn, clusdata, celldataX, celldataY);
       pmdcont->Add(pmdcl);
     }
@@ -355,10 +355,12 @@ int AliPMDClusteringV1::CrClust(double ave, double cutoff, int nmx1)
       }
     }
   }
+
   //  for(icell=0; icell<=cellcount; icell++){
   //    ofl0 << fInfcl[0][icell] << " " << fInfcl[1][icell] << " " <<
   //      fInfcl[2][icell] << endl;
-  //  }
+  //}
+
   return cellcount;
 }
 // ------------------------------------------------------------------------ //
@@ -373,8 +375,8 @@ void AliPMDClusteringV1::RefClust(int incr)
   int ig, nsupcl, lev1[20], lev2[20];
   double x[4500], y[4500], z[4500], x1, y1, z1, x2, y2, z2, dist;
   double xc[4500], yc[4500], zc[4500], cells[4500], sum, rc[4500], rr;
-
-
+  
+  
   //asso
   Int_t t[4500],cellCount[4500];
   for(i=0; i<4500; i++)
@@ -382,8 +384,8 @@ void AliPMDClusteringV1::RefClust(int incr)
       t[i]=-1;
       cellCount[i]=0;
     }
-
-
+  
+  
   // fClno counts the final clusters
   // nsupcl =  # of superclusters; ncl[i]= # of cells in supercluster i
   // x, y and z store (x,y) coordinates of and energy deposited in a cell
@@ -391,26 +393,26 @@ void AliPMDClusteringV1::RefClust(int incr)
   // zc stores the energy deposited in a cluster
   // rc is cluster radius
   // finally the cluster information is put in 2-dimensional array clusters
-  // ofstream ofl1("checking.5",ios::app);
+  //ofstream ofl1("checking.5",ios::app);
   fClno  = -1;
   nsupcl = -1;
   for(i=0; i<4500; i++){ncl[i]=-1;}
-  for(i=0; i<incr; i++){
+  for(i=0; i<incr; i++){
     if(fInfcl[0][i] != nsupcl){ nsupcl=nsupcl+1; }
     if (nsupcl > 4500) {
       AliWarning("RefClust: Too many superclusters!");
       nsupcl = 4500;
       break;
     }
+     
     ncl[nsupcl]=ncl[nsupcl]+1;
   }
 
   AliDebug(1,Form("Number of cells = %d Number of Superclusters = %d",
                  incr+1,nsupcl+1));
-
   id=-1;
   icl=-1;
-  for(i=0; i<nsupcl; i++){
+  for(i=0; i<=nsupcl; i++) {
     if(ncl[i] == 0){
       id=id+1;
       icl=icl+1;
@@ -430,8 +432,7 @@ void AliPMDClusteringV1::RefClust(int incr)
       fClusters[3][fClno] = 1.;
       fClusters[4][fClno] = 0.5;
 
-
-      //asso
+      //association
 
       fClTr[0][fClno]=fCellTrNo[i1][i2];
       for(Int_t icltr=1;icltr<14;icltr++)
@@ -441,6 +442,7 @@ void AliPMDClusteringV1::RefClust(int incr)
       
       //ofl1 << icl << " " << fCoord[0][i1][i2] << " " << fCoord[1][i1][i2] <<
       //" " << fEdepCell[i1][i2] << " " << fClusters[3][fClno] <<endl;
+      
     }else if(ncl[i] == 1){
       // two cell super-cluster --> single cluster
       // cluster center is at ener. dep.-weighted mean of two cells
@@ -486,7 +488,7 @@ void AliPMDClusteringV1::RefClust(int incr)
 
 
       //ofl1 << icl << " " << fClusters[0][fClno] << " " << fClusters[1][fClno]
-      //   << " " << fClusters[2][fClno] << " " <<fClusters[3][fClno] <<endl;
+      //  << " " << fClusters[2][fClno] << " " <<fClusters[3][fClno] <<endl;
     }
     else{
       
diff --git a/PMD/AliPMDSDigitsRead.C b/PMD/AliPMDSDigitsRead.C
new file mode 100644 (file)
index 0000000..7719a1d
--- /dev/null
@@ -0,0 +1,106 @@
+// ----------------------------------------------------//
+//                                                     //
+//    This is a macro  to read PMD.SDigits.root         //
+//                                                     //
+// ----------------------------------------------------//
+#include <Riostream.h>
+void AliPMDSDigitsRead(Int_t nevt = 1) 
+{
+  TStopwatch timer;
+  timer.Start();
+
+
+  TH2F *h2 = new TH2F("h2","Y vs. X",200,-100.,100.,200,-100.,100.);
+  AliRunLoader *fRunLoader = AliRunLoader::Open("galice.root");
+  if (!fRunLoader)
+    { 
+      Error("Open","Can not open session for file %s.",file);
+    }
+  
+  //  fRunLoader->LoadgAlice();
+  fRunLoader->LoadHeader();
+  gAlice = fRunLoader->GetAliRun();
+  
+  if (gAlice)
+    {
+      printf("AliRun object found on file.\n");
+    }
+  else
+    {
+      printf("Could not find AliRun object.\n");
+    }
+  fPMD  = (AliPMD*)gAlice->GetDetector("PMD");
+  fPMDLoader = fRunLoader->GetLoader("PMDLoader");
+  if (fPMDLoader == 0x0)
+    {
+      cerr<<"OpengAlice : Can not find PMD or PMDLoader\n";
+    }
+  
+  fPMDLoader->LoadSDigits("READ");
+  TClonesArray *fSDigits; 
+  
+  // -------------------------------------------------------------- //
+  
+  Int_t    det = 0,smn = 0;
+  Int_t    xpos, ypos;
+  Int_t    xpad, ypad;
+  Float_t  edep;
+  Float_t  xx,yy;
+
+  AliPMDUtility cc;  
+  
+  for (Int_t ievt = 0; ievt <nevt; ievt++)
+    {
+      fRunLoader->GetEvent(ievt);
+      fTreeS = fPMDLoader->TreeS();
+      if (fTreeS == 0x0)
+       {
+         cout << " Can not get TreeD" << endl;
+       }
+      AliPMDsdigit  *pmdsdigit;
+      TBranch *branch = fTreeS->GetBranch("PMDSDigit");
+      branch->SetAddress(&fSDigits);
+      
+      Int_t nmodules = (Int_t) fTreeS->GetEntries();
+      
+      for (Int_t imodule = 0; imodule < nmodules; imodule++)
+       {
+         fTreeS->GetEntry(imodule); 
+         Int_t nentries = fSDigits->GetLast();
+         for (Int_t ient = 0; ient < nentries+1; ient++)
+           {
+             pmdsdigit = (AliPMDsdigit*)fSDigits->UncheckedAt(ient);
+             
+             det    = pmdsdigit->GetDetector();
+             smn    = pmdsdigit->GetSMNumber();
+             xpos   = pmdsdigit->GetRow();
+             ypos   = pmdsdigit->GetColumn();
+             edep   = pmdsdigit->GetCellEdep();
+             Int_t trno   = pmdsdigit->GetTrackNumber();
+             
+             if(smn <12)
+               {
+                 xpad = ypos;
+                 ypad = xpos;
+               }
+             else if(smn >=12 && smn < 24)
+               {
+                 xpad = xpos;
+                 ypad = ypos;
+               }
+             
+             if(det == 1)
+               {
+                 cc.RectGeomCellPos(smn,xpad,ypad,xx,yy);
+                 h2->Fill(xx,yy);  
+               }
+           }
+       } // modules
+
+    }   
+  h2->Draw();
+
+  timer.Stop();
+  timer.Print();
+}
+