Line data Source code
1 : // $Id$
2 :
3 : //**************************************************************************
4 : //* This file is property of and copyright by the ALICE HLT Project *
5 : //* ALICE Experiment at CERN, All rights reserved. *
6 : //* *
7 : //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no> *
8 : //* for The ALICE HLT Project. *
9 : //* *
10 : //* Permission to use, copy, modify and distribute this software and its *
11 : //* documentation strictly for non-commercial purposes is hereby granted *
12 : //* without fee, provided that the above copyright notice appears in all *
13 : //* copies and that both the copyright notice and this permission notice *
14 : //* appear in the supporting documentation. The authors make no claims *
15 : //* about the suitability of this software for any purpose. It is *
16 : //* provided "as is" without express or implied warranty. *
17 : //**************************************************************************
18 :
19 : /// @file AliHLTTPCSpacePointContainer.h
20 : /// @author Matthias Richter
21 : /// @date 2011-04-29
22 : /// @brief Helper class for handling of HLT TPC space point data blocks
23 : ///
24 :
25 : #include "AliHLTTPCSpacePointContainer.h"
26 : #include "AliHLTTPCSpacePointData.h"
27 : #include "AliHLTTPCDefinitions.h"
28 : #include "AliHLTTPCGeometry.h"
29 : #include "AliHLTComponent.h"
30 : #include "AliHLTTemplates.h"
31 : #include "TMath.h"
32 : #include <memory>
33 : #include <algorithm>
34 :
35 : /** ROOT macro for the implementation of ROOT specific class methods */
36 6 : ClassImp(AliHLTTPCSpacePointContainer)
37 :
38 : AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointContainer()
39 0 : : AliHLTSpacePointContainer()
40 0 : , fClusters()
41 0 : , fSelections()
42 0 : {
43 : // see header file for class documentation
44 : // or
45 : // refer to README to build package
46 : // or
47 : // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
48 0 : }
49 :
50 : AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointContainer(const AliHLTTPCSpacePointContainer& c)
51 0 : : AliHLTSpacePointContainer(c)
52 0 : , fClusters(c.fClusters.begin(), c.fClusters.end())
53 0 : , fSelections()
54 0 : {
55 : /// copy constructor
56 0 : }
57 :
58 : AliHLTTPCSpacePointContainer& AliHLTTPCSpacePointContainer::operator=(const AliHLTTPCSpacePointContainer& c)
59 : {
60 : /// assignment operator
61 0 : if (&c==this) return *this;
62 0 : AliHLTSpacePointContainer::operator=(c);
63 0 : fClusters=c.fClusters;
64 :
65 0 : return *this;
66 0 : }
67 :
68 0 : AliHLTTPCSpacePointContainer::~AliHLTTPCSpacePointContainer()
69 0 : {
70 : // destructor
71 :
72 0 : Clear();
73 0 : }
74 :
75 : int AliHLTTPCSpacePointContainer::AddInputBlock(const AliHLTComponentBlockData* pDesc)
76 : {
77 : // add input block to the collection
78 0 : if (!pDesc) return -EINVAL;
79 : int count=0;
80 0 : if (pDesc->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType) {
81 0 : HLTWarning("ignoring data block of type %s", AliHLTComponent::DataType2Text(pDesc->fDataType).c_str());
82 0 : return 0;
83 : }
84 0 : if (!pDesc->fPtr || pDesc->fSize<sizeof(AliHLTTPCClusterData)) return -ENODATA;
85 :
86 : // consistency check of the input block
87 0 : const AliHLTTPCClusterData* pClusterData=reinterpret_cast<AliHLTTPCClusterData*>(pDesc->fPtr);
88 0 : if (pClusterData->fSpacePointCnt*sizeof(AliHLTTPCSpacePointData)+sizeof(AliHLTTPCClusterData)!=pDesc->fSize) {
89 0 : HLTError("data block of type %s corrupted: number of entries %d is not consistent with block size %d",
90 : AliHLTComponent::DataType2Text(pDesc->fDataType).c_str(), pClusterData->fSpacePointCnt, pDesc->fSize);
91 0 : return -EBADMSG;
92 : }
93 :
94 0 : AliHLTUInt8_t minslice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification );
95 0 : AliHLTUInt8_t maxslice = AliHLTTPCDefinitions::GetMaxSliceNr( pDesc->fSpecification );
96 0 : AliHLTUInt8_t minpart = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification );
97 0 : AliHLTUInt8_t maxpart = AliHLTTPCDefinitions::GetMaxPatchNr( pDesc->fSpecification );
98 0 : bool bIsSinglePartition=(pDesc->fSpecification==kAliHLTVoidDataSpec?false:minslice==maxslice && minpart==maxpart);
99 :
100 0 : for (UInt_t i=0; i<pClusterData->fSpacePointCnt; i++) {
101 0 : AliHLTUInt32_t clusterID=~(AliHLTUInt32_t)0;
102 0 : if (bIsSinglePartition) {
103 : // cluster ID has to match ID encoded from slice, partition and index
104 0 : clusterID=AliHLTTPCSpacePointData::GetID(minslice, minpart, i);
105 0 : if (clusterID!=pClusterData->fSpacePoints[i].fID) {
106 0 : HLTWarning("cluster ID 0x%08x does not match slice %d partition %d index %d",
107 : pClusterData->fSpacePoints[i].fID, minslice, minpart, i);
108 : }
109 : } else {
110 : // check the cluster ID for correct bounds
111 0 : clusterID=pClusterData->fSpacePoints[i].fID;
112 0 : UInt_t clusterSlice =AliHLTTPCSpacePointData::GetSlice(clusterID);
113 0 : UInt_t clusterPart =AliHLTTPCSpacePointData::GetPatch(clusterID);
114 0 : UInt_t clusterNumber=AliHLTTPCSpacePointData::GetNumber(clusterID);
115 0 : if (pDesc->fSpecification!=kAliHLTVoidDataSpec) {
116 0 : if (clusterSlice<minslice || clusterSlice>maxslice ||
117 0 : clusterPart<minpart || clusterPart>maxpart) {
118 0 : HLTWarning("cluster ID 0x%08d out of range indicated by data specification 0x%08x", clusterID, pDesc->fSpecification);
119 0 : } else if (clusterNumber!=i) {
120 0 : HLTWarning("cluster ID 0x%08d does not match index %d in block 0x%08x", clusterID, i, pDesc->fSpecification);
121 : }
122 : }
123 : }
124 : {
125 : // consistency check for x and row number
126 : // UInt_t clusterSlice =AliHLTTPCSpacePointData::GetSlice(clusterID);
127 : // UInt_t clusterPart =AliHLTTPCSpacePointData::GetPatch(clusterID);
128 : // int row=AliHLTTPCGeometry::GetPadRow(pClusterData->fSpacePoints[i].fX);
129 : // if (row<AliHLTTPCGeometry::GetFirstRow(clusterPart) || row>AliHLTTPCGeometry::GetLastRow(clusterPart)) {
130 : // HLTError("row number %d calculated from x value %f is outside slice %d partition %d, expected row %d"
131 : // , row, pClusterData->fSpacePoints[i].fX, clusterSlice, clusterPart, pClusterData->fSpacePoints[i].fPadRow);
132 : // }
133 : }
134 :
135 0 : if (fClusters.find(clusterID)==fClusters.end()) {
136 : // new cluster
137 0 : fClusters[clusterID]=AliHLTTPCSpacePointProperties(&pClusterData->fSpacePoints[i]);
138 0 : count++;
139 0 : } else {
140 0 : HLTError("cluster with ID 0x%08x already existing, skipping cluster %d of data block 0x%08x",
141 : clusterID, i, pDesc->fSpecification);
142 : }
143 0 : }
144 :
145 : return count;
146 0 : }
147 :
148 : int AliHLTTPCSpacePointContainer::GetClusterIDs(vector<AliHLTUInt32_t>& tgt) const
149 : {
150 : // get array of cluster IDs
151 0 : tgt.clear();
152 0 : transform(fClusters.begin(), fClusters.end(), back_inserter(tgt), HLT::AliGetKey());
153 0 : return tgt.size();
154 : }
155 :
156 : bool AliHLTTPCSpacePointContainer::Check(AliHLTUInt32_t clusterID) const
157 : {
158 : // check if the cluster is available
159 0 : return fClusters.find(clusterID)!=fClusters.end();
160 : }
161 :
162 : const vector<AliHLTUInt32_t>* AliHLTTPCSpacePointContainer::GetClusterIDs(AliHLTUInt32_t mask)
163 : {
164 : // get array of cluster IDs filtered by mask
165 0 : if (fSelections.find(mask)!=fSelections.end()) {
166 : // return existing selection
167 0 : return fSelections.find(mask)->second;
168 : }
169 : // create new collection
170 0 : vector<AliHLTUInt32_t>* selected=new vector<AliHLTUInt32_t>;
171 0 : if (!selected) return NULL;
172 0 : UInt_t slice=AliHLTTPCSpacePointData::GetSlice(mask);
173 0 : UInt_t partition=AliHLTTPCSpacePointData::GetPatch(mask);
174 0 : HLTInfo("creating collection 0x%08x", mask);
175 0 : for (std::map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator cl=fClusters.begin();
176 0 : cl!=fClusters.end(); cl++) {
177 0 : UInt_t s=AliHLTTPCSpacePointData::GetSlice(cl->first);
178 0 : UInt_t p=AliHLTTPCSpacePointData::GetPatch(cl->first);
179 0 : if ((slice>=(unsigned)AliHLTTPCGeometry::GetNSlice() || s==slice) &&
180 0 : (partition>=(unsigned)AliHLTTPCGeometry::GetNumberOfPatches() || p==partition)) {
181 0 : selected->push_back(cl->first);
182 0 : }
183 : }
184 0 : HLTInfo("collection 0x%08x with %d spacepoints", mask, selected->size());
185 0 : fSelections[mask]=selected;
186 : return selected;
187 0 : }
188 :
189 : float AliHLTTPCSpacePointContainer::GetX(AliHLTUInt32_t clusterID) const
190 : {
191 : // get X coordinate
192 0 : if (fClusters.find(clusterID)==fClusters.end() ||
193 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
194 : // FIXME: understand deviation from the nominal x value
195 : // there is a small deviation in the x coordinate - padrow number correlation
196 : // in principle, the clusterfinder only uses the mapping to set the x parameter.
197 : // now extracting the x value from the padrow no.
198 : //return fClusters.find(clusterID)->second.Data()->fX;
199 0 : return AliHLTTPCGeometry::Row2X(fClusters.find(clusterID)->second.Data()->fPadRow);
200 0 : }
201 :
202 : float AliHLTTPCSpacePointContainer::GetXWidth(AliHLTUInt32_t clusterID) const
203 : {
204 : // get error for X coordinate
205 0 : if (fClusters.find(clusterID)==fClusters.end() ||
206 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
207 0 : return 0.0; // fixed in padrow number
208 0 : }
209 :
210 : float AliHLTTPCSpacePointContainer::GetY(AliHLTUInt32_t clusterID) const
211 : {
212 : // get Y coordinate
213 0 : if (fClusters.find(clusterID)==fClusters.end() ||
214 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
215 0 : return fClusters.find(clusterID)->second.Data()->fY;
216 0 : }
217 :
218 : float AliHLTTPCSpacePointContainer::GetYWidth(AliHLTUInt32_t clusterID) const
219 : {
220 : // get error for Y coordinate
221 0 : if (fClusters.find(clusterID)==fClusters.end() ||
222 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
223 0 : return fClusters.find(clusterID)->second.Data()->fSigmaY2;
224 0 : }
225 :
226 : float AliHLTTPCSpacePointContainer::GetZ(AliHLTUInt32_t clusterID) const
227 : {
228 : // get Z coordinate
229 0 : if (fClusters.find(clusterID)==fClusters.end() ||
230 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
231 0 : return fClusters.find(clusterID)->second.Data()->fZ;
232 0 : }
233 :
234 : float AliHLTTPCSpacePointContainer::GetZWidth(AliHLTUInt32_t clusterID) const
235 : {
236 : // get error for Z coordinate
237 0 : if (fClusters.find(clusterID)==fClusters.end() ||
238 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
239 0 : return fClusters.find(clusterID)->second.Data()->fSigmaZ2;
240 0 : }
241 :
242 : float AliHLTTPCSpacePointContainer::GetCharge(AliHLTUInt32_t clusterID) const
243 : {
244 : // get charge
245 0 : if (fClusters.find(clusterID)==fClusters.end() ||
246 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
247 0 : return fClusters.find(clusterID)->second.Data()->fCharge;
248 0 : }
249 :
250 : float AliHLTTPCSpacePointContainer::GetMaxSignal(AliHLTUInt32_t clusterID) const
251 : {
252 : // get charge
253 0 : if (fClusters.find(clusterID)==fClusters.end() ||
254 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0.0;
255 0 : return fClusters.find(clusterID)->second.Data()->fQMax;
256 0 : }
257 :
258 : float AliHLTTPCSpacePointContainer::GetPhi(AliHLTUInt32_t clusterID) const
259 : {
260 : // get charge
261 :
262 : // phi can be derived directly from the id, no need to search
263 : // for existing cluster
264 0 : int slice=AliHLTTPCSpacePointData::GetSlice(clusterID);
265 0 : return ( slice + 0.5 ) * TMath::Pi() / 9.0;
266 : }
267 :
268 : UChar_t AliHLTTPCSpacePointContainer::GetPadRow(AliHLTUInt32_t clusterID) const
269 : {
270 : // get pad row
271 0 : if (fClusters.find(clusterID)==fClusters.end() ||
272 0 : fClusters.find(clusterID)->second.Data()==NULL) return 0;
273 0 : return fClusters.find(clusterID)->second.Data()->fPadRow;
274 0 : }
275 :
276 : void AliHLTTPCSpacePointContainer::Clear(Option_t * option)
277 : {
278 : // clear the object and reset pointer references
279 0 : fClusters.clear();
280 :
281 0 : for (std::map<AliHLTUInt32_t, vector<AliHLTUInt32_t>*>::iterator selection=fSelections.begin();
282 0 : selection!=fSelections.end(); selection++) {
283 0 : if (selection->second) delete selection->second;
284 : }
285 0 : fSelections.clear();
286 :
287 0 : AliHLTSpacePointContainer::Clear(option);
288 0 : }
289 :
290 : void AliHLTTPCSpacePointContainer::Print(ostream& out, Option_t */*option*/) const
291 : {
292 : // print to stream
293 0 : out << "AliHLTTPCSpacePointContainer::Print" << endl;
294 0 : out << "n clusters: " << fClusters.size() << endl;
295 0 : for (std::map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator cl=fClusters.begin();
296 0 : cl!=fClusters.end(); cl++) {
297 0 : out << " " << cl->first << cl->second << endl;
298 : }
299 0 : }
300 :
301 : AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::SelectByMask(AliHLTUInt32_t mask, bool /*bAlloc*/) const
302 : {
303 : /// create a collection of clusters for a space point mask
304 0 : std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer);
305 0 : if (!c.get()) return NULL;
306 :
307 0 : UInt_t slice=AliHLTTPCSpacePointData::GetSlice(mask);
308 0 : UInt_t partition=AliHLTTPCSpacePointData::GetPatch(mask);
309 0 : for (std::map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator cl=fClusters.begin();
310 0 : cl!=fClusters.end(); cl++) {
311 0 : UInt_t s=AliHLTTPCSpacePointData::GetSlice(cl->first);
312 0 : UInt_t p=AliHLTTPCSpacePointData::GetPatch(cl->first);
313 0 : if ((slice>=(unsigned)AliHLTTPCGeometry::GetNSlice() || s==slice) &&
314 0 : (partition>=(unsigned)AliHLTTPCGeometry::GetNumberOfPatches() || p==partition)) {
315 0 : c->fClusters[cl->first]=cl->second;
316 0 : }
317 : }
318 0 : return c.release();
319 0 : }
320 :
321 : AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::SelectByTrack(int trackId, bool /*bAlloc*/) const
322 : {
323 : /// create a collection of clusters for a specific track
324 0 : std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer);
325 0 : if (!c.get()) return NULL;
326 :
327 0 : HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, int>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::GetTrackId,trackId));
328 0 : return c.release();
329 0 : }
330 :
331 : AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::SelectByMC(int mcId, bool /*bAlloc*/) const
332 : {
333 : /// create a collection of clusters for a specific MC track
334 0 : std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer);
335 0 : if (!c.get()) return NULL;
336 :
337 0 : HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, int>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::GetMCId,mcId));
338 0 : return c.release();
339 0 : }
340 :
341 : AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::UsedClusters(bool /*bAlloc*/) const
342 : {
343 : /// create a collection of all used clusters
344 0 : std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer);
345 0 : if (!c.get()) return NULL;
346 :
347 0 : HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, bool>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::IsUsed,true));
348 0 : return c.release();
349 0 : }
350 :
351 : AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::UnusedClusters(bool /*bAlloc*/) const
352 : {
353 : /// create a collection of all unused clusters
354 0 : std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer);
355 0 : if (!c.get()) return NULL;
356 :
357 0 : HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, bool>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::IsUsed,false));
358 0 : return c.release();
359 0 : }
360 :
361 : int AliHLTTPCSpacePointContainer::MarkUsed(const AliHLTUInt32_t* clusterIDs, int arraySize)
362 : {
363 : /// mark the clusters with specified IDs as used
364 0 : if (!clusterIDs) return -EINVAL;
365 : int iCount=0;
366 0 : for (int i=0; i<arraySize; i++) {
367 0 : if (fClusters.find(clusterIDs[i])==fClusters.end()) continue;
368 0 : fClusters[clusterIDs[i]].MarkUsed();
369 0 : iCount++;
370 0 : }
371 : return iCount;
372 0 : }
373 :
374 : int AliHLTTPCSpacePointContainer::SetTrackID(int trackID, const AliHLTUInt32_t* clusterIDs, int arraySize)
375 : {
376 : /// set track id for specified clusters
377 0 : if (!clusterIDs) return -EINVAL;
378 : int iCount=0;
379 0 : for (int i=0; i<arraySize; i++) {
380 0 : if (fClusters.find(clusterIDs[i])==fClusters.end()) continue;
381 0 : fClusters[clusterIDs[i]].SetTrackId(trackID);
382 0 : iCount++;
383 0 : }
384 : return iCount;
385 0 : }
386 :
387 : int AliHLTTPCSpacePointContainer::GetTrackID(AliHLTUInt32_t clusterID) const
388 : {
389 : /// get track id for specified cluster
390 0 : map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator element=fClusters.find(clusterID);
391 0 : if (element==fClusters.end()) return -1;
392 0 : return element->second.GetTrackId();
393 0 : }
394 :
395 : int AliHLTTPCSpacePointContainer::SetMCID(int mcID, const AliHLTUInt32_t* clusterIDs, int arraySize)
396 : {
397 : /// set mc id for specified clusters
398 0 : if (!clusterIDs) return -EINVAL;
399 : int iCount=0;
400 0 : for (int i=0; i<arraySize; i++) {
401 0 : if (fClusters.find(clusterIDs[i])==fClusters.end()) continue;
402 0 : fClusters[clusterIDs[i]].SetMCId(mcID);
403 0 : iCount++;
404 0 : }
405 : return iCount;
406 0 : }
407 :
408 : AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::AliHLTTPCSpacePointProperties()
409 0 : : fCluster(NULL)
410 0 : , fUsed(false)
411 0 : , fTrackId(-1)
412 0 : , fMCId(-1)
413 0 : {
414 : // constructor
415 0 : }
416 :
417 : AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::AliHLTTPCSpacePointProperties(const AliHLTTPCSpacePointData* pCluster)
418 0 : : fCluster(pCluster)
419 0 : , fUsed(false)
420 0 : , fTrackId(-1)
421 0 : , fMCId(-1)
422 0 : {
423 : // constructor
424 0 : }
425 :
426 : AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::AliHLTTPCSpacePointProperties(const AliHLTTPCSpacePointProperties& src)
427 0 : : fCluster(src.fCluster)
428 0 : , fUsed(src.fUsed)
429 0 : , fTrackId(src.fTrackId)
430 0 : , fMCId(src.fMCId)
431 0 : {
432 : // copy constructor
433 0 : }
434 :
435 : AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties& AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::operator=(const AliHLTTPCSpacePointProperties& src)
436 : {
437 : // assignment operator
438 0 : if (&src==this) return *this;
439 0 : fCluster=src.fCluster;
440 0 : fUsed=src.fUsed;
441 0 : fTrackId=src.fTrackId;
442 0 : fMCId=src.fMCId;
443 :
444 0 : return *this;
445 0 : }
446 :
447 : void AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::Print(ostream& out, Option_t */*option*/) const
448 : {
449 : // print info
450 0 : if (!Data()) {
451 0 : out << "no data";
452 0 : return;
453 : }
454 0 : const AliHLTTPCSpacePointData* data=Data();
455 0 : out << " " << data->fX << " " << data->fY << " " << data->fZ << " " << (UInt_t)data->fPadRow
456 0 : << " " << data->GetSigmaY2() << " " << data->GetSigmaZ2()
457 0 : << " " << data->GetCharge() << " " << data->GetQMax()
458 0 : << " " << fTrackId << " " << fMCId << " " << fUsed;
459 0 : }
460 :
461 : ostream& operator<<(ostream &out, const AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties& p)
462 : {
463 0 : p.Print(out);
464 0 : return out;
465 : }
|