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>
54 //_____________________________________________________________________________
55 class AliGDCInterruptHandler : public TSignalHandler {
57 AliGDCInterruptHandler();
58 Bool_t Notify() {fStop = kTRUE; return kTRUE;};
59 Bool_t Stop() const {return fStop;};
61 Bool_t fStop; // CTRL-C pressed
64 //_____________________________________________________________________________
65 AliGDCInterruptHandler::AliGDCInterruptHandler() :
66 TSignalHandler(kSigInterrupt, kFALSE)
72 //_____________________________________________________________________________
74 int main(int argc, char** argv)
76 // set ROOT in batch mode
80 FILE* file = fopen("monitorGDC.log", "w");
84 // get data from a file or online from this node
87 status = monitorSetDataSource(argv[1]);
89 status = monitorSetDataSource(":");
91 if (status) ::Fatal("monitorSetDataSource", monitorDecodeError(status));
93 // monitor only a sample of physics events
94 char* table[] = {"Physics event", "yes", NULL};
95 status = monitorDeclareTable(table);
96 if (status) ::Fatal("monitorDeclareTable", monitorDecodeError(status));
98 // declare this monitoring program to DATE
99 status = monitorDeclareMp("GDC physics monitoring");
100 if (status) ::Fatal("monitorDeclareMp", monitorDecodeError(status));
102 // initialize HLT transformations
103 if (!AliL3Transform::Init("./", kFALSE)) {
104 ::Fatal("AliL3Transform::Init", "HLT initialization failed");
106 AliESD *esd = new AliESD;
107 AliKalmanTrack::SetConvConst(
108 1000/0.299792458/AliL3Transform::GetSolenoidField()
110 AliITSgeom *geom = new AliITSgeom();
111 geom->ReadNewFile("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymmFMD.det");
113 AliMagF* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., AliMagFMaps::k5kG);
114 AliTracker::SetFieldMap(field);
116 // create the signal handler
117 AliGDCInterruptHandler* handler = new AliGDCInterruptHandler;
118 gSystem->AddSignalHandler(handler);
122 while (!handler->Stop()) {
123 // get the next event
124 status = monitorGetEventDynamic(&ptr);
125 if (status == (Int_t)MON_ERR_EOF) break;
126 if (status) ::Fatal("monitorGetEventDynamic", monitorDecodeError(status));
130 gSystem->Sleep(1000); // sleep for 1 second
134 AliRawReaderDate rawReader(ptr);
135 // if ((rawReader.GetAttributes()[0] & 0x02) != 0) {
137 // Int_t errorCode = rawReader.CheckData();
139 if (errorCode && (errorCode != AliRawReader::kErrSize)) {
141 if (file) fprintf(file, "%s\n", time.AsString());
142 if (file) fprintf(file, "run: %d event: %d %d\n",
143 rawReader.GetRunNumber(),
144 rawReader.GetEventId()[0],
145 rawReader.GetEventId()[1]);
146 fprintf(file, "ERROR: %d\n\n", errorCode);
150 AliL3Benchmark *fBenchmark = new AliL3Benchmark();
151 fBenchmark->Start("Overall timing");
153 // ITS clusterer and vertexer
154 fBenchmark->Start("ITS Clusterer");
155 AliL3ITSclusterer clusterer(geom);
156 AliRawReader *itsrawreader=new AliRawReaderDate(ptr);
157 TTree* treeClusters = new TTree("TreeL3ITSclusters"," "); //make a tree
158 clusterer.Digits2Clusters(itsrawreader,treeClusters);
159 fBenchmark->Stop("ITS Clusterer");
161 AliL3ITSVertexerZ vertexer;
162 AliESDVertex *vertex = vertexer.FindVertexForCurrentEvent(geom,treeClusters);
164 Double_t vtxErr[3]={0.005,0.005,0.010};
165 vertex->GetXYZ(vtxPos);
166 // vertex->GetSigmaXYZ(vtxErr);
167 esd->SetVertex(vertex);
169 // TPC Hough reconstruction
170 Float_t ptmin = 0.1*AliL3Transform::GetSolenoidField();
171 Float_t zvertex = vtxPos[2];
173 // Run the Hough Transformer
174 fBenchmark->Start("Init");
175 AliL3Hough *hough1 = new AliL3Hough();
177 hough1->SetThreshold(4);
178 hough1->CalcTransformerParams(ptmin);
179 hough1->SetPeakThreshold(70,-1);
180 // printf("Pointer is %x\n",ptr);
181 hough1->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
182 hough1->SetAddHistograms();
183 fBenchmark->Stop("Init");
185 fBenchmark->Start("Init");
186 AliL3Hough *hough2 = new AliL3Hough();
188 hough2->SetThreshold(4);
189 hough2->CalcTransformerParams(ptmin);
190 hough2->SetPeakThreshold(70,-1);
191 // printf("Pointer is %x\n",ptr);
192 hough2->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
193 hough2->SetAddHistograms();
194 fBenchmark->Stop("Init");
196 Int_t nglobaltracks = 0;
198 hough1->StartProcessInThread(0,17);
199 hough2->StartProcessInThread(18,35);
201 // gSystem->Sleep(20000);
202 if(hough1->WaitForThreadFinish())
203 ::Fatal("AliL3Hough::WaitForThreadFinish"," Can not join the required thread! ");
204 if(hough2->WaitForThreadFinish())
205 ::Fatal("AliL3Hough::WaitForThreadFinish"," Can not join the required thread! ");
207 gSystem->MakeDirectory("hough1");
208 hough1->WriteTracks("./hough1");
209 gSystem->MakeDirectory("hough2");
210 hough2->WriteTracks("./hough2");
213 for(int slice=0; slice<=17; slice++)
215 // cout<<"Processing slice "<<slice<<endl;
216 fBenchmark->Start("ReadData");
217 hough1->ReadData(slice,0);
218 fBenchmark->Stop("ReadData");
219 fBenchmark->Start("Transform");
221 fBenchmark->Stop("Transform");
222 hough1->AddAllHistogramsRows();
223 hough1->FindTrackCandidatesRow();
224 fBenchmark->Start("AddTracks");
226 fBenchmark->Stop("AddTracks");
228 // AliL3TrackArray* tracks = (AliL3TrackArray*)hough1->GetTracks(0);
229 // nglobaltracks += tracks->GetNTracks();
231 for(int slice=18; slice<=35; slice++)
233 // cout<<"Processing slice "<<slice<<endl;
234 fBenchmark->Start("ReadData");
235 hough2->ReadData(slice,0);
236 fBenchmark->Stop("ReadData");
237 fBenchmark->Start("Transform");
239 fBenchmark->Stop("Transform");
240 hough2->AddAllHistogramsRows();
241 hough2->FindTrackCandidatesRow();
242 fBenchmark->Start("AddTracks");
244 fBenchmark->Stop("AddTracks");
246 // AliL3TrackArray* tracks = (AliL3TrackArray*)hough2->GetTracks(0);
247 // nglobaltracks += tracks->GetNTracks();
250 nglobaltracks += hough1->FillESD(esd);
251 nglobaltracks += hough2->FillESD(esd);
254 AliL3ITStracker itsTracker(geom);
255 itsTracker.SetVertex(vtxPos,vtxErr);
257 itsTracker.LoadClusters(treeClusters);
258 itsTracker.Clusters2Tracks(esd);
259 itsTracker.UnloadClusters();
261 fBenchmark->Stop("Overall timing");
263 if (file) fprintf(file, "%s\n", time.AsString());
264 if (file) fprintf(file, "run: %d event: %d %d\n",
265 rawReader.GetRunNumber(),
266 rawReader.GetEventId()[0],
267 rawReader.GetEventId()[1]);
268 if (errorCode) fprintf(file, "ERROR: %d\n", errorCode);
270 if (file) fprintf(file, "Hough Transformer found %d tracks\n",nglobaltracks);
272 hough1->DoBench("hough1");
273 hough2->DoBench("hough2");
274 fBenchmark->Analyze("overall");
276 FILE* bench = fopen("hough1.dat", "r");
277 while (bench && !feof(bench)) {
279 if (!fgets(buffer, 256, bench)) break;
280 fprintf(file, "%s", buffer);
285 FILE* bench = fopen("hough2.dat", "r");
286 while (bench && !feof(bench)) {
288 if (!fgets(buffer, 256, bench)) break;
289 fprintf(file, "%s", buffer);
294 FILE* bench = fopen("overall.dat", "r");
295 while (bench && !feof(bench)) {
297 if (!fgets(buffer, 256, bench)) break;
298 fprintf(file, "%s", buffer);
301 fprintf(file, "\n\n");
312 // read run, event, detector, DDL numbers and data size
313 AliRawReaderDate rawReader(ptr);
315 printf("\n%s\n", time.AsString());
316 printf("run: %d event: %d %d\n", rawReader.GetRunNumber(),
317 rawReader.GetEventId()[0], rawReader.GetEventId()[1]);
318 while (rawReader.ReadMiniHeader()) {
319 printf(" detector: %d DDL: %d size: %d\n",
320 rawReader.GetDetectorID(), rawReader.GetDDLID(),
321 rawReader.GetDataSize());
326 gSystem->Sleep(100); // sleep for 0.1 second
329 gSystem->ProcessEvents();
330 if (file) fflush(file);
334 gSystem->RemoveSignalHandler(handler);
335 if (file) fclose(file);
341 int main(int /*argc*/, char** /*argv*/)
343 ::Fatal("main", "this program was compiled without DATE");