]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
New classes that updates histograms contniously during run.
authorphille <phille@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 11 Mar 2007 20:09:21 +0000 (20:09 +0000)
committerphille <phille@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sun, 11 Mar 2007 20:09:21 +0000 (20:09 +0000)
TH1F Histograms for each crystal are written to file
at the end of run. There is one root file for each
readoutpaartition.

HLT/PHOS/AliHLTPHOSRcuHistogramProducer.cxx [new file with mode: 0644]
HLT/PHOS/AliHLTPHOSRcuHistogramProducer.h [new file with mode: 0644]
HLT/PHOS/AliHLTPHOSRcuHistogramProducerComponent.cxx [new file with mode: 0644]
HLT/PHOS/AliHLTPHOSRcuHistogramProducerComponent.h [new file with mode: 0644]

diff --git a/HLT/PHOS/AliHLTPHOSRcuHistogramProducer.cxx b/HLT/PHOS/AliHLTPHOSRcuHistogramProducer.cxx
new file mode 100644 (file)
index 0000000..bfad40f
--- /dev/null
@@ -0,0 +1,212 @@
+/**************************************************************************
+ * Copyright(c) 2006, ALICE Experiment at CERN, All rights reserved.      *
+ *                                                                        *
+ * Authors: Boris Polichtchouk & Per Thomas Hille for the ALICE           *
+ * offline/HLT Project. Contributors are mentioned in the code where      *
+ * appropriate.                                                           *
+ *                                                                        *
+ * 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.                  *
+ **************************************************************************/
+
+#include "AliHLTPHOSRcuHistogramProducer.h"
+#include <iostream>
+#include "stdio.h"
+#include <cstdlib>
+#include "AliHLTPHOSRcuCellEnergyDataStruct.h"
+#include "TFile.h"
+
+/*************************************************************************
+* Class AliHLTPHOSRcuHistogramProducer accumulating histograms     *
+* with amplitudes per PHOS channel                                       *
+* It is intended to run at the HLT farm                                  *
+* and it fills the histograms with amplitudes per channel.               * 
+* Usage example see in PHOS/macros/Shuttle/AliPHOSCalibHistoProducer.C   *
+**************************************************************************/
+AliHLTPHOSRcuHistogramProducer:: AliHLTPHOSRcuHistogramProducer(): fEventCount(0), fEquippmentID(0), fModuleID(0), fRcuX(0), fRcuZ(0)
+{
+  cout << "WARNING: You cannot invoke the AliHLTPHOSRcuHistogramProducer without arguments" << endl;
+  cout << "Usage AliHLTPHOSRcuHistogramProducer(ModuleID, X. Z)" << endl;
+  //  Reset();
+  //  Init();
+} 
+
+AliHLTPHOSRcuHistogramProducer::AliHLTPHOSRcuHistogramProducer(AliHLTUInt8_t moduleID, AliHLTUInt8_t rcuX, AliHLTUInt8_t rcuZ)
+{
+  SetModuleID(moduleID);
+  SetRcuX(rcuX); 
+  SetRcuZ(rcuZ); 
+  Init();
+}
+
+AliHLTPHOSRcuHistogramProducer::~ AliHLTPHOSRcuHistogramProducer()
+{
+
+}
+
+
+void
+AliHLTPHOSRcuHistogramProducer::Init()
+{
+  char tmpHistoName[256];
+  int geomx;
+  int geomz;
+    
+
+  for(int x = 0; x < N_XCOLUMNS_RCU; x ++)
+    {
+      for(int z = 0; z < N_ZROWS_RCU; z ++)
+       {
+         for(int gain = 0; gain < N_GAINS; gain ++)
+           {
+             geomx = x + N_XCOLUMNS_RCU*fRcuX;
+             geomz = z + N_ZROWS_RCU*fRcuZ;
+
+             fEnergyAverageValues[x][z][gain] = 0; 
+             fAccumulatedValues[x][z][gain]   = 0;
+             fTimingAverageValues[x][z][gain] = 0; 
+             fHits[x][z][gain]                = 0;
+             sprintf(tmpHistoName, "Edistribution_%d_%d_%d_%d",(int)fModuleID,  geomx, geomz, gain);
+             fEnergyHistogramPtrs[x][z][gain] = 0;
+
+             fEnergyHistogramPtrs[x][z][gain] = new TH1F( tmpHistoName, tmpHistoName, N_BINS, XBIN_LOW, XBIN_UP);
+
+             sprintf(tmpHistoName, "TOFdistribution_%d_%d_%d_%d",(int)fModuleID,  geomx, geomz, gain);
+       
+             fTimingHistogramPtrs[x][z][gain] = 0;
+
+             fTimingHistogramPtrs[x][z][gain] = new TH1F(tmpHistoName , tmpHistoName, N_BINS, XBIN_LOW, XBIN_UP);
+
+             fCellAccEnergy.fAccumulatedEnergies[x][z][gain] = 0;
+             fCellAccEnergy.fHits[x][z][gain] = 0;
+             fCellAccEnergy.fModuleID = 0;
+             fCellAccEnergy.fRcuX = 0;
+             fCellAccEnergy.fRcuZ = 0; 
+           }
+       } 
+    } 
+}
+
+int 
+AliHLTPHOSRcuHistogramProducer::GetEquippmentId()
+{
+  return  fEquippmentID;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducer::SetRcuX(AliHLTUInt8_t X)
+{
+  fRcuX = X; 
+  fCellAccEnergy.fRcuX = X;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducer::SetRcuZ(AliHLTUInt8_t Z)
+{
+  fRcuZ = Z; 
+  fCellAccEnergy.fRcuZ = Z;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducer::SetModuleID(AliHLTUInt8_t moduleID)
+{
+ fModuleID = moduleID;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducer::SetEquippmentId(int id)
+{
+  fEquippmentID = id;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducer::FillEnergy(AliHLTUInt8_t x, AliHLTUInt8_t z,  AliHLTUInt8_t gain, float energy)
+{
+  fCellAccEnergy.fAccumulatedEnergies[x][z][gain] += energy;
+  fCellAccEnergy.fHits[x][z][gain] ++; 
+  fEnergyHistogramPtrs[x][z][gain]->Fill(energy);
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducer::FillTime(AliHLTUInt8_t x, AliHLTUInt8_t z, AliHLTUInt8_t gain, float time)
+{
+  fTimingHistogramPtrs[x][z][gain]->Fill(time);
+}
+
+const AliHLTPHOSRcuCellAccumulatedEnergyDataStruct& 
+AliHLTPHOSRcuHistogramProducer::GetCellAccumulatedEnergies()
+{
+  //  return &fCellAccEnergy ;
+  return fCellAccEnergy ;
+}
+
+int 
+AliHLTPHOSRcuHistogramProducer::IncrementEventCounter()
+{
+   fEventCount ++;
+   return  fEventCount;
+}
+
+
+void
+AliHLTPHOSRcuHistogramProducer::Reset()
+{
+  for(int x = 0; x < N_XCOLUMNS_RCU; x ++)
+    {
+      for(int z = 0; z < N_ZROWS_RCU; z ++)
+       {
+         for(int gain = 0; gain < N_GAINS; gain ++)
+           {
+             fEnergyAverageValues[x][z][gain] = 0; 
+             fAccumulatedValues[x][z][gain]   = 0;
+             fTimingAverageValues[x][z][gain] = 0; 
+             fHits[x][z][gain]                = 0;
+           }
+       } 
+    }
+  
+  for(int i = 0; i <ALTRO_MAX_SAMPLES;  i++)
+    {
+      fTmpChannelData[i] = 0;
+    }
+}
+
+void 
+AliHLTPHOSRcuHistogramProducer::WriteEnergyHistograms()
+{
+  char tmpFileName[256];
+  sprintf(tmpFileName,"/home/aliphoshlt/rundir/outdata/calibHisto_%d_%d_%d.root", (int)fModuleID, (int)fRcuX, (int)fRcuZ);
+  TFile *histoFile =  new TFile(tmpFileName,"update");
+  char hname[128];
+  if(!histoFile) return;
+  if(!histoFile->IsOpen()) return;
+
+  cout <<"printing histograms"<< endl;
+  cout <<"histofile-Getname() =" << histoFile->GetName() << endl;
+
+    for(int x = 0; x <  N_XCOLUMNS_RCU; x ++)
+    {
+      for(int z = 0; z < N_ZROWS_RCU; z ++)
+       {
+         for(int gain = 0; gain < N_GAINS; gain ++)
+           {
+             fEnergyHistogramPtrs[x][z][gain]->Write();
+           }
+       } 
+    }
+
+    cout << "printing histograms, finished"<< endl;
+    histoFile->Close();
+
+}
diff --git a/HLT/PHOS/AliHLTPHOSRcuHistogramProducer.h b/HLT/PHOS/AliHLTPHOSRcuHistogramProducer.h
new file mode 100644 (file)
index 0000000..ccbf605
--- /dev/null
@@ -0,0 +1,65 @@
+#ifndef ALIHLTPHOSRCUHISTOGRAMPRODUCER_H
+#define ALIHLTPHOSRCUHISTOGRAMPRODUCER_H 
+
+/* Copyright(c) 2006, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice  */ 
+
+#include "AliHLTPHOSDefinitions.h"
+#include "AliHLTPHOSCommonDefs.h"
+#include "TH1.h"
+//#include "AliHLTPHOSModuleCellAccumulatedEnergyDataStruct.h"
+#include "AliHLTPHOSRcuCellAccumulatedEnergyDataStruct.h"
+
+#define XBIN_LOW  0
+#define XBIN_UP   1023
+#define N_BINS    1023
+
+class AliHLTPHOSRcuHistogramProducer
+{
+ public:
+  AliHLTPHOSRcuHistogramProducer();
+  AliHLTPHOSRcuHistogramProducer(AliHLTUInt8_t moduleID, AliHLTUInt8_t rcuX, AliHLTUInt8_t rcuZ);
+
+  virtual ~AliHLTPHOSRcuHistogramProducer();
+  int GetEquippmentId();
+  const AliHLTPHOSRcuCellAccumulatedEnergyDataStruct& GetCellAccumulatedEnergies(); 
+  int IncrementEventCounter();
+  void Init();
+  void SetEquippmentId(int id);
+  void SetRcuX(AliHLTUInt8_t X);
+  void SetRcuZ(AliHLTUInt8_t Z);
+  void SetModuleID(AliHLTUInt8_t moduleID); 
+  void FillEnergy(AliHLTUInt8_t x, AliHLTUInt8_t z,  AliHLTUInt8_t gain, float energy);
+  void FillTime(AliHLTUInt8_t x,   AliHLTUInt8_t z,  AliHLTUInt8_t gain, float time); 
+  void Reset();
+  void WriteEnergyHistograms();
+
+
+ protected:
+ private:
+  AliHLTPHOSRcuHistogramProducer(const AliHLTPHOSRcuHistogramProducer & );
+  AliHLTPHOSRcuHistogramProducer & operator = (const AliHLTPHOSRcuHistogramProducer &)
+   {
+      return *this;
+   };
+
+  TH1F *fEnergyHistogramPtrs[N_XCOLUMNS_RCU][N_ZROWS_RCU][N_GAINS];
+  TH1F *fTimingHistogramPtrs[N_XCOLUMNS_RCU][N_ZROWS_RCU][N_GAINS];
+
+  Float_t fEnergyAverageValues[N_XCOLUMNS_RCU][N_ZROWS_RCU][N_GAINS];  
+  Double_t fAccumulatedValues[N_XCOLUMNS_RCU][N_ZROWS_RCU][N_GAINS];
+  Float_t fTimingAverageValues[N_XCOLUMNS_RCU][N_ZROWS_RCU][N_GAINS]; 
+  AliHLTUInt32_t fHits[N_XCOLUMNS_RCU][N_ZROWS_RCU][N_GAINS];
+  int fEventCount;
+  AliHLTUInt32_t fEquippmentID;
+  Double_t fTmpChannelData[ALTRO_MAX_SAMPLES];
+  //  AliHLTPHOSModuleCellAccumulatedEnergyDataStruct*  fOutPtr;
+  AliHLTPHOSRcuCellAccumulatedEnergyDataStruct fCellAccEnergy;
+  //  AliHLTPHOSRcuCellAccumulatedEnergyDataStruct*  fOutPtr;
+  AliHLTUInt8_t fModuleID;
+  AliHLTUInt8_t fRcuX;
+  AliHLTUInt8_t fRcuZ;
+};
+
+#endif
diff --git a/HLT/PHOS/AliHLTPHOSRcuHistogramProducerComponent.cxx b/HLT/PHOS/AliHLTPHOSRcuHistogramProducerComponent.cxx
new file mode 100644 (file)
index 0000000..705f3d6
--- /dev/null
@@ -0,0 +1,318 @@
+/**************************************************************************
+ * Copyright(c) 2006, ALICE Experiment at CERN, All rights reserved.      *
+ *                                                                        *
+ * Authors: Boris Polichtchouk & Per Thomas Hille for the ALICE           *
+ * offline/HLT Project. Contributors are mentioned in the code where      *
+ * appropriate.                                                           *
+ *                                                                        *
+ * 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.                  *
+ **************************************************************************/
+
+#include <iostream>
+#include "stdio.h"
+#include <cstdlib>
+#include "AliHLTPHOSRcuCellEnergyDataStruct.h"
+//#include "AliHLTPHOSModuleCellAccumulatedEnergyDataStruct.h"
+#include "AliHLTPHOSRcuHistogramProducer.h"
+#include "AliHLTPHOSRcuHistogramProducerComponent.h"
+
+
+
+const AliHLTComponentDataType  AliHLTPHOSRcuHistogramProducerComponent::inputDataTypes[]={kAliHLTVoidDataType,{0,"",""}}; //'zero' terminated array
+const AliHLTComponentDataType  AliHLTPHOSRcuHistogramProducerComponent::outputDataType=kAliHLTVoidDataType;
+AliHLTPHOSRcuHistogramProducerComponent gAliHLTPHOSRcuHistogramProducerComponent;
+
+//AliHLTPHOSHistogramProducerComponent gAliHLTPHOSHistogramProducerComponent;
+/*************************************************************************
+* Class AliHLTPHOSRcuHistogramProducerComponent accumulating histograms  *
+* with amplitudes per PHOS channel                                       *
+* It is intended to run at the HLT farm                                  *
+* and it fills the histograms with amplitudes per channel.               * 
+* Usage example see in PHOS/macros/Shuttle/AliPHOSCalibHistoProducer.C   *
+**************************************************************************/
+AliHLTPHOSRcuHistogramProducerComponent:: AliHLTPHOSRcuHistogramProducerComponent():AliHLTProcessor(), fEventCount(0), fRcuHistoProducerPtr(0)
+{
+  //  Reset();
+} 
+
+
+
+AliHLTPHOSRcuHistogramProducerComponent::~ AliHLTPHOSRcuHistogramProducerComponent()
+{
+
+}
+
+
+AliHLTPHOSRcuHistogramProducerComponent::AliHLTPHOSRcuHistogramProducerComponent(const  AliHLTPHOSRcuHistogramProducerComponent & ) : AliHLTProcessor(), fEventCount(0), fRcuHistoProducerPtr(0)
+{
+
+}
+
+
+int 
+AliHLTPHOSRcuHistogramProducerComponent::Deinit()
+{
+  cout << "AliHLTPHOSRcuHistogramProducerComponent::Deinit()" << endl;
+  fRcuHistoProducerPtr->WriteEnergyHistograms();
+  return 0;
+}
+
+
+int 
+AliHLTPHOSRcuHistogramProducerComponent::DoDeinit()
+{
+  Logging(kHLTLogInfo, "HLT", "PHOS", ",AliHLTPHOSRcuHistogramProducer DoDeinit");
+  return 0;
+}
+
+
+const char* 
+AliHLTPHOSRcuHistogramProducerComponent::GetComponentID()
+{
+  return "RcuHistogramProducer";
+}
+
+
+void
+ AliHLTPHOSRcuHistogramProducerComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
+{
+  const AliHLTComponentDataType* pType=inputDataTypes;
+  while (pType->fID!=0) 
+    {
+      list.push_back(*pType);
+      pType++;
+    }
+}
+
+
+AliHLTComponentDataType 
+AliHLTPHOSRcuHistogramProducerComponent::GetOutputDataType()
+{
+  return AliHLTPHOSDefinitions::gkCellEnergyDataType;
+}
+
+
+void
+AliHLTPHOSRcuHistogramProducerComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier )
+{
+  constBase = 30;
+  inputMultiplier = 1;
+}
+
+
+int  AliHLTPHOSRcuHistogramProducerComponent::DoEvent( const AliHLTComponentEventData& evtData, const AliHLTComponentBlockData* blocks, 
+                                             AliHLTComponentTriggerData& trigData, AliHLTUInt8_t* outputPtr, 
+                                             AliHLTUInt32_t& size, vector<AliHLTComponentBlockData>& outputBlocks )
+{
+  unsigned long ndx       = 0;
+  UInt_t offset           = 0; 
+  UInt_t mysize           = 0;
+  UInt_t tSize            = 0;
+  const AliHLTComponentBlockData* iter = NULL;   
+  AliHLTPHOSRcuCellEnergyDataStruct *cellDataPtr;
+  AliHLTUInt8_t* outBPtr;
+  // outBPtr = outputPtr;
+  // fOutPtr =  (AliHLTPHOSRcuCellAccumulatedEnergyDataStruct*)outBPtr;
+  int tmpCnt;
+
+  for( ndx = 0; ndx < evtData.fBlockCnt; ndx++ )
+    {
+      iter = blocks+ndx;
+      cellDataPtr = (AliHLTPHOSRcuCellEnergyDataStruct*)( iter->fPtr);
+      tmpCnt =  cellDataPtr->fCnt;
+
+      for(int i= 0; i <= tmpCnt; i ++)
+       {
+         fRcuHistoProducerPtr->FillEnergy(cellDataPtr->fValidData[i].fX,
+                                          cellDataPtr->fValidData[i].fZ, 
+                                          cellDataPtr->fValidData[i].fGain, 
+                                          cellDataPtr->fValidData[i].fEnergy);
+       }
+    }
+  
+  outBPtr = outputPtr;
+  fOutPtr =  (AliHLTPHOSRcuCellAccumulatedEnergyDataStruct*)outBPtr;
+  const AliHLTPHOSRcuCellAccumulatedEnergyDataStruct  &innPtr = fRcuHistoProducerPtr->GetCellAccumulatedEnergies();
+
+  fOutPtr->fModuleID = fModuleID;
+  fOutPtr->fRcuX     = fRcuX;
+  fOutPtr->fRcuZ     = fRcuZ;
+
+
+  for(int x=0; x < N_XCOLUMNS_RCU; x ++)
+    {
+      for(int z=0; z < N_XCOLUMNS_RCU; z ++)
+       {
+         for(int gain =0;  gain < N_GAINS; gain ++)
+           {
+             fOutPtr->fAccumulatedEnergies[x][z][gain] = innPtr.fAccumulatedEnergies[x][z][gain];
+             fOutPtr->fHits[x][z][gain] = innPtr.fHits[x][z][gain];
+           }
+       }
+    }
+
+
+  //pushing data to shared output memory
+  mysize += sizeof(AliHLTPHOSRcuCellAccumulatedEnergyDataStruct);
+  AliHLTComponentBlockData bd;
+  FillBlockData( bd );
+  bd.fOffset = offset;
+  bd.fSize = mysize;
+  bd.fDataType = AliHLTPHOSDefinitions::gkCellAccumulatedEnergyDataType;
+  bd.fSpecification = 0xFFFFFFFF;
+  outputBlocks.push_back( bd );
+  tSize += mysize;
+  outBPtr += mysize;
+
+  if( tSize > size )
+    {
+      Logging( kHLTLogFatal, "HLT::AliHLTRcuHistogramProducerComponent::DoEvent", "Too much data",
+              "Data written over allowed buffer. Amount written: %lu, allowed amount: %lu."
+              , tSize, size );
+      return EMSGSIZE;
+    }
+
+  fEventCount++; 
+  return 0;
+}//end DoEvent
+
+
+int
+AliHLTPHOSRcuHistogramProducerComponent::DoInit( int argc, const char** argv )
+{
+  int iResult=0;
+  TString argument="";
+  //  fRcuHistoProducerPtr = new AliHLTPHOSRcuHistogramProducer();
+  AliHLTUInt8_t tmpRcuX;
+  AliHLTUInt8_t tmpRcuZ; 
+  AliHLTUInt8_t tmpModuleID;
+  Bool_t isSetEquippmentID = kFALSE;
+
+  for(int i=0; i<argc && iResult>=0; i++) 
+    {
+      argument=argv[i];
+      
+      if (argument.IsNull()) 
+       {
+         continue;
+       }
+      if (argument.CompareTo("-equipmentID") == 0) 
+       {
+         if(i+1 <= argc)
+           {
+             fEquippmentID = atoi(argv[i+1]);
+             isSetEquippmentID = kTRUE;
+           }
+         else
+           {
+              iResult= -1;
+              Logging( kHLTLogFatal, "HLT::AliHLTPHOSRcuHistogramProducerComponent::DoInt( int argc, const char** argv )", "Missing argument",
+                       "The argument -equippmentID expects a number");
+              return  iResult;   
+           }
+       }
+
+      int rcuIndex =  (fEquippmentID - 1792)%N_RCUS_PER_MODULE;
+      //     fModuleID = (fEquippmentID  -1792 -rcuIndex)/N_RCUS_PER_MODULE;
+      tmpModuleID = ((fEquippmentID -1792 -rcuIndex)/N_RCUS_PER_MODULE);
+      SetModuleID(tmpModuleID);
+
+      if(rcuIndex == 0)
+       {
+         tmpRcuX = 0; 
+         tmpRcuZ = 0;
+       }
+      
+      if(rcuIndex == 1)
+       {
+        tmpRcuX = 0; 
+        tmpRcuZ = 1;
+       }
+      
+      if(rcuIndex == 2)
+       {
+         tmpRcuX = 1; 
+         tmpRcuZ = 0;
+       }
+      
+      if(rcuIndex == 3)
+       {
+         tmpRcuX = 1; 
+         tmpRcuZ = 1;
+       }
+
+      SetRcuX(tmpRcuX);
+      SetRcuZ(tmpRcuZ); 
+      cout <<"********InitInfo************"<< endl;
+      cout <<"AliHLTPHOSRcuHistogramProducerComponent::SetCoordinate"<< endl;
+      cout <<"Equpippment ID =\t"<< fEquippmentID <<endl;
+      cout <<"Module ID =\t"<<  (int) tmpModuleID<<endl;
+      cout <<"RCUX =\t\t" << (int)tmpRcuX << endl;
+      cout <<"RCUZ =\t\t" << (int)tmpRcuZ << endl;
+     }
+
+  if(isSetEquippmentID == kFALSE)
+    {
+       Logging( kHLTLogFatal, "HLT::AliHLTPHOSRcuHistogramProducerComponent::DoInt( int argc, const char** argv )", "Missing argument",
+                       "The argument equippmentID is not set: set it with a component argumet like this: -equippmentID  <number>");
+       iResult = -2; 
+    }
+
+  
+  fRcuHistoProducerPtr = new AliHLTPHOSRcuHistogramProducer( tmpModuleID, tmpRcuX, tmpRcuZ);
+
+
+  return  iResult;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducerComponent::SetRcuX(AliHLTUInt8_t X)
+{
+  fRcuX = X;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducerComponent::SetRcuZ(AliHLTUInt8_t Z)
+{
+  fRcuZ = Z;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducerComponent::SetModuleID(AliHLTUInt8_t moduleID)
+{
+  fModuleID = moduleID;
+}
+
+
+void 
+AliHLTPHOSRcuHistogramProducerComponent::SetEquippmentId(int id)
+{
+  fEquippmentID = id;
+  fRcuHistoProducerPtr->SetEquippmentId(id);
+}
+
+
+int 
+AliHLTPHOSRcuHistogramProducerComponent::GetEquippmentId()
+{
+  return  fEquippmentID;
+}
+
+
+AliHLTComponent*
+AliHLTPHOSRcuHistogramProducerComponent::Spawn()
+{
+  return new AliHLTPHOSRcuHistogramProducerComponent;
+}
+
+
diff --git a/HLT/PHOS/AliHLTPHOSRcuHistogramProducerComponent.h b/HLT/PHOS/AliHLTPHOSRcuHistogramProducerComponent.h
new file mode 100644 (file)
index 0000000..3efba03
--- /dev/null
@@ -0,0 +1,61 @@
+#ifndef ALIHLTPHOSRCUHISTOGRAMPRODUCERCOMPONENT_H
+#define ALIHLTPHOSRCUHISTOGRAMPRODUCERCOMPONENT_H 
+
+/* Copyright(c) 2006, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice  */ 
+
+#include "AliHLTProcessor.h"
+#include "AliHLTPHOSDefinitions.h"
+#include "AliHLTPHOSCommonDefs.h"
+//#include "AliHLTPHOSModuleCellAccumulatedEnergyDataStruct.h"
+#include "AliHLTPHOSRcuCellAccumulatedEnergyDataStruct.h"
+#include "Rtypes.h"
+
+class AliHLTPHOSRcuHistogramProducer;
+
+class AliHLTPHOSRcuHistogramProducerComponent:public AliHLTProcessor
+{
+ public:
+  AliHLTPHOSRcuHistogramProducerComponent();
+  virtual ~AliHLTPHOSRcuHistogramProducerComponent();
+  virtual int DoInit( int argc, const char** argv );
+  virtual int Deinit();
+  virtual int DoDeinit();
+  virtual int DoEvent(const AliHLTComponentEventData&, const AliHLTComponentBlockData*, AliHLTComponentTriggerData&, AliHLTUInt8_t*,
+                     AliHLTUInt32_t&, std::vector<AliHLTComponentBlockData, std::allocator<AliHLTComponentBlockData> >&);
+  virtual void GetInputDataTypes(std::vector<AliHLTComponentDataType, std::allocator<AliHLTComponentDataType> >&);
+  virtual AliHLTComponentDataType GetOutputDataType();
+  virtual void GetOutputDataSize(unsigned long& constBase, double& inputMultiplier);
+  virtual AliHLTComponent* Spawn();
+
+  virtual const char* GetComponentID();
+  int GetEquippmentId();
+  void SetRcuX(AliHLTUInt8_t X);
+  void SetRcuZ(AliHLTUInt8_t Z);
+  void SetModuleID(AliHLTUInt8_t moduleID);
+  void SetEquippmentId(int id);
+
+
+ private:
+  AliHLTPHOSRcuHistogramProducerComponent(const AliHLTPHOSRcuHistogramProducerComponent & );
+  AliHLTPHOSRcuHistogramProducerComponent & operator = (const AliHLTPHOSRcuHistogramProducerComponent &)
+   {
+      return *this;
+   };
+  int fEventCount;
+  AliHLTUInt8_t fRcuX;
+  AliHLTUInt8_t fRcuZ; 
+  AliHLTUInt8_t fModuleID;
+  AliHLTPHOSRcuHistogramProducer* fRcuHistoProducerPtr;
+  AliHLTPHOSRcuCellAccumulatedEnergyDataStruct*  fOutPtr;
+
+  //  AliHLTPHOSRcuCellAccumulatedEnergyDataStruct*  fInnPtr;
+  //  AliHLTPHOSRcuCellAccumulatedEnergyDataStruct&  fInnPtr;
+
+  static const AliHLTComponentDataType inputDataTypes[];
+  static const AliHLTComponentDataType outputDataType;
+
+  int fEquippmentID; 
+};
+
+#endif