]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Monitoring component for EMCAL electron trigger
authorfronchet <fronchet@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 18 Aug 2011 13:47:29 +0000 (13:47 +0000)
committerfronchet <fronchet@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 18 Aug 2011 13:47:29 +0000 (13:47 +0000)
HLT/trigger/AliHLTEmcalElectronMonitor.cxx [new file with mode: 0644]
HLT/trigger/AliHLTEmcalElectronMonitor.h [new file with mode: 0644]
HLT/trigger/AliHLTEmcalElectronMonitorComponent.cxx [new file with mode: 0644]
HLT/trigger/AliHLTEmcalElectronMonitorComponent.h [new file with mode: 0644]

diff --git a/HLT/trigger/AliHLTEmcalElectronMonitor.cxx b/HLT/trigger/AliHLTEmcalElectronMonitor.cxx
new file mode 100644 (file)
index 0000000..e51997f
--- /dev/null
@@ -0,0 +1,82 @@
+#include "AliHLTEmcalElectronMonitor.h"
+#include "AliHLTScalars.h"
+#include "TString.h"
+#include "TMath.h"
+
+
+ClassImp(AliHLTEmcalElectronMonitor);
+
+AliHLTEmcalElectronMonitor::AliHLTEmcalElectronMonitor():
+  hList(NULL),
+  hTracksPt(NULL),
+  hClusterEn(NULL),
+  hdEta(NULL),
+  hdPhi(NULL),
+  hdR(NULL),
+  hEoverP(NULL)
+{
+
+  // book histograms
+  hList = new TObjArray;
+
+  hTracksPt   = new TH1F("hTracksPt","Tracks pT (GeV/c)", 500, 0, 100);
+  hList->Add(hTracksPt);
+
+  hClusterEn  = new TH1F("hClusterEn","Cluster Energy (GeV)", 500, 0, 100);
+  hList->Add(hClusterEn);
+
+  hdEta = new TH1F("hdEta", "#Delta #eta (Cluster-Track)", 200,-.1,.1);
+  hList->Add(hdEta);
+
+  hdPhi = new TH1F("hdPhi", "#Delta #phi (Cluster-Track)", 200,-.1,.1);
+  hList->Add(hdPhi);
+  
+  hdR = new TH1F("hdR","#Delta R (Track-Cluster);#Delta R(#sqrt{#Delta #eta ^{2} +#Delta #Phi ^{2}})",200,0.,.1);
+  hList->Add(hdR);
+  
+  hEoverP=new TH1F("hEoverP","E/P for matched tracks;E/P",200,0.,10.);
+  hList->Add(hEoverP);
+  
+}
+//___________________________________________________________________________________________________________________________________________________
+
+AliHLTEmcalElectronMonitor::~AliHLTEmcalElectronMonitor()
+{
+
+  // default destructor
+
+}
+//___________________________________________________________________________________________________________________________________________________
+
+TObjArray* AliHLTEmcalElectronMonitor::GetHistograms()
+{
+
+  // pointer to histogram objects
+  
+  return hList;
+
+}
+//___________________________________________________________________________________________________________________________________________________
+
+Int_t AliHLTEmcalElectronMonitor::MakeHisto(AliHLTScalars *scalar)
+{
+
+  // make the histograms
+  
+  hTracksPt->Fill( scalar->GetScalar("TracksPt").Value() );
+  
+  hClusterEn->Fill( scalar->GetScalar("ClusterEn").Value() );
+  
+  hdEta->Fill( scalar->GetScalar("dEta").Value() );
+  
+  hdPhi->Fill( scalar->GetScalar("dPhi").Value() );
+  
+  hdR->Fill( scalar->GetScalar("dR").Value() );
+  
+  hEoverP->Fill( scalar->GetScalar("EoverP").Value() );
+  
+
+  return 0;
+
+}
+  
diff --git a/HLT/trigger/AliHLTEmcalElectronMonitor.h b/HLT/trigger/AliHLTEmcalElectronMonitor.h
new file mode 100644 (file)
index 0000000..561da28
--- /dev/null
@@ -0,0 +1,43 @@
+#ifndef ALIHLTEMCALELECTRONMONITOR_H
+#define ALIHLTEMCALELECTRONMONITOR_H
+
+#include "TH1F.h"
+#include "TObjArray.h"
+#include "TString.h"
+#include "AliHLTScalars.h"
+
+class AliHLTEmcalElectronMonitor : public TObject
+{
+
+ public:
+  // constructor
+  AliHLTEmcalElectronMonitor();
+
+  // destructor
+  virtual ~AliHLTEmcalElectronMonitor();
+
+  // make histos
+  Int_t MakeHisto(AliHLTScalars *scalar);
+
+  // retrieve histograms
+  TObjArray* GetHistograms();
+
+ private:
+  TObjArray *hList;
+  TH1F      *hTracksPt;
+  TH1F      *hClusterEn;
+  TH1F      *hdEta;
+  TH1F      *hdPhi;
+  TH1F      *hdR;
+  TH1F      *hEoverP;
+  
+  AliHLTEmcalElectronMonitor(const AliHLTEmcalElectronMonitor &);
+  AliHLTEmcalElectronMonitor & operator = (const AliHLTEmcalElectronMonitor &);
+
+  ClassDef(AliHLTEmcalElectronMonitor, 0);
+
+};
+
+#endif
+
+
diff --git a/HLT/trigger/AliHLTEmcalElectronMonitorComponent.cxx b/HLT/trigger/AliHLTEmcalElectronMonitorComponent.cxx
new file mode 100644 (file)
index 0000000..48bdc68
--- /dev/null
@@ -0,0 +1,147 @@
+#include "AliHLTEmcalElectronMonitorComponent.h"
+#include "AliHLTEmcalElectronMonitor.h"
+#include "AliHLTScalars.h"
+
+#include "TFile.h"
+#include "TString.h"
+
+AliHLTEmcalElectronMonitorComponent::AliHLTEmcalElectronMonitorComponent() :
+  fRootFileName("EmcalElectrontrigger_histos.root"),
+  fPushFraction(10),
+  fLocalEventCount(0),
+  fVerbose(0),
+  fHistoPtr(NULL)
+{
+
+  // default constructor
+
+}
+//____________________________________________________________________________________________________________
+
+AliHLTEmcalElectronMonitorComponent::~AliHLTEmcalElectronMonitorComponent()
+{
+
+  // default destructor
+
+}
+//____________________________________________________________________________________________________________
+
+int AliHLTEmcalElectronMonitorComponent::DoInit(int argc, const char **argv)
+{
+  // initialize
+
+  fHistoPtr = new AliHLTEmcalElectronMonitor();
+  for (int i = 0; i < argc; i++) {
+    if (!strcmp("-roothistofilename", argv[i]))
+      fRootFileName = argv[i+1];
+    if (!strcmp("-pushfraction", argv[i]))
+      fPushFraction = atoi(argv[i+1]);
+    if (!strcmp("-verbose", argv[i]))
+      fVerbose = atoi(argv[i+1]);
+  }
+
+  return 0;
+
+}
+//____________________________________________________________________________________________________________
+
+int AliHLTEmcalElectronMonitorComponent::Deinit()
+{
+  // de-initialize
+
+  if (fHistoPtr) {
+    delete fHistoPtr;
+    fHistoPtr = NULL;
+  }
+  return 0;
+
+}
+//____________________________________________________________________________________________________________
+const char* AliHLTEmcalElectronMonitorComponent::GetComponentID()
+{
+  // component id
+  
+  return "EmcalElectronMonitor";
+
+}
+//____________________________________________________________________________________________________________
+
+void AliHLTEmcalElectronMonitorComponent::GetInputDataTypes(vector<AliHLTComponentDataType> &list)
+{
+  // define input data types
+  
+  list.clear();
+  list.push_back(kAliHLTDataTypeEventStatistics|kAliHLTDataOriginHLT);
+
+}
+//____________________________________________________________________________________________________________
+
+AliHLTComponentDataType AliHLTEmcalElectronMonitorComponent::GetOutputDataType()
+{
+
+  // return output data types
+  
+  return kAliHLTDataTypeHistogram | kAliHLTDataOriginEMCAL;
+
+}
+//____________________________________________________________________________________________________________
+
+void AliHLTEmcalElectronMonitorComponent::GetOutputDataSize(unsigned long &constBase, double &inputMultiplier)
+{
+
+  // calculate output data size
+  
+  constBase = 0;
+  inputMultiplier = 100;
+
+}
+//____________________________________________________________________________________________________________
+
+int AliHLTEmcalElectronMonitorComponent::DoEvent(const AliHLTComponentEventData &evtData, const AliHLTComponentBlockData *blocks,
+                                          AliHLTComponentTriggerData &/*trigData*/, AliHLTUInt8_t */*outputPtr*/, AliHLTUInt32_t &/*size*/,
+                                          std::vector<AliHLTComponentBlockData> &/*outputBlocks*/)
+{
+
+  // do event
+
+  const AliHLTComponentBlockData *iter = NULL;
+  UInt_t specification = 0;
+  
+  for (unsigned long ij = 0; ij < evtData.fBlockCnt; ij++) {
+    AliHLTScalars *scalarPtr = NULL;
+    iter = blocks + ij;
+    if (fVerbose) PrintComponentDataTypeInfo(iter->fDataType);
+
+    if (iter->fDataType == kAliHLTDataTypeEventStatistics)  { 
+      scalarPtr = reinterpret_cast<AliHLTScalars*>(iter->fPtr);
+    }
+    else {
+      if (fVerbose)  HLTWarning("FastJet Monitor: Data block does not contain event stats - check if flag is set for histograming for AliHLTTriggerFastJet \n");
+    }
+    
+    specification |= iter->fSpecification;
+    fHistoPtr->MakeHisto(scalarPtr);
+  }
+
+  fLocalEventCount++;
+    
+  TFile rootHistFile(fRootFileName, "RECREATE");
+  
+  fHistoPtr->GetHistograms()->Write();
+
+  if (fLocalEventCount%fPushFraction == 0) {
+    if (fVerbose) cout << "Emcal Electron Monitor: pushback done at " << fLocalEventCount << " evens " << endl;
+    PushBack(fHistoPtr->GetHistograms(), kAliHLTDataTypeTObjArray | kAliHLTDataOriginEMCAL, specification);
+  }
+  
+  return 0;
+}
+//____________________________________________________________________________________________________________
+
+AliHLTComponent* AliHLTEmcalElectronMonitorComponent::Spawn()
+{
+  // spawn
+
+  return new AliHLTEmcalElectronMonitorComponent();
+
+}
diff --git a/HLT/trigger/AliHLTEmcalElectronMonitorComponent.h b/HLT/trigger/AliHLTEmcalElectronMonitorComponent.h
new file mode 100644 (file)
index 0000000..cf0da33
--- /dev/null
@@ -0,0 +1,63 @@
+#ifndef ALIHLTEMCALELECTRONMONITORCOMPONENT_H
+#define ALIHLTEMCALELECTRONMONITORCOMPONENT_H
+
+#include "AliHLTProcessor.h"
+
+class AliHLTEmcalElectronMonitor;
+
+class AliHLTEmcalElectronMonitorComponent : public AliHLTProcessor
+{
+ public:
+  // Constructor
+  AliHLTEmcalElectronMonitorComponent();
+
+  // Destructor
+  virtual ~AliHLTEmcalElectronMonitorComponent();
+
+  // interface function
+  const char* GetComponentID();
+
+  // interface function
+  void GetInputDataTypes(std::vector<AliHLTComponentDataType> &list);
+
+  // interface function
+  AliHLTComponentDataType GetOutputDataType();
+
+  // interface function
+  void GetOutputDataSize(unsigned long &constBase, double &inputMultiplier);
+
+  // interface function
+  int DoEvent(const AliHLTComponentEventData &evtData, const AliHLTComponentBlockData *blocks,
+             AliHLTComponentTriggerData &/*trigData*/, AliHLTUInt8_t */*outputPtr*/, AliHLTUInt32_t &/*size*/,
+             std::vector<AliHLTComponentBlockData> &/*outputBlocks*/);
+  
+  // interface function
+  AliHLTComponent* Spawn();
+
+ protected:
+  
+  // interface function
+  int DoInit(int argc, const char **argv);
+  int DoDeinit() {return 0;};
+  
+  using AliHLTProcessor::DoEvent;
+
+  // interface function
+  virtual int Deinit();
+
+
+ private:
+  TString fRootFileName;
+  int     fPushFraction;
+  int     fLocalEventCount;
+  int     fVerbose;
+
+  // pointer to the histo maker itself
+  AliHLTEmcalElectronMonitor *fHistoPtr;
+
+  AliHLTEmcalElectronMonitorComponent(const AliHLTEmcalElectronMonitorComponent &);
+  AliHLTEmcalElectronMonitorComponent & operator = (const AliHLTEmcalElectronMonitorComponent &);
+
+};
+
+#endif