]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
AliHLTFileWriter moved to libAliHLTUtil
authorrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 27 Apr 2007 11:24:47 +0000 (11:24 +0000)
committerrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 27 Apr 2007 11:24:47 +0000 (11:24 +0000)
HLT/BASE/util/AliHLTFileWriter.cxx [new file with mode: 0644]
HLT/BASE/util/AliHLTFileWriter.h [new file with mode: 0644]

diff --git a/HLT/BASE/util/AliHLTFileWriter.cxx b/HLT/BASE/util/AliHLTFileWriter.cxx
new file mode 100644 (file)
index 0000000..e0cfd10
--- /dev/null
@@ -0,0 +1,304 @@
+// $Id$
+
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Authors: Matthias Richter <Matthias.Richter@ift.uib.no>                *
+ *          for The ALICE Off-line 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   AliHLTFileWriter.cxx
+    @author Matthias Richter
+    @date   
+    @brief  HLT file writer component implementation. */
+
+#if __GNUC__>= 3
+using namespace std;
+#endif
+
+#include "AliHLTFileWriter.h"
+#include <TObjArray.h>
+#include <TObjString.h>
+#include <TMath.h>
+#include <TFile.h>
+
+/** the global object for component registration */
+AliHLTFileWriter gAliHLTFileWriter;
+
+/** ROOT macro for the implementation of ROOT specific class methods */
+ClassImp(AliHLTFileWriter)
+
+AliHLTFileWriter::AliHLTFileWriter()
+  :
+  AliHLTDataSink(),
+  fBaseName(""),
+  fExtension(""),
+  fDirectory(""),
+  fCurrentFileName(""),
+  fMode(0)
+{
+  // see header file for class documentation
+  // or
+  // refer to README to build package
+  // or
+  // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
+}
+
+AliHLTFileWriter::AliHLTFileWriter(const AliHLTFileWriter&)
+  :
+  AliHLTDataSink(),
+  fBaseName(""),
+  fExtension(""),
+  fDirectory(""),
+  fCurrentFileName(""),
+  fMode(0)
+{
+  // see header file for class documentation
+  HLTFatal("copy constructor untested");
+}
+
+AliHLTFileWriter& AliHLTFileWriter::operator=(const AliHLTFileWriter&)
+{ 
+  // see header file for class documentation
+  HLTFatal("assignment operator untested");
+  return *this;
+}
+
+AliHLTFileWriter::~AliHLTFileWriter()
+{
+  // see header file for class documentation
+
+  // file list and file name list are owner of their objects and
+  // delete all the objects
+}
+
+const char* AliHLTFileWriter::GetComponentID()
+{
+  // see header file for class documentation
+  return "FileWriter";
+}
+
+void AliHLTFileWriter::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
+{
+  // see header file for class documentation
+  list.clear();
+  list.push_back(kAliHLTAnyDataType);
+}
+
+AliHLTComponent* AliHLTFileWriter::Spawn()
+{
+  // see header file for class documentation
+  return new AliHLTFileWriter;
+}
+
+int AliHLTFileWriter::DoInit( int argc, const char** argv )
+{
+  // see header file for class documentation
+  int iResult=0;
+  TString argument="";
+  int bMissingParam=0;
+  for (int i=0; i<argc && iResult>=0; i++) {
+    argument=argv[i];
+    if (argument.IsNull()) continue;
+
+    // -basename
+    if (argument.CompareTo("-datafile")==0) {
+      if ((bMissingParam=(++i>=argc))) break;
+      fBaseName=argv[i];
+      TObjArray* pTokens=fBaseName.Tokenize(".");
+      if (pTokens) {
+       int iEntries=pTokens->GetEntries();
+       if (iEntries>1) {
+         int i=0;
+         fBaseName=((TObjString*)pTokens->At(i++))->GetString();
+         while (i<iEntries-1) {
+           fBaseName+="." + ((TObjString*)pTokens->At(i++))->GetString();
+         }
+         fExtension=((TObjString*)pTokens->At(i))->GetString();
+       }
+       delete pTokens;
+      }
+
+      // -directory
+    } else if (argument.CompareTo("-directory")==0) {
+      if ((bMissingParam=(++i>=argc))) break;
+      fDirectory=argv[i];
+
+      // -enumeration
+    } else if (argument.CompareTo("-enumerate")==0) {
+      SetMode(kEnumerate);
+
+      // -concatenate-blocks
+    } else if (argument.CompareTo("-concatenate-blocks")==0) {
+      SetMode(kConcatenateBlocks);
+
+      // -concatenate-events
+    } else if (argument.CompareTo("-concatenate-events")==0) {
+      SetMode(kConcatenateEvents);
+
+    } else {
+      if ((iResult=ScanArgument(argc-i, &argv[i]))==-EINVAL) {
+       HLTError("unknown argument %s", argument.Data());
+       break;
+      } else if (iResult==-EPROTO) {
+       bMissingParam=1;
+       break;
+      } else if (iResult>=0) {
+       i+=iResult;
+       iResult=0;
+      }
+    }
+  }
+  if (bMissingParam) {
+    HLTError("missing parameter for argument %s", argument.Data());
+    iResult=-EINVAL;
+  }
+  if (iResult>=0) {
+    iResult=InitWriter();
+  }
+
+  return iResult;
+}
+
+int AliHLTFileWriter::InitWriter()
+{
+  // see header file for class documentation
+  return 0; // note: this doesn't mean 'error'
+}
+
+int AliHLTFileWriter::ScanArgument(int argc, const char** argv)
+{
+  // see header file for class documentation
+
+  // there are no other arguments than the standard ones
+  if (argc==0 && argv==NULL) {
+    // this is just to get rid of the warning "unused parameter"
+  }
+  return -EINVAL;
+}
+
+int AliHLTFileWriter::DoDeinit()
+{
+  // see header file for class documentation
+  int iResult=CloseWriter();
+  ClearMode(kEnumerate);
+  return iResult;
+}
+
+int AliHLTFileWriter::CloseWriter()
+{
+  // see header file for class documentation
+  return 0; // note: this doesn't mean 'error'
+}
+
+int AliHLTFileWriter::DumpEvent( const AliHLTComponentEventData& evtData,
+                        const AliHLTComponentBlockData* blocks, 
+                        AliHLTComponentTriggerData& trigData )
+{
+  // see header file for class documentation
+  int iResult=0;
+  if (CheckMode(kConcatenateEvents)==0) {
+    // reset the current file name in order to open a new file
+    // for the first block. If events are concatenated, the current
+    // file name stays in order to be opended in append mode.
+    fCurrentFileName="";
+  }
+  for (int n=0; n<(int)evtData.fBlockCnt; n++ ) {
+    //HLTDebug("block %d out of %d", n, evtData.fBlockCnt);
+    TString filename;
+    HLTDebug("dataspec 0x%x", blocks[n].fSpecification);
+    iResult=BuildFileName(evtData.fEventID, n, blocks[n].fDataType, filename);
+    ios::openmode filemode=(ios::openmode)0;
+    if (fCurrentFileName.CompareTo(filename)==0) {
+      // append to the file
+      filemode=ios::app;
+    } else {
+      // store the file for the next block
+      fCurrentFileName=filename;
+    }
+    if (iResult>=0) {
+      ofstream dump(filename.Data(), filemode);
+      if (dump.good()) {
+       dump.write((const char*)blocks[n].fPtr, blocks[n].fSize);
+       HLTDebug("wrote %d byte(s) to file %s", blocks[n].fSize, filename.Data());
+      } else {
+       HLTError("can not open file %s for writing", filename.Data());
+       iResult=-EBADF;
+      }
+      dump.close();
+    }
+  }
+  if (trigData.fStructSize==0) {
+    // this is just to get rid of the warning "unused parameter"
+  }
+  return iResult;
+}
+
+int AliHLTFileWriter::BuildFileName(const AliHLTEventID_t eventID, const int blockID, const AliHLTComponentDataType& dataType, TString& filename)
+{
+  // see header file for class documentation
+  int iResult=0;
+  //HLTDebug("build file name for event %d block %d", eventID, blockID);
+  filename="";
+  if (!fDirectory.IsNull()) {
+    filename+=fDirectory;
+    if (!fDirectory.EndsWith("/"))
+      filename+="/";
+  }
+  if (!fBaseName.IsNull())
+    filename+=fBaseName;
+  else
+    filename+="event";
+  if (!CheckMode(kConcatenateEvents)) {
+    if (!CheckMode(kEnumerate)) {
+      if (eventID!=kAliHLTVoidEventID) {
+       filename+=Form("_0x%08x", eventID);
+      }
+    } else {
+      filename+=Form("_%d", GetEventCount());
+    }
+  }
+  if (blockID>=0 && !CheckMode(kConcatenateBlocks)) {
+    filename+=Form("_0x%x", blockID);
+    if (dataType!=kAliHLTVoidDataType) {
+      filename+="_";
+      filename+=AliHLTComponent::DataType2Text(dataType).data();
+    }
+  }
+  if (!fExtension.IsNull())
+    filename+="." + fExtension;
+  filename.ReplaceAll(" ", "");
+  return iResult;
+}
+
+int AliHLTFileWriter::SetMode(Short_t mode) 
+{
+  // see header file for class documentation
+  fMode|=mode;
+  //HLTDebug("mode set to 0x%x", fMode);
+  return fMode;
+}
+
+int AliHLTFileWriter::ClearMode(Short_t mode)
+{
+  // see header file for class documentation
+  fMode&=~mode;
+  //HLTDebug("mode set to 0x%x", fMode);
+  return fMode;
+}
+
+int AliHLTFileWriter::CheckMode(Short_t mode)
+{
+  // see header file for class documentation
+
+  //HLTDebug("check mode 0x%x for flag 0x%x: %d", fMode, mode, (fMode&mode)!=0);
+  return (fMode&mode)!=0;
+}
diff --git a/HLT/BASE/util/AliHLTFileWriter.h b/HLT/BASE/util/AliHLTFileWriter.h
new file mode 100644 (file)
index 0000000..bd994c7
--- /dev/null
@@ -0,0 +1,190 @@
+// @(#) $Id$
+
+#ifndef ALIHLTFILEWRITER_H
+#define ALIHLTFILEWRITER_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/** @file   AliHLTFileWriter.h
+    @author Matthias Richter
+    @date   
+    @brief  An HLT file dump (data sink) component.
+*/
+
+#include "AliHLTDataSink.h"
+#include <TList.h>
+
+/**
+ * @class AliHLTFileWriter
+ * An HLT data sink component which writes data to file(s).
+ *
+ * Component ID: \b FileWriter <br>
+ * Library: \b libHLTBase (in order to use the component from the external
+ * interface, it might be necessary to specify a dummy library with the
+ * \em -componentlibrary argument).
+ *
+ * Mandatory arguments: <br>
+ * <!-- NOTE: ignore the \li. <i> and </i>: it's just doxygen formating -->
+ *
+ * Optional arguments: <br>
+ * <!-- NOTE: ignore the \li. <i> and </i>: it's just doxygen formating -->
+ * \li -datafile     <i> filename   </i> <br>
+ *      file name base
+ * \li -directory    <i> directory  </i> <br>
+ *      target directory
+ * \li -enumerate <br>
+ *      don't use the event number but an event counter beginning from 0
+ * \li -concatenate-blocks <br>
+ *      concatenate all blocks of one event into one file, this skips
+ *      the block no, and the block data type in the file name
+ * \li -concatenate-events <br>
+ *      concatenate all events into one file, this skips the event no,
+ *      the block no, and the block data type in the file name. Currently,
+ *      this implies the -concatenate-blocks option.
+ *
+ * The file name is built from the basename, the event number, the block
+ * number and the data type in the format:
+ * <pre>
+ * basename_eventno_blockno_dt
+ * </pre>
+ * If the basename was not given, \em 'event' ist used instead. A file
+ * extension after the last dot is separated from the basename and appended
+ * to the final name.
+ *
+ * The class can be used as a base class for file writers. Additional
+ * argument scan can be implemented in @ref ScanArgument which is called
+ * for each unknown argument.
+ * @ingroup alihlt_component
+ */
+class AliHLTFileWriter : public AliHLTDataSink  {
+ public:
+  /** standard constructor */
+  AliHLTFileWriter();
+  /** not a valid copy constructor, defined according to effective C++ style */
+  AliHLTFileWriter(const AliHLTFileWriter&);
+  /** not a valid assignment op, but defined according to effective C++ style */
+  AliHLTFileWriter& operator=(const AliHLTFileWriter&);
+  /** destructor */
+  virtual ~AliHLTFileWriter();
+
+  virtual const char* GetComponentID();
+  virtual void GetInputDataTypes( vector<AliHLTComponentDataType>& list);
+  virtual AliHLTComponent* Spawn();
+
+ protected:
+  /**
+   * Init method.
+   */
+  int DoInit( int argc, const char** argv );
+
+  /**
+   * Deinit method.
+   */
+  int DoDeinit();
+
+  /**
+   * Init the writer.
+   * The DoInit function is not available for child classes. InitWriter is the
+   * corresponding function for classes derived from AliHLTFileWriter.
+   */
+  virtual int InitWriter();
+
+  /**
+   * Init the writer.
+   * The DoDeinit function is not available for child classes. CloseWriter is the
+   * corresponding function for classes derived from AliHLTFileWriter.
+   */
+  virtual int CloseWriter();
+
+  /**
+   * Data processing method for the component.
+   * The function can be overloaded by other file writer components.
+   * @param evtData       event data structure
+   * @param blocks        input data block descriptors
+   * @param trigData     trigger data structure
+   */
+  virtual int DumpEvent( const AliHLTComponentEventData& evtData,
+                        const AliHLTComponentBlockData* blocks, 
+                        AliHLTComponentTriggerData& trigData );
+
+  /**
+   * Scan one argument and adjacent parameters.
+   * Can be overloaded by child classes in order to add additional arguments
+   * beyond the standard arguments of the file publisher. The method is called
+   * whenever a non-standard argument is recognized. Make sure to return 
+   * <tt> -EPROTO </tt> if the argument is not recognized be the child.
+   * @param argc           size of the argument array
+   * @param argv           agument array for component initialization
+   * @return number of processed members of the argv <br>
+   *         -EINVAL unknown argument <br>
+   *         -EPROTO parameter for argument missing <br>
+   */
+  virtual int ScanArgument(int argc, const char** argv);
+
+  /**
+   * Build file name from eventID data type and the specified directory and basename.
+   * @param eventID [in]   the ID of the event
+   * @param blockID [in]   the ID of the current block
+   *                       no block string appended if -1
+   * @param dataType [in]  the data type of the data block
+   *                       no type string appanded if @ref kAliHLTVoidDataType
+   * @param filename [out] string to receive the file name
+   */
+  int BuildFileName(const AliHLTEventID_t eventID, const int blockID, const AliHLTComponentDataType& dataType, TString& filename);
+
+  /**
+   * Set a mode flag.
+   * @return current mode flags
+   */
+  int SetMode(Short_t mode);
+    
+  /**
+   * Clear a mode flag.
+   * @return current mode flags
+   */
+  int ClearMode(Short_t mode);
+
+  /**
+   * Check a mode flag.
+   * @return 1 if flag is set, 0 if not
+   */
+  int CheckMode(Short_t mode);
+
+  /**
+   * Working modes of the writer
+   * @internal
+   */
+  enum TWriterMode {
+    /**
+     * flag to indicate whether to write each incoming block to separate files
+     * or all blocks of one event to one file. set = concatenate (one file).
+     */
+    kConcatenateBlocks = 0x1,
+
+    /**
+     * flag to indicate whether to concatenate incoming blocks of the same type
+     * for all events to one file. If also @ref kConcatenateBlocks is set,
+     * or all blocks of all events are written to the same file.
+     */
+    kConcatenateEvents = 0x2,
+
+    /** event enumeration flag */
+    kEnumerate = 0x4
+  };
+
+ private:
+  /** the basename of the output file */
+  TString    fBaseName;                                            // see above
+  /** the extension of the output file */
+  TString    fExtension;                                           // see above
+  /** target directory */
+  TString    fDirectory;                                           // see above
+  /** enumeration format string */
+  TString    fCurrentFileName;                                     // see above
+
+  /** mode specifier, see @ref TWriterMode */
+  Short_t    fMode;                                                // see above
+
+  ClassDef(AliHLTFileWriter, 1)
+};
+#endif