]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MONITOR/monitorGDC.cxx
- Dipole rotated wr to ALICE coordinate system
[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   AliESD *esd = new AliESD;
108   //  AliKalmanTrack::SetConvConst(
109   //     1000/0.299792458/AliHLTTransform::GetSolenoidField()
110   //  );
111   AliITSgeom *geom = new AliITSgeom();
112   geom->ReadNewFile("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymmFMD.det");
113   if (!geom) return 1;
114   Int_t sfield = 0;
115   switch ((Int_t)(AliHLTTransform::GetSolenoidField()+0.5)) {
116   case 2:
117     sfield = AliMagFMaps::k2kG;
118     break;
119   case 4:
120     sfield = AliMagFMaps::k4kG;
121     break;
122   case 5:
123     sfield = AliMagFMaps::k5kG;
124     break;
125   default:
126     ::Fatal("AliHLTTransform::GetSolenoidField", "Incorrect magnetic field");
127   }
128   AliMagF* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., sfield);
129   AliTracker::SetFieldMap(field,kTRUE);
130
131   // Init PID
132   AliPID pid;
133
134   // create the signal handler
135   AliGDCInterruptHandler* handler = new AliGDCInterruptHandler;
136   gSystem->AddSignalHandler(handler);
137
138   // endless loop
139   void* ptr = NULL;
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));
145
146     // if no new event
147     if (!ptr) {
148       gSystem->Sleep(1000);   // sleep for 1 second
149       continue;
150     }
151
152     AliRawReaderDate rawReader(ptr);
153     //    if ((rawReader.GetAttributes()[0] & 0x02) != 0) {
154
155     //      Int_t errorCode = rawReader.CheckData();
156     Int_t errorCode = 0;
157       if (errorCode && (errorCode != AliRawReader::kErrSize)) {
158         time.Set();
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);
165
166       } else {
167
168         AliHLTBenchmark *fBenchmark = new AliHLTBenchmark();
169         fBenchmark->Start("Overall timing");
170
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");
178         
179         AliHLTITSVertexerZ vertexer;
180         AliESDVertex *vertex = vertexer.FindVertexForCurrentEvent(geom,treeClusters);
181         Double_t vtxPos[3];
182         Double_t vtxErr[3]={0.005,0.005,0.010};
183         vertex->GetXYZ(vtxPos);
184         //      vertex->GetSigmaXYZ(vtxErr);
185         esd->SetVertex(vertex);
186
187         // TPC Hough reconstruction
188         Float_t ptmin = 0.1*AliHLTTransform::GetSolenoidField();
189         Float_t zvertex = vtxPos[2];
190
191         // Run the Hough Transformer
192         fBenchmark->Start("Init");
193         AliHLTHough *hough1 = new AliHLTHough();
194
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");
202
203         fBenchmark->Start("Init");
204         AliHLTHough *hough2 = new AliHLTHough();
205
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");
213
214         Int_t nglobaltracks = 0;
215         /*
216         hough1->StartProcessInThread(0,17);
217         hough2->StartProcessInThread(18,35);
218
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! ");
224
225         gSystem->MakeDirectory("hough1");
226         hough1->WriteTracks("./hough1");
227         gSystem->MakeDirectory("hough2");
228         hough2->WriteTracks("./hough2");
229         */
230
231         for(int slice=0; slice<=17; slice++)
232         {
233           //      cout<<"Processing slice "<<slice<<endl;
234           fBenchmark->Start("ReadData");
235           hough1->ReadData(slice,0);
236           fBenchmark->Stop("ReadData");
237           fBenchmark->Start("Transform");
238           hough1->Transform();
239           fBenchmark->Stop("Transform");
240           hough1->AddAllHistogramsRows();
241           hough1->FindTrackCandidatesRow();
242           fBenchmark->Start("AddTracks");
243           hough1->AddTracks();
244           fBenchmark->Stop("AddTracks");
245
246           //      AliHLTTrackArray* tracks = (AliHLTTrackArray*)hough1->GetTracks(0);
247           //      nglobaltracks += tracks->GetNTracks();
248         }
249         for(int slice=18; slice<=35; slice++)
250         {
251           //      cout<<"Processing slice "<<slice<<endl;
252           fBenchmark->Start("ReadData");
253           hough2->ReadData(slice,0);
254           fBenchmark->Stop("ReadData");
255           fBenchmark->Start("Transform");
256           hough2->Transform();
257           fBenchmark->Stop("Transform");
258           hough2->AddAllHistogramsRows();
259           hough2->FindTrackCandidatesRow();
260           fBenchmark->Start("AddTracks");
261           hough2->AddTracks();
262           fBenchmark->Stop("AddTracks");
263
264           //      AliHLTTrackArray* tracks = (AliHLTTrackArray*)hough2->GetTracks(0);
265           //      nglobaltracks += tracks->GetNTracks();
266         }
267
268         nglobaltracks += hough1->FillESD(esd);
269         nglobaltracks += hough2->FillESD(esd);
270
271         // ITS tracker
272         AliHLTITStracker itsTracker(0);
273         itsTracker.SetVertex(vtxPos,vtxErr);
274
275         itsTracker.LoadClusters(treeClusters);
276         itsTracker.Clusters2Tracks(esd);
277         itsTracker.UnloadClusters();
278
279         fBenchmark->Stop("Overall timing");
280         time.Set();
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);
287
288         if (file) fprintf(file, "Hough Transformer found %d tracks\n",nglobaltracks);
289
290         hough1->DoBench("hough1");
291         hough2->DoBench("hough2");
292         fBenchmark->Analyze("overall");
293         if (file) {
294           FILE* bench = fopen("hough1.dat", "r");
295           while (bench && !feof(bench)) {
296             char buffer[256];
297             if (!fgets(buffer, 256, bench)) break;
298             fprintf(file, "%s", buffer);
299           }
300           fclose(bench);
301         }
302         if (file) {
303           FILE* bench = fopen("hough2.dat", "r");
304           while (bench && !feof(bench)) {
305             char buffer[256];
306             if (!fgets(buffer, 256, bench)) break;
307             fprintf(file, "%s", buffer);
308           }
309           fclose(bench);
310         }
311         if (file) {
312           FILE* bench = fopen("overall.dat", "r");
313           while (bench && !feof(bench)) {
314             char buffer[256];
315             if (!fgets(buffer, 256, bench)) break;
316             fprintf(file, "%s", buffer);
317           }
318           fclose(bench);
319           fprintf(file, "\n\n");
320         }
321
322         delete hough1;
323         delete hough2;
324
325         esd->Reset();
326       }
327     //    }
328
329     /*
330     // read run, event, detector, DDL numbers and data size
331     AliRawReaderDate rawReader(ptr);
332     time.Set();
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());
340     }
341
342     */
343
344     gSystem->Sleep(100);   // sleep for 0.1 second
345     free(ptr);
346
347     gSystem->ProcessEvents();
348     if (file) fflush(file);
349
350   }
351
352   gSystem->RemoveSignalHandler(handler);
353   if (file) fclose(file);
354
355   return 0;
356 }
357
358 #else
359 int main(int /*argc*/, char** /*argv*/)
360 {
361   ::Fatal("main", "this program was compiled without DATE");
362
363   return 1;
364 }
365 #endif