]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
style modifications (Markus)
authorabercuci <abercuci@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 2 Jul 2010 07:34:46 +0000 (07:34 +0000)
committerabercuci <abercuci@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 2 Jul 2010 07:34:46 +0000 (07:34 +0000)
PWG1/TRD/AliTRDcheckDET.cxx
PWG1/TRD/AliTRDcheckDET.h

index abdb5b2626ea9f4f754deaeb67b5e1945c84df5b..2b370c82cf1e82d70e87488e93b875a900f36937 100644 (file)
@@ -37,6 +37,7 @@
 #include <TLinearFitter.h>
 #include <TMath.h>
 #include <TMap.h>
+#include <TProfile2D.h>
 #include <TObjArray.h>
 #include <TObject.h>
 #include <TObjString.h>
@@ -257,6 +258,13 @@ void AliTRDcheckDET::MakeSummary(){
   cOut->cd(); GetRefFigure(kFigPH);
   cOut->SaveAs(Form("TRDsummary%s2.gif", GetName())); 
   delete cOut;
+
+  // Third Plot: Mean Number of clusters as function of eta, phi and layer
+   cOut = new TCanvas(Form("summary%s3", GetName()), Form("Summary 3 for task %s", GetName()), 1024, 768);
+  cOut->cd(); MakePlotMeanClustersLayer();
+  cOut->SaveAs(Form("TRDsummary%s3.gif", GetName())); 
+  delete cOut;
+
 }
 
 //_______________________________________________________
@@ -394,6 +402,18 @@ TObjArray *AliTRDcheckDET::Histos(){
   } else h->Reset();
   fContainer->AddAt(h, kNclustersTrack);
 
+  TObjArray *arr = new TObjArray(AliTRDgeometry::kNlayer);
+  arr->SetOwner(kTRUE);  arr->SetName("clusters");
+  fContainer->AddAt(arr, kNclustersLayer);
+  for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
+    if(!(h = (TProfile2D *)gROOT->FindObject(Form("hNcl%d", ily)))){
+      h = new TProfile2D(Form("hNcl%d", ily), Form("Mean Number of clusters in Layer %d", ily), 100, -1.0, 1.0, 50, -1.1*TMath::Pi(), 1.1*TMath::Pi());
+      h->GetXaxis()->SetTitle("#eta");
+      h->GetYaxis()->SetTitle("#phi");
+    } else h->Reset();
+    arr->AddAt(h, ily);
+  }
+
   if(!(h = (TH1F *)gROOT->FindObject("hNclTls"))){
     h = new TH1F("hNclTls","N_{clusters} / tracklet", 51, -0.5, 50.5);
     h->GetXaxis()->SetTitle("N_{clusters}");
@@ -468,7 +488,7 @@ TObjArray *AliTRDcheckDET::Histos(){
   }
   fContainer->AddAt(h, kTrackStatus);
 
-  TObjArray *arr = new TObjArray(AliTRDgeometry::kNlayer);
+  arr = new TObjArray(AliTRDgeometry::kNlayer);
   arr->SetOwner(kTRUE);  arr->SetName("TrackletStatus");
   fContainer->AddAt(arr, kTrackletStatus);
   for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
@@ -636,15 +656,28 @@ TH1 *AliTRDcheckDET::PlotNClustersTracklet(const AliTRDtrackV1 *track){
     AliDebug(4, "No Track defined.");
     return NULL;
   }
+  AliExternalTrackParam *par = fkTrack->GetTrackIn() ? fkTrack->GetTrackIn() : fkTrack->GetTrackOut();
   TH1 *h = NULL;
+  TProfile2D *hlayer = NULL;
+  Double_t eta = 0., phi = 0.;
   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNclustersTracklet)))){
     AliWarning("No Histogram defined.");
     return NULL;
   }
   AliTRDseedV1 *tracklet = NULL;
