]>
Commit | Line | Data |
---|---|---|
7e7b2c34 | 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" | |
c7585a2a | 28 | #include "AliHLTTPCTransform.h" |
7e7b2c34 | 29 | #include "AliHLTComponent.h" |
30 | #include "AliHLTTemplates.h" | |
b39a860d | 31 | #include "TMath.h" |
7e7b2c34 | 32 | #include <memory> |
33 | #include <algorithm> | |
34 | ||
35 | /** ROOT macro for the implementation of ROOT specific class methods */ | |
36 | ClassImp(AliHLTTPCSpacePointContainer) | |
37 | ||
38 | AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointContainer() | |
39 | : AliHLTSpacePointContainer() | |
40 | , fClusters() | |
9003c201 | 41 | , fSelections() |
7e7b2c34 | 42 | { |
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 | } | |
49 | ||
50 | AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointContainer(const AliHLTTPCSpacePointContainer& c) | |
51 | : AliHLTSpacePointContainer(c) | |
52 | , fClusters(c.fClusters.begin(), c.fClusters.end()) | |
9003c201 | 53 | , fSelections() |
7e7b2c34 | 54 | { |
55 | /// copy constructor | |
56 | } | |
57 | ||
58 | AliHLTTPCSpacePointContainer& AliHLTTPCSpacePointContainer::operator=(const AliHLTTPCSpacePointContainer& c) | |
59 | { | |
60 | /// assignment operator | |
61 | if (&c==this) return *this; | |
62 | AliHLTSpacePointContainer::operator=(c); | |
63 | fClusters=c.fClusters; | |
64 | ||
65 | return *this; | |
66 | } | |
67 | ||
68 | AliHLTTPCSpacePointContainer::~AliHLTTPCSpacePointContainer() | |
69 | { | |
70 | // destructor | |
9003c201 | 71 | |
5ecd199d | 72 | Clear(); |
7e7b2c34 | 73 | } |
74 | ||
75 | int AliHLTTPCSpacePointContainer::AddInputBlock(const AliHLTComponentBlockData* pDesc) | |
76 | { | |
77 | // add input block to the collection | |
78 | if (!pDesc) return -EINVAL; | |
79 | int count=0; | |
80 | if (pDesc->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType) { | |
81 | HLTWarning("ignoring data block of type %s", AliHLTComponent::DataType2Text(pDesc->fDataType).c_str()); | |
82 | return 0; | |
83 | } | |
84 | if (!pDesc->fPtr || pDesc->fSize<sizeof(AliHLTTPCClusterData)) return -ENODATA; | |
85 | ||
86 | // consistency check of the input block | |
87 | const AliHLTTPCClusterData* pClusterData=reinterpret_cast<AliHLTTPCClusterData*>(pDesc->fPtr); | |
88 | if (pClusterData->fSpacePointCnt*sizeof(AliHLTTPCSpacePointData)+sizeof(AliHLTTPCClusterData)!=pDesc->fSize) { | |
89 | 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 | return -EBADMSG; | |
92 | } | |
93 | ||
94 | AliHLTUInt8_t minslice = AliHLTTPCDefinitions::GetMinSliceNr( pDesc->fSpecification ); | |
95 | AliHLTUInt8_t maxslice = AliHLTTPCDefinitions::GetMaxSliceNr( pDesc->fSpecification ); | |
96 | AliHLTUInt8_t minpart = AliHLTTPCDefinitions::GetMinPatchNr( pDesc->fSpecification ); | |
97 | AliHLTUInt8_t maxpart = AliHLTTPCDefinitions::GetMaxPatchNr( pDesc->fSpecification ); | |
98 | bool bIsSinglePartition=(pDesc->fSpecification==kAliHLTVoidDataSpec?false:minslice==maxslice && minpart==maxpart); | |
99 | ||
100 | for (UInt_t i=0; i<pClusterData->fSpacePointCnt; i++) { | |
101 | AliHLTUInt32_t clusterID=~(AliHLTUInt32_t)0; | |
102 | if (bIsSinglePartition) { | |
103 | // cluster ID has to match ID encoded from slice, partition and index | |
104 | clusterID=AliHLTTPCSpacePointData::GetID(minslice, minpart, i); | |
105 | if (clusterID!=pClusterData->fSpacePoints[i].fID) { | |
106 | 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 | clusterID=pClusterData->fSpacePoints[i].fID; | |
112 | UInt_t clusterSlice =AliHLTTPCSpacePointData::GetSlice(clusterID); | |
113 | UInt_t clusterPart =AliHLTTPCSpacePointData::GetPatch(clusterID); | |
114 | UInt_t clusterNumber=AliHLTTPCSpacePointData::GetNumber(clusterID); | |
115 | if (pDesc->fSpecification!=kAliHLTVoidDataSpec) { | |
116 | if (clusterSlice<minslice || clusterSlice>maxslice || | |
117 | clusterPart<minpart || clusterPart>maxpart) { | |
118 | HLTWarning("cluster ID 0x%08d out of range indicated by data specification 0x%08x", clusterID, pDesc->fSpecification); | |
119 | } else if (clusterNumber!=i) { | |
120 | HLTWarning("cluster ID 0x%08d does not match index %d in block 0x%08x", clusterID, i, pDesc->fSpecification); | |
121 | } | |
122 | } | |
123 | } | |
c7585a2a | 124 | { |
9003c201 | 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=AliHLTTPCTransform::GetPadRow(pClusterData->fSpacePoints[i].fX); | |
129 | // if (row<AliHLTTPCTransform::GetFirstRow(clusterPart) || row>AliHLTTPCTransform::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 | // } | |
c7585a2a | 133 | } |
134 | ||
7e7b2c34 | 135 | if (fClusters.find(clusterID)==fClusters.end()) { |
136 | // new cluster | |
137 | fClusters[clusterID]=AliHLTTPCSpacePointProperties(&pClusterData->fSpacePoints[i]); | |
138 | count++; | |
139 | } else { | |
140 | HLTError("cluster with ID 0x%08x already existing, skipping cluster %d of data block 0x%08x", | |
141 | clusterID, i, pDesc->fSpecification); | |
142 | } | |
143 | } | |
144 | ||
145 | return count; | |
146 | } | |
147 | ||
148 | int AliHLTTPCSpacePointContainer::GetClusterIDs(vector<AliHLTUInt32_t>& tgt) const | |
149 | { | |
150 | // get array of cluster IDs | |
151 | tgt.clear(); | |
152 | transform(fClusters.begin(), fClusters.end(), back_inserter(tgt), HLT::AliGetKey()); | |
153 | return tgt.size(); | |
154 | } | |
155 | ||
9003c201 | 156 | bool AliHLTTPCSpacePointContainer::Check(AliHLTUInt32_t clusterID) const |
157 | { | |
158 | // check if the cluster is available | |
159 | 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 | if (fSelections.find(mask)!=fSelections.end()) { | |
166 | // return existing selection | |
167 | return fSelections.find(mask)->second; | |
168 | } | |
169 | // create new collection | |
170 | vector<AliHLTUInt32_t>* selected=new vector<AliHLTUInt32_t>; | |
171 | if (!selected) return NULL; | |
172 | UInt_t slice=AliHLTTPCSpacePointData::GetSlice(mask); | |
173 | UInt_t partition=AliHLTTPCSpacePointData::GetPatch(mask); | |
174 | HLTInfo("creating collection 0x%08x", mask); | |
175 | for (std::map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator cl=fClusters.begin(); | |
176 | cl!=fClusters.end(); cl++) { | |
177 | UInt_t s=AliHLTTPCSpacePointData::GetSlice(cl->first); | |
178 | UInt_t p=AliHLTTPCSpacePointData::GetPatch(cl->first); | |
179 | if ((slice>=(unsigned)AliHLTTPCTransform::GetNSlice() || s==slice) && | |
180 | (partition>=(unsigned)AliHLTTPCTransform::GetNumberOfPatches() || p==partition)) { | |
181 | selected->push_back(cl->first); | |
182 | } | |
183 | } | |
184 | HLTInfo("collection 0x%08x with %d spacepoints", mask, selected->size()); | |
185 | fSelections[mask]=selected; | |
186 | return selected; | |
187 | } | |
188 | ||
7e7b2c34 | 189 | float AliHLTTPCSpacePointContainer::GetX(AliHLTUInt32_t clusterID) const |
190 | { | |
191 | // get X coordinate | |
192 | if (fClusters.find(clusterID)==fClusters.end() || | |
193 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
c7585a2a | 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 | return AliHLTTPCTransform::Row2X(fClusters.find(clusterID)->second.Data()->fPadRow); | |
7e7b2c34 | 200 | } |
201 | ||
202 | float AliHLTTPCSpacePointContainer::GetXWidth(AliHLTUInt32_t clusterID) const | |
203 | { | |
204 | // get error for X coordinate | |
205 | if (fClusters.find(clusterID)==fClusters.end() || | |
206 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
207 | return 0.0; // fixed in padrow number | |
208 | } | |
209 | ||
210 | float AliHLTTPCSpacePointContainer::GetY(AliHLTUInt32_t clusterID) const | |
211 | { | |
212 | // get Y coordinate | |
213 | if (fClusters.find(clusterID)==fClusters.end() || | |
214 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
215 | return fClusters.find(clusterID)->second.Data()->fY; | |
216 | } | |
217 | ||
218 | float AliHLTTPCSpacePointContainer::GetYWidth(AliHLTUInt32_t clusterID) const | |
219 | { | |
220 | // get error for Y coordinate | |
221 | if (fClusters.find(clusterID)==fClusters.end() || | |
222 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
223 | return fClusters.find(clusterID)->second.Data()->fSigmaY2; | |
224 | } | |
225 | ||
226 | float AliHLTTPCSpacePointContainer::GetZ(AliHLTUInt32_t clusterID) const | |
227 | { | |
228 | // get Z coordinate | |
229 | if (fClusters.find(clusterID)==fClusters.end() || | |
230 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
231 | return fClusters.find(clusterID)->second.Data()->fZ; | |
232 | } | |
233 | ||
234 | float AliHLTTPCSpacePointContainer::GetZWidth(AliHLTUInt32_t clusterID) const | |
235 | { | |
236 | // get error for Z coordinate | |
237 | if (fClusters.find(clusterID)==fClusters.end() || | |
238 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
239 | return fClusters.find(clusterID)->second.Data()->fSigmaZ2; | |
240 | } | |
241 | ||
242 | float AliHLTTPCSpacePointContainer::GetCharge(AliHLTUInt32_t clusterID) const | |
243 | { | |
244 | // get charge | |
245 | if (fClusters.find(clusterID)==fClusters.end() || | |
246 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
247 | return fClusters.find(clusterID)->second.Data()->fCharge; | |
248 | } | |
249 | ||
5ecd199d | 250 | float AliHLTTPCSpacePointContainer::GetMaxSignal(AliHLTUInt32_t clusterID) const |
251 | { | |
252 | // get charge | |
253 | if (fClusters.find(clusterID)==fClusters.end() || | |
254 | fClusters.find(clusterID)->second.Data()==NULL) return 0.0; | |
255 | return fClusters.find(clusterID)->second.Data()->fQMax; | |
256 | } | |
257 | ||
b39a860d | 258 | float AliHLTTPCSpacePointContainer::GetPhi(AliHLTUInt32_t clusterID) const |
259 | { | |
260 | // get charge | |
c7585a2a | 261 | |
262 | // phi can be derived directly from the id, no need to search | |
263 | // for existing cluster | |
264 | int slice=AliHLTTPCSpacePointData::GetSlice(clusterID); | |
b39a860d | 265 | return ( slice + 0.5 ) * TMath::Pi() / 9.0; |
266 | } | |
267 | ||
5ecd199d | 268 | void AliHLTTPCSpacePointContainer::Clear(Option_t * option) |
7e7b2c34 | 269 | { |
270 | // clear the object and reset pointer references | |
5ecd199d | 271 | fClusters.clear(); |
272 | ||
273 | for (std::map<AliHLTUInt32_t, vector<AliHLTUInt32_t>*>::iterator selection=fSelections.begin(); | |
274 | selection!=fSelections.end(); selection++) { | |
275 | if (selection->second) delete selection->second; | |
276 | } | |
277 | fSelections.clear(); | |
278 | ||
279 | AliHLTSpacePointContainer::Clear(option); | |
7e7b2c34 | 280 | } |
281 | ||
282 | void AliHLTTPCSpacePointContainer::Print(ostream& out, Option_t */*option*/) const | |
283 | { | |
284 | // print to stream | |
285 | out << "AliHLTTPCSpacePointContainer::Print" << endl; | |
286 | out << "n clusters: " << fClusters.size() << endl; | |
287 | for (std::map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator cl=fClusters.begin(); | |
288 | cl!=fClusters.end(); cl++) { | |
289 | out << " " << cl->first << cl->second << endl; | |
290 | } | |
291 | } | |
292 | ||
9003c201 | 293 | AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::SelectByMask(AliHLTUInt32_t mask, bool /*bAlloc*/) const |
294 | { | |
295 | /// create a collection of clusters for a space point mask | |
296 | std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer); | |
297 | if (!c.get()) return NULL; | |
298 | ||
299 | UInt_t slice=AliHLTTPCSpacePointData::GetSlice(mask); | |
300 | UInt_t partition=AliHLTTPCSpacePointData::GetPatch(mask); | |
301 | for (std::map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator cl=fClusters.begin(); | |
302 | cl!=fClusters.end(); cl++) { | |
303 | UInt_t s=AliHLTTPCSpacePointData::GetSlice(cl->first); | |
304 | UInt_t p=AliHLTTPCSpacePointData::GetPatch(cl->first); | |
305 | if ((slice>=(unsigned)AliHLTTPCTransform::GetNSlice() || s==slice) && | |
306 | (partition>=(unsigned)AliHLTTPCTransform::GetNumberOfPatches() || p==partition)) { | |
307 | c->fClusters[cl->first]=cl->second; | |
308 | } | |
309 | } | |
310 | return c.release(); | |
311 | } | |
312 | ||
7e7b2c34 | 313 | AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::SelectByTrack(int trackId, bool /*bAlloc*/) const |
314 | { | |
315 | /// create a collection of clusters for a specific track | |
316 | std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer); | |
317 | if (!c.get()) return NULL; | |
318 | ||
319 | HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, int>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::GetTrackId,trackId)); | |
320 | return c.release(); | |
321 | } | |
322 | ||
323 | AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::SelectByMC(int mcId, bool /*bAlloc*/) const | |
324 | { | |
325 | /// create a collection of clusters for a specific MC track | |
326 | std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer); | |
327 | if (!c.get()) return NULL; | |
328 | ||
329 | HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, int>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::GetMCId,mcId)); | |
330 | return c.release(); | |
331 | } | |
332 | ||
333 | AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::UsedClusters(bool /*bAlloc*/) const | |
334 | { | |
335 | /// create a collection of all used clusters | |
336 | std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer); | |
337 | if (!c.get()) return NULL; | |
338 | ||
339 | HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, bool>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::IsUsed,true)); | |
340 | return c.release(); | |
341 | } | |
342 | ||
343 | AliHLTSpacePointContainer* AliHLTTPCSpacePointContainer::UnusedClusters(bool /*bAlloc*/) const | |
344 | { | |
345 | /// create a collection of all unused clusters | |
346 | std::auto_ptr<AliHLTTPCSpacePointContainer> c(new AliHLTTPCSpacePointContainer); | |
347 | if (!c.get()) return NULL; | |
348 | ||
349 | HLT::copy_map_if(fClusters.begin(), fClusters.end(), c->fClusters, HLT::AliHLTUnaryPredicate<AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties, bool>(&AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::IsUsed,false)); | |
350 | return c.release(); | |
351 | } | |
352 | ||
353 | int AliHLTTPCSpacePointContainer::MarkUsed(const AliHLTUInt32_t* clusterIDs, int arraySize) | |
354 | { | |
355 | /// mark the clusters with specified IDs as used | |
356 | if (!clusterIDs) return -EINVAL; | |
357 | int iCount=0; | |
358 | for (int i=0; i<arraySize; i++) { | |
359 | if (fClusters.find(clusterIDs[i])==fClusters.end()) continue; | |
360 | fClusters[clusterIDs[i]].MarkUsed(); | |
361 | iCount++; | |
362 | } | |
363 | return iCount; | |
364 | } | |
365 | ||
366 | int AliHLTTPCSpacePointContainer::SetTrackID(int trackID, const AliHLTUInt32_t* clusterIDs, int arraySize) | |
367 | { | |
368 | /// set track id for specified clusters | |
369 | if (!clusterIDs) return -EINVAL; | |
370 | int iCount=0; | |
371 | for (int i=0; i<arraySize; i++) { | |
372 | if (fClusters.find(clusterIDs[i])==fClusters.end()) continue; | |
373 | fClusters[clusterIDs[i]].SetTrackId(trackID); | |
374 | iCount++; | |
375 | } | |
376 | return iCount; | |
377 | } | |
378 | ||
c7585a2a | 379 | int AliHLTTPCSpacePointContainer::GetTrackID(AliHLTUInt32_t clusterID) const |
380 | { | |
381 | /// get track id for specified cluster | |
382 | map<AliHLTUInt32_t, AliHLTTPCSpacePointProperties>::const_iterator element=fClusters.find(clusterID); | |
383 | if (element==fClusters.end()) return -1; | |
384 | return element->second.GetTrackId(); | |
385 | } | |
386 | ||
7e7b2c34 | 387 | int AliHLTTPCSpacePointContainer::SetMCID(int mcID, const AliHLTUInt32_t* clusterIDs, int arraySize) |
388 | { | |
389 | /// set mc id for specified clusters | |
390 | if (!clusterIDs) return -EINVAL; | |
391 | int iCount=0; | |
392 | for (int i=0; i<arraySize; i++) { | |
393 | if (fClusters.find(clusterIDs[i])==fClusters.end()) continue; | |
394 | fClusters[clusterIDs[i]].SetMCId(mcID); | |
395 | iCount++; | |
396 | } | |
397 | return iCount; | |
398 | } | |
399 | ||
400 | AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::AliHLTTPCSpacePointProperties() | |
401 | : fCluster(NULL) | |
402 | , fUsed(false) | |
403 | , fTrackId(-1) | |
404 | , fMCId(-1) | |
405 | { | |
406 | // constructor | |
407 | } | |
408 | ||
409 | AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::AliHLTTPCSpacePointProperties(const AliHLTTPCSpacePointData* pCluster) | |
410 | : fCluster(pCluster) | |
411 | , fUsed(false) | |
412 | , fTrackId(-1) | |
413 | , fMCId(-1) | |
414 | { | |
415 | // constructor | |
416 | } | |
417 | ||
418 | AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::AliHLTTPCSpacePointProperties(const AliHLTTPCSpacePointProperties& src) | |
419 | : fCluster(src.fCluster) | |
420 | , fUsed(src.fUsed) | |
421 | , fTrackId(src.fTrackId) | |
422 | , fMCId(src.fMCId) | |
423 | { | |
424 | // copy constructor | |
425 | } | |
426 | ||
5795ff9b | 427 | AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties& AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::operator=(const AliHLTTPCSpacePointProperties& src) |
7e7b2c34 | 428 | { |
429 | // assignment operator | |
430 | if (&src==this) return *this; | |
431 | fCluster=src.fCluster; | |
432 | fUsed=src.fUsed; | |
433 | fTrackId=src.fTrackId; | |
434 | fMCId=src.fMCId; | |
435 | ||
436 | return *this; | |
437 | } | |
438 | ||
439 | void AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties::Print(ostream& out, Option_t */*option*/) const | |
440 | { | |
441 | // print info | |
442 | if (!Data()) { | |
443 | out << "no data"; | |
444 | return; | |
445 | } | |
446 | const AliHLTTPCSpacePointData* data=Data(); | |
447 | out << " " << data->fX << " " << data->fY << " " << data->fZ << " " << (UInt_t)data->fPadRow | |
5ecd199d | 448 | << " " << data->GetSigmaY2() << " " << data->GetSigmaZ2() |
449 | << " " << data->GetCharge() << " " << data->GetQMax() | |
7e7b2c34 | 450 | << " " << fTrackId << " " << fMCId << " " << fUsed; |
451 | } | |
452 | ||
453 | ostream& operator<<(ostream &out, const AliHLTTPCSpacePointContainer::AliHLTTPCSpacePointProperties& p) | |
454 | { | |
455 | p.Print(out); | |
456 | return out; | |
457 | } |