]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/AliHLTTPCClusterDumpComponent.cxx
Added a new component (AliHLTTPCHWCFDataReverterComponent) which reads the data,...
[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     fDirectory(""),
39     fSlice(-1)
40 {
41   // see header file for class documentation
42   // or
43   // refer to README to build package
44   // or
45   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
46 }
47
48 AliHLTTPCClusterDumpComponent::~AliHLTTPCClusterDumpComponent()
49 {
50   // see header file for class documentation
51 }
52
53 const char* AliHLTTPCClusterDumpComponent::GetComponentID()
54 {
55   // see header file for class documentation
56   return "TPCClusterDump";
57 }
58
59 void AliHLTTPCClusterDumpComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
60 {
61   // see header file for class documentation
62   list.clear();
63   list.push_back(AliHLTTPCDefinitions::fgkClustersDataType);
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
82   TString argument="";
83   bool bMissingParam=0;
84   int i=0;
85   do {
86     if (i>=argc || (argument=argv[i]).IsNull()) continue;
87     
88     // -directory
89     if (argument.CompareTo("-directory-clusterdump")==0) {
90       if ((bMissingParam=(++i>=argc))) break;
91       fDirectory=argv[i];
92       break;
93     }
94     //-slice
95     if (argument.CompareTo("-slice")==0) {
96       if ((bMissingParam=(++i>=argc))) break;
97       TString str= argv[i];
98       fSlice=str.Atoi();
99       break;
100     }
101   }while(0);
102
103   HLTWarning("AliHLTTPCClusterDumpComponent does not have any arguments at this time");
104   return 0;
105 }
106
107 int AliHLTTPCClusterDumpComponent::CloseWriter()
108 {
109   // see header file for class documentation
110   return 0;
111 }
112
113 int AliHLTTPCClusterDumpComponent::DumpEvent( const AliHLTComponentEventData& /*evtData*/,
114                                               const AliHLTComponentBlockData* /*blocks*/, 
115                                               AliHLTComponentTriggerData& /*trigData*/ )
116 {
117   // see header file for class documentation
118
119   HLTDebug("Entering DumpEvent");
120
121   int iResult=0;
122   int blockno=0;
123   const AliHLTComponentBlockData* pDesc=NULL;
124
125   Int_t spacePointCounter=0;
126
127   //building the filename
128   fCurrentFileName="";
129   if (!fDirectory.IsNull()) {
130     fCurrentFileName+=fDirectory;
131   }
132   fCurrentFileName+="ClusterDump";
133   fCurrentFileName+=Form("_RunNo-%d",GetRunNo());
134   if(fSlice!=-1){
135     fCurrentFileName+=Form("_Slice-%d", fSlice);
136   }
137   fCurrentFileName+=Form("_Event-%d", GetEventCount());
138
139   ofstream dump;
140   dump.open(fCurrentFileName.Data());
141
142   for (pDesc=GetFirstInputBlock(AliHLTTPCDefinitions::fgkClustersDataType); pDesc!=NULL; pDesc=GetNextInputBlock(), blockno++) {
143     HLTDebug("event %Lu block %d: %s 0x%08x size %d", GetEventId(), blockno, DataType2Text(pDesc->fDataType).c_str(), pDesc->fSpecification, pDesc->fSize);
144
145     if(pDesc->fDataType!=AliHLTTPCDefinitions::fgkClustersDataType){continue;}
146  
147     if (dump.good() || 1) {//the || 1 is there since dump.good() will return false( EOF )
148        iResult=1;
149        const AliHLTTPCClusterData* clusterData = (const AliHLTTPCClusterData*) pDesc->fPtr;
150        Int_t nSpacepoints = (Int_t) clusterData->fSpacePointCnt;
151        AliHLTTPCSpacePointData *clusters = (AliHLTTPCSpacePointData*) &clusterData->fSpacePoints;
152        
153        for(int i=0;i<nSpacepoints;i++){
154          dump << "" << endl;
155          dump << "ClusterNumber: " << spacePointCounter << endl;
156          dump << "Slice:         " << (Int_t)(clusters[i].fID/10) << endl;//quick fix to get the partiion and slice numbers to the clusterdump
157          dump << "Partition:     " << (Int_t)(clusters[i].fID%10) << endl;//quick fix to get the partiion and slice numbers to the clusterdump
158          dump << "[X,Y,Z]:       [" << clusters[i].fX<<" , "<<clusters[i].fY<<" , "<<clusters[i].fZ <<"]"<< endl;
159          Float_t xyz[3]={clusters[i].fX,clusters[i].fY,clusters[i].fZ};
160          AliHLTTPCTransform::LocHLT2Raw(xyz,(Int_t)(clusters[i].fID/10),(Int_t)(clusters[i].fID%10));
161          dump << "[R,P,T]:       [" << xyz[0]<<" , "<<xyz[1]<<" , "<<xyz[2] <<"]"<< endl;
162          dump << "Total Charge:  " << clusters[i].fCharge         << endl;
163          dump << "Q Max:         " << clusters[i].fQMax           << endl;
164          spacePointCounter++;
165        }
166        
167      }
168      else {
169        HLTError("can not open file %s for writing", fCurrentFileName.Data());
170        iResult=-EBADF;
171      }
172
173   }
174   dump.close();
175   return iResult;
176 }