2 /**************************************************************************
3 * This file is property of and copyright by the ALICE HLT Project *
4 * ALICE Experiment at CERN, All rights reserved. *
6 * Primary Authors: Kalliopi Kanaki <Kalliopi.Kanaki@ift.uib.no> *
7 * for The ALICE HLT Project. *
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 **************************************************************************/
18 /** @file AliHLTTPCCalibrationComponent.cxx
19 @author Kalliopi Kanaki
24 // see header file for class documentation
26 // refer to README to build package
28 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
34 #include "AliHLTTPCCalibrationComponent.h"
35 #include "AliHLTTPCDefinitions.h"
36 #include "AliHLTTPCAnalysisTaskcalib.h"
37 #include "AliHLTReadoutList.h"
39 #include "AliAnalysisManager.h"
40 #include "AliESDEvent.h"
41 //#include "AliESDInputHandler.h"
43 #include "AliTPCcalibTime.h"
44 #include "AliTPCcalibTimeGain.h"
45 #include "AliTPCseed.h"
48 #include "TObjArray.h"
49 #include "TTimeStamp.h"
54 ClassImp(AliHLTTPCCalibrationComponent) // ROOT macro for the implementation of ROOT specific class methods
56 AliHLTTPCCalibrationComponent::AliHLTTPCCalibrationComponent()
68 fEnableAnalysis(kTRUE)
70 // see header file for class documentation
72 // refer to README to build package
74 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
77 AliHLTTPCCalibrationComponent::~AliHLTTPCCalibrationComponent() {
78 // see header file for class documentation
81 const char* AliHLTTPCCalibrationComponent::GetComponentID() {
82 // see header file for class documentation
84 return "TPCCalibration";
87 void AliHLTTPCCalibrationComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list) {
88 // see header file for class documentation
91 list.push_back( kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC ); // output of TPCCalibSeedMaker
92 list.push_back( kAliHLTDataTypeESDObject|kAliHLTDataOriginOut ); // output of global esd converter
96 AliHLTComponentDataType AliHLTTPCCalibrationComponent::GetOutputDataType() {
97 // see header file for class documentation
99 return AliHLTTPCDefinitions::fgkCalibCEDataType;
102 void AliHLTTPCCalibrationComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier ) {
103 // see header file for class documentation
106 inputMultiplier = (2.0); // to be estimated
109 AliHLTComponent* AliHLTTPCCalibrationComponent::Spawn() {
110 // see header file for class documentation
112 return new AliHLTTPCCalibrationComponent();
116 Int_t AliHLTTPCCalibrationComponent::ScanArgument( Int_t argc, const char** argv ) {
117 // see header file for class documentation
120 TString argument = "";
121 TString parameter = "";
123 if(!argc) return -EINVAL;
125 argument = argv[iResult];
127 if(argument.IsNull()) return -EINVAL;
129 if( argument.CompareTo("-enable-analysis") == 0 ){
130 HLTInfo( "Analysis before shipping data to FXS enabled." );
131 fEnableAnalysis = kTRUE;
140 Int_t AliHLTTPCCalibrationComponent::InitCalibration() {
141 // see header file for class documentation
143 if(fCalibTask) return EINPROGRESS;
144 fCalibTask = new AliHLTTPCAnalysisTaskcalib("TPC Calibration Task");
146 if(fCalibTime) return EINPROGRESS;
147 //fCalibTime = new AliTPCcalibTime();
148 fCalibTime = new AliTPCcalibTime("calibTime","time dependent Vdrift calibration",-2, 2, 1);
150 fCalibTime->SetDebugLevel(20);
151 fCalibTime->SetStreamLevel(10);
152 fCalibTime->SetTriggerMask(-1,-1,kFALSE); //accept everything
154 if(fCalibTimeGain) return EINPROGRESS;
155 fCalibTimeGain = new AliTPCcalibTimeGain("calibTimeGain","time dependent gain calibration",-2, 2, 1);
157 fCalibTask->AddJob(fCalibTime);
158 fCalibTask->AddJob(fCalibTimeGain);
159 fCalibTask->GetJobs();
164 Int_t AliHLTTPCCalibrationComponent::DeinitCalibration() {
165 // see header file for class documentation
167 if(fCalibTask) delete fCalibTask; fCalibTask = NULL;
168 if(fCalibTime) delete fCalibTime; fCalibTime = NULL;
169 if(fCalibTimeGain) delete fCalibTimeGain; fCalibTimeGain = NULL;
174 Int_t AliHLTTPCCalibrationComponent::ProcessCalibration( const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/ ){
175 // see header file for class documentation
177 if(GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR )) return 0;
179 TObject *iter = NULL;
181 //--------------- output over TObjArray of AliTPCseed objects (output of TPCSeedMaker) -------------------//
183 for(iter = (TObject*)GetFirstInputObject(kAliHLTDataTypeTObjArray|kAliHLTDataOriginTPC); iter != NULL; iter = (TObject*)GetNextInputObject()){
185 if(GetDataType(iter) != (kAliHLTDataTypeTObjArray | kAliHLTDataOriginTPC)) continue;
186 fSeedArray = dynamic_cast<TObjArray*>(iter);
189 //----------- loop over output of global esd converter ----------------//
191 for(iter = (TObject*)GetFirstInputObject(kAliHLTDataTypeESDObject | kAliHLTDataOriginOut); iter != NULL; iter = (TObject*)GetNextInputObject()){
193 if(GetDataType(iter) != (kAliHLTDataTypeESDObject | kAliHLTDataOriginOut)) continue;
195 AliHLTUInt8_t slice = AliHLTTPCDefinitions::GetMinSliceNr(GetSpecification(iter));
196 AliHLTUInt8_t partition = AliHLTTPCDefinitions::GetMinPatchNr(GetSpecification(iter));
198 if( partition < fMinPartition ) fMinPartition = partition;
199 if( partition > fMaxPartition ) fMaxPartition = partition;
200 if( slice < fMinSlice ) fMinSlice = slice;
201 if( slice > fMaxSlice ) fMaxSlice = slice;
203 fESDEvent = dynamic_cast<AliESDEvent*>(iter);
204 fESDEvent->CreateStdContent();
206 HLTDebug("# Seeds: %i\n", fSeedArray->GetEntriesFast());
208 for(Int_t i=0; i<fSeedArray->GetEntriesFast(); i++){
210 AliTPCseed *seed = (AliTPCseed*)fSeedArray->UncheckedAt(i);
212 AliESDtrack *esd = fESDEvent->GetTrack(i);
213 AliTPCseed *seedCopy = new AliTPCseed(*seed, kTRUE);
214 esd->AddCalibObject(seedCopy);
218 fCalibTask->Process(fESDEvent);
220 fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification( fMinSlice, fMaxSlice, fMinPartition, fMaxPartition );
221 PushBack( (TObject*) fCalibTask, AliHLTTPCDefinitions::fgkCalibCEDataType | kAliHLTDataOriginOut, fSpecification);
226 Int_t AliHLTTPCCalibrationComponent::ShipDataToFXS( const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/ ){
227 // see header file for class documentation
229 if(fEnableAnalysis) fCalibTask->Analyze();
230 static AliHLTReadoutList rdList(AliHLTReadoutList::kTPC);
231 PushToFXS( (TObject*) fCalibTask, "TPC", "CALIB", &rdList ) ;