]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDPattern.cxx
Clean-up compile-time warnings
[u/mrichter/AliRoot.git] / FMD / AliFMDPattern.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, 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 /* $Id$ */
16 /** @file    AliFMDPattern.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:39:09 2006
19     @brief   FMD 2D Event display 
20 */
21 //___________________________________________________________________
22 //
23 // The classes defined here, are utility classes for reading in data
24 // for the FMD.  They are  put in a seperate library to not polute the
25 // normal libraries.  The classes are intended to be used as base
26 // classes for customized class that do some sort of analysis on the
27 // various types of data produced by the FMD. 
28 //
29 // Latest changes by Christian Holm Christensen
30 //
31
32 // #include <iostream>
33
34 // #include <TApplication.h>
35 // #include <TButton.h>
36 #include <TCanvas.h>
37 #include <TH2F.h>
38 #include <TMath.h>
39 #include <TPad.h>
40 #include <TRandom.h>
41 // #include <TSlider.h>
42 #include <TStyle.h>
43 // #include <TSystem.h>
44 #include <TVector2.h>
45 // #include <TView.h>
46 #include <TGraph.h>
47 #include "AliFMDPattern.h"      // ALIFMDDISPLAY_H
48 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
49 //#include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
50 #include "AliFMDRing.h"
51 // #include "AliFMDDetector.h"
52 #include "AliFMDHit.h"
53 // #include <AliLog.h>
54 #include "AliFMDDebug.h" // Better debug macros
55 class AliFMDDetector;
56
57 //____________________________________________________________________
58 ClassImp(AliFMDPattern)
59 #if 0
60   ; // This is here to keep Emacs for indenting the next line
61 #endif
62
63 //____________________________________________________________________
64 AliFMDPattern::AliFMDPatternDetector::AliFMDPatternDetector(UShort_t id) 
65   : fId(id),
66     fCounts(0), 
67     fGraphs(0), 
68     fFrame(0)
69 {
70   // CTOR 
71   // 
72   // Parameters: 
73   // 
74   //   ID       Identifier 
75 }
76
77 //____________________________________________________________________
78 AliFMDPattern::AliFMDPatternDetector::~AliFMDPatternDetector()
79 {
80   // DTOR 
81   // Destructor - 
82   // deletes mother frame 
83   if (fFrame) delete fFrame;
84 }
85
86 //____________________________________________________________________
87 void
88 AliFMDPattern::AliFMDPatternDetector::DrawShape(TObjArray& a) 
89 {
90   // Draw all shapes. 
91   // 
92   // Paramters 
93   // 
94   //    a       Array of shapes 
95   //
96   TIter next(&a);
97   TGraph* g = 0;
98   while ((g = static_cast<TGraph*>(next()))) {
99     g->DrawClone("f same");
100     g->DrawClone("l same");
101   }
102 }
103
104 //____________________________________________________________________
105 void
106 AliFMDPattern::AliFMDPatternDetector::Begin(Int_t nlevel, Double_t r, 
107                                TObjArray& inners, TObjArray& outers)
108 {
109   // Start of a run. 
110   // 
111   // Parameters 
112   // 
113   //    nlevel          Number of levels 
114   //    r               Radius 
115   //    inners          Array of inner shapes 
116   //    outers          Array of outer shapes 
117   //
118
119   // To make code-checker shut up
120   TStyle* style = gStyle;  
121   if (nlevel < 1) nlevel = style->GetNumberOfColors();
122   fCounts.Set(nlevel);
123   if (!fFrame) {
124     // The code-checker thinks this is not using the declaration of
125     // TH2F - what a morron!   
126     fFrame = new TH2F(Form("fmd%dFrame", fId), Form("FMD%d", fId), 
127                       10, -r, r, 10, -r, r);
128     fFrame->SetStats(kFALSE);
129     fFrame->Draw();
130   }
131   DrawShape(inners);
132   if (fId != 1) DrawShape(outers);
133   for (Int_t i = 0; i < nlevel; i++) { 
134     TGraph* g = new TGraph;
135     Int_t idx = Int_t(Float_t(i) / nlevel * style->GetNumberOfColors());
136     Int_t col = style->GetColorPalette(idx);
137     g->SetName(Form("FMD%d_L%02d", fId, i));
138     g->SetMarkerColor(col);
139     g->SetLineColor(col);
140     g->SetFillColor(col);
141     g->SetMarkerSize(i * .2 + .2);
142     g->SetMarkerStyle(2);
143     g->Draw("same p");
144     fGraphs.AddAtAndExpand(g, i);
145   }
146   TIter   next(&fGraphs);
147 }
148
149 //____________________________________________________________________
150 void
151 AliFMDPattern::AliFMDPatternDetector::Clear() 
152 {
153   // Clear this display.  
154   // Simply reset counters to zero. 
155   // Avoid deleting memory. 
156   fCounts.Reset(0);
157 }
158
159 //____________________________________________________________________
160 void
161 AliFMDPattern::AliFMDPatternDetector::End()
162 {
163   // Called when displaying the data.  
164   // Simply resets number of points at each level to 
165   // the seen number of hits at that level. 
166   // Avoid deleting memory. 
167   TIter   next(&fGraphs);
168   TGraph* g = 0;
169   Int_t   i = 0;
170   while ((g = static_cast<TGraph*>(next()))) g->Set(fCounts[i++]);
171 }
172 //____________________________________________________________________
173 void
174 AliFMDPattern::AliFMDPatternDetector::AddMarker(Double_t x, Double_t y, Float_t s, 
175                                    Float_t max)
176 {
177   // Add a marker at (X,Y,Z).  The marker color and size is chosen
178   // relative to the MAX argument. 
179   // 
180   // Parameters 
181   // 
182   //    X,Y,Z           Coordiantes 
183   //    MAX             Maximum value. 
184   // 
185   /** Sigh, for some odd reason, the code-checker does not recognise
186       this a usage of the TMath namespace declaration! Idiot */
187   Int_t i = TMath::Min(Int_t(fCounts.fN * s / max),  
188                        Int_t(fGraphs.GetEntries()-1));
189   TGraph* g = static_cast<TGraph*>(fGraphs.At(i));
190   if (!g) return;
191   g->SetPoint(fCounts[i]++, x, y);
192 }
193
194
195 //____________________________________________________________________
196 AliFMDPattern::AliFMDPattern(const char* gAliceFile)
197   : AliFMDDisplay(kTRUE, gAliceFile),
198     fInners(0), 
199     fOuters(0),
200     fInnerMax(0), 
201     fOuterMax(0),
202     fFMD1Pad(0),
203     fFMD1(1),
204     fFMD2Pad(0),
205     fFMD2(2),
206     fFMD3Pad(0),
207     fFMD3(3),
208     fSummary(0),
209     fEvent(.1, .8, "Event #"),
210     fFMD1Sum(.2, .7, "# in FMD1: "),
211     fFMD2Sum(.2, .6, "# in FMD2: "),
212     fFMD3Sum(.2, .5, "# in FMD3: "),
213     fLine(.15, .47, .85, .47),
214     fTotal(.2, .35, "Total:   ")
215 {
216   // Constructor. 
217   // 
218   // Parameters 
219   // 
220   //   gAliceFile       The galice.root file to use - if any. 
221   // 
222
223   // RemoveLoad(kGeometry);
224   fEvent.SetBit(TLatex::kTextNDC);
225   fFMD1Sum.SetBit(TLatex::kTextNDC);
226   fFMD2Sum.SetBit(TLatex::kTextNDC);
227   fFMD3Sum.SetBit(TLatex::kTextNDC);
228   fLine.SetBit(TLine::kLineNDC);
229   fTotal.SetBit(TLatex::kTextNDC);
230 }
231
232 //____________________________________________________________________
233 AliFMDPattern::~AliFMDPattern()
234 {
235   // DTOR 
236   // Free all allocated shapes. 
237   // note, that most members are real objects, so we do not need to
238   // deal with them here. 
239   fInners.Delete();
240   fOuters.Delete();
241 }
242
243     
244 //____________________________________________________________________
245 Bool_t 
246 AliFMDPattern::Init()
247 {
248   // Initialize.  Get transforms and such, 
249   if (!AliFMDInput::Init()) return kFALSE;
250   AliFMDGeometry* geom = AliFMDGeometry::Instance();
251   if (!geom) return kFALSE;
252   geom->Init();
253   geom->InitTransformations();
254   
255   Char_t rs[] = { 'I' , 'O', '\0' };
256   Char_t *r   = rs;
257   do {
258     AliFMDRing* ring = geom->GetRing(*r);
259     if (!ring) continue;
260     const TObjArray& vs = ring->GetVerticies();
261     TObjArray&       gs = (*r == 'I' ? fInners   : fOuters);
262     Float_t&         mr = (*r == 'I' ? fInnerMax : fOuterMax);
263     Int_t            nm = ring->GetNModules();
264     AliInfo(Form("Making %d modules for %c", nm, *r));
265     for (Int_t m = 0; m < nm; m++) {
266       Int_t          nv = vs.GetEntries();
267       Double_t       a  = TMath::Pi() / 180 * (m * 2 + 1) * ring->GetTheta();
268       TGraph*        g  = new TGraph(nv+1);
269       Double_t       x0 = 0, y0 = 0;
270       gs.AddAtAndExpand(g, m);
271       for (Int_t c = 0; c < nv; c++) {
272         TVector2* v = static_cast<TVector2*>(vs.At(c));
273         mr          = TMath::Max(mr, Float_t(v->Mod()));
274         TVector2  w(v->Rotate(a));
275         if (c == 0) { x0 = w.X(); y0 = w.Y(); }
276         g->SetPoint(c, w.X(), w.Y());
277       }
278       g->SetName(Form("FMDX%c_%02d", *r, m));
279       g->SetPoint(nv, x0, y0);
280       g->SetFillColor((*rs == 'I' ? 
281                        (m % 2 == 0 ? 18 : 17) :
282                        (m % 2 == 0 ? 20 : 23)));
283       g->SetFillStyle(3001);
284       g->SetLineColor(1);
285       g->SetLineWidth(1);
286       g->SetLineStyle(2);
287     }
288   } while (*(++r));
289     
290   return kTRUE;
291 }
292
293 //____________________________________________________________________
294 Bool_t 
295 AliFMDPattern::Begin(Int_t event) 
296 {
297   // Called at the begining of an event. 
298   // 
299   // Parameters 
300   //
301   //    EVENT           The event number 
302   // 
303   MakeAux();
304   if (!fCanvas) {
305     const char* which[] = { "Continue", "Redisplay", 0 };
306     MakeCanvas(which);
307     
308     AliFMDGeometry* geom = AliFMDGeometry::Instance();
309     // AliFMDDetector* det;
310     if ((/* det = */ geom->GetDetector(1))) {
311       fPad->cd();
312       fFMD1Pad = new TPad("FMD1", "FMD1", 0.0, 0.50, 0.5, 1.0, 0, 0);
313       fFMD1Pad->Draw();
314       fFMD1Pad->cd();
315       fFMD1.Begin(-1, fInnerMax, fInners, fOuters);
316     }
317     if ((/* det = */ geom->GetDetector(2))) {
318       fPad->cd();
319       fFMD2Pad = new TPad("FMD2", "FMD2", 0.5, 0.50, 1.0, 1.0, 0, 0);
320       fFMD2Pad->Draw();
321       fFMD2Pad->cd();
322       fFMD2.Begin(-1, fOuterMax, fInners, fOuters);
323     }
324     if ((/* det = */ geom->GetDetector(3))) {
325       fPad->cd();
326       fFMD3Pad = new TPad("FMD3", "FMD3", 0.0, 0.0, .5, .5, 0, 0);
327       fFMD3Pad->Draw();
328       fFMD3Pad->cd();
329       fFMD3.Begin(-1, fOuterMax, fInners, fOuters);
330     }
331     fPad->cd();
332     fSummary = new TPad("display", "Display", 0.5, 0.0, 1.0, 0.5, 0, 0);
333     fSummary->Draw();
334     fSummary->cd();
335     fEvent.Draw();
336     fFMD1Sum.Draw();
337     fFMD2Sum.Draw();
338     fFMD3Sum.Draw();
339     fLine.Draw();
340     fTotal.Draw();
341   }
342   fEvent.SetTitle(Form("Event # %6d", event));
343
344   fCanvas->Modified();
345   fCanvas->Update();
346   fCanvas->cd();
347   fFMD1.Clear();
348   fFMD2.Clear();
349   fFMD3.Clear();
350   return AliFMDInput::Begin(event);
351 }
352
353 //____________________________________________________________________
354 void
355 AliFMDPattern::Redisplay()
356 {
357   // Redraw the displayu 
358   fFMD1.Clear();
359   fFMD2.Clear();
360   fFMD3.Clear();
361   AliFMDDisplay::Redisplay();
362 }
363
364 //____________________________________________________________________
365 void
366 AliFMDPattern::AtEnd()
367 {
368   // Called at the end of an event. 
369   DrawAux();
370   
371   Int_t total = 0;
372   
373   fFMD1.End();
374   fFMD1Pad->Modified();
375   fFMD1Sum.SetTitle(Form("# hits in FMD1: %5d", fFMD1.Total()));
376   total += fFMD1.Total();
377
378   fFMD2.End();
379   fFMD2Pad->Modified();
380   fFMD2Sum.SetTitle(Form("# hits in FMD2: %5d", fFMD2.Total()));
381   total += fFMD2.Total();
382
383   fFMD3.End();
384   fFMD3Pad->Modified();
385   fFMD3Sum.SetTitle(Form("# hits in FMD3: %5d", fFMD3.Total()));
386   total += fFMD3.Total();
387
388   fTotal.SetTitle(Form("Total:    %5d/51200 (%3d%%)", 
389                       total, Int_t(100. / 51200 * total)));
390   fSummary->Modified();
391   fCanvas->Modified();
392   fCanvas->Update();
393   fCanvas->cd();
394 }
395
396 //____________________________________________________________________
397 Bool_t 
398 AliFMDPattern::ProcessHit(AliFMDHit* hit, TParticle*) 
399 {
400   // Process a hit. 
401   // 
402   // Parameters 
403   // 
404   //    HIT             The hit to process. 
405   // 
406   // The TParticle argument is never used. 
407   switch (hit->Detector()) {
408   case 1: fFMD1.AddMarker(hit->X(), hit->Y(), hit->Edep(), 1); break;
409   case 2: fFMD2.AddMarker(hit->X(), hit->Y(), hit->Edep(), 1); break;
410   case 3: fFMD3.AddMarker(hit->X(), hit->Y(), hit->Edep(), 1); break;
411   }
412   return kTRUE;
413 }
414
415
416 //____________________________________________________________________
417 void
418 AliFMDPattern::AddMarker(UShort_t det, Char_t rng, 
419                          UShort_t sec, UShort_t str,
420                          TObject*, Float_t s, Float_t max)
421 {
422   // Add a marker to the display
423   //
424   //    det     Detector
425   //    rng     Ring
426   //    sec     Sector 
427   //    str     Strip
428   //    o       Object to refer to
429   //    s       Signal 
430   //    max     Maximum of signal 
431   //
432   AliFMDPatternDetector* d = 0;
433   switch (det) {
434   case 1: d = &fFMD1; break;
435   case 2: d = &fFMD2; break;
436   case 3: d = &fFMD3; break;
437   }
438   if (!d) return;
439   AliFMDGeometry*   geom = AliFMDGeometry::Instance();
440   Double_t x, y, z;
441   geom->Detector2XYZ(det, rng, sec, str, x, y, z);
442   // Make code-checker shut the f**k up 
443   TRandom* rand = gRandom;
444   if (true) {
445     AliFMDRing* r  = geom->GetRing(rng);
446     Double_t    t  = .9 * r->GetTheta() / 2;
447     Double_t    a  = rand->Uniform(-t,t) * TMath::Pi() / 180;
448     Double_t    x1 = x * TMath::Cos(a) - y * TMath::Sin(a);
449     Double_t    y1 = x * TMath::Sin(a) + y * TMath::Cos(a);
450     x = x1;
451     y = y1;
452   }
453   d->AddMarker(x, y, s, max);
454 }
455
456 //____________________________________________________________________
457 //
458 // EOF
459 //