]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/calibration/AliHLTTPCCalibrationComponent.cxx
correcting memory leaks on error conditions
[u/mrichter/AliRoot.git] / HLT / TPCLib / calibration / AliHLTTPCCalibrationComponent.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: Kalliopi Kanaki <Kalliopi.Kanaki@ift.uib.no>          *
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   AliHLTTPCCalibrationComponent.cxx
19     @author Kalliopi Kanaki
20     @date   
21     @brief  
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 "AliHLTTPCCalibrationComponent.h"
35 #include "AliHLTTPCDefinitions.h"
36 #include "AliHLTTPCAnalysisTaskcalib.h"
37 #include "AliHLTReadoutList.h"
38
39 #include "AliAnalysisManager.h"
40 #include "AliESDEvent.h"
41 //#include "AliESDInputHandler.h"
42
43 #include "AliTPCcalibTime.h"
44 #include "AliTPCcalibTimeGain.h"
45 #include "AliTPCseed.h"
46
47 #include "TString.h"
48 #include "TObjArray.h"
49 #include "TTimeStamp.h"
50
51 #include <cstdlib>
52 #include <cerrno>
53
54 ClassImp(AliHLTTPCCalibrationComponent) // ROOT macro for the implementation of ROOT specific class methods
55
56 AliHLTTPCCalibrationComponent::AliHLTTPCCalibrationComponent()
57   :
58   fCalibTask(NULL),
59   fCalibTime(NULL),
60   fCalibTimeGain(NULL),
61   fESDEvent(NULL),
62   fSeedArray(NULL),
63   fMinPartition(5),
64   fMaxPartition(0),
65   fMinSlice(35),
66   fMaxSlice(0),
67   fSpecification(0) ,
68   fEnableAnalysis(kTRUE)
69 {
70   // see header file for class documentation
71   // or
72   // refer to README to build package
73   // or
74   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
75 }
76
77 AliHLTTPCCalibrationComponent::~AliHLTTPCCalibrationComponent() {
78 // see header file for class documentation
79 }
80
81 const char* AliHLTTPCCalibrationComponent::GetComponentID() {
82 // see header file for class documentation
83
84   return "TPCCalibration";
85 }
86
87 void AliHLTTPCCalibrationComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
88 // see header file for class documentation
89
90   list.clear();     
91   list.push_back( kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC );  // output of TPCCalibSeedMaker
92   list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOut ); // output of global esd converter 
93
94 }
95
96 AliHLTComponentDataType AliHLTTPCCalibrationComponent::GetOutputDataType() {
97 // see header file for class documentation
98
99   return AliHLTTPCDefinitions::fgkCalibCEDataType;
100 }
101
102 void AliHLTTPCCalibrationComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
103 // see header file for class documentation
104
105   constBase = 0;
106   inputMultiplier = (2.0); // to be estimated
107 }
108
109 AliHLTComponent* AliHLTTPCCalibrationComponent::Spawn() {
110 // see header file for class documentation
111
112   return new AliHLTTPCCalibrationComponent();
113 }  
114
115
116 Int_t AliHLTTPCCalibrationComponent::ScanArgument( Int_t argc, const char** argv ) {
117 // see header file for class documentation
118
119   Int_t iResult = 0;
120   TString argument = "";
121   TString parameter = "";
122
123   if(!argc) return -EINVAL;
124
125   argument = argv[iResult];
126   
127   if(argument.IsNull()) return -EINVAL;
128
129   if( argument.CompareTo("-enable-analysis") == 0 ){
130     HLTInfo( "Analysis before shipping data to FXS enabled." );
131     fEnableAnalysis = kTRUE;
132   }
133   else {
134     iResult = -EINVAL;
135   }
136       
137   return iResult;
138 }
139
140 Int_t AliHLTTPCCalibrationComponent::InitCalibration() {
141 // see header file for class documentation
142   
143   if(fCalibTask) return EINPROGRESS;
144   fCalibTask = new AliHLTTPCAnalysisTaskcalib("TPC Calibration Task");
145   
146   if(fCalibTime) return EINPROGRESS;
147   //fCalibTime = new AliTPCcalibTime();
148   fCalibTime = new AliTPCcalibTime("calibTime","time dependent Vdrift calibration",-2, 2, 1);
149   
150   fCalibTime->SetDebugLevel(20);
151   fCalibTime->SetStreamLevel(10);
152   fCalibTime->SetTriggerMask(-1,-1,kFALSE); //accept everything 
153
154   if(fCalibTimeGain) return EINPROGRESS;
155   fCalibTimeGain = new AliTPCcalibTimeGain("calibTimeGain","time dependent gain calibration",-2, 2, 1);
156   
157   fCalibTask->AddJob(fCalibTime);
158   fCalibTask->AddJob(fCalibTimeGain);
159   fCalibTask->GetJobs();
160  
161   return 0;
162 }
163
164 Int_t AliHLTTPCCalibrationComponent::DeinitCalibration() {
165 // see header file for class documentation
166
167   if(fCalibTask)     delete fCalibTask;     fCalibTask     = NULL;
168   if(fCalibTime)     delete fCalibTime;     fCalibTime     = NULL;
169   if(fCalibTimeGain) delete fCalibTimeGain; fCalibTimeGain = NULL;
170
171   return 0;
172 }
173
174 Int_t AliHLTTPCCalibrationComponent::ProcessCalibration( const AliHLTComponentEventData& /*evtData*/,  AliHLTComponentTriggerData& /*trigData*/ ){
175 // see header file for class documentation
176   
177   if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
178
179   TObject *iter = NULL;
180
181   //--------------- output over TObjArray of AliTPCseed objects (output of TPCSeedMaker) -------------------//
182   
183   for(iter = (TObject*)GetFirstInputObject(kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC); iter != NULL; iter = (TObject*)GetNextInputObject()){  
184               
185       if(GetDataType(iter) != (kAliHLTDataTypeTObjArray | kAliHLTDataOriginTPC)) continue;
186       fSeedArray = dynamic_cast<TObjArray*>(iter); 
187   }
188   
189   //----------- loop over output of global esd converter ----------------//
190   
191   for(iter = (TObject*)GetFirstInputObject(kAliHLTDataTypeESDObject | kAliHLTDataOriginOut); iter != NULL; iter = (TObject*)GetNextInputObject()){   
192       
193       if(GetDataType(iter) != (kAliHLTDataTypeESDObject | kAliHLTDataOriginOut)) continue;
194   
195       AliHLTUInt8_t slice     = AliHLTTPCDefinitions::GetMinSliceNr(GetSpecification(iter)); 
196       AliHLTUInt8_t partition = AliHLTTPCDefinitions::GetMinPatchNr(GetSpecification(iter));
197       
198       if( partition < fMinPartition ) fMinPartition = partition;
199       if( partition > fMaxPartition ) fMaxPartition = partition;
200       if( slice < fMinSlice ) fMinSlice = slice;
201       if( slice > fMaxSlice ) fMaxSlice = slice;
202       
203       fESDEvent = dynamic_cast<AliESDEvent*>(iter);
204       fESDEvent->CreateStdContent();
205      
206       HLTDebug("# Seeds: %i\n", fSeedArray->GetEntriesFast());
207       
208       for(Int_t i=0; i<fSeedArray->GetEntriesFast(); i++){
209           
210           AliTPCseed *seed = (AliTPCseed*)fSeedArray->UncheckedAt(i);
211           if(!seed) continue;
212           AliESDtrack *esd = fESDEvent->GetTrack(i);    
213           AliTPCseed *seedCopy = new AliTPCseed(*seed, kTRUE); 
214           esd->AddCalibObject(seedCopy);
215       }      
216  }   
217   
218   fCalibTask->Process(fESDEvent);
219
220   fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification( fMinSlice, fMaxSlice, fMinPartition, fMaxPartition );
221   PushBack( (TObject*) fCalibTask, AliHLTTPCDefinitions::fgkCalibCEDataType | kAliHLTDataOriginOut, fSpecification);
222   
223   return 0;
224 }
225
226 Int_t AliHLTTPCCalibrationComponent::ShipDataToFXS( const AliHLTComponentEventData& /*evtData*/,  AliHLTComponentTriggerData& /*trigData*/ ){
227   // see header file for class documentation
228     
229   if(fEnableAnalysis) fCalibTask->Analyze();  
230   static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
231   PushToFXS( (TObject*) fCalibTask, "TPC", "CALIB", &rdList ) ;
232   
233   return 0;
234
235
236