]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
added monitoring example component and macro
authorrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 19 Feb 2008 11:45:26 +0000 (11:45 +0000)
committerrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 19 Feb 2008 11:45:26 +0000 (11:45 +0000)
HLT/SampleLib/AliHLTAgentSample.cxx
HLT/SampleLib/AliHLTSampleMonitoringComponent.cxx [new file with mode: 0644]
HLT/SampleLib/AliHLTSampleMonitoringComponent.h [new file with mode: 0644]
HLT/exa/monitoring.C [new file with mode: 0644]
HLT/libAliHLTSample.pkg

index 8d6978c88cf323c462ab777d55eac18334c989a1..5f97332dd38a5aebcf354a864c67fa3c68a99f21 100644 (file)
@@ -31,6 +31,7 @@
 #include "AliHLTDummyComponent.h"
 #include "AliHLTSampleComponent1.h"
 #include "AliHLTSampleComponent2.h"
 #include "AliHLTDummyComponent.h"
 #include "AliHLTSampleComponent1.h"
 #include "AliHLTSampleComponent2.h"
+#include "AliHLTSampleMonitoringComponent.h"
 
 // header file of the module preprocessor
 #include "AliHLTSamplePreprocessor.h"
 
 // header file of the module preprocessor
 #include "AliHLTSamplePreprocessor.h"
@@ -124,6 +125,7 @@ int AliHLTAgentSample::RegisterComponents(AliHLTComponentHandler* pHandler) cons
   pHandler->AddComponent(new AliHLTDummyComponent);
   pHandler->AddComponent(new AliHLTSampleComponent1);
   pHandler->AddComponent(new AliHLTSampleComponent2);
   pHandler->AddComponent(new AliHLTDummyComponent);
   pHandler->AddComponent(new AliHLTSampleComponent1);
   pHandler->AddComponent(new AliHLTSampleComponent2);
+  pHandler->AddComponent(new AliHLTSampleMonitoringComponent);
 
   return 0;
 }
 
   return 0;
 }
