2 //**************************************************************************
3 //* This file is property of and copyright by the *
4 //* ALICE Experiment at CERN, All rights reserved. *
6 //* Primary Authors: Matthias Richter <Matthias.Richter@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 AliHLTTPCDataCompressionComponent.cxx
19 /// @author Matthias Richter
21 /// @brief TPC component for data compression
24 #include "AliHLTTPCDataCompressionComponent.h"
25 #include "AliHLTTPCDefinitions.h"
26 #include "AliHLTTPCTrackGeometry.h"
27 #include "AliHLTTPCSpacePointContainer.h"
28 #include "AliHLTTPCHWCFSpacePointContainer.h"
29 #include "AliHLTGlobalBarrelTrack.h"
30 #include "AliHLTComponentBenchmark.h"
31 #include "AliHLTDataDeflaterSimple.h"
32 #include "AliHLTDataDeflaterHuffman.h"
33 #include "AliHLTTPCTransform.h"
34 #include "AliHLTTPCClusterMCData.h"
35 #include "AliHLTTPCClusterTransformation.h"
36 #include "AliHLTErrorGuard.h"
37 #include "AliRawDataHeader.h"
38 #include "AliCDBManager.h"
39 #include "AliCDBPath.h"
41 #include "AliCDBMetaData.h"
42 #include "AliCDBEntry.h"
47 ClassImp(AliHLTTPCDataCompressionComponent)
49 AliHLTTPCDataCompressionComponent::AliHLTTPCDataCompressionComponent()
53 , fVerificationMode(0)
54 , fMaxDeltaPad(AliHLTTPCDefinitions::GetMaxClusterDeltaPad())
55 , fMaxDeltaTime(AliHLTTPCDefinitions::GetMaxClusterDeltaTime())
56 , fRawInputClusters(NULL)
57 , fInputClusters(NULL)
59 , fSpacePointGrid(NULL)
60 , fpDataDeflater(NULL)
61 , fHistoCompFactor(NULL)
62 , fHistoResidualPad(NULL)
63 , fHistoResidualTime(NULL)
64 , fHistoClustersOnTracks(NULL)
65 , fHistoClusterRatio(NULL)
66 , fHistoTrackClusterRatio(NULL)
68 , fTrainingTableOutput()
70 , fpWrittenAssociatedClusterIds(NULL)
71 , fDriftTimeFactorA(1.)
72 , fDriftTimeOffsetA(0.)
73 , fDriftTimeFactorC(1.)
74 , fDriftTimeOffsetC(0.)
79 AliHLTTPCDataCompressionComponent::~AliHLTTPCDataCompressionComponent()
82 if (fpWrittenAssociatedClusterIds) delete fpWrittenAssociatedClusterIds;
86 const char* AliHLTTPCDataCompressionComponent::GetComponentID()
88 /// inherited from AliHLTComponent: id of the component
89 return "TPCDataCompressor";
93 void AliHLTTPCDataCompressionComponent::GetInputDataTypes( AliHLTComponentDataTypeList& tgtList)
95 /// inherited from AliHLTComponent: list of data types in the vector reference
97 tgtList.push_back(AliHLTTPCDefinitions::HWClustersDataType());
98 tgtList.push_back(AliHLTTPCDefinitions::ClustersDataType());
99 tgtList.push_back(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
102 AliHLTComponentDataType AliHLTTPCDataCompressionComponent::GetOutputDataType()
104 /// inherited from AliHLTComponent: output data type of the component.
105 return kAliHLTMultipleDataType;
108 int AliHLTTPCDataCompressionComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
110 /// inherited from AliHLTComponent: multiple output data types of the component.
112 tgtList.push_back(AliHLTTPCDefinitions::RawClustersDataType());
113 tgtList.push_back(AliHLTTPCDefinitions::RemainingClustersCompressedDataType());
114 tgtList.push_back(AliHLTTPCDefinitions::RemainingClusterIdsDataType());
115 tgtList.push_back(AliHLTTPCDefinitions::ClusterTracksCompressedDataType());
116 tgtList.push_back(AliHLTTPCDefinitions::ClusterIdTracksDataType());
117 return tgtList.size();
120 void AliHLTTPCDataCompressionComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
122 /// inherited from AliHLTComponent: output data size estimator
124 inputMultiplier=1.; // there should not be more data than input
125 inputMultiplier+=.3; // slightly more data when using the old HWCF data with 20 Byte and raw clusters 22 Byte
126 if (fpWrittenAssociatedClusterIds) inputMultiplier+=.3; // space for optional cluster id array
129 AliHLTComponent* AliHLTTPCDataCompressionComponent::Spawn()
131 /// inherited from AliHLTComponent: spawn function.
132 return new AliHLTTPCDataCompressionComponent;
135 void AliHLTTPCDataCompressionComponent::GetOCDBObjectDescription(TMap* const targetMap)
137 /// Get a list of OCDB object needed for the particular component
138 if (!targetMap) return;
140 targetMap->Add(new TObjString("HLT/ConfigTPC/TPCDataCompressor"),
141 new TObjString("component arguments"));
142 if (fDeflaterMode==2) {
143 targetMap->Add(new TObjString("HLT/ConfigTPC/TPCDataCompressorHuffmanTables"),
144 new TObjString("huffman tables for deflater mode 2"));
148 int AliHLTTPCDataCompressionComponent::DoEvent( const AliHLTComponentEventData& /*evtData*/,
149 const AliHLTComponentBlockData* /*inputBlocks*/,
150 AliHLTComponentTriggerData& /*trigData*/,
151 AliHLTUInt8_t* outputPtr,
152 AliHLTUInt32_t& size,
153 AliHLTComponentBlockDataList& outputBlocks )
155 /// inherited from AliHLTProcessor: data processing
157 AliHLTUInt32_t capacity=size;
160 if (!IsDataEvent()) return 0;
162 if (!fRawInputClusters) {
166 if (GetBenchmarkInstance()) {
167 GetBenchmarkInstance()->StartNewEvent();
168 GetBenchmarkInstance()->Start(0);
172 // Loop over all input blocks in the event
173 bool bHaveMC=(GetFirstInputBlock(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo | kAliHLTDataOriginTPC))!=NULL;
174 if ((bHaveMC || fVerificationMode>0) && fpWrittenAssociatedClusterIds==NULL) {
175 fpWrittenAssociatedClusterIds=new vector<AliHLTUInt32_t>;
178 const AliHLTComponentBlockData* pDesc=NULL;
180 AliHLTUInt8_t minSlice=0xFF, maxSlice=0xFF, minPatch=0xFF, maxPatch=0xFF;
181 AliHLTUInt32_t inputRawClusterSize=0;
182 AliHLTUInt32_t outputDataSize=0;
184 int associatedClusters=0;
187 /// input track array
188 vector<AliHLTGlobalBarrelTrack> inputTrackArray;
190 if (GetBenchmarkInstance()) {
191 GetBenchmarkInstance()->Start(2);
194 // transformed clusters
195 if (fMode==10) { // FIXME: condition to be adjusted
196 for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType);
197 pDesc!=NULL; pDesc=GetNextInputBlock()) {
198 if (GetBenchmarkInstance()) {
199 GetBenchmarkInstance()->AddInput(pDesc->fSize);
201 AliHLTUInt8_t slice = 0;
202 AliHLTUInt8_t patch = 0;
203 slice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification );
204 patch = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification );
205 if ( minSlice==0xFF || slice<minSlice ) minSlice = slice;
206 if ( maxSlice==0xFF || slice>maxSlice ) maxSlice = slice;
207 if ( minPatch==0xFF || patch<minPatch ) minPatch = patch;
208 if ( maxPatch==0xFF || patch>maxPatch ) maxPatch = patch;
209 if (fInputClusters) {
210 fInputClusters->AddInputBlock(pDesc);
213 if (GetBenchmarkInstance()) {
214 GetBenchmarkInstance()->Stop(2);
215 GetBenchmarkInstance()->Start(3);
219 vector<int> trackindexmap; // stores index for every track id
222 if (fMode==2 || fMode==4) {
223 for (pDesc=GetFirstInputBlock(kAliHLTDataTypeTrack|kAliHLTDataOriginTPC);
224 pDesc!=NULL; pDesc=GetNextInputBlock()) {
225 if (GetBenchmarkInstance()) {
226 GetBenchmarkInstance()->AddInput(pDesc->fSize);
228 AliHLTUInt8_t slice = 0;
229 AliHLTUInt8_t patch = 0;
230 slice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification );
231 patch = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification );
232 if ( minSlice==0xFF || slice<minSlice ) minSlice = slice;
233 if ( maxSlice==0xFF || slice>maxSlice ) maxSlice = slice;
234 if ( minPatch==0xFF || patch<minPatch ) minPatch = patch;
235 if ( maxPatch==0xFF || patch>maxPatch ) maxPatch = patch;
236 const AliHLTTracksData* pTracks=reinterpret_cast<const AliHLTTracksData*>(pDesc->fPtr);
237 if ((iResult=AliHLTGlobalBarrelTrack::ConvertTrackDataArray(pTracks, pDesc->fSize, inputTrackArray))<0) {
240 trackindexmap.resize(inputTrackArray.size(), -1);
244 if (GetBenchmarkInstance()) {
245 GetBenchmarkInstance()->Stop(3);
246 GetBenchmarkInstance()->Start(4);
251 for (vector<AliHLTGlobalBarrelTrack>::iterator track=inputTrackArray.begin();
252 track!=inputTrackArray.end();
253 track++, trackindex++) {
254 int trackID=track->GetID();
256 // FIXME: error guard
257 HLTError("invalid track ID");
260 if (trackID>=(int)trackindexmap.size())
261 trackindexmap.resize(trackID+1, -1);
262 trackindexmap[trackID]=trackindex;
265 UInt_t nofPoints=track->GetNumberOfPoints();
266 const UInt_t* points=track->GetPoints();
267 for (unsigned i=0; i<nofPoints; i++) {
268 int slice=AliHLTTPCSpacePointData::GetSlice(points[i]);
269 int partition=AliHLTTPCSpacePointData::GetPatch(points[i]);
270 int number=AliHLTTPCSpacePointData::GetNumber(points[i]);
271 HLTInfo("track %d point %d id 0x%08x slice %d partition %d number %d", track->GetID(), i, points[i], slice, partition, number);
275 AliHLTTPCTrackGeometry* trackpoints=new AliHLTTPCTrackGeometry;
276 if (!trackpoints) continue;
277 HLTDebug("track %d id %d:", trackindex, trackID);
279 // in order to avoid rounding errors the track points are
280 // calculated in exactly the same way as in the decoding
281 // Thats why the track instance can not be used directly
282 // but a new instance is created from the values in the
284 // think about moving that to some common code used by
285 // both compression and decoding
286 AliHLTExternalTrackParam param;
287 memset(¶m, 0, sizeof(param));
288 float alpha=track->GetAlpha();
289 while (alpha<0.) alpha+=TMath::TwoPi();
290 while (alpha>TMath::TwoPi()) alpha-=TMath::TwoPi();
291 AliHLTUInt8_t tSlice=AliHLTUInt8_t(9*alpha/TMath::Pi());
292 param.fAlpha =( tSlice + 0.5 ) * TMath::Pi() / 9.0;
293 if (param.fAlpha>TMath::TwoPi()) param.fAlpha-=TMath::TwoPi();
294 param.fX = track->GetX();
295 param.fY = track->GetY();
296 param.fZ = track->GetZ();
297 param.fSinPsi = track->GetSnp();
298 param.fTgl = track->GetTgl();
299 param.fq1Pt = track->GetSigned1Pt();
300 AliHLTGlobalBarrelTrack ctrack(param);
301 ctrack.CalculateHelixParams(bz);
302 trackpoints->InitDriftTimeTransformation(fDriftTimeFactorA, fDriftTimeOffsetA, fDriftTimeFactorC, fDriftTimeOffsetC);
303 trackpoints->SetTrackId(trackID);
304 trackpoints->CalculateTrackPoints(ctrack);
305 trackpoints->RegisterTrackPoints(fTrackGrid);
306 track->SetTrackGeometry(trackpoints);
309 for (vector<AliHLTGlobalBarrelTrack>::const_iterator track=inputTrackArray.begin();
310 track!=inputTrackArray.end();
312 AliHLTTrackGeometry* trackpoints=track->GetTrackGeometry();
313 if (!trackpoints) continue;
314 trackpoints->FillTrackPoints(fTrackGrid);
320 if (GetBenchmarkInstance()) {
321 GetBenchmarkInstance()->Stop(4);
322 GetBenchmarkInstance()->Start(5);
325 // loop over raw cluster blocks, assign to tracks and write
326 // unassigned clusters
327 for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkHWClustersDataType);
328 pDesc!=NULL; pDesc=GetNextInputBlock()) {
329 if (pDesc->fSize<=sizeof(AliRawDataHeader)) continue;
330 if (GetBenchmarkInstance()) {
331 GetBenchmarkInstance()->Start(1);
332 GetBenchmarkInstance()->AddInput(pDesc->fSize);
334 AliHLTUInt8_t slice = 0;
335 AliHLTUInt8_t patch = 0;
336 slice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification );
337 patch = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification );
338 if ( minSlice==0xFF || slice<minSlice ) minSlice = slice;
339 if ( maxSlice==0xFF || slice>maxSlice ) maxSlice = slice;
340 if ( minPatch==0xFF || patch<minPatch ) minPatch = patch;
341 if ( maxPatch==0xFF || patch>maxPatch ) maxPatch = patch;
342 inputRawClusterSize+=pDesc->fSize;
344 // add the data and populate the index grid
345 fRawInputClusters->AddInputBlock(pDesc);
346 fRawInputClusters->PopulateAccessGrid(fSpacePointGrid, pDesc->fSpecification);
347 if (fVerbosity>0 && fSpacePointGrid->GetNumberOfSpacePoints()>0) {
348 HLTInfo("index grid slice %d partition %d", slice, patch);
349 fSpacePointGrid->Print();
350 for (AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid::iterator& cl=fSpacePointGrid->begin();
351 cl!=fSpacePointGrid->end(); cl++) {
352 AliHLTUInt32_t id=cl.Data().fId;
353 float row=fRawInputClusters->GetX(id);
354 float pad=fRawInputClusters->GetY(id);
355 float time=fRawInputClusters->GetZ(id);
356 HLTInfo(" cluster id 0x%08x: row %f pad %f time %f", id, row, pad, time);
359 if (GetBenchmarkInstance()) {
360 GetBenchmarkInstance()->Stop(1);
361 GetBenchmarkInstance()->Start(4);
364 // process the clusters per padrow and check the track grid
365 // for tracks crossing that particular padrow
366 if (GetBenchmarkInstance()) {
367 GetBenchmarkInstance()->Stop(4);
368 GetBenchmarkInstance()->Start(5);
370 allClusters+=fSpacePointGrid->GetNumberOfSpacePoints();
371 iResult=ProcessTrackClusters(&inputTrackArray[0], inputTrackArray.size(), fTrackGrid, trackindexmap, fSpacePointGrid, fRawInputClusters, slice, patch);
372 int assignedInThisPartition=0;
374 assignedInThisPartition=iResult;
375 associatedClusters+=iResult;
377 iResult=ProcessRemainingClusters(&inputTrackArray[0], inputTrackArray.size(), fTrackGrid, trackindexmap, fSpacePointGrid, fRawInputClusters, slice, patch);
379 if (fSpacePointGrid->GetNumberOfSpacePoints()>0) {
380 if (fVerbosity>0) HLTInfo("associated %d (%d) of %d clusters in slice %d partition %d", iResult+assignedInThisPartition, assignedInThisPartition, fSpacePointGrid->GetNumberOfSpacePoints(), slice, patch);
382 associatedClusters+=iResult;
385 // write all remaining clusters not yet assigned to tracks
386 // the index grid is used to write sorted in padrow
387 // FIXME: decoder index instead of data specification to be used
388 // use an external access grid to reduce allocated memory
389 // set to NULL after writing the clusters
390 const char* writeoptions="";
391 if (fpWrittenAssociatedClusterIds) {
392 writeoptions="write-cluster-ids";
394 fRawInputClusters->SetSpacePointPropertyGrid(pDesc->fSpecification, fSpacePointGrid);
395 iResult=fRawInputClusters->Write(outputPtr+size, capacity-size, outputBlocks, fpDataDeflater, writeoptions);
396 fRawInputClusters->SetSpacePointPropertyGrid(pDesc->fSpecification, NULL);
399 outputDataSize+=iResult;
400 // the size of the optional cluster id array must be subtracted
401 if (fpWrittenAssociatedClusterIds && outputBlocks.size()>0 &&
402 outputBlocks.back().fDataType==AliHLTTPCDefinitions::RemainingClusterIdsDataType()) {
403 outputDataSize-=outputBlocks.back().fSize;
405 if (GetBenchmarkInstance()) GetBenchmarkInstance()->AddOutput(iResult);
407 if (GetBenchmarkInstance()) {
408 GetBenchmarkInstance()->Stop(5);
411 fSpacePointGrid->Clear();
413 if (fHistoClusterRatio && allClusters>0) {
414 if (fVerbosity>0) HLTInfo("associated %d of %d clusters to tracks", associatedClusters, allClusters);
415 float ratio=associatedClusters; ratio/=allClusters;
416 fHistoClusterRatio->Fill(ratio);
419 // output of track model clusters
421 AliHLTUInt32_t tracksBufferOffset=sizeof(AliHLTTPCTrackModelBlock);
422 if (capacity-size<tracksBufferOffset) {
426 if (fpWrittenAssociatedClusterIds) fpWrittenAssociatedClusterIds->clear();
427 AliHLTTPCTrackModelBlock* trackModelBlock=reinterpret_cast<AliHLTTPCTrackModelBlock*>(outputPtr+size);
428 trackModelBlock->fVersion=1;
429 trackModelBlock->fDeflaterMode=fpDataDeflater?fpDataDeflater->GetDeflaterVersion():0;
430 trackModelBlock->fTrackCount=inputTrackArray.size();
431 trackModelBlock->fClusterCount=0;
432 trackModelBlock->fGlobalParameterCnt=5;
433 tracksBufferOffset+=trackModelBlock->fGlobalParameterCnt*sizeof(trackModelBlock->fGlobalParameters);
434 if (capacity-size<tracksBufferOffset) {
439 AliHLTUInt32_t parameterIndex=0;
440 trackModelBlock->fGlobalParameters[parameterIndex++]=bz;
441 trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeFactorA;
442 trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeOffsetA;
443 trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeFactorC;
444 trackModelBlock->fGlobalParameters[parameterIndex++]=fDriftTimeOffsetC;
445 if (parameterIndex!=trackModelBlock->fGlobalParameterCnt) {
446 HLTError("internal error, size of parameter array has changed without providing all values");
451 if (trackindexmap.size()>0) {// condition for track model compression
452 iResult=WriteTrackClusters(inputTrackArray, fRawInputClusters, fpDataDeflater, outputPtr+size+tracksBufferOffset, capacity-size-tracksBufferOffset);
454 AliHLTComponent_BlockData bd;
457 bd.fSize = tracksBufferOffset+iResult;
458 bd.fDataType = AliHLTTPCDefinitions::ClusterTracksCompressedDataType();
459 bd.fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(minSlice, maxSlice, minPatch, maxPatch);
460 outputBlocks.push_back(bd);
462 outputDataSize+=bd.fSize;
463 HLTBenchmark("track data block of %d tracks: size %d", inputTrackArray.size(), bd.fSize);
465 if (fpWrittenAssociatedClusterIds && fpWrittenAssociatedClusterIds->size()>0) {
466 AliHLTComponent::FillBlockData(bd);
468 bd.fSize = fpWrittenAssociatedClusterIds->size()*sizeof(vector<AliHLTUInt32_t>::value_type);
469 if (capacity-size>bd.fSize) {
470 memcpy(outputPtr+bd.fOffset, &(*fpWrittenAssociatedClusterIds)[0], bd.fSize);
471 bd.fDataType = AliHLTTPCDefinitions::ClusterIdTracksDataType();
472 bd.fSpecification = AliHLTTPCDefinitions::EncodeDataSpecification(minSlice, maxSlice, minPatch, maxPatch);
473 outputBlocks.push_back(bd);
479 fpWrittenAssociatedClusterIds->clear();
485 fRawInputClusters->Clear();
487 float compressionFactor=(float)inputRawClusterSize;
488 if ((outputDataSize)>0) compressionFactor/=outputDataSize;
489 else compressionFactor=0.;
490 if (fHistoCompFactor) fHistoCompFactor->Fill(compressionFactor);
492 if (GetBenchmarkInstance() && allClusters>0) {
493 GetBenchmarkInstance()->Stop(0);
494 if (fDeflaterMode!=3) {
495 HLTBenchmark("%s - compression factor %.2f", GetBenchmarkInstance()->GetStatistics(), compressionFactor);
497 HLTBenchmark("%s", GetBenchmarkInstance()->GetStatistics());
501 if (fInputClusters) {
502 fInputClusters->Clear();
504 if (fRawInputClusters) {
505 fRawInputClusters->Clear();
512 for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkAliHLTDataTypeClusterMCInfo | kAliHLTDataOriginTPC);
513 pDesc!=NULL; pDesc=GetNextInputBlock()) {
514 outputBlocks.push_back(*pDesc);
520 int AliHLTTPCDataCompressionComponent::ProcessTrackClusters(AliHLTGlobalBarrelTrack* pTracks, unsigned nofTracks,
521 AliHLTTrackGeometry::AliHLTTrackGrid* pTrackIndex,
522 const vector<int>& trackIndexMap,
523 AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid* pClusterIndex,
524 AliHLTSpacePointContainer* pClusters,
525 int slice, int partition) const
527 // process to assigned track clusters
529 int assignedClusters=0;
530 if (!pTracks || nofTracks==0) return 0;
532 vector<int> processedTracks(nofTracks, -1);
533 for (AliHLTTrackGeometry::AliHLTTrackGrid::iterator& trackId=pTrackIndex->begin(slice, partition, -1);
534 trackId!=pTrackIndex->end(); trackId++) {
535 if (trackId.Data()>=trackIndexMap.size()) {
536 HLTError("can not find track id %d in index map of size %d", trackId.Data(), trackIndexMap.size());
539 int trackindex=trackIndexMap[trackId.Data()];
540 if (trackindex<0 || trackindex>=(int)nofTracks) {
541 HLTError("invalid index %d found for track id %d", trackindex, trackId.Data());
544 if (processedTracks[trackindex]>0) continue;
545 processedTracks[trackindex]=1;
546 AliHLTGlobalBarrelTrack& track=pTracks[trackindex];
547 if (!track.GetTrackGeometry()) {
548 HLTError("can not find track geometry for track %d", trackId.Data());
551 AliHLTTPCTrackGeometry* pTrackPoints=dynamic_cast<AliHLTTPCTrackGeometry*>(track.GetTrackGeometry());
553 HLTError("invalid track geometry type for track %d, expecting AliHLTTPCTrackGeometry", trackId.Data());
557 UInt_t nofTrackPoints=track.GetNumberOfPoints();
558 const UInt_t* trackPoints=track.GetPoints();
559 for (unsigned i=0; i<nofTrackPoints; i++) {
560 const AliHLTUInt32_t& clusterId=trackPoints[i];
561 if (AliHLTTPCSpacePointData::GetSlice(clusterId)!=(unsigned)slice ||
562 AliHLTTPCSpacePointData::GetPatch(clusterId)!=(unsigned)partition) {
563 // not in the current partition;
567 int clusterrow=(int)pClusters->GetX(clusterId);
568 AliHLTUInt32_t pointId=AliHLTTPCSpacePointData::GetID(slice, partition, clusterrow);
569 AliHLTTrackGeometry::AliHLTTrackPoint* point=pTrackPoints->GetRawTrackPoint(pointId);
571 //HLTError("can not find track point slice %d partition %d padrow %d (0x%08x) of track %d", slice, partition, clusterrow, pointId, trackId.Data());
574 float pad=point->GetU();
575 float time=point->GetV();
577 iResult=FindCellClusters(trackId.Data(), clusterrow, pad, time, pClusterIndex, pClusters, point, clusterId);
578 if (iResult>0) assignedClusters+=iResult;
581 return assignedClusters;
584 int AliHLTTPCDataCompressionComponent::ProcessRemainingClusters(AliHLTGlobalBarrelTrack* pTracks, unsigned nofTracks,
585 AliHLTTrackGeometry::AliHLTTrackGrid* pTrackIndex,
586 const vector<int>& trackIndexMap,
587 AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid* pClusterIndex,
588 AliHLTSpacePointContainer* pClusters,
589 int slice, int partition) const
591 // assign remaining clusters to tracks
593 int associatedClusters=0;
594 if (!pTracks || nofTracks==0) return 0;
596 for (int padrow=0; padrow<AliHLTTPCTransform::GetNRows(partition); padrow++) {
597 for (AliHLTTrackGeometry::AliHLTTrackGrid::iterator& trackId=pTrackIndex->begin(slice, partition, padrow);
598 trackId!=pTrackIndex->end(); trackId++) {
599 if (trackId.Data()>=trackIndexMap.size()) {
600 HLTError("can not find track id %d in index map of size %d", trackId.Data(), trackIndexMap.size());
603 int trackindex=trackIndexMap[trackId.Data()];
604 if (trackindex<0 || trackindex>=(int)nofTracks) {
605 HLTError("invalid index %d found for track id %d", trackindex, trackId.Data());
608 AliHLTGlobalBarrelTrack& track=pTracks[trackindex];
609 if (!track.GetTrackGeometry()) {
610 HLTError("can not find track geometry for track %d", trackId.Data());
613 AliHLTTPCTrackGeometry* pTrackPoints=dynamic_cast<AliHLTTPCTrackGeometry*>(track.GetTrackGeometry());
615 HLTError("invalid track geometry type for track %d, expecting AliHLTTPCTrackGeometry", trackId.Data());
618 AliHLTUInt32_t pointId=AliHLTTPCSpacePointData::GetID(slice, partition, padrow);
619 AliHLTTrackGeometry::AliHLTTrackPoint* point=pTrackPoints->GetRawTrackPoint(pointId);
621 //HLTError("can not find track point slice %d partition %d padrow %d (0x%08x) of track %d", slice, partition, padrow, pointId, trackId.Data());
624 float pad=point->GetU();
625 float time=point->GetV();
627 iResult=FindCellClusters(trackId.Data(), padrow, pad, time, pClusterIndex, pClusters, point);
628 if (iResult>0) associatedClusters+=iResult;
630 HLTInfo("trackpoint track %d slice %d partition %d padrow %d: %.3f \t%.3f - associated %d", track.GetID(), slice, partition, padrow, pad, time, iResult);
634 if (iResult<0) return iResult;
635 return associatedClusters;
638 int AliHLTTPCDataCompressionComponent::FindCellClusters(int trackId, int padrow, float pad, float time,
639 AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid* pClusterIndex,
640 AliHLTSpacePointContainer* pClusters,
641 AliHLTTrackGeometry::AliHLTTrackPoint* pTrackPoint,
642 AliHLTUInt32_t clusterId) const
644 // check index cell for entries and assign to track
646 // search a 4x4 matrix out of the 9x9 matrix around the cell addressed by
648 int rowindex=pClusterIndex->GetXIndex((float)padrow);
649 int padstartindex=pClusterIndex->GetYIndex(pad);
650 int timestartindex=pClusterIndex->GetZIndex(time);
651 int cellindex=pClusterIndex->Index(rowindex, padstartindex, timestartindex);
652 float centerpad=pClusterIndex->GetCenterY(cellindex);
653 float centertime=pClusterIndex->GetCenterZ(cellindex);
654 if ((TMath::Abs(centerpad-pad)>fMaxDeltaPad && pad>0.) ||
655 (TMath::Abs(centertime-time)>fMaxDeltaTime && time>0.)) {
656 ALIHLTERRORGUARD(20, "invalid pad center calculation, please check dimensions if dimensions of index grid match the maximum possible deviation");
661 if (centerpad>pad) paddirection=-1;
662 if (centertime>time) timedirection=-1;
663 for (int padcount=0, padindex=padstartindex; padcount<2; padcount++, padindex+=paddirection) {
664 if (padindex<0) continue;
665 if (padindex>=pClusterIndex->GetDimensionY()) break;
666 for (int timecount=0, timeindex=timestartindex; timecount<2; timecount++, timeindex+=timedirection) {
667 if (timeindex<0) continue;
668 if (timeindex>=pClusterIndex->GetDimensionZ()) break;
669 cellindex=pClusterIndex->Index(rowindex, padindex, timeindex);
670 float cellpad=pClusterIndex->GetCenterY(cellindex);
671 float celltime=pClusterIndex->GetCenterZ(cellindex);
672 for (AliHLTSpacePointContainer::AliHLTSpacePointPropertyGrid::iterator& cl=pClusterIndex->begin((float)padrow, cellpad, celltime);
673 cl!=pClusterIndex->end(); cl++) {
674 if (cl.Data().fTrackId>=0) continue;
675 if (clusterId!=(~(AliHLTUInt32_t)0) && clusterId!=cl.Data().fId) continue;
676 if (TMath::Abs(padrow-pClusters->GetX(cl.Data().fId))>=1.) {
677 HLTError("cluster 0x%08x: mismatch on padrow: trackpoint %d cluster %f", cl.Data().fId, padrow, pClusters->GetX(cl.Data().fId));
680 float clusterpad=pClusters->GetY(cl.Data().fId);
681 float clustertime=pClusters->GetZ(cl.Data().fId);
682 if (TMath::Abs(clusterpad-pad)<fMaxDeltaPad &&
683 TMath::Abs(clustertime-time)<fMaxDeltaTime) {
684 // add this cluster to the track point and mark in the index grid
685 cl.Data().fTrackId=trackId;
686 pTrackPoint->AddAssociatedSpacePoint(cl.Data().fId, clusterpad-pad, clustertime-time);
689 if (clusterId!=(~(AliHLTUInt32_t)0)) break;
696 int AliHLTTPCDataCompressionComponent::WriteTrackClusters(const vector<AliHLTGlobalBarrelTrack>& tracks,
697 AliHLTSpacePointContainer* pSpacePoints,
698 AliHLTDataDeflater* pDeflater,
699 AliHLTUInt8_t* outputPtr,
700 AliHLTUInt32_t capacity) const
702 // write the track data block including all associated clusters
703 AliHLTUInt32_t size=0;
704 for (vector<AliHLTGlobalBarrelTrack>::const_iterator track=tracks.begin();
707 if (!track->GetTrackGeometry()) {
708 HLTError("can not find track geometry for track %d", track->GetID());
711 AliHLTTPCTrackGeometry* pTrackPoints=dynamic_cast<AliHLTTPCTrackGeometry*>(track->GetTrackGeometry());
713 HLTError("invalid track geometry type for track %d, expecting AliHLTTPCTrackGeometry", track->GetID());
717 int result=pTrackPoints->Write(*track, pSpacePoints, pDeflater, outputPtr+size, capacity-size, fpWrittenAssociatedClusterIds);
718 if (result<0) return result;
721 UInt_t nofTrackPoints=track->GetNumberOfPoints();
722 const UInt_t* trackPoints=track->GetPoints();
724 int assignedPoints=0;
725 int assignedTrackPoints=0;
726 const vector<AliHLTTrackGeometry::AliHLTTrackPoint>& rawPoints=pTrackPoints->GetRawPoints();
727 for (vector<AliHLTTrackGeometry::AliHLTTrackPoint>::const_iterator point=rawPoints.begin();
728 point!=rawPoints.end(); point++) {
729 const vector<AliHLTTrackGeometry::AliHLTTrackSpacepoint>& spacePoints=point->GetSpacepoints();
730 for (vector<AliHLTTrackGeometry::AliHLTTrackSpacepoint>::const_iterator spacePoint=spacePoints.begin();
731 spacePoint!=spacePoints.end(); spacePoint++) {
732 float dpad=spacePoint->GetResidual(0);
733 float dtime=spacePoint->GetResidual(1);
734 if (dpad>-1000 && dtime>-1000 && fHistoResidualPad && fHistoResidualTime) {
735 fHistoResidualPad->Fill(dpad);
736 fHistoResidualTime->Fill(dtime);
739 for (unsigned i=0; i<nofTrackPoints; i++) {
740 if (trackPoints[i]==spacePoint->fId) {
741 assignedTrackPoints++;
747 if (fHistoClustersOnTracks) {
748 fHistoClustersOnTracks->Fill(assignedPoints);
750 if (fHistoTrackClusterRatio && nofTrackPoints>0) {
751 float ratio=assignedTrackPoints; ratio/=nofTrackPoints;
752 fHistoTrackClusterRatio->Fill(ratio);
758 int AliHLTTPCDataCompressionComponent::DoInit( int argc, const char** argv )
760 /// inherited from AliHLTComponent: component initialisation and argument scan.
763 // component configuration
764 //Stage 1: default initialization.
768 TString cdbPath("HLT/ConfigTPC/");
769 cdbPath += GetComponentID();
771 iResult = ConfigureFromCDBTObjString(cdbPath);
775 //Stage 3: command line arguments.
776 if (argc && (iResult = ConfigureFromArgumentString(argc, argv)) < 0)
779 std::auto_ptr<AliHLTComponentBenchmark> benchmark(new AliHLTComponentBenchmark);
780 if (benchmark.get()) {
781 benchmark->SetTimer(0,"total");
782 benchmark->SetTimer(1,"rawclusterinput");
783 benchmark->SetTimer(2,"clusterinput");
784 benchmark->SetTimer(3,"trackinput");
785 benchmark->SetTimer(4,"processing");
786 benchmark->SetTimer(5,"output");
791 unsigned spacePointContainerMode=0;
792 if (fMode==2 || fMode==4) {
793 // initialize map data for cluster access in the track association loop
794 spacePointContainerMode|=AliHLTTPCHWCFSpacePointContainer::kModeCreateMap;
796 if (fMode==3 || fMode==4) {
797 // optimized storage format: differential pad and time storage
798 spacePointContainerMode|=AliHLTTPCHWCFSpacePointContainer::kModeDifferentialPadTime;
800 std::auto_ptr<AliHLTTPCHWCFSpacePointContainer> rawInputClusters(new AliHLTTPCHWCFSpacePointContainer(spacePointContainerMode));
801 std::auto_ptr<AliHLTTPCSpacePointContainer> inputClusters(new AliHLTTPCSpacePointContainer);
803 std::auto_ptr<TH1F> histoCompFactor(new TH1F("CompressionFactor",
804 "HLT TPC data compression factor",
806 std::auto_ptr<TH1F> histoResidualPad(new TH1F("PadResidual",
807 "HLT TPC pad residual",
808 100, -fMaxDeltaPad, fMaxDeltaPad));
809 std::auto_ptr<TH1F> histoResidualTime(new TH1F("TimeResidual",
810 "HLT TPC time residual",
811 100, -fMaxDeltaTime, fMaxDeltaTime));
812 std::auto_ptr<TH1F> histoClustersOnTracks(new TH1F("ClustersOnTracks",
813 "Clusters in track model compression",
815 std::auto_ptr<TH1F> histoClusterRatio(new TH1F("ClusterRatio",
816 "Fraction of clusters in track model compression",
818 std::auto_ptr<TH1F> histoTrackClusterRatio(new TH1F("UsedTrackClusters",
819 "Fraction of track clusters in track model compression",
822 // track grid: 36 slices, each 6 partitions with max 33 rows
823 fTrackGrid=new AliHLTTrackGeometry::AliHLTTrackGrid(36, 1, 6, 1, 33, 1, 20000);
824 fSpacePointGrid=AliHLTTPCHWCFSpacePointContainer::AllocateIndexGrid();
826 if (!rawInputClusters.get() ||
827 !inputClusters.get() ||
830 if (fTrackGrid) delete fTrackGrid; fTrackGrid=NULL;
831 if (fSpacePointGrid) delete fSpacePointGrid; fSpacePointGrid=NULL;
835 if (fDeflaterMode>0 && (iResult=InitDeflater(fDeflaterMode))<0)
838 fpBenchmark=benchmark.release();
839 fRawInputClusters=rawInputClusters.release();
840 fInputClusters=inputClusters.release();
842 // initialize the histograms if stored at the end
843 // condition might be extended
844 if (!fHistogramFile.IsNull()) {
845 fHistoCompFactor=histoCompFactor.release();
846 fHistoResidualPad=histoResidualPad.release();
847 fHistoResidualTime=histoResidualTime.release();
848 fHistoClustersOnTracks=histoClustersOnTracks.release();
849 fHistoClusterRatio=histoClusterRatio.release();
850 fHistoTrackClusterRatio=histoTrackClusterRatio.release();
853 if (iResult>=0 && (iResult=InitDriftTimeTransformation())<0) return iResult;
858 int AliHLTTPCDataCompressionComponent::InitDeflater(int mode)
860 /// init the data deflater
862 if (mode==2 || mode==3) {
864 std::auto_ptr<AliHLTDataDeflaterHuffman> deflater(new AliHLTDataDeflaterHuffman(mode==3));
865 if (!deflater.get()) return -ENOMEM;
867 if (!deflater->IsTrainingMode()) {
868 TString cdbPath("HLT/ConfigTPC/");
869 cdbPath += GetComponentID();
870 cdbPath += "HuffmanTables";
871 TObject* pConf=LoadAndExtractOCDBObject(cdbPath);
872 if (!pConf) return -ENOENT;
873 if (dynamic_cast<TList*>(pConf)==NULL) {
874 HLTError("huffman table configuration object of inconsistent type");
877 iResult=deflater->InitDecoders(dynamic_cast<TList*>(pConf));
878 if (iResult<0) return iResult;
881 if (!fHistogramFile.IsNull())
882 deflater->EnableStatistics();
884 unsigned nofParameters=AliHLTTPCDefinitions::GetNumberOfClusterParameterDefinitions();
886 for (; p<nofParameters; p++) {
887 const AliHLTTPCDefinitions::AliClusterParameter& parameter=AliHLTTPCDefinitions::fgkClusterParameterDefinitions[p];
888 // use the pad/time length as reference for the calculation of ratio for residuals
889 unsigned refLength=0;
890 unsigned refLengthPad=0;
891 unsigned refLengthTime=0;
892 if (parameter.fId==AliHLTTPCDefinitions::kPad) refLengthPad=parameter.fBitLength;
893 else if (parameter.fId==AliHLTTPCDefinitions::kTime) refLengthTime=parameter.fBitLength;
894 else if (parameter.fId==AliHLTTPCDefinitions::kResidualPad) refLength=refLengthPad;
895 else if (parameter.fId==AliHLTTPCDefinitions::kResidualTime) refLength=refLengthTime;
897 if (deflater->AddParameterDefinition(parameter.fName,
898 parameter.fBitLength,
899 refLength)!=(int)parameter.fId) {
900 // for performance reason the parameter id is simply used as index in the array of
901 // definitions, the position must match the id
902 HLTFatal("mismatch between parameter id and position in array for parameter %s, rearrange definitions!", parameter.fName);
906 fpDataDeflater=deflater.release();
910 std::auto_ptr<AliHLTDataDeflaterSimple> deflater(new AliHLTDataDeflaterSimple);
911 if (!deflater.get()) return -ENOMEM;
913 if (!fHistogramFile.IsNull())
914 deflater->EnableStatistics();
916 unsigned nofParameters=AliHLTTPCDefinitions::GetNumberOfClusterParameterDefinitions();
918 for (; p<nofParameters; p++) {
919 const AliHLTTPCDefinitions::AliClusterParameter& parameter=AliHLTTPCDefinitions::fgkClusterParameterDefinitions[p];
920 if (deflater->AddParameterDefinition(parameter.fName,
921 parameter.fBitLength,
922 parameter.fOptional)!=(int)parameter.fId) {
923 // for performance reason the parameter id is simply used as index in the array of
924 // definitions, the position must match the id
925 HLTFatal("mismatch between parameter id and position in array for parameter %s, rearrange definitions!", parameter.fName);
929 fpDataDeflater=deflater.release();
932 HLTError("invalid deflater mode %d, allowed 1=simple 2=huffman", mode);
936 int AliHLTTPCDataCompressionComponent::DoDeinit()
938 /// inherited from AliHLTComponent: component cleanup
940 if (fpBenchmark) delete fpBenchmark; fpBenchmark=NULL;
941 if (fRawInputClusters) delete fRawInputClusters; fRawInputClusters=NULL;
942 if (fInputClusters) delete fInputClusters; fInputClusters=NULL;
943 if (!fHistogramFile.IsNull()) {
944 TFile out(fHistogramFile, "RECREATE");
945 if (!out.IsZombie()) {
947 if (fHistoCompFactor) fHistoCompFactor->Write();
948 if (fHistoResidualPad) fHistoResidualPad->Write();
949 if (fHistoResidualTime) fHistoResidualTime->Write();
950 if (fHistoClusterRatio) fHistoClusterRatio->Write();
951 if (fHistoClustersOnTracks) fHistoClustersOnTracks->Write();
952 if (fHistoTrackClusterRatio) fHistoTrackClusterRatio->Write();
956 if (fHistoCompFactor) delete fHistoCompFactor;
957 fHistoCompFactor=NULL;
958 if (fHistoResidualPad) delete fHistoResidualPad;
959 fHistoResidualPad=NULL;
960 if (fHistoResidualTime) delete fHistoResidualTime;
961 fHistoResidualTime=NULL;
962 if (fHistoClustersOnTracks) delete fHistoClustersOnTracks;
963 fHistoClustersOnTracks=NULL;
964 if (fHistoClusterRatio) delete fHistoClusterRatio;
965 fHistoClusterRatio=NULL;
966 if (fHistoTrackClusterRatio) delete fHistoTrackClusterRatio;
967 fHistoTrackClusterRatio=NULL;
969 if (fpDataDeflater) {
970 if (!fHistogramFile.IsNull()) {
971 TString filename=fHistogramFile;
972 filename.ReplaceAll(".root", "-deflater.root");
973 fpDataDeflater->SaveAs(filename);
975 if (fDeflaterMode==3) {
976 if (fTrainingTableOutput.IsNull()) {
977 fTrainingTableOutput=GetComponentID();
978 fTrainingTableOutput+="-huffman.root";
980 // TODO: currently, the code tables are also calculated in FindObject
981 // check if a different function is more appropriate
982 TObject* pConf=fpDataDeflater->FindObject("DeflaterConfiguration");
984 TString cdbEntryPath("HLT/ConfigTPC/");
985 cdbEntryPath += GetComponentID();
986 cdbEntryPath += "HuffmanTables";
987 AliCDBPath cdbPath(cdbEntryPath);
988 AliCDBId cdbId(cdbPath, AliCDBManager::Instance()->GetRun(), AliCDBRunRange::Infinity(), 0, 0);
989 AliCDBMetaData* cdbMetaData=new AliCDBMetaData;
990 cdbMetaData->SetResponsible("ALICE HLT Matthias.Richter@cern.ch");
991 cdbMetaData->SetComment("Huffman encoder configuration");
992 AliCDBEntry* entry=new AliCDBEntry(pConf, cdbId, cdbMetaData, kTRUE);
994 entry->SaveAs(fTrainingTableOutput);
995 // this is a small memory leak
996 // seg fault in ROOT object handling if the two objects are deleted
999 //delete cdbMetaData;
1002 delete fpDataDeflater;
1004 fpDataDeflater=NULL;
1007 if (fTrackGrid) delete fTrackGrid; fTrackGrid=NULL;
1008 if (fSpacePointGrid) delete fSpacePointGrid; fSpacePointGrid=NULL;
1013 int AliHLTTPCDataCompressionComponent::ScanConfigurationArgument(int argc, const char** argv)
1015 /// inherited from AliHLTComponent: argument scan
1017 if (argc<1) return 0;
1018 int bMissingParam=0;
1020 TString argument=argv[i];
1024 if (argument.CompareTo("-mode")==0) {
1025 if ((bMissingParam=(++i>=argc))) break;
1026 TString parameter=argv[i];
1027 if (parameter.IsDigit()) {
1028 fMode=parameter.Atoi();
1031 HLTError("invalid parameter for argument %s, expecting number instead of %s", argument.Data(), parameter.Data());
1037 if (argument.CompareTo("-deflater-mode")==0) {
1038 if ((bMissingParam=(++i>=argc))) break;
1039 TString parameter=argv[i];
1040 if (parameter.IsDigit()) {
1041 fDeflaterMode=parameter.Atoi();
1044 HLTError("invalid parameter for argument %s, expecting number instead of %s", argument.Data(), parameter.Data());
1050 if (argument.CompareTo("-histogram-file")==0) {
1051 if ((bMissingParam=(++i>=argc))) break;
1052 fHistogramFile=argv[i++];
1055 // -save-histogram-table
1056 if (argument.CompareTo("-save-huffman-table")==0) {
1057 if ((bMissingParam=(++i>=argc))) break;
1058 fTrainingTableOutput=argv[i++];
1061 // -cluster-verification
1062 if (argument.CompareTo("-cluster-verification")==0) {
1063 if ((bMissingParam=(++i>=argc))) break;
1064 TString parameter=argv[i];
1065 if (parameter.IsDigit()) {
1066 fVerificationMode=parameter.Atoi();
1069 HLTError("invalid parameter for argument %s, expecting number instead of %s", argument.Data(), parameter.Data());
1073 } while (0); // using do-while only to have break available
1075 if (bMissingParam) {
1076 HLTError("missing parameter for argument %s", argument.Data());
1083 int AliHLTTPCDataCompressionComponent::InitDriftTimeTransformation()
1085 /// calculate correction factor and offset for a linear approximation of the
1086 /// drift time transformation, separately for A and C side
1088 AliHLTTPCClusterTransformation transform;
1089 if ((iResult=transform.Init( GetBz(), GetTimeStamp()))<0) {
1090 HLTError("failed to init AliHLTTPCClusterTransformation: %d", iResult);
1094 if ((iResult=CalculateDriftTimeTransformation(transform, 0, 0, fDriftTimeFactorA, fDriftTimeOffsetA))<0) return iResult;
1095 if (fVerbosity>0) HLTInfo("drift time transformation A side: m=%f n=%f", fDriftTimeFactorA, fDriftTimeOffsetA);
1096 if ((iResult=CalculateDriftTimeTransformation(transform, 18, 0, fDriftTimeFactorC, fDriftTimeOffsetC))<0) return iResult;
1097 if (fVerbosity>0) HLTInfo("drift time transformation C side: m=%f n=%f", fDriftTimeFactorC, fDriftTimeOffsetC);
1102 int AliHLTTPCDataCompressionComponent::CalculateDriftTimeTransformation(AliHLTTPCClusterTransformation& transform,
1103 int slice, int padrow,
1104 float& m, float& n) const
1106 /// calculate correction factor and offset for a linear approximation of the
1107 /// drift time transformation by just probing the range of timebins with
1108 /// AliHLTTPCClusterTransformation
1109 const int nofSteps=100;
1110 vector<float> zvalues;
1112 int nofTimebins=AliHLTTPCTransform::GetNTimeBins();
1113 int stepWidth=nofTimebins/nofSteps;
1118 for (time=0; time<nofTimebins; time+=stepWidth, count++) {
1120 transform.Transform(slice, padrow, 0, time, xyz);
1121 zvalues.push_back(xyz[2]);
1125 if (count==0) count=1;
1131 for (vector<float>::const_iterator z=zvalues.begin();
1132 z!=zvalues.end(); z++, time+=stepWidth) {
1133 sumTZ+=(time-meanT)*((*z)-meanZ);
1134 sumT2+=(time-meanT)*(time-meanT);