]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCCalibCEComponent.cxx
added new helper components to libAliHLTUtil (EsdCollector and AliHLTOUTPublisher...
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCCalibCEComponent.cxx
1 // $Id$
2 /**************************************************************************
3  * This file is property of and copyright by the ALICE HLT Project        * 
4  * ALICE Experiment at CERN, All rights reserved.                         *
5  *                                                                        *
6  * Primary Authors: Jochen Thaeder <thaeder@kip.uni-heidelberg.de>        *
7  *                  for The ALICE HLT Project.                            *
8  *                                                                        *
9  * Permission to use, copy, modify and distribute this software and its   *
10  * documentation strictly for non-commercial purposes is hereby granted   *
11  * without fee, provided that the above copyright notice appears in all   *
12  * copies and that both the copyright notice and this permission notice   *
13  * appear in the supporting documentation. The authors make no claims     *
14  * about the suitability of this software for any purpose. It is          *
15  * provided "as is" without express or implied warranty.                  *
16  **************************************************************************/
17
18 /** @file   AliHLTTPCCalibCEComponent.cxx
19     @author Jochen Thaeder
20     @date   
21     @brief  A pedestal calibration component for the TPC.
22 */
23
24 // see header file for class documentation
25 // or
26 // refer to README to build package
27 // or
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt   
29
30 #if __GNUC__>= 3
31 using namespace std;
32 #endif
33
34 //#include "AliHLTTPCLogging.h"
35 #include "AliHLTTPCDefinitions.h"
36 //#include "AliHLTTPCTransform.h"
37
38 #include "AliHLTTPCCalibCEComponent.h"
39
40 //#include "AliRawDataHeader.h"
41 #include "AliRawReaderMemory.h"
42 #include "AliTPCRawStream.h"
43
44 #include "AliTPCCalibCE.h"
45
46 #include <cstdlib>
47 #include <cerrno>
48 #include "TString.h"
49
50 /** ROOT macro for the implementation of ROOT specific class methods */
51 ClassImp(AliHLTTPCCalibCEComponent)
52
53 AliHLTTPCCalibCEComponent::AliHLTTPCCalibCEComponent()
54   :
55   fRawReader(NULL),
56   fRawStream(NULL),
57   fCalibCE(NULL),
58   fMinPatch(5),
59   fMaxPatch(0),
60   fSpecification(0) ,
61   fEnableAnalysis(kFALSE) {
62   // see header file for class documentation
63   // or
64   // refer to README to build package
65   // or
66   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
67 }
68
69 AliHLTTPCCalibCEComponent::~AliHLTTPCCalibCEComponent() {
70   // see header file for class documentation
71 }
72
73 // Public functions to implement AliHLTComponent's interface.
74 // These functions are required for the registration process
75
76 const char* AliHLTTPCCalibCEComponent::GetComponentID() {
77   // see header file for class documentation
78
79   return "TPCCalibCE";
80 }
81
82 void AliHLTTPCCalibCEComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
83   // see header file for class documentation
84
85   list.clear(); 
86   list.push_back( AliHLTTPCDefinitions::fgkDDLPackedRawDataType );
87 }
88
89 AliHLTComponentDataType AliHLTTPCCalibCEComponent::GetOutputDataType() {
90   // see header file for class documentation
91
92   return AliHLTTPCDefinitions::fgkCalibCEDataType;
93 }
94
95 void AliHLTTPCCalibCEComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
96   // see header file for class documentation
97
98   // XXX TODO: Find more realistic values.  
99   constBase = 0;
100   inputMultiplier = (2.0);
101 }
102
103 AliHLTComponent* AliHLTTPCCalibCEComponent::Spawn() {
104   // see header file for class documentation
105
106   return new AliHLTTPCCalibCEComponent();
107 }  
108
109
110 Int_t AliHLTTPCCalibCEComponent::ScanArgument( Int_t argc, const char** argv ) {
111   // see header file for class documentation
112
113   Int_t iResult = 0;
114   TString argument = "";
115   TString parameter = "";
116
117   if ( !argc ) 
118     return -EINVAL;
119
120   argument = argv[iResult];
121   
122   if ( argument.IsNull() ) 
123     return -EINVAL;
124
125   // -rcuformat
126   if ( argument.CompareTo("-rcuformat") == 0 ) {
127
128     if ( ++iResult >= argc  ) {
129       iResult = -EPROTO;
130     }
131     else {
132       parameter = argv[1];
133       if ( parameter.CompareTo("old") == 0 ) {
134         HLTInfo( "RCU Format is set to old." );
135       }
136       else if ( parameter.CompareTo("new") == 0 ) {
137         HLTInfo( "RCU Format is set to new." );
138       }
139       else {
140         HLTError( "Cannot convert rcu format specifier '%s'.", argv[1] );
141         iResult = -EPROTO;
142       }
143     } 
144   }
145   else if ( argument.CompareTo("-enableanalysis") == 0 ) {
146     HLTInfo( "Analysis before shipping data to FXS enabled." );
147     fEnableAnalysis = kTRUE;
148   }
149   else {
150     iResult = -EINVAL;
151   }
152       
153   return iResult;
154 }
155
156 Int_t AliHLTTPCCalibCEComponent::InitCalibration() {
157   // see header file for class documentation
158     
159   // ** Create pedestal calibration
160   if ( fCalibCE )
161     return EINPROGRESS;
162   
163   fCalibCE = new AliTPCCalibCE();
164
165   // **  Create AliRoot Memory Reader
166   if (fRawReader)
167     return EINPROGRESS;
168
169   fRawReader = new AliRawReaderMemory();
170
171   return 0;
172 }
173
174 Int_t AliHLTTPCCalibCEComponent::DeinitCalibration() {
175   // see header file for class documentation
176
177   if ( fRawReader )
178     delete fRawReader;
179   fRawReader = NULL;
180
181   if ( fCalibCE )
182     delete fCalibCE;
183   fCalibCE = NULL;
184
185   return 0;
186 }
187
188 /*
189  * --- setter for rcuformat need in AliTPCCalibCE class
190  */
191 Int_t AliHLTTPCCalibCEComponent::ProcessCalibration( const AliHLTComponentEventData& /*evtData*/, 
192                                                            AliHLTComponentTriggerData& /*trigData*/ ) {
193   // see header file for class documentation
194   
195   const AliHLTComponentBlockData* iter = NULL;
196
197   AliHLTUInt8_t slice=0, patch=0;
198   Int_t ddlId = 0;
199     
200   // ** Loop over all input blocks and specify which data format should be read - only select Raw Data
201   iter = GetFirstInputBlock( kAliHLTDataTypeDDLRaw | kAliHLTDataOriginTPC);
202   
203   while ( iter != NULL ) {
204     
205     // ** Print Debug output which data format was received
206     //    char tmp1[14], tmp2[14];
207     //    DataType2Text( iter->fDataType, tmp1 );
208     //    DataType2Text( AliHLTTPCDefinitions::fgkDDLPackedRawDataType, tmp2 );
209     //    HLTDebug ( "Event received - Event 0x%08LX (%Lu) received datatype: %s - required datatype: %s", 
210     //    evtData.fEventID, evtData.fEventID, tmp1, tmp2 );
211
212     // ** Get DDL ID in order to tell the memory reader which slice/patch to use
213     slice = AliHLTTPCDefinitions::GetMinSliceNr( *iter );
214     patch = AliHLTTPCDefinitions::GetMinPatchNr( *iter );
215
216     if (patch < 2) ddlId = 768 + (2 * slice) + patch;
217     else ddlId = 838 + (4*slice) + patch;
218
219     HLTDebug ( "Input Raw Data - Slice/Patch: %d/%d - EquipmentID : %d.", slice, patch, ddlId );
220
221     // ** Get min and max patch, used for output specification
222     if ( patch < fMinPatch ) fMinPatch =  patch;
223     if ( patch > fMaxPatch ) fMaxPatch =  patch;
224
225     // ** Init TPCRawStream
226     fRawReader->SetMemory( reinterpret_cast<UChar_t*>( iter->fPtr ), iter->fSize );
227     fRawReader->SetEquipmentID(ddlId);
228
229     fRawStream = new AliTPCRawStream( fRawReader );
230
231     // ** Process actual Pedestal Calibration - Fill histograms
232     fCalibCE->ProcessEvent( fRawStream );
233   
234     // ** Delete TPCRawStream
235     if ( fRawStream )
236       delete fRawStream;
237     fRawStream = NULL;    
238
239     // ** Get next input block, with the same specification as defined in GetFirstInputBlock()
240     iter = GetNextInputBlock();
241
242   } //  while ( iter != NULL ) {
243     
244   // ** Get output specification
245   fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification( slice, slice, fMinPatch, fMaxPatch );
246
247   // ** PushBack data to shared memory ... 
248   PushBack( (TObject*) fCalibCE, AliHLTTPCDefinitions::fgkCalibCEDataType, fSpecification);
249   
250   return 0;
251 } // Int_t AliHLTTPCCalibCEComponent::ProcessCalibration( const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& trigData ) {
252
253
254 Int_t AliHLTTPCCalibCEComponent::ShipDataToFXS( const AliHLTComponentEventData& /*evtData*/, 
255                                                       AliHLTComponentTriggerData& /*trigData*/ ) {
256   // see header file for class documentation
257     
258   if ( fEnableAnalysis )
259     fCalibCE->Analyse();
260   
261   // ** PushBack data to FXS ...
262   PushToFXS( (TObject*) fCalibCE, "TPC", "CE" ) ;
263   
264   return 0;
265 } // Int_t AliHLTTPCCalibCEComponent::ShipDataToFXS( const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& trigData ) {
266
267