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 <AliL3StandardIncludes.h>
38 #include <AliL3Transform.h>
39 #include <AliL3MemHandler.h>
40 #include <AliL3TrackArray.h>
41 #include <AliL3HoughMaxFinder.h>
42 #include <AliL3HoughBaseTransformer.h>
43 #include <AliL3Hough.h>
44 #include <AliL3Benchmark.h>
45 #include <AliKalmanTrack.h>
46 #include "AliITSgeom.h"
48 #include "AliMagFMaps.h"
49 #include <AliL3ITSclusterer.h>
50 #include <AliL3ITSVertexerZ.h>
51 #include <AliL3ITStracker.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 (!AliL3Transform::Init("./", kFALSE)) {
105 ::Fatal("AliL3Transform::Init", "HLT initialization failed");
107 AliESD *esd = new AliESD;
108 AliKalmanTrack::SetConvConst(
109 1000/0.299792458/AliL3Transform::GetSolenoidField()
111 AliITSgeom *geom = new AliITSgeom();
112 geom->ReadNewFile("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymmFMD.det");
114 AliMagF* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., AliMagFMaps::k5kG);
115 AliTracker::SetFieldMap(field);
117 // create the signal handler
118 AliGDCInterruptHandler* handler = new AliGDCInterruptHandler;
119 gSystem->AddSignalHandler(handler);
123 while (!handler->Stop()) {
124 // get the next event
125 status = monitorGetEventDynamic(&ptr);
126 if (status == (Int_t)MON_ERR_EOF) break;
127 if (status) ::Fatal("monitorGetEventDynamic", monitorDecodeError(status));
131 gSystem->Sleep(1000); // sleep for 1 second
135 AliRawReaderDate rawReader(ptr);
136 // if ((rawReader.GetAttributes()[0] & 0x02) != 0) {
138 // Int_t errorCode = rawReader.CheckData();
140 if (errorCode && (errorCode != AliRawReader::kErrSize)) {
142 if (file) fprintf(file, "%s\n", time.AsString());
143 if (file) fprintf(file, "run: %d event: %d %d\n",
144 rawReader.GetRunNumber(),
145 rawReader.GetEventId()[0],
146 rawReader.GetEventId()[1]);
147 fprintf(file, "ERROR: %d\n\n", errorCode);
151 AliL3Benchmark *fBenchmark = new AliL3Benchmark();
152 fBenchmark->Start("Overall timing");
154 // ITS clusterer and vertexer
155 fBenchmark->Start("ITS Clusterer");
156 AliL3ITSclusterer clusterer(geom);
157 AliRawReader *itsrawreader=new AliRawReaderDate(ptr);
158 TTree* treeClusters = new TTree("TreeL3ITSclusters"," "); //make a tree
159 clusterer.Digits2Clusters(itsrawreader,treeClusters);
160 fBenchmark->Stop("ITS Clusterer");
162 AliL3ITSVertexerZ vertexer;
163 AliESDVertex *vertex = vertexer.FindVertexForCurrentEvent(geom,treeClusters);
165 Double_t vtxErr[3]={0.005,0.005,0.010};
166 vertex->GetXYZ(vtxPos);
167 // vertex->GetSigmaXYZ(vtxErr);
168 esd->SetVertex(vertex);
170 // TPC Hough reconstruction
171 Float_t ptmin = 0.1*AliL3Transform::GetSolenoidField();
172 Float_t zvertex = vtxPos[2];
174 // Run the Hough Transformer
175 fBenchmark->Start("Init");
176 AliL3Hough *hough1 = new AliL3Hough();
178 hough1->SetThreshold(4);
179 hough1->CalcTransformerParams(ptmin);
180 hough1->SetPeakThreshold(70,-1);
181 // printf("Pointer is %x\n",ptr);
182 hough1->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
183 hough1->SetAddHistograms();
184 fBenchmark->Stop("Init");
186 fBenchmark->Start("Init");
187 AliL3Hough *hough2 = new AliL3Hough();
189 hough2->SetThreshold(4);
190 hough2->CalcTransformerParams(ptmin);
191 hough2->SetPeakThreshold(70,-1);
192 // printf("Pointer is %x\n",ptr);
193 hough2->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
194 hough2->SetAddHistograms();
195 fBenchmark->Stop("Init");
197 Int_t nglobaltracks = 0;
199 hough1->StartProcessInThread(0,17);
200 hough2->StartProcessInThread(18,35);
202 // gSystem->Sleep(20000);
203 if(hough1->WaitForThreadFinish())
204 ::Fatal("AliL3Hough::WaitForThreadFinish"," Can not join the required thread! ");
205 if(hough2->WaitForThreadFinish())
206 ::Fatal("AliL3Hough::WaitForThreadFinish"," Can not join the required thread! ");
208 gSystem->MakeDirectory("hough1");
209 hough1->WriteTracks("./hough1");
210 gSystem->MakeDirectory("hough2");
211 hough2->WriteTracks("./hough2");
214 for(int slice=0; slice<=17; slice++)
216 // cout<<"Processing slice "<<slice<<endl;
217 fBenchmark->Start("ReadData");
218 hough1->ReadData(slice,0);
219 fBenchmark->Stop("ReadData");
220 fBenchmark->Start("Transform");
222 fBenchmark->Stop("Transform");
223 hough1->AddAllHistogramsRows();
224 hough1->FindTrackCandidatesRow();
225 fBenchmark->Start("AddTracks");
227 fBenchmark->Stop("AddTracks");
229 // AliL3TrackArray* tracks = (AliL3TrackArray*)hough1->GetTracks(0);
230 // nglobaltracks += tracks->GetNTracks();
232 for(int slice=18; slice<=35; slice++)
234 // cout<<"Processing slice "<<slice<<endl;
235 fBenchmark->Start("ReadData");
236 hough2->ReadData(slice,0);
237 fBenchmark->Stop("ReadData");
238 fBenchmark->Start("Transform");
240 fBenchmark->Stop("Transform");
241 hough2->AddAllHistogramsRows();
242 hough2->FindTrackCandidatesRow();
243 fBenchmark->Start("AddTracks");
245 fBenchmark->Stop("AddTracks");
247 // AliL3TrackArray* tracks = (AliL3TrackArray*)hough2->GetTracks(0);
248 // nglobaltracks += tracks->GetNTracks();
251 nglobaltracks += hough1->FillESD(esd);
252 nglobaltracks += hough2->FillESD(esd);
255 AliL3ITStracker itsTracker(geom);
256 itsTracker.SetVertex(vtxPos,vtxErr);
258 itsTracker.LoadClusters(treeClusters);
259 itsTracker.Clusters2Tracks(esd);
260 itsTracker.UnloadClusters();
262 fBenchmark->Stop("Overall timing");
264 if (file) fprintf(file, "%s\n", time.AsString());
265 if (file) fprintf(file, "run: %d event: %d %d\n",
266 rawReader.GetRunNumber(),
267 rawReader.GetEventId()[0],
268 rawReader.GetEventId()[1]);
269 if (errorCode) fprintf(file, "ERROR: %d\n", errorCode);
271 if (file) fprintf(file, "Hough Transformer found %d tracks\n",nglobaltracks);
273 hough1->DoBench("hough1");
274 hough2->DoBench("hough2");
275 fBenchmark->Analyze("overall");
277 FILE* bench = fopen("hough1.dat", "r");
278 while (bench && !feof(bench)) {
280 if (!fgets(buffer, 256, bench)) break;
281 fprintf(file, "%s", buffer);
286 FILE* bench = fopen("hough2.dat", "r");
287 while (bench && !feof(bench)) {
289 if (!fgets(buffer, 256, bench)) break;
290 fprintf(file, "%s", buffer);
295 FILE* bench = fopen("overall.dat", "r");
296 while (bench && !feof(bench)) {
298 if (!fgets(buffer, 256, bench)) break;
299 fprintf(file, "%s", buffer);
302 fprintf(file, "\n\n");
313 // read run, event, detector, DDL numbers and data size
314 AliRawReaderDate rawReader(ptr);
316 printf("\n%s\n", time.AsString());
317 printf("run: %d event: %d %d\n", rawReader.GetRunNumber(),
318 rawReader.GetEventId()[0], rawReader.GetEventId()[1]);
319 while (rawReader.ReadMiniHeader()) {
320 printf(" detector: %d DDL: %d size: %d\n",
321 rawReader.GetDetectorID(), rawReader.GetDDLID(),
322 rawReader.GetDataSize());
327 gSystem->Sleep(100); // sleep for 0.1 second
330 gSystem->ProcessEvents();
331 if (file) fflush(file);
335 gSystem->RemoveSignalHandler(handler);
336 if (file) fclose(file);
342 int main(int /*argc*/, char** /*argv*/)
344 ::Fatal("main", "this program was compiled without DATE");