]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
added generation of trigger map
authorjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 4 Dec 2007 18:36:47 +0000 (18:36 +0000)
committerjgrosseo <jgrosseo@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 4 Dec 2007 18:36:47 +0000 (18:36 +0000)
PWG0/trigger/AliGenTriggerMapSelector.cxx
PWG0/trigger/AliGenTriggerMapSelector.h
PWG0/trigger/runGenTriggerMapSelector.C

index 0614058cccee895e92b76d12a8414f99ee7c2ae7..19c913069714199fa4c2c9202f0857cede66517d 100644 (file)
@@ -230,6 +230,7 @@ Bool_t AliGenTriggerMapSelector::Process(Long64_t entry)
         if (chip[currentLayer] == -1)
         {
           chip[currentLayer] = moduleIndex * 5 + isChip;
+          //Printf("%d --> %d %d", label, moduleIndex, isChip);
           nClusters[currentLayer]++;
         }
         else
@@ -335,3 +336,35 @@ void AliGenTriggerMapSelector::ReadHistograms(const char* filename)
   fTracklets  = dynamic_cast<TNtuple*> (file->Get("fTracklets"));
   fChipsFired  = dynamic_cast<TH2F*> (file->Get("fChipsFired"));
 }
+
+void AliGenTriggerMapSelector::GenerateTriggerMap(Bool_t clean)
+{
+  Int_t nEntries = fTracklets->GetEntries();
+  //nEntries = 10;
+
+  TH2F* hist = new TH2F("hist", ";layer1;layer2", 400, -0.5, 399.5, 800, 399.5, 1199.5);
+
+  Long64_t nSkipped = 0;
+  for (Int_t i = 0; i < nEntries; ++i)
+  {
+    fTracklets->GetEntry(i);
+    Float_t* vars = fTracklets->GetArgs();
+
+    if (clean && (Bool_t) vars[3] != kFALSE)
+    {
+      nSkipped++;
+      continue;
+    }
+
+    hist->Fill(vars[1], vars[2]);
+  }
+  if (clean)
+    Printf("Skipped %lld dirty entries", nSkipped);
+
+  hist->Draw("COLZ");
+
+  for (Int_t i = 1; i <= hist->GetNbinsX(); ++i)
+    for (Int_t j = 1; j <= hist->GetNbinsY(); ++j)
+      if (hist->GetBinContent(i, j) > 0)
+        Printf("(%d,%d) ", TMath::Nint(hist->GetXaxis()->GetBinCenter(i)), TMath::Nint(hist->GetYaxis()->GetBinCenter(j)));
+}
index 433db6896940ba056491b83f29318cdc8edfb813..b29ba86a9a3563a95fc8e644e0912f0c2d61e3a7 100644 (file)
@@ -25,6 +25,7 @@ class AliGenTriggerMapSelector : public AliSelectorRL {
 
     void WriteHistograms(const char* filename = "triggerMap.root");
     void ReadHistograms(const char* filename = "triggerMap.root");
+    void GenerateTriggerMap(Bool_t clean = kFALSE);
 
  protected:
     TH2F* fChipsFired;
index 2daae22575be62d5cb7afa9669d4d077a2897039..34a2ca659dd9d2aedfa2b2e4eb2c5abc3359cfe0 100644 (file)
@@ -36,3 +36,11 @@ void runGenTriggerMapSelector(Char_t* data, Int_t nRuns=20, Int_t offset=0, Bool
   executeQuery(chain, &inputList, selectorName, option);
 }
 
+void generateTriggerMap(Bool_t clean = kFALSE)
+{
+  gSystem->Load("libPWG0base");
+  gROOT->ProcessLine(".L AliGenTriggerMapSelector.cxx+");
+  AliGenTriggerMapSelector selector;
+  selector.ReadHistograms();
+  selector.GenerateTriggerMap(clean);
+}