]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCClusterDumpComponent.cxx
Cleanup: hardcoded bit operations with TPC cluster ID's removed. Instead static...
[u/mrichter/AliRoot.git] / HLT / TPCLib / AliHLTTPCClusterDumpComponent.cxx
1 // $Id$
2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project        * 
4 //* ALICE Experiment at CERN, All rights reserved.                         *
5 //*                                                                        *
6 //* Primary Authors: Kenneth Aamodt <Kenneth.Aamodt@cern.ch>               *
7 //*                  for The ALICE HLT Project.                            *
8 //*                                                                        *
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 //**************************************************************************
17
18 /** @file   AliHLTTPCClusterDumpComponent.cxx
19     @author Kenneth Aamodt
20     @date   
21     @brief  Special file writer converting TPC clusters input to readable
22             ASCII format. 
23 */
24
25 #include <cassert>
26 #include "AliHLTTPCClusterDumpComponent.h"
27 #include "AliHLTTPCDefinitions.h"
28 #include "AliHLTTPCSpacePointData.h"
29 #include "AliHLTTPCClusterDataFormat.h"
30 #include "AliHLTTPCTransform.h"
31
32 /** ROOT macro for the implementation of ROOT specific class methods */
33 ClassImp(AliHLTTPCClusterDumpComponent)
34
35   AliHLTTPCClusterDumpComponent::AliHLTTPCClusterDumpComponent()
36     :
37     AliHLTFileWriter()
38 //fSlice(-1)
39 {
40   // see header file for class documentation
41   // or
42   // refer to README to build package
43   // or
44   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
45 }
46
47 AliHLTTPCClusterDumpComponent::~AliHLTTPCClusterDumpComponent()
48 {
49   // see header file for class documentation
50 }
51
52 const char* AliHLTTPCClusterDumpComponent::GetComponentID()
53 {
54   // see header file for class documentation
55   return "TPCClusterDump";
56 }
57
58 void AliHLTTPCClusterDumpComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
59 {
60   // see header file for class documentation
61   list.clear();
62   list.push_back(AliHLTTPCDefinitions::fgkClustersDataType);
63   list.push_back(AliHLTTPCDefinitions::fgkAlterClustersDataType);
64 }
65
66 AliHLTComponent* AliHLTTPCClusterDumpComponent::Spawn()
67 {
68   // see header file for class documentation
69   return new AliHLTTPCClusterDumpComponent;
70 }
71
72 int AliHLTTPCClusterDumpComponent::InitWriter()
73 {
74   // see header file for class documentation
75   return 0;
76 }
77
78 int AliHLTTPCClusterDumpComponent::ScanArgument(int /*argc*/, const char** /*argv*/)
79 {
80   // see header file for class documentation
81   int iResult=0;
82   TString argument="";
83   bool bMissingParam=0;
84   int i=0;
85   
86   if (bMissingParam) iResult=-EPROTO;
87   else if (iResult>=0) iResult=i;
88   
89   return iResult;
90 }
91
92 int AliHLTTPCClusterDumpComponent::CloseWriter()
93 {
94   // see header file for class documentation
95   return 0;
96 }
97
98 int AliHLTTPCClusterDumpComponent::DumpEvent( const AliHLTComponentEventData& evtData,
99                                               AliHLTComponentTriggerData& /*trigData*/ )
100 {
101   // see header file for class documentation
102
103   HLTDebug("Entering DumpEvent");
104
105   int iResult=0;
106   int blockno=0;
107   const AliHLTComponentBlockData* pDesc=NULL;
108
109   if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) )
110     return 0;
111
112   Int_t spacePointCounter=0;
113   
114   for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); pDesc!=NULL; pDesc=GetNextInputBlock(), blockno++) {
115     TString filename;
116     iResult=BuildFileName(evtData.fEventID, 0, pDesc->fDataType, 0, filename);
117     ios::openmode filemode=(ios::openmode)0;
118     if (fCurrentFileName.CompareTo(filename)==0) {
119       filemode=ios::app;
120     } else {
121       fCurrentFileName=filename;
122     }
123     
124     if (iResult>=0) {
125       ofstream dump(fCurrentFileName.Data(), filemode);
126       if (dump.good()) {
127         
128         if(pDesc->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType){continue;}
129         
130         //if (dump.good() || 1) {//the || 1 is there since dump.good() will return false( EOF )
131         iResult=1;
132         const AliHLTTPCClusterData* clusterData = (const AliHLTTPCClusterData*) pDesc->fPtr;
133         Int_t nSpacepoints = (Int_t) clusterData->fSpacePointCnt;
134         AliHLTTPCSpacePointData *clusters = (AliHLTTPCSpacePointData*) &clusterData->fSpacePoints;
135         
136         for(int i=0;i<nSpacepoints;i++){
137           UInt_t idCluster = clusters[i].fID;
138           Int_t slice = AliHLTTPCSpacePointData::GetSlice(idCluster);
139           Int_t patch = AliHLTTPCSpacePointData::GetPatch(idCluster);
140           
141           dump << "" << endl;
142           dump << "ClusterNumber: " << spacePointCounter << endl;
143           dump << "Slice:         " << slice << endl;
144           dump << "Partition:     " << patch << endl;
145           dump << "[X,Y,Z]:       [" << clusters[i].fX<<" , "<<clusters[i].fY<<" , "<<clusters[i].fZ <<"]"<< endl;
146           //Float_t xyz[3]={clusters[i].fX,clusters[i].fY,clusters[i].fZ};
147           //AliHLTTPCTransform::LocHLT2Raw(xyz,(Int_t)(clusters[i].fID/10),(Int_t)(clusters[i].fID%10));
148           //dump << "[R,P,T]:       [" << xyz[0]<<" , "<<xyz[1]<<" , "<<xyz[2] <<"]"<< endl;
149           dump << "Total Charge:  " << clusters[i].fCharge         << endl;
150           dump << "Q Max:         " << clusters[i].fQMax           << endl;
151           spacePointCounter++;
152         }
153       }
154       else {
155         HLTError("can not open file %s for writing", fCurrentFileName.Data());
156         iResult=-EBADF;
157       }
158       dump.close();
159     }
160   }
161   
162   for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkAlterClustersDataType); pDesc!=NULL; pDesc=GetNextInputBlock(), blockno++) {
163     TString filename;
164     iResult=BuildFileName(evtData.fEventID, 0, pDesc->fDataType, 0, filename);
165     ios::openmode filemode=(ios::openmode)0;
166     if (fCurrentFileName.CompareTo(filename)==0) {
167       filemode=ios::app;
168     } else {
169       fCurrentFileName=filename;
170     }
171     
172     if (iResult>=0) {
173       ofstream dump(fCurrentFileName.Data(), filemode);
174       if (dump.good()) {
175         
176         if(pDesc->fDataType!=AliHLTTPCDefinitions::fgkAlterClustersDataType){continue;}
177         
178         //if (dump.good() || 1) {//the || 1 is there since dump.good() will return false( EOF )
179         iResult=1;
180         const AliHLTTPCClusterData* clusterData = (const AliHLTTPCClusterData*) pDesc->fPtr;
181         Int_t nSpacepoints = (Int_t) clusterData->fSpacePointCnt;
182         AliHLTTPCSpacePointData *clusters = (AliHLTTPCSpacePointData*) &clusterData->fSpacePoints;
183         
184         for(int i=0;i<nSpacepoints;i++){
185           UInt_t idCluster = clusters[i].fID;
186           Int_t slice = AliHLTTPCSpacePointData::GetSlice(idCluster);
187           Int_t patch = AliHLTTPCSpacePointData::GetPatch(idCluster);
188           
189           dump << "" << endl;
190           dump << "ClusterNumber: " << spacePointCounter << endl;
191           dump << "Slice:         " << slice << endl;
192           dump << "Partition:     " << patch << endl;
193           dump << "[X,Y,Z]:       [" << clusters[i].fX<<" , "<<clusters[i].fY<<" , "<<clusters[i].fZ <<"]"<< endl;
194           //Float_t xyz[3]={clusters[i].fX,clusters[i].fY,clusters[i].fZ};
195           //AliHLTTPCTransform::LocHLT2Raw(xyz,(Int_t)(clusters[i].fID/10),(Int_t)(clusters[i].fID%10));
196           //dump << "[R,P,T]:       [" << xyz[0]<<" , "<<xyz[1]<<" , "<<xyz[2] <<"]"<< endl;
197           dump << "Total Charge:  " << clusters[i].fCharge         << endl;
198           dump << "Q Max:         " << clusters[i].fQMax           << endl;
199           spacePointCounter++;
200         }
201       }
202       else {
203         HLTError("can not open file %s for writing", fCurrentFileName.Data());
204         iResult=-EBADF;
205       }
206       dump.close();
207     }
208   }
209   return iResult;
210 }