1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // this program performs local monitoring on a GDC by running the HLT code //
22 // If an argument is given, this is taken as the name of a date file which //
23 // is used instead of the local node. //
24 // The program can be stopped by pressing CTRL-C. //
26 ///////////////////////////////////////////////////////////////////////////////
29 #include <TSysEvtHandler.h>
34 #include "AliRawReaderDate.h"
37 #include <AliHLTStandardIncludes.h>
38 #include <AliHLTTransform.h>
39 #include <AliHLTMemHandler.h>
40 #include <AliHLTTrackArray.h>
41 #include <AliHLTHoughMaxFinder.h>
42 #include <AliHLTHoughBaseTransformer.h>
43 #include <AliHLTHough.h>
44 #include <AliHLTBenchmark.h>
45 #include <AliKalmanTrack.h>
46 #include "AliITSgeomTGeo.h"
47 #include "AliITSgeom.h"
49 #include "AliMagFMaps.h"
50 #include <AliHLTITSclusterer.h>
51 #include <AliHLTITSVertexerZ.h>
52 #include <AliHLTITStracker.h>
55 //_____________________________________________________________________________
56 class AliGDCInterruptHandler : public TSignalHandler {
58 AliGDCInterruptHandler();
59 Bool_t Notify() {fStop = kTRUE; return kTRUE;};
60 Bool_t Stop() const {return fStop;};
62 Bool_t fStop; // CTRL-C pressed
65 //_____________________________________________________________________________
66 AliGDCInterruptHandler::AliGDCInterruptHandler() :
67 TSignalHandler(kSigInterrupt, kFALSE)
73 //_____________________________________________________________________________
75 int main(int argc, char** argv)
77 // set ROOT in batch mode
81 FILE* file = fopen("monitorGDC.log", "w");
85 // get data from a file or online from this node
88 status = monitorSetDataSource(argv[1]);
90 status = monitorSetDataSource(":");
92 if (status) ::Fatal("monitorSetDataSource", monitorDecodeError(status));
94 // monitor only a sample of physics events
95 char* table[] = {"Physics event", "yes", NULL};
96 status = monitorDeclareTable(table);
97 if (status) ::Fatal("monitorDeclareTable", monitorDecodeError(status));
99 // declare this monitoring program to DATE
100 status = monitorDeclareMp("GDC physics monitoring");
101 if (status) ::Fatal("monitorDeclareMp", monitorDecodeError(status));
103 // initialize HLT transformations
104 if (!AliHLTTransform::Init("./", kFALSE)) {
105 ::Fatal("AliHLTTransform::Init", "HLT initialization failed");
107 AliESD *esd = new AliESD;
108 // AliKalmanTrack::SetConvConst(
109 // 1000/0.299792458/AliHLTTransform::GetSolenoidField()
111 AliITSgeom *geom = new AliITSgeom();
112 geom->ReadNewFile("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymmFMD.det");
115 switch ((Int_t)(AliHLTTransform::GetSolenoidField()+0.5)) {
117 sfield = AliMagFMaps::k2kG;
120 sfield = AliMagFMaps::k4kG;
123 sfield = AliMagFMaps::k5kG;
126 ::Fatal("AliHLTTransform::GetSolenoidField", "Incorrect magnetic field");
128 AliMagF* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., sfield);
129 AliTracker::SetFieldMap(field,kTRUE);
134 // create the signal handler
135 AliGDCInterruptHandler* handler = new AliGDCInterruptHandler;
136 gSystem->AddSignalHandler(handler);
140 while (!handler->Stop()) {
141 // get the next event
142 status = monitorGetEventDynamic(&ptr);
143 if (status == (Int_t)MON_ERR_EOF) break;
144 if (status) ::Fatal("monitorGetEventDynamic", monitorDecodeError(status));
148 gSystem->Sleep(1000); // sleep for 1 second
152 AliRawReaderDate rawReader(ptr);
153 // if ((rawReader.GetAttributes()[0] & 0x02) != 0) {
155 // Int_t errorCode = rawReader.CheckData();
157 if (errorCode && (errorCode != AliRawReader::kErrSize)) {
159 if (file) fprintf(file, "%s\n", time.AsString());
160 if (file) fprintf(file, "run: %d event: %d %d\n",
161 rawReader.GetRunNumber(),
162 rawReader.GetEventId()[0],
163 rawReader.GetEventId()[1]);
164 fprintf(file, "ERROR: %d\n\n", errorCode);
168 AliHLTBenchmark *fBenchmark = new AliHLTBenchmark();
169 fBenchmark->Start("Overall timing");
171 // ITS clusterer and vertexer
172 fBenchmark->Start("ITS Clusterer");
173 AliHLTITSclusterer clusterer(0);
174 AliRawReader *itsrawreader=new AliRawReaderDate(ptr);
175 TTree* treeClusters = new TTree("TreeL3ITSclusters"," "); //make a tree
176 clusterer.Digits2Clusters(itsrawreader,treeClusters);
177 fBenchmark->Stop("ITS Clusterer");
179 AliHLTITSVertexerZ vertexer;
180 AliESDVertex *vertex = vertexer.FindVertexForCurrentEvent(geom,treeClusters);
182 Double_t vtxErr[3]={0.005,0.005,0.010};
183 vertex->GetXYZ(vtxPos);
184 // vertex->GetSigmaXYZ(vtxErr);
185 esd->SetVertex(vertex);
187 // TPC Hough reconstruction
188 Float_t ptmin = 0.1*AliHLTTransform::GetSolenoidField();
189 Float_t zvertex = vtxPos[2];
191 // Run the Hough Transformer
192 fBenchmark->Start("Init");
193 AliHLTHough *hough1 = new AliHLTHough();
195 hough1->SetThreshold(4);
196 hough1->CalcTransformerParams(ptmin);
197 hough1->SetPeakThreshold(70,-1);
198 // printf("Pointer is %x\n",ptr);
199 hough1->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
200 hough1->SetAddHistograms();
201 fBenchmark->Stop("Init");
203 fBenchmark->Start("Init");
204 AliHLTHough *hough2 = new AliHLTHough();
206 hough2->SetThreshold(4);
207 hough2->CalcTransformerParams(ptmin);
208 hough2->SetPeakThreshold(70,-1);
209 // printf("Pointer is %x\n",ptr);
210 hough2->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
211 hough2->SetAddHistograms();
212 fBenchmark->Stop("Init");
214 Int_t nglobaltracks = 0;
216 hough1->StartProcessInThread(0,17);
217 hough2->StartProcessInThread(18,35);
219 // gSystem->Sleep(20000);
220 if(hough1->WaitForThreadFinish())
221 ::Fatal("AliHLTHough::WaitForThreadFinish"," Can not join the required thread! ");
222 if(hough2->WaitForThreadFinish())
223 ::Fatal("AliHLTHough::WaitForThreadFinish"," Can not join the required thread! ");
225 gSystem->MakeDirectory("hough1");
226 hough1->WriteTracks("./hough1");
227 gSystem->MakeDirectory("hough2");
228 hough2->WriteTracks("./hough2");
231 for(int slice=0; slice<=17; slice++)
233 // cout<<"Processing slice "<<slice<<endl;
234 fBenchmark->Start("ReadData");
235 hough1->ReadData(slice,0);
236 fBenchmark->Stop("ReadData");
237 fBenchmark->Start("Transform");
239 fBenchmark->Stop("Transform");
240 hough1->AddAllHistogramsRows();
241 hough1->FindTrackCandidatesRow();
242 fBenchmark->Start("AddTracks");
244 fBenchmark->Stop("AddTracks");
246 // AliHLTTrackArray* tracks = (AliHLTTrackArray*)hough1->GetTracks(0);
247 // nglobaltracks += tracks->GetNTracks();
249 for(int slice=18; slice<=35; slice++)
251 // cout<<"Processing slice "<<slice<<endl;
252 fBenchmark->Start("ReadData");
253 hough2->ReadData(slice,0);
254 fBenchmark->Stop("ReadData");
255 fBenchmark->Start("Transform");
257 fBenchmark->Stop("Transform");
258 hough2->AddAllHistogramsRows();
259 hough2->FindTrackCandidatesRow();
260 fBenchmark->Start("AddTracks");
262 fBenchmark->Stop("AddTracks");
264 // AliHLTTrackArray* tracks = (AliHLTTrackArray*)hough2->GetTracks(0);
265 // nglobaltracks += tracks->GetNTracks();
268 nglobaltracks += hough1->FillESD(esd);
269 nglobaltracks += hough2->FillESD(esd);
272 AliHLTITStracker itsTracker(0);
273 itsTracker.SetVertex(vtxPos,vtxErr);
275 itsTracker.LoadClusters(treeClusters);
276 itsTracker.Clusters2Tracks(esd);
277 itsTracker.UnloadClusters();
279 fBenchmark->Stop("Overall timing");
281 if (file) fprintf(file, "%s\n", time.AsString());
282 if (file) fprintf(file, "run: %d event: %d %d\n",
283 rawReader.GetRunNumber(),
284 rawReader.GetEventId()[0],
285 rawReader.GetEventId()[1]);
286 if (errorCode) fprintf(file, "ERROR: %d\n", errorCode);
288 if (file) fprintf(file, "Hough Transformer found %d tracks\n",nglobaltracks);
290 hough1->DoBench("hough1");
291 hough2->DoBench("hough2");
292 fBenchmark->Analyze("overall");
294 FILE* bench = fopen("hough1.dat", "r");
295 while (bench && !feof(bench)) {
297 if (!fgets(buffer, 256, bench)) break;
298 fprintf(file, "%s", buffer);
303 FILE* bench = fopen("hough2.dat", "r");
304 while (bench && !feof(bench)) {
306 if (!fgets(buffer, 256, bench)) break;
307 fprintf(file, "%s", buffer);
312 FILE* bench = fopen("overall.dat", "r");
313 while (bench && !feof(bench)) {
315 if (!fgets(buffer, 256, bench)) break;
316 fprintf(file, "%s", buffer);
319 fprintf(file, "\n\n");
330 // read run, event, detector, DDL numbers and data size
331 AliRawReaderDate rawReader(ptr);
333 printf("\n%s\n", time.AsString());
334 printf("run: %d event: %d %d\n", rawReader.GetRunNumber(),
335 rawReader.GetEventId()[0], rawReader.GetEventId()[1]);
336 while (rawReader.ReadMiniHeader()) {
337 printf(" detector: %d DDL: %d size: %d\n",
338 rawReader.GetDetectorID(), rawReader.GetDDLID(),
339 rawReader.GetDataSize());
344 gSystem->Sleep(100); // sleep for 0.1 second
347 gSystem->ProcessEvents();
348 if (file) fflush(file);
352 gSystem->RemoveSignalHandler(handler);
353 if (file) fclose(file);
359 int main(int /*argc*/, char** /*argv*/)
361 ::Fatal("main", "this program was compiled without DATE");