]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTSpacePointContainer.cxx
create 2D projections of space points, adding interface method to retrieve phi of...
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTSpacePointContainer.cxx
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   AliHLTSpacePointContainer.cxx
20 /// @author Matthias Richter
21 /// @date   2011-04-29
22 /// @brief  Base helper class for handling of HLT space point data blocks
23 ///
24
25 #include "AliHLTSpacePointContainer.h"
26 #include "AliHLTComponent.h"
27 #include "TFile.h"
28 #include "TString.h"
29 #include "TArrayC.h"
30 #include "TObjArray.h"
31 #include "TH2F.h"
32 #include "TMath.h"
33 #include <memory>
34 #include <iostream>
35
36 /** ROOT macro for the implementation of ROOT specific class methods */
37 ClassImp(AliHLTSpacePointContainer)
38
39 AliHLTSpacePointContainer::AliHLTSpacePointContainer()
40   : TObject(), AliHLTLogging()
41   , fBuffers()
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 AliHLTSpacePointContainer::AliHLTSpacePointContainer(const AliHLTSpacePointContainer& c)
51   : TObject(c), AliHLTLogging()
52   , fBuffers()
53 {
54   /// copy constructor
55
56   // currently no need to copy or allocate, the fBuffers of the source object are expected
57   // to be persistent throughout the lifetime of the new object. Only pointers are copied.
58 }
59
60 AliHLTSpacePointContainer& AliHLTSpacePointContainer::operator=(const AliHLTSpacePointContainer& c)
61 {
62   /// assignment operator
63   if (&c==this) return *this;
64
65   // currently no need to copy or allocate
66   return *this;
67 }
68
69 AliHLTSpacePointContainer::~AliHLTSpacePointContainer()
70 {
71   // destructor
72   Clear();
73 }
74
75 void AliHLTSpacePointContainer::Clear(Option_t * /*option*/)
76 {
77   // internal cleanup
78   vector<TArrayC*>::iterator element=fBuffers.begin();
79   while (element!=fBuffers.end()) {
80     TArrayC* pBuffer=*element;
81     element=fBuffers.erase(element);
82     delete pBuffer;
83   }
84 }
85
86 void AliHLTSpacePointContainer::Print(Option_t *option) const
87 {
88   // print info
89   Print(cout, option);
90 }
91
92 void AliHLTSpacePointContainer::Print(ostream& out, Option_t */*option*/) const
93 {
94   // print to stream
95   out << "AliHLTSpacePointContainer::Print" << endl;
96 }
97
98 int AliHLTSpacePointContainer::AddInputBlock(const char* filename, AliHLTComponentDataType dt, unsigned specification)
99 {
100   // open block from file and add to collection
101   if (!filename) return -EINVAL;
102   
103   TString input=filename;
104   input+="?filetype=raw";
105   std::auto_ptr<TFile> pFile(new TFile(input));
106   if (!pFile.get()) return -ENOMEM;
107   if (pFile->IsZombie()) return -ENOENT;
108
109   int iResult=0;
110   pFile->Seek(0);
111   std::auto_ptr<TArrayC> buffer(new TArrayC);
112   if (!buffer.get()) return -ENOMEM;
113
114   buffer->Set(pFile->GetSize());
115   if (pFile->ReadBuffer(buffer->GetArray(), buffer->GetSize())==0) {
116     AliHLTComponentBlockData bd;
117     AliHLTComponent::FillBlockData(bd);
118     bd.fPtr=buffer->GetArray();
119     bd.fSize=buffer->GetSize();
120     bd.fDataType=dt;
121     bd.fSpecification=specification;
122     HLTDebug("adding data block %d byte(s) from file %s", pFile->GetSize(), filename);
123     iResult=AddInputBlock(&bd);
124   } else {
125     HLTError("failed reading %d byte(s) from file %s", pFile->GetSize(), filename);
126     iResult=-ENODATA;
127   }
128
129   fBuffers.push_back(buffer.release());
130   return iResult;
131 }
132
133 int AliHLTSpacePointContainer::AddInputBlocks(const char* listfile, AliHLTComponentDataType dt)
134 {
135   // open blank separated list of files and add data
136   ifstream list(listfile);
137   if (!list.good()) return -ENOENT;
138
139   int count=0;
140   TString file;
141   while (file.ReadLine(list)) {
142     HLTInfo("adding data from file %s", file.Data());
143     int iResult=AddInputBlock(file.Data(), dt, kAliHLTVoidDataSpec);
144     if (iResult<0) {
145       HLTInfo("failed to add data from file %s: error %d", file.Data(), iResult);
146       return iResult;
147     }
148     count+=iResult;
149   }
150
151   // std::auto_ptr<TObjArray> tokens(files.Tokenize(" "));
152   // if (!tokens.get()) return 0;
153   // for (int i=0; i<tokens->GetEntriesFast(); i++) {
154   //   if (!tokens->At(i)) continue;
155   //   cout << "adding data from file " << tokens->At(i)->GetName() << endl;
156   //   AddInputBlock(tokens->At(i)->GetName(), dt, kAliHLTVoidDataSpec);
157   // }
158
159   return count;
160 }
161
162 AliHLTSpacePointContainer* AliHLTSpacePointContainer::SelectByTrack(int /*trackId*/, bool /*bAlloc*/) const
163 {
164   /// create a collection of clusters for a specific track
165   /// default implementation, nothing to do
166   return NULL;
167 }
168
169 AliHLTSpacePointContainer* AliHLTSpacePointContainer::SelectByMC(int /*mcId*/, bool /*bAlloc*/) const
170 {
171   /// create a collection of clusters for a specific MC track
172   /// default implementation, nothing to do
173   return NULL;
174 }
175
176 AliHLTSpacePointContainer* AliHLTSpacePointContainer::UsedClusters(bool /*bAlloc*/) const
177 {
178   /// create a collection of all used clusters
179   /// default implementation, nothing to do
180   return NULL;
181 }
182
183 AliHLTSpacePointContainer* AliHLTSpacePointContainer::UnusedClusters(bool /*bAlloc*/) const
184 {
185   /// create a collection of all unused clusters
186   /// default implementation, nothing to do
187   return NULL;
188 }
189
190 int AliHLTSpacePointContainer::MarkUsed(const AliHLTUInt32_t* /*clusterIDs*/, int /*arraySize*/)
191 {
192   /// default implementation, nothing to do
193   return -ENOSYS;
194 }
195
196 int AliHLTSpacePointContainer::SetTrackID(int /*trackID*/, const AliHLTUInt32_t* /*clusterIDs*/, int /*arraySize*/)
197 {
198   /// default implementation, nothing to do
199   return -ENOSYS;
200 }
201
202 int AliHLTSpacePointContainer::SetMCID(int /*mcID*/, const AliHLTUInt32_t* /*clusterIDs*/, int /*arraySize*/)
203 {
204   /// default implementation, nothing to do
205   return -ENOSYS;
206 }
207
208 TH1* AliHLTSpacePointContainer::DrawProjection(const char* plane, const vector<AliHLTUInt32_t>& selection) const
209 {
210   // draw the projection of space points in a specified plane
211   // "xy", "xz", "yz"
212   
213   int mode=0;
214   if (strcmp(plane, "xy")==0) mode=0;
215   else if (strcmp(plane, "xz")==0) mode=1;
216   else if (strcmp(plane, "yz")==0) mode=2;
217   else {
218     HLTError("invalid plane specification %s, allowed 'xy', 'xz', 'yz'", plane);
219     return NULL;
220   }
221   
222   float maxXAxis=100.0;
223   float maxYAxis=100.0;
224   vector<AliHLTUInt32_t> clusterIDs;
225   GetClusterIDs(clusterIDs);
226   vector<AliHLTUInt32_t>::iterator clusterID=clusterIDs.begin();
227   while (clusterID!=clusterIDs.end()) {
228     vector<AliHLTUInt32_t>::const_iterator element=selection.begin();
229     for (; element!=selection.end(); element++) {
230       if (((*element)&0xffff0000)==((*clusterID)&0xffff0000)) break;
231     }
232     if (element==selection.end() && selection.size()>0) {
233       clusterID=clusterIDs.erase(clusterID);
234       continue;
235     }
236     float XValue=0.0;
237     float YValue=0.0;
238     if (mode==0) {
239       XValue=GetX(*clusterID);
240       YValue=GetY(*clusterID);
241     } else if (mode==1) {
242       XValue=GetX(*clusterID);
243       YValue=GetZ(*clusterID);
244     } else if (mode==2) {
245       XValue=GetY(*clusterID);
246       YValue=GetZ(*clusterID);
247     }
248     if (maxXAxis<XValue) maxXAxis=XValue;
249     if (maxYAxis<YValue) maxYAxis=YValue;
250     clusterID++;
251   }
252   if (maxXAxis<maxYAxis) maxXAxis=maxYAxis;
253   else maxYAxis=maxXAxis;
254
255   TH2F* histogram=new TH2F(plane, plane, 1000, -maxXAxis, maxXAxis, 1000, -maxYAxis, maxYAxis);
256   if (!histogram) return NULL;
257     if (mode==0) {
258       histogram->GetXaxis()->SetTitle("X [cm]");
259       histogram->GetYaxis()->SetTitle("Y [cm]");
260     } else if (mode==1) {
261       histogram->GetXaxis()->SetTitle("Z [cm]");
262       histogram->GetYaxis()->SetTitle("X [cm]");
263     } else if (mode==2) {
264       histogram->GetXaxis()->SetTitle("Z [cm]");
265       histogram->GetYaxis()->SetTitle("Y [cm]");
266     }
267
268   for (clusterID=clusterIDs.begin(); clusterID!=clusterIDs.end(); clusterID++) {
269     float XValue=GetX(*clusterID);
270     float YValue=GetY(*clusterID);
271     float ZValue=GetZ(*clusterID);
272     Float_t phi     = GetPhi(*clusterID);
273     Float_t cos     = TMath::Cos( phi );
274     Float_t sin     = TMath::Sin( phi );
275     if (mode==0) {
276       histogram->SetTitle("XY projection");
277       histogram->Fill(cos*XValue - sin*YValue, sin*XValue + cos*YValue);
278     } else if (mode==1) {
279       // should be maybe 'ZX' but somehow 'XZ' is more 'lingual convenient'
280       histogram->SetTitle("XZ projection");
281       histogram->Fill(ZValue, cos*XValue - sin*YValue);
282     } else if (mode==2) {
283       // same comment
284       histogram->SetTitle("YZ projection");
285       histogram->Fill(ZValue, sin*XValue + cos*YValue);
286     }
287   }
288
289   return histogram;
290 }
291
292 ostream& operator<<(ostream &out, const AliHLTSpacePointContainer& c)
293 {
294   c.Print(out);
295   return out;
296 }