]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MONITOR/monitorGDC.cxx
New class AliESDEvent, backward compatibility with the old AliESD (Christian)
[u/mrichter/AliRoot.git] / MONITOR / monitorGDC.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // this program performs local monitoring on a GDC by running the HLT code   //
21 //                                                                           //
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.                            //
25 //                                                                           //
26 ///////////////////////////////////////////////////////////////////////////////
27
28 #include <TError.h>
29 #include <TSysEvtHandler.h>
30 #ifdef ALI_DATE
31 #include <TROOT.h>
32 #include <TSystem.h>
33 #include <TDatime.h>
34 #include "AliRawReaderDate.h"
35 #include "event.h"
36 #include "monitor.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"
48 #include "AliMagF.h"
49 #include "AliMagFMaps.h"
50 #include <AliHLTITSclusterer.h>
51 #include <AliHLTITSVertexerZ.h>
52 #include <AliHLTITStracker.h>
53 #endif
54
55 //_____________________________________________________________________________
56 class AliGDCInterruptHandler : public TSignalHandler {
57 public:
58   AliGDCInterruptHandler();
59   Bool_t Notify() {fStop = kTRUE; return kTRUE;};
60   Bool_t Stop() const {return fStop;};
61 private:
62   Bool_t fStop;  // CTRL-C pressed
63 };
64
65 //_____________________________________________________________________________
66 AliGDCInterruptHandler::AliGDCInterruptHandler() : 
67   TSignalHandler(kSigInterrupt, kFALSE) 
68 {
69   fStop = kFALSE;
70 };
71
72
73 //_____________________________________________________________________________
74 #ifdef ALI_DATE
75 int main(int argc, char** argv)
76 {
77   // set ROOT in batch mode
78   gROOT->SetBatch();   
79
80   // open a log file
81   FILE* file = fopen("monitorGDC.log", "w");
82
83   TDatime time;
84
85   // get data from a file or online from this node
86   Int_t status = 0;
87   if (argc > 1) {
88     status = monitorSetDataSource(argv[1]);
89   } else {
90     status = monitorSetDataSource(":");
91   }
92   if (status) ::Fatal("monitorSetDataSource", monitorDecodeError(status));
93
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));
98
99   // declare this monitoring program to DATE
100   status = monitorDeclareMp("GDC physics monitoring");
101   if (status) ::Fatal("monitorDeclareMp", monitorDecodeError(status));
102
103   // initialize HLT transformations
104   if (!AliHLTTransform::Init("./", kFALSE)) {
105     ::Fatal("AliHLTTransform::Init", "HLT initialization failed");
106   }
107   AliESDEvent *esd = new AliESDEvent;
108   esd->CreateStdContent();
109   //  AliKalmanTrack::SetConvConst(
110   //     1000/0.299792458/AliHLTTransform::GetSolenoidField()
111   //  );
112   AliITSgeom *geom = new AliITSgeom();
113   geom->ReadNewFile("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymmFMD.det");
114   if (!geom) return 1;
115   Int_t sfield = 0;
116   switch ((Int_t)(AliHLTTransform::GetSolenoidField()+0.5)) {
117   case 2:
118     sfield = AliMagFMaps::k2kG;
119     break;
120   case 4:
121     sfield = AliMagFMaps::k4kG;
122     break;
123   case 5:
124     sfield = AliMagFMaps::k5kG;
125     break;
126   default:
127     ::Fatal("AliHLTTransform::GetSolenoidField", "Incorrect magnetic field");
128   }
129   AliMagF* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., sfield);
130   AliTracker::SetFieldMap(field,kTRUE);
131
132   // Init PID
133   AliPID pid;
134
135   // create the signal handler
136   AliGDCInterruptHandler* handler = new AliGDCInterruptHandler;
137   gSystem->AddSignalHandler(handler);
138
139   // endless loop
140   void* ptr = NULL;
141   while (!handler->Stop()) {
142     // get the next event
143     status = monitorGetEventDynamic(&ptr);
144     if (status == (Int_t)MON_ERR_EOF) break;
145     if (status) ::Fatal("monitorGetEventDynamic", monitorDecodeError(status));
146
147     // if no new event
148     if (!ptr) {
149       gSystem->Sleep(1000);   // sleep for 1 second
150       continue;
151     }
152
153     AliRawReaderDate rawReader(ptr);
154     //    if ((rawReader.GetAttributes()[0] & 0x02) != 0) {
155
156     //      Int_t errorCode = rawReader.CheckData();
157     Int_t errorCode = 0;
158       if (errorCode && (errorCode != AliRawReader::kErrSize)) {
159         time.Set();
160         if (file) fprintf(file, "%s\n", time.AsString());
161         if (file) fprintf(file, "run: %d  event: %d %d\n", 
162                           rawReader.GetRunNumber(), 
163                           rawReader.GetEventId()[0], 
164                           rawReader.GetEventId()[1]);
165         fprintf(file, "ERROR: %d\n\n", errorCode);
166
167       } else {
168
169         AliHLTBenchmark *fBenchmark = new AliHLTBenchmark();
170         fBenchmark->Start("Overall timing");
171
172         // ITS clusterer and vertexer
173         fBenchmark->Start("ITS Clusterer");
174         AliHLTITSclusterer clusterer(0);
175         AliRawReader *itsrawreader=new AliRawReaderDate(ptr);
176         TTree* treeClusters = new TTree("TreeL3ITSclusters"," "); //make a tree
177         clusterer.Digits2Clusters(itsrawreader,treeClusters);
178         fBenchmark->Stop("ITS Clusterer");
179         
180         AliHLTITSVertexerZ vertexer;
181         AliESDVertex *vertex = vertexer.FindVertexForCurrentEvent(geom,treeClusters);
182         Double_t vtxPos[3];
183         Double_t vtxErr[3]={0.005,0.005,0.010};
184         vertex->GetXYZ(vtxPos);
185         //      vertex->GetSigmaXYZ(vtxErr);
186         esd->SetVertex(vertex);
187
188         // TPC Hough reconstruction
189         Float_t ptmin = 0.1*AliHLTTransform::GetSolenoidField();
190         Float_t zvertex = vtxPos[2];
191
192         // Run the Hough Transformer
193         fBenchmark->Start("Init");
194         AliHLTHough *hough1 = new AliHLTHough();
195
196         hough1->SetThreshold(4);
197         hough1->CalcTransformerParams(ptmin);
198         hough1->SetPeakThreshold(70,-1);
199         //      printf("Pointer is %x\n",ptr);
200         hough1->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
201         hough1->SetAddHistograms();
202         fBenchmark->Stop("Init");
203
204         fBenchmark->Start("Init");
205         AliHLTHough *hough2 = new AliHLTHough();
206
207         hough2->SetThreshold(4);
208         hough2->CalcTransformerParams(ptmin);
209         hough2->SetPeakThreshold(70,-1);
210         //      printf("Pointer is %x\n",ptr);
211         hough2->Init("./", kFALSE, 100, kFALSE,4,0,(Char_t*)ptr,zvertex);
212         hough2->SetAddHistograms();
213         fBenchmark->Stop("Init");
214
215         Int_t nglobaltracks = 0;
216         /*
217         hough1->StartProcessInThread(0,17);
218         hough2->StartProcessInThread(18,35);
219
220         //      gSystem->Sleep(20000);
221         if(hough1->WaitForThreadFinish())
222           ::Fatal("AliHLTHough::WaitForThreadFinish"," Can not join the required thread! ");
223         if(hough2->WaitForThreadFinish())
224           ::Fatal("AliHLTHough::WaitForThreadFinish"," Can not join the required thread! ");
225
226         gSystem->MakeDirectory("hough1");
227         hough1->WriteTracks("./hough1");
228         gSystem->MakeDirectory("hough2");
229         hough2->WriteTracks("./hough2");
230         */
231
232         for(int slice=0; slice<=17; slice++)
233         {
234           //      cout<<"Processing slice "<<slice<<endl;
235           fBenchmark->Start("ReadData");
236           hough1->ReadData(slice,0);
237           fBenchmark->Stop("ReadData");
238           fBenchmark->Start("Transform");
239           hough1->Transform();
240           fBenchmark->Stop("Transform");
241           hough1->AddAllHistogramsRows();
242           hough1->FindTrackCandidatesRow();
243           fBenchmark->Start("AddTracks");
244           hough1->AddTracks();
245           fBenchmark->Stop("AddTracks");
246
247           //      AliHLTTrackArray* tracks = (AliHLTTrackArray*)hough1->GetTracks(0);
248           //      nglobaltracks += tracks->GetNTracks();
249         }
250         for(int slice=18; slice<=35; slice++)
251         {
252           //      cout<<"Processing slice "<<slice<<endl;
253           fBenchmark->Start("ReadData");
254           hough2->ReadData(slice,0);
255           fBenchmark->Stop("ReadData");
256           fBenchmark->Start("Transform");
257           hough2->Transform();
258           fBenchmark->Stop("Transform");
259           hough2->AddAllHistogramsRows();
260           hough2->FindTrackCandidatesRow();
261           fBenchmark->Start("AddTracks");
262           hough2->AddTracks();
263           fBenchmark->Stop("AddTracks");
264
265           //      AliHLTTrackArray* tracks = (AliHLTTrackArray*)hough2->GetTracks(0);
266           //      nglobaltracks += tracks->GetNTracks();
267         }
268
269         nglobaltracks += hough1->FillESD(esd);
270         nglobaltracks += hough2->FillESD(esd);
271
272         // ITS tracker
273         AliHLTITStracker itsTracker(0);
274         itsTracker.SetVertex(vtxPos,vtxErr);
275
276         itsTracker.LoadClusters(treeClusters);
277         itsTracker.Clusters2Tracks(esd);
278         itsTracker.UnloadClusters();
279
280         fBenchmark->Stop("Overall timing");
281         time.Set();
282         if (file) fprintf(file, "%s\n", time.AsString());
283         if (file) fprintf(file, "run: %d  event: %d %d\n", 
284                           rawReader.GetRunNumber(), 
285                           rawReader.GetEventId()[0], 
286                           rawReader.GetEventId()[1]);
287         if (errorCode) fprintf(file, "ERROR: %d\n", errorCode);
288
289         if (file) fprintf(file, "Hough Transformer found %d tracks\n",nglobaltracks);
290
291         hough1->DoBench("hough1");
292         hough2->DoBench("hough2");
293         fBenchmark->Analyze("overall");
294         if (file) {
295           FILE* bench = fopen("hough1.dat", "r");
296           while (bench && !feof(bench)) {
297             char buffer[256];
298             if (!fgets(buffer, 256, bench)) break;
299             fprintf(file, "%s", buffer);
300           }
301           fclose(bench);
302         }
303         if (file) {
304           FILE* bench = fopen("hough2.dat", "r");
305           while (bench && !feof(bench)) {
306             char buffer[256];
307             if (!fgets(buffer, 256, bench)) break;
308             fprintf(file, "%s", buffer);
309           }
310           fclose(bench);
311         }
312         if (file) {
313           FILE* bench = fopen("overall.dat", "r");
314           while (bench && !feof(bench)) {
315             char buffer[256];
316             if (!fgets(buffer, 256, bench)) break;
317             fprintf(file, "%s", buffer);
318           }
319           fclose(bench);
320           fprintf(file, "\n\n");
321         }
322
323         delete hough1;
324         delete hough2;
325
326         esd->Reset();
327       }
328     //    }
329
330     /*
331     // read run, event, detector, DDL numbers and data size
332     AliRawReaderDate rawReader(ptr);
333     time.Set();
334     printf("\n%s\n", time.AsString());
335     printf("run: %d  event: %d %d\n", rawReader.GetRunNumber(), 
336            rawReader.GetEventId()[0], rawReader.GetEventId()[1]);
337     while (rawReader.ReadMiniHeader()) {
338       printf(" detector: %d   DDL: %d  size: %d\n", 
339              rawReader.GetDetectorID(), rawReader.GetDDLID(), 
340              rawReader.GetDataSize());
341     }
342
343     */
344
345     gSystem->Sleep(100);   // sleep for 0.1 second
346     free(ptr);
347
348     gSystem->ProcessEvents();
349     if (file) fflush(file);
350
351   }
352
353   gSystem->RemoveSignalHandler(handler);
354   if (file) fclose(file);
355
356   return 0;
357 }
358
359 #else
360 int main(int /*argc*/, char** /*argv*/)
361 {
362   ::Fatal("main", "this program was compiled without DATE");
363
364   return 1;
365 }
366 #endif