]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
deactivating generation of the default detector/ddl list in the FXS header
authorrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 4 Nov 2009 01:45:07 +0000 (01:45 +0000)
committerrichterm <richterm@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 4 Nov 2009 01:45:07 +0000 (01:45 +0000)
the calibration component has to generate the list explicitely.
only ddls for one detector can be set

HLT/BASE/AliHLTCalibrationProcessor.cxx
HLT/BASE/AliHLTCalibrationProcessor.h
HLT/PHOS/AliHLTPHOSRcuDAComponent.cxx
HLT/TPCLib/calibration/AliHLTTPCCalibTimeComponent.cxx
HLT/TPCLib/calibration/AliHLTTPCCalibTimeGainComponent.cxx
HLT/TPCLib/calibration/AliHLTTPCCalibrationComponent.cxx
HLT/TPCLib/offline/AliHLTTPCOfflineCalibrationComponent.cxx
HLT/comp/AliHLTCOMPHuffmanAltroCalibComponent.cxx

index ea9fcc25ce7d15b19cf5201e5144cc49cc7bb959..a5caab2fd2a1f9672e5dccc51f626fd46c94f7ad 100644 (file)
@@ -210,7 +210,8 @@ Int_t AliHLTCalibrationProcessor::DoEvent( const AliHLTComponentEventData& evtDa
 
   // ** if event Type is not SOR or EOR -> fill DDLNumber list and process data
   if ( ! blkEOR  && !blkSOR ) {
-    
+
+#if 0    
     // ** Set DDLNumber List
     if ( trigData.fData != NULL) {
       AliHLTEventTriggerData* trg = ( AliHLTEventTriggerData* ) trigData.fData;
@@ -240,6 +241,7 @@ Int_t AliHLTCalibrationProcessor::DoEvent( const AliHLTComponentEventData& evtDa
        } // if ( wordNdx > 0 ) {
       } // if ( trg != NULL) {
     } // if ( trigData.fData != NULL) {
+#endif
     
     // ** ProcessData
     iResult = ProcessCalibration( evtData, blocks, trigData, outputPtr, size, outputBlocks );
index 6a1fa60c4140e141c45194ac8f35da31e75b6b40..5987012dc3b42f9ee332b58f2880176942fc3ec5 100644 (file)
@@ -65,7 +65,7 @@ class AliHLTCalibrationProcessor : public AliHLTProcessor {
    *                   Will be filled automatically if not supplied by the component.
    * @return neg. error code if failed 
    */
-  Int_t PushToFXS(TObject* pObject, const char* pDetector, const char* pFileID, AliHLTEventDDL* pDDLList = NULL);
+  Int_t PushToFXS(TObject* pObject, const char* pDetector, const char* pFileID, AliHLTEventDDL* pDDLList);
 
   /**
    * Insert an object into the output. FXS header will be inserted before the root object.
@@ -79,7 +79,7 @@ class AliHLTCalibrationProcessor : public AliHLTProcessor {
    *                   Will be filled automatically if not supplied by the component.
    * @return neg. error code if failed 
    */
-   Int_t PushToFXS(void* pBuffer, int iSize, const char* pDetector, const char* pFileID, AliHLTEventDDL* pDDLList = NULL);
+   Int_t PushToFXS(void* pBuffer, int iSize, const char* pDetector, const char* pFileID, AliHLTEventDDL* pDDLList);
 
   /** Constants  */ 
   static const AliHLTUInt32_t fgkFXSProtocolHeaderSize;
index 56ab03c1300362ab06d78add525828dc883a3010..4d66208a36c2589bcd643d5cddbd81b1c14c6273 100644 (file)
@@ -24,6 +24,7 @@
 //#include "AliHLTPHOSRcuCellEnergyDataStruct.h"
 #include "TObjArray.h"
 #include "AliHLTPHOSUtilities.h"
+#include "AliHLTReadoutList.h"
 
 //#include <iostream>
 
@@ -217,7 +218,8 @@ AliHLTPHOSRcuDAComponent::ShipDataToFXS( const AliHLTComponentEventData& /*evtDa
   TFile *outFile =  new TFile(filename, "recreate");
   calibPtr->Write(); 
   outFile->Close();
-  PushToFXS( (TObject*)fPHOSDAPtr->GetHistoContainer(), "PHOS",  filename);
+  static AliHLTReadoutList rdList(AliHLTReadoutList::kPHOS);
+  PushToFXS( (TObject*)fPHOSDAPtr->GetHistoContainer(), "PHOS",  filename, rdList.Buffer());
   cout << "Finnished pushing data to HLT FXS" << endl;
   return 0;
 }  
index 4469ecb8c9900f5ceec59c5e1c4d7bec676e3df5..4f0c8dd8e96dfa63a8e092eedcb368d277052ea2 100644 (file)
@@ -33,6 +33,7 @@ using namespace std;
 
 #include "AliHLTTPCCalibTimeComponent.h"
 #include "AliHLTTPCDefinitions.h"
+#include "AliHLTReadoutList.h"
 
 #include "AliESDEvent.h"
 #include "AliESDtrack.h"
@@ -204,7 +205,8 @@ Int_t AliHLTTPCCalibTimeComponent::ShipDataToFXS( const AliHLTComponentEventData
 // see header file for class documentation
     
   if(fEnableAnalysis) fCalibTime->Analyze();  
-  PushToFXS( (TObject*)fCalibTime, "TPC", "Time" ) ;
+  static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
+  PushToFXS( (TObject*)fCalibTime, "TPC", "Time", rdList.Buffer() ) ;
   
   return 0;
 } 
index 82172d308fc4fb776c69d9b76a9b953d633163b8..a570bdfa2924f9edd8e78291c630ea3264a2f311 100644 (file)
@@ -33,6 +33,7 @@ using namespace std;
 
 #include "AliHLTTPCCalibTimeGainComponent.h"
 #include "AliHLTTPCDefinitions.h"
+#include "AliHLTReadoutList.h"
 
 #include "AliESDEvent.h"
 #include "AliESDtrack.h"
@@ -203,7 +204,8 @@ Int_t AliHLTTPCCalibTimeGainComponent::ShipDataToFXS( const AliHLTComponentEvent
   // see header file for class documentation
     
   if(fEnableAnalysis) fCalibTimeGain->Analyze();  
-  PushToFXS( (TObject*) fCalibTimeGain, "TPC", "TimeGain" ) ;
+  static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
+  PushToFXS( (TObject*) fCalibTimeGain, "TPC", "TimeGain", rdList.Buffer() ) ;
   
   return 0;
 } 
index 42e10c62b575e382934e6849acb9246cb4ced52c..e03072c25992a6a01bbd4f2f392a3db4a82f5efd 100644 (file)
@@ -34,6 +34,7 @@ using namespace std;
 #include "AliHLTTPCCalibrationComponent.h"
 #include "AliHLTTPCDefinitions.h"
 #include "AliHLTTPCAnalysisTaskcalib.h"
+#include "AliHLTReadoutList.h"
 
 #include "AliAnalysisManager.h"
 #include "AliESDEvent.h"
@@ -226,7 +227,8 @@ Int_t AliHLTTPCCalibrationComponent::ShipDataToFXS( const AliHLTComponentEventDa
   // see header file for class documentation
     
   if(fEnableAnalysis) fCalibTask->Analyze();  
-  PushToFXS( (TObject*) fCalibTask, "TPC", "CALIB" ) ;
+  static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
+  PushToFXS( (TObject*) fCalibTask, "TPC", "CALIB", rdList.Buffer() ) ;
   
   return 0;
 } 
index f993f529a7f05508deffa6cc8e3054f04c1f62cc..ef317324d077dfdd0fae720a26e065de16c4e78a 100644 (file)
@@ -38,6 +38,7 @@
 #include "AliTPCcalibTracksCuts.h"
 #include "AliTPCClusterParam.h"
 #include "AliHLTTPCDefinitions.h"
+#include "AliHLTReadoutList.h"
 
 /** ROOT macro for the implementation of ROOT specific class methods */
 ClassImp(AliHLTTPCOfflineCalibrationComponent)
@@ -254,9 +255,10 @@ Int_t AliHLTTPCOfflineCalibrationComponent::ShipDataToFXS( const AliHLTComponent
       fTPCcalibTracksGain->Analyze();
       fTPCcalibTracks->Analyze();
    }
-   PushToFXS((TObject*)fTPCcalibAlign, "TPC", "TPCcalibAlign") ;
-   PushToFXS((TObject*)fTPCcalibTracksGain, "TPC", "TPCcalibTracksGain") ;
-   PushToFXS((TObject*)fTPCcalibTracks, "TPC", "TPCcalibTracks") ;
+   static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
+   PushToFXS((TObject*)fTPCcalibAlign, "TPC", "TPCcalibAlign", rdList.Buffer()) ;
+   PushToFXS((TObject*)fTPCcalibTracksGain, "TPC", "TPCcalibTracksGain", rdList.Buffer()) ;
+   PushToFXS((TObject*)fTPCcalibTracks, "TPC", "TPCcalibTracks", rdList.Buffer()) ;
 
 return 0;
 }
index 8b90937b6a599e799abd3989b4dee9f79f48f3bb..bbed3c65d0476a7ec1496827b20781627ebd57c0 100644 (file)
@@ -27,6 +27,7 @@ using namespace std;
 #include "AliHLTCOMPHuffmanAltro.h"
 #include "AliHLTCompDefinitions.h"
 #include "AliHLTStdIncludes.h"
+#include "AliHLTReadoutList.h"
 #include "TFile.h" // necessary for HuffmanData writing
 
 ClassImp(AliHLTCOMPHuffmanAltroCalibComponent)
@@ -365,7 +366,8 @@ Int_t AliHLTCOMPHuffmanAltroCalibComponent::ShipDataToFXS( const AliHLTComponent
   Int_t dataspec = (Int_t) fSpecification;
 
   fHuffmanData->SetOCDBSpecifications(fOrigin, dataspec);
-  PushToFXS( (TObject*) fHuffmanData, "TPC", "HuffmanData" ) ;
+  static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
+  PushToFXS( (TObject*) fHuffmanData, "TPC", "HuffmanData", rdList.Buffer() ) ;
   
   return 0;
 } // Int_t AliHLTCOMPHuffmanAltroCalibComponent::ShipDataToFXS( const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& trigData ) {