+  TObjArray *histosLayer = dynamic_cast<TObjArray *>(fContainer->At(kNclustersLayer));
+  if(!histosLayer){
+    AliWarning("No Histograms for single layer defined");
+  }
   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
     h->Fill(tracklet->GetN2());
+    if(histosLayer && par){
+      if((hlayer = dynamic_cast<TProfile2D *>(histosLayer->At(itl)))){
+        GetEtaPhiAt(par, tracklet->GetX0(), eta, phi);
+        hlayer->Fill(eta, phi, tracklet->GetN2());
+      }
+    }
   }
   return h;
 }
@@ -1212,6 +1245,20 @@ void AliTRDcheckDET::GetDistanceToTracklet(Double_t *dist, AliTRDseedV1 * const
   dist[1] = c->GetZ() - tracklet->GetZat(x);
 }
 
+//________________________________________________________
+void AliTRDcheckDET::GetEtaPhiAt(AliExternalTrackParam *track, Double_t x, Double_t &eta, Double_t &phi){
+  //
+  // Get phi and eta at a given radial position
+  // 
+  AliExternalTrackParam workpar(*track);
+
+  Double_t posLocal[3];
+  Bool_t sucPos = workpar.GetXYZAt(x, fEventInfo->GetRunInfo()->GetMagneticField(), posLocal);
+  Double_t sagPhi = sucPos ? TMath::ATan2(posLocal[1], posLocal[0]) : 0.;
+  phi = sagPhi;
+  eta = workpar.Eta();
+}
+
 
 //_______________________________________________________
 TH1* AliTRDcheckDET::MakePlotChi2()
@@ -1416,6 +1463,27 @@ Bool_t AliTRDcheckDET::MakePlotPulseHeight(){
   return kTRUE;
 }
 
+//________________________________________________________
+void AliTRDcheckDET::MakePlotMeanClustersLayer(){
+  //
+  // Create Summary plot for the mean number of clusters per layer
+  //
+  TCanvas *output = gPad->GetCanvas();
+  output->Divide(3,2);
+  TObjArray *histos = (TObjArray *)fContainer->At(kNclustersLayer);
+  if(!histos){
+    AliWarning("Histos for each layer not found");
+    return;
+  }
+  TProfile2D *hlayer = NULL;
+  for(Int_t ily = 0; ily < AliTRDgeometry::kNlayer; ily++){
+    hlayer = dynamic_cast<TProfile2D *>(histos->At(ily));
+    output->cd(ily + 1);
+    gPad->SetGrid(0,0);
+    hlayer->Draw("colz");
+  }
+}
+
 //________________________________________________________
 Bool_t AliTRDcheckDET::MakeBarPlot(TH1 *histo, Int_t color){
   //
index 577979100483ce571e274798c2e0f24874130fa3..1722524a12862c626d326b6d7dda03c91fc11071 100644 (file)
@@ -18,6 +18,7 @@ class TObjArray;
 class TH1;
 class TMap;
 class AliESDHeader;
+class AliExternalTrackParam;
 class AliTRDcluster;
 class AliTRDseedV1;
 class AliTRDgeometry;
@@ -46,7 +47,8 @@ public:
     kTriggerPurity      = 15,
     kTrackStatus        = 16,
     kTrackletStatus     = 17,
-    kNTrackletsP        = 18
+    kNTrackletsP        = 18,
+    kNclustersLayer     = 19
   };
   enum FigureType_t{
     kFigNclustersTrack,
@@ -110,7 +112,9 @@ private:
   TH1* MakePlotNTracklets();
   Bool_t MakePlotPulseHeight();
   void MakePlotnTrackletsVsP();
+  void MakePlotMeanClustersLayer();
   Bool_t MakeBarPlot(TH1 *histo, Int_t Color);
+  void GetEtaPhiAt(AliExternalTrackParam *track, Double_t x, Double_t &eta, Double_t &phi);
 
   AliTRDeventInfo *fEventInfo;         //! ESD Header
   TMap *fTriggerNames;                 //! Containing trigger class names