diff --git a/HLT/SampleLib/AliHLTSampleMonitoringComponent.cxx b/HLT/SampleLib/AliHLTSampleMonitoringComponent.cxx
new file mode 100644 (file)
index 0000000..812ce4c
--- /dev/null
@@ -0,0 +1,253 @@
+// $Id$
+
+/**************************************************************************
+ * This file is property of and copyright by the ALICE HLT Project        * 
+ * ALICE Experiment at CERN, All rights reserved.                         *
+ *                                                                        *
+ * Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
+ *                  for The ALICE HLT Project.                            *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+/** @file   AliHLTSampleMonitoringComponent.cxx
+    @author Matthias Richter
+    @date   
+    @brief  A sample monitoring component for the HLT.
+*/
+
+#include "AliHLTSampleMonitoringComponent.h"
+#include "TH1F.h"
+#include "TTree.h"
+#include "TString.h"
+#include "TObjString.h"
+#include "TObjArray.h"
+#include "AliCDBEntry.h"
+#include "AliCDBManager.h"
+
+/** ROOT macro for the implementation of ROOT specific class methods */
+ClassImp(AliHLTSampleMonitoringComponent)
+
+AliHLTSampleMonitoringComponent::AliHLTSampleMonitoringComponent()
+  :
+  fPushHistograms(false),
+  fPushTTree(false),
+  fPushTObjArray(false),
+  fHpx(NULL),
+  fHpy(NULL)
+{
+  // see header file for class documentation
+  // or
+  // refer to README to build package
+  // or
+  // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
+}
+
+AliHLTSampleMonitoringComponent::~AliHLTSampleMonitoringComponent()
+{
+  // see header file for class documentation
+}
+
+const char* AliHLTSampleMonitoringComponent::GetComponentID()
+{
+  // see header file for class documentation
+  return "Sample-MonitoringComponent";
+}
+
+void AliHLTSampleMonitoringComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
+{
+  // see header file for class documentation
+  list.push_back(kAliHLTAnyDataType);
+}
+
+AliHLTComponentDataType AliHLTSampleMonitoringComponent::GetOutputDataType()
+{
+  // see header file for class documentation
+  return kAliHLTVoidDataType;
+}
+
+void AliHLTSampleMonitoringComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
+{
+  // see header file for class documentation
+  constBase = 100000;
+  inputMultiplier = 0;
+}
+
+AliHLTComponent* AliHLTSampleMonitoringComponent::Spawn()
+{
+  // see header file for class documentation
+  return new AliHLTSampleMonitoringComponent;
+}
+
+int AliHLTSampleMonitoringComponent::DoInit( int argc, const char** argv )
+{
+  // see header file for class documentation
+  int iResult=0;
+
+  TString argument="";
+  TString configuration=""; 
+  int bMissingParam=0;
+
+  for (int i=0; i<argc && iResult>=0; i++) {
+    argument=argv[i];
+    if (argument.IsNull()) continue;
+
+    // -push-histograms
+    if (argument.CompareTo("-push-histograms")==0) {
+      fPushHistograms=true;
+
+    // -push-ttree
+    } else if (argument.CompareTo("-push-ttree")==0) {
+      fPushTTree=true;
+
+    // -push-array
+    } else if (argument.CompareTo("-push-array")==0) {
+      fPushTObjArray=true;
+
+    } else {
+      // the remaining arguments are treated as configuration
+      if (!configuration.IsNull()) configuration+=" ";
+      configuration+=argument;
+    }
+  }
+  if (bMissingParam) {
+    HLTError("missing parameter for argument %s", argument.Data());
+    iResult=-EINVAL;
+  }
+
+  // choose fPushTTree as default if none is set
+  if (!(fPushTTree || fPushTObjArray || fPushHistograms)) fPushTTree=true;
+
+  // strictly speaking I would prefer to use local or dynamic variables
+  // locally in DoEvent, but there is a ROOT bug or feature (related to
+  // garbage collection) which causes seg faults after a while. 
+  fHpx = new TH1F("hpx","px distribution",100,-4,4);
+  fHpy = new TH1F("hpy","py distribution",100,-10,10);
+
+  return iResult;
+}
+
+int AliHLTSampleMonitoringComponent::DoDeinit()
+{
+  // see header file for class documentation
+  if (fHpx) delete fHpx;
+  fHpx=NULL;
+  if (fHpy) delete fHpy;
+  fHpy=NULL;
+  return 0;
+}
+
+int AliHLTSampleMonitoringComponent::DoEvent( const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& /*trigData*/)
+{
+  // see header file for class documentation
+
+  // the function ignores all input blocks and fakes some monitoring histogram
+  int iResult=0;
+
+  fHpx->Reset();
+  fHpx->FillRandom("gaus",100*(GetEventCount()+1));
+
+  fHpy->Reset();
+  fHpy->FillRandom("gaus",500*(GetEventCount()+1));
+
+  if (fPushHistograms) {
+    PushBack(fHpx, "ROOTTH1F", "EXPL", 0);
+    PushBack(fHpy, "ROOTTH1F", "EXPL", 1);
+  }
+
+  if (fPushTTree) {
+    TString event;
+    TTree *pTree = new TTree("T","A Root Tree");
+    if (pTree) {
+      pTree->SetDirectory(0);
+      event.Form("event_%d_hpx", GetEventCount());
+      pTree->Branch(event, "TH1F", &fHpx, 32000, 0);
+      event.Form("event_%d_hpy", GetEventCount());
+      pTree->Branch(event, "TH1F", &fHpy, 32000, 0);
+
+      PushBack(pTree, "ROOTTREE", "EXPL");
+      delete pTree;
+    } else {
+      iResult=-ENOMEM;
+    }
+  }
+
+  if (fPushTObjArray) {
+    TObjArray* pArray=new TObjArray;
+    if (pArray) {
+      pArray->Add(fHpx);
+      pArray->Add(fHpy);
+      
+      PushBack(pArray, "ROOTOBJA", "EXPL");
+      delete pArray;
+    }
+  }
+
+  return iResult;
+}
+
+int AliHLTSampleMonitoringComponent::Configure(const char* arguments)
+{
+  // see header file for class documentation
+  int iResult=0;
+  if (!arguments) return iResult;
+
+  TString allArgs=arguments;
+  TString argument;
+  int bMissingParam=0;
+
+  TObjArray* pTokens=allArgs.Tokenize(" ");
+  if (pTokens) {
+    for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
+      argument=((TObjString*)pTokens->At(i))->GetString();
+      if (argument.IsNull()) continue;
+
+      // -reset
+      if (argument.CompareTo("-reset")==0) {
+
+      } else {
+       HLTError("unknown argument %s", argument.Data());
+       iResult=-EINVAL;
+       break;
+      }
+    }
+    delete pTokens;
+  }
+  if (bMissingParam) {
+    HLTError("missing parameter for argument %s", argument.Data());
+    iResult=-EINVAL;
+  }
+  return iResult;
+}
+
+int AliHLTSampleMonitoringComponent::Reconfigure(const char* cdbEntry, const char* chainId)
+{
+  // see header file for class documentation
+  int iResult=0;
+  const char* path=NULL;
+  const char* defaultNotify="";
+  if (cdbEntry) {
+    path=cdbEntry;
+    defaultNotify=" (default)";
+  }
+  if (path) {
+    AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
+    if (pEntry) {
+      TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
+      if (pString) {
+       iResult=Configure(pString->GetString().Data());
+      } else {
+       HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
+      }
+    } else {
+      HLTError("can not fetch object \"%s\" from CDB", path);
+    }
+  }
+  return iResult;
+}
diff --git a/HLT/SampleLib/AliHLTSampleMonitoringComponent.h b/HLT/SampleLib/AliHLTSampleMonitoringComponent.h
new file mode 100644 (file)
index 0000000..3528f1c
--- /dev/null
@@ -0,0 +1,146 @@
+//-*- Mode: C++ -*-
+// $Id$
+#ifndef ALIHLTSAMPLEMONITORINGCOMPONENT_H
+#define ALIHLTSAMPLEMONITORINGCOMPONENT_H
+
+//* This file is property of and copyright by the ALICE HLT Project        * 
+//* ALICE Experiment at CERN, All rights reserved.                         *
+//* See cxx source for full Copyright notice                               */
+
+/** @file   AliHLTSampleMonitoringComponent.h
+    @author Matthias Richter
+    @date   
+    @brief  A sample monitoring component for the HLT.
+*/
+
+#include "AliHLTProcessor.h"
+
+class TH1F;
+
+/**
+ * @class AliHLTSampleMonitoringComponent
+ * An HLT sample component.
+ * This component does not any data processing at all. It just
+ * illustrates the existence of several components in ine library and
+ * allows to set up a very simple chain with different components.
+ *
+ * The component generates two histograms and accumulates random data
+ * as events are coming in. The two histograms are sent out via tree
+ * different approaches, see output data types.
+ *
+ * <h2>General properties:</h2>
+ *
+ * Component ID: \b Sample-MonitoringComponent <br>
+ * Library: \b libAliHLTSample.so     <br>
+ * Input Data Types: @ref kAliHLTAnyDataType <br>
+ *    in fact, the component ignores all incoming data blocks
+ *
+ * Output Data Types:
+ * \it {ROOT_TH1,EXPL}                                              <br>
+ *     one histogram per data block, specification identifies the
+ *     specific histogram.
+ * \it {ROOTOBJA,EXPL}                                              <br>
+ *     the two histograms are added to a TOBjArray which is pushed
+ *     to the output stream
+ * \it {ROOTTREE,EXPL}                                              <br>
+ *     the two histograms are added to a TTree which is pushed
+ *     to the output stream
+ *
+ * <h2>Mandatory arguments:</h2>
+ * <!-- NOTE: ignore the \li. <i> and </i>: it's just doxygen formating -->
+ *
+ * <h2>Optional arguments:</h2>
+ * <!-- NOTE: ignore the \li. <i> and </i>: it's just doxygen formating -->
+ * \li -push-histograms                              <br>
+ *      push histograms idividually
+ * \li -push-ttree (default)                         <br>
+ *      push histograms embedded in TTree
+ * \li -push-array                                   <br>
+ *      push histograms embedded in TObjArray
+ *
+ * <h2>Configuration:</h2>
+ * <!-- NOTE: ignore the \li. <i> and </i>: it's just doxygen formating -->
+ *
+ * <h2>Default CDB entries:</h2>
+ * The component has no default CDB entry.
+ * It does not load any configuration from the global <tt>ConfigHLT</tt>
+ * folder.
+ * 
+ * For re-configuration/steering the component expects a TObjString object
+ * holding a string with the configuration parameters explained above.
+ *
+ * <h2>Performance:</h2>
+ * The component does not any event data processing.
+ *
+ * <h2>Memory consumption:</h2>
+ * The component does not any event data processing.
+ *
+ * <h2>Output size:</h2>
+ * The component has no output data.
+ *
+ * The component implements the @ref alihltcomponent-low-level-interface.
+ * for data processing.
+ *
+ * Using the latter case, a component can also be reconfigured. Special
+ * events are propageted through the chain in order to trigger the re-
+ * configuration. The component needs to implement the function
+ * @ref Reconfigure(). The simplest version of a configuration object is
+ * a string object (TObjString) containing configuration arguments.
+ *
+ * @ingroup alihlt_tutorial
+ */
+class AliHLTSampleMonitoringComponent : public AliHLTProcessor {
+public:
+  AliHLTSampleMonitoringComponent();
+  virtual ~AliHLTSampleMonitoringComponent();
+
+  // AliHLTComponent interface functions
+  const char* GetComponentID();
+  void GetInputDataTypes( vector<AliHLTComponentDataType>& list);
+  AliHLTComponentDataType GetOutputDataType();
+  void GetOutputDataSize( unsigned long& constBase, double& inputMultiplier );
+  AliHLTComponent* Spawn();
+
+ protected:
+  // AliHLTComponent interface functions
+  int DoInit( int argc, const char** argv );
+  int DoDeinit();
+  int DoEvent( const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& trigData);
+  int Reconfigure(const char* cdbEntry, const char* chainId);
+
+  using AliHLTProcessor::DoEvent;
+
+private:
+  /** copy constructor prohibited */
+  AliHLTSampleMonitoringComponent(const AliHLTSampleMonitoringComponent&);
+  /** assignment operator prohibited */
+  AliHLTSampleMonitoringComponent& operator=(const AliHLTSampleMonitoringComponent&);
+
+  /**
+   * Configure the component.
+   * Parse a string for the configuration arguments and set the component
+   * properties.
+   *
+   * This function illustrates the scanning of an argument string. The string
+   * was presumably fetched from the CDB.
+   */
+  int Configure(const char* arguments);
+
+  /** send histograms individually as data blocks */
+  bool fPushHistograms; // !transient
+
+  /** send histograms TTree embedded */
+  bool fPushTTree; // !transient
+
+  /** send histograms Object Array embedded */
+  bool fPushTObjArray; // !transient
+
+  /** test histogram */
+  TH1F* fHpx; //!transient
+
+  /** test histogram */
+  TH1F* fHpy; //!transient
+
+  ClassDef(AliHLTSampleMonitoringComponent, 0)
+};
+#endif
diff --git a/HLT/exa/monitoring.C b/HLT/exa/monitoring.C
new file mode 100644 (file)
index 0000000..e394521
--- /dev/null
@@ -0,0 +1,62 @@
+// $Id$
+/**
+ * @file monitoring.C
+ * @brief Sample macro for the a monitoring component;
+ *
+ * Usage:
+ * <pre>
+ *   aliroot -b -q monitoring.C | tee monitoring.log
+ * </pre>
+ *
+ * This macro illustrates the creation of an HLT monitoring component.
+ * The \b Sample-MonitoringComponent component ignores all input data
+ * and just fakes some histograms. The histograms can be sent out via
+ * different ways. A real component has to have input data and a parent
+ * component providing the data in the chain.
+ *
+ * The ROOTFileWriter component (AliHLTRootFileWriterComponent) provides
+ * a simple means to write ROOT objects to a ROOT file.
+ *
+ * See AliHLTSampleMonitoringComponent for detailed description.
+ *
+ * @author Matthias.Richter@ift.uib.no
+ * @ingroup alihlt_tutorial
+ */
+{
+  /////////////////////////////////////////////////////////////////////////
+  /////////////////////////////////////////////////////////////////////////
+  /////////////////////////////////////////////////////////////////////////
+  // global initialization of the HLT
+
+  // this is just a tool to switch the logging systems
+  AliHLTLogging log;
+  //log.SwitchAliLog(0);
+
+  AliHLTSystem gHLT;
+  gHLT.SetGlobalLoggingLevel(0x3c);
+
+  // load the component library
+  gHLT.LoadComponentLibraries("libAliHLTUtil.so");
+  gHLT.LoadComponentLibraries("libAliHLTSample.so");
+
+  /////////////////////////////////////////////////////////////////////////
+  /////////////////////////////////////////////////////////////////////////
+  /////////////////////////////////////////////////////////////////////////
+  // now we build up a small chain
+
+  // publisher for the reconfigure event
+  TString arg;
+  AliHLTConfiguration monitoring("monitoring", "Sample-MonitoringComponent", NULL , "-push-histograms");
+
+  AliHLTConfiguration writer("writer", "ROOTFileWriter", "monitoring" , "");
+
+  // run the chain
+  gHLT.BuildTaskList("writer");
+  gHLT.Run(1);
+
+  /////////////////////////////////////////////////////////////////////////
+  /////////////////////////////////////////////////////////////////////////
+  /////////////////////////////////////////////////////////////////////////
+  // cleanup
+
+}
index b375f8f53be1fb821adcf571c2041e9db40239ff..957f37e7690c9618c3b7dbc4c7d7a6c25b0dda9f 100644 (file)
@@ -6,23 +6,19 @@
 # and linking process. For further information refer to the 
 # README.
 
 # and linking process. For further information refer to the 
 # README.
 
