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