-# library sources
-MODULE_SRCS=   AliHLTSampleComponent1.cxx \
-               AliHLTSampleComponent2.cxx \
-               AliHLTAgentSample.cxx \
-               AliHLTSamplePreprocessor.cxx \
-               AliHLTSampleOfflineSinkComponent.cxx \
-               AliHLTDummyComponent.cxx
-
 # class header files, the link definition for the root dictionary
 # will be created from the names of the header files
 # class header files, the link definition for the root dictionary
 # will be created from the names of the header files
-CLASS_HDRS:=   AliHLTSampleComponent1.h \
+CLASS_HDRS:=   AliHLTSampleComponent1.h \
                AliHLTSampleComponent2.h \
                AliHLTSampleComponent2.h \
+               AliHLTSampleMonitoringComponent.h \
                AliHLTAgentSample.h \
                AliHLTSamplePreprocessor.h \
                AliHLTSampleOfflineSinkComponent.h \
                AliHLTDummyComponent.h
 
                AliHLTAgentSample.h \
                AliHLTSamplePreprocessor.h \
                AliHLTSampleOfflineSinkComponent.h \
                AliHLTDummyComponent.h
 
+# library sources, generated from the class headers
+MODULE_SRCS=   $(CLASS_HDRS:.h=.cxx)
+
 # library headers
 # in most cases you might have already added all the header files to
 # the CLASS_HDRS variable. So we just use the content of this. You
 # library headers
 # in most cases you might have already added all the header files to
 # the CLASS_HDRS variable. So we just use the content of this. You