]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDPattern.cxx
Compilation of AliAnalysisGoodies.cxx only if Root was compiled with XML support
[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
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
55 //____________________________________________________________________
56 ClassImp(AliFMDPattern)
57 #if 0
58   ; // This is here to keep Emacs for indenting the next line
59 #endif
60
61 //____________________________________________________________________
62 AliFMDPattern::Detector::Detector(UShort_t id) 
63   : fId(id),
64     fCounts(0), 
65     fGraphs(0), 
66     fFrame(0)
67 {}
68
69 //____________________________________________________________________
70 AliFMDPattern::Detector::~Detector()
71 {
72   if (fFrame) delete fFrame;
73 }
74
75 //____________________________________________________________________
76 void
77 AliFMDPattern::Detector::DrawShape(TObjArray& a) 
78 {
79   TIter next(&a);
80   TGraph* g = 0;
81   while ((g = static_cast<TGraph*>(next()))) {
82     g->DrawClone("f same");
83     g->DrawClone("l same");
84   }
85 }
86
87 //____________________________________________________________________
88 void
89 AliFMDPattern::Detector::Begin(Int_t nlevel, Double_t r, 
90                                TObjArray& inners, TObjArray& outers)
91 {
92   if (nlevel < 1) nlevel = gStyle->GetNumberOfColors();
93   fCounts.Set(nlevel);
94   if (!fFrame) {
95     fFrame = new TH2F(Form("fmd%dFrame", fId), Form("FMD%d", fId), 
96                       10, -r, r, 10, -r, r);
97     fFrame->SetStats(kFALSE);
98     fFrame->Draw();
99   }
100   DrawShape(inners);
101   if (fId != 1) DrawShape(outers);
102   for (Int_t i = 0; i < nlevel; i++) { 
103     TGraph* g = new TGraph;
104     Int_t idx = Int_t(Float_t(i) / nlevel * gStyle->GetNumberOfColors());
105     Int_t col = gStyle->GetColorPalette(idx);
106     g->SetName(Form("FMD%d_L%02d", fId, i));
107     g->SetMarkerColor(col);
108     g->SetLineColor(col);
109     g->SetFillColor(col);
110     g->SetMarkerSize(i * .2 + .2);
111     g->SetMarkerStyle(2);
112     g->Draw("same p");
113     fGraphs.AddAtAndExpand(g, i);
114   }
115   TIter   next(&fGraphs);
116 }
117
118 //____________________________________________________________________
119 void
120 AliFMDPattern::Detector::Clear() 
121 {
122   fCounts.Reset(0);
123 }
124
125 //____________________________________________________________________
126 void
127 AliFMDPattern::Detector::End()
128 {
129   TIter   next(&fGraphs);
130   TGraph* g = 0;
131   Int_t   i = 0;
132   while ((g = static_cast<TGraph*>(next()))) g->Set(fCounts[i++]);
133 }
134 //____________________________________________________________________
135 void
136 AliFMDPattern::Detector::AddMarker(Double_t x, Double_t y, Float_t s, 
137                                    Float_t max)
138 {
139   Int_t i = TMath::Min(Int_t(fCounts.fN * s / max), 
140                        Int_t(fGraphs.GetEntries()-1));
141   TGraph* g = static_cast<TGraph*>(fGraphs.At(i));
142   if (!g) return;
143   g->SetPoint(fCounts[i]++, x, y);
144 }
145
146
147 //____________________________________________________________________
148 AliFMDPattern::AliFMDPattern(const char* gAliceFile)
149   : AliFMDDisplay(kTRUE, gAliceFile),
150     fInnerMax(0), 
151     fOuterMax(0),
152     fFMD1Pad(0),
153     fFMD1(1),
154     fFMD2Pad(0),
155     fFMD2(2),
156     fFMD3Pad(0),
157     fFMD3(3),
158     fEvent(.1, .8, "Event #"),
159     fFMD1Sum(.2, .7, "# in FMD1: "),
160     fFMD2Sum(.2, .6, "# in FMD2: "),
161     fFMD3Sum(.2, .5, "# in FMD3: "),
162     fLine(.15, .47, .85, .47),
163     fTotal(.2, .35, "Total:   ")
164 {
165   // RemoveLoad(kGeometry);
166   fEvent.SetBit(TLatex::kTextNDC);
167   fFMD1Sum.SetBit(TLatex::kTextNDC);
168   fFMD2Sum.SetBit(TLatex::kTextNDC);
169   fFMD3Sum.SetBit(TLatex::kTextNDC);
170   fLine.SetBit(TLine::kLineNDC);
171   fTotal.SetBit(TLatex::kTextNDC);
172 }
173
174 //____________________________________________________________________
175 AliFMDPattern::~AliFMDPattern()
176 {
177   fInners.Delete();
178   fOuters.Delete();
179 }
180
181     
182 //____________________________________________________________________
183 Bool_t 
184 AliFMDPattern::Init()
185 {
186   // Initialize.  GEt transforms and such, 
187   if (!AliFMDInput::Init()) return kFALSE;
188   AliFMDGeometry* geom = AliFMDGeometry::Instance();
189   geom->Init();
190   geom->InitTransformations();
191   
192   Char_t rs[] = { 'I' , 'O', '\0' };
193   Char_t *r   = rs;
194   do {
195     AliFMDRing* ring = geom->GetRing(*r);
196     if (!ring) continue;
197     const TObjArray& vs = ring->GetVerticies();
198     TObjArray&       gs = (*r == 'I' ? fInners   : fOuters);
199     Float_t&         mr = (*r == 'I' ? fInnerMax : fOuterMax);
200     Int_t            nm = ring->GetNModules();
201     AliInfo(Form("Making %d modules for %c", nm, *r));
202     for (Int_t m = 0; m < nm; m++) {
203       Int_t          nv = vs.GetEntries();
204       Double_t       a  = TMath::Pi() / 180 * (m * 2 + 1) * ring->GetTheta();
205       TGraph*        g  = new TGraph(nv+1);
206       Double_t       x0 = 0, y0 = 0;
207       gs.AddAtAndExpand(g, m);
208       for (Int_t c = 0; c < nv; c++) {
209         TVector2* v = static_cast<TVector2*>(vs.At(c));
210         mr          = TMath::Max(mr, Float_t(v->Mod()));
211         TVector2  w(v->Rotate(a));
212         if (c == 0) { x0 = w.X(); y0 = w.Y(); }
213         g->SetPoint(c, w.X(), w.Y());
214       }
215       g->SetName(Form("FMDX%c_%02d", *r, m));
216       g->SetPoint(nv, x0, y0);
217       g->SetFillColor((*rs == 'I' ? 
218                        (m % 2 == 0 ? 18 : 17) :
219                        (m % 2 == 0 ? 20 : 23)));
220       g->SetFillStyle(3001);
221       g->SetLineColor(1);
222       g->SetLineWidth(1);
223       g->SetLineStyle(2);
224     }
225   } while (*(++r));
226     
227   return kTRUE;
228 }
229
230 //____________________________________________________________________
231 Bool_t 
232 AliFMDPattern::Begin(Int_t event) 
233 {
234   MakeAux();
235   if (!fCanvas) {
236     const char* which[] = { "Continue", "Redisplay", 0 };
237     MakeCanvas(which);
238     
239     AliFMDGeometry* geom = AliFMDGeometry::Instance();
240     AliFMDDetector* det;
241     if ((det = geom->GetDetector(1))) {
242       fPad->cd();
243       fFMD1Pad = new TPad("FMD1", "FMD1", 0.0, 0.50, 0.5, 1.0, 0, 0);
244       fFMD1Pad->Draw();
245       fFMD1Pad->cd();
246       fFMD1.Begin(-1, fInnerMax, fInners, fOuters);
247     }
248     if ((det = geom->GetDetector(2))) {
249       fPad->cd();
250       fFMD2Pad = new TPad("FMD2", "FMD2", 0.5, 0.50, 1.0, 1.0, 0, 0);
251       fFMD2Pad->Draw();
252       fFMD2Pad->cd();
253       fFMD2.Begin(-1, fOuterMax, fInners, fOuters);
254     }
255     if ((det = geom->GetDetector(3))) {
256       fPad->cd();
257       fFMD3Pad = new TPad("FMD3", "FMD3", 0.0, 0.0, .5, .5, 0, 0);
258       fFMD3Pad->Draw();
259       fFMD3Pad->cd();
260       fFMD3.Begin(-1, fOuterMax, fInners, fOuters);
261     }
262     fPad->cd();
263     fSummary = new TPad("display", "Display", 0.5, 0.0, 1.0, 0.5, 0, 0);
264     fSummary->Draw();
265     fSummary->cd();
266     fEvent.Draw();
267     fFMD1Sum.Draw();
268     fFMD2Sum.Draw();
269     fFMD3Sum.Draw();
270     fLine.Draw();
271     fTotal.Draw();
272   }
273   fEvent.SetTitle(Form("Event # %6d", event));
274
275   fCanvas->Modified();
276   fCanvas->Update();
277   fCanvas->cd();
278   fFMD1.Clear();
279   fFMD2.Clear();
280   fFMD3.Clear();
281   return AliFMDInput::Begin(event);
282 }
283
284 //____________________________________________________________________
285 void
286 AliFMDPattern::Redisplay()
287 {
288   fFMD1.Clear();
289   fFMD2.Clear();
290   fFMD3.Clear();
291   AliFMDDisplay::Redisplay();
292 }
293
294 //____________________________________________________________________
295 void
296 AliFMDPattern::AtEnd()
297 {
298   DrawAux();
299   
300   Int_t total = 0;
301   
302   fFMD1.End();
303   fFMD1Pad->Modified();
304   fFMD1Sum.SetTitle(Form("# hits in FMD1: %5d", fFMD1.Total()));
305   total += fFMD1.Total();
306
307   fFMD2.End();
308   fFMD2Pad->Modified();
309   fFMD2Sum.SetTitle(Form("# hits in FMD2: %5d", fFMD2.Total()));
310   total += fFMD2.Total();
311
312   fFMD3.End();
313   fFMD3Pad->Modified();
314   fFMD3Sum.SetTitle(Form("# hits in FMD3: %5d", fFMD3.Total()));
315   total += fFMD3.Total();
316
317   fTotal.SetTitle(Form("Total:    %5d/51200 (%3d%%)", 
318                       total, Int_t(100. / 51200 * total)));
319   fSummary->Modified();
320   fCanvas->Modified();
321   fCanvas->Update();
322   fCanvas->cd();
323 }
324
325 //____________________________________________________________________
326 Bool_t 
327 AliFMDPattern::ProcessHit(AliFMDHit* hit, TParticle*) 
328 {
329   switch (hit->Detector()) {
330   case 1: fFMD1.AddMarker(hit->X(), hit->Y(), hit->Edep(), 1); break;
331   case 2: fFMD2.AddMarker(hit->X(), hit->Y(), hit->Edep(), 1); break;
332   case 3: fFMD3.AddMarker(hit->X(), hit->Y(), hit->Edep(), 1); break;
333   }
334   return kTRUE;
335 }
336
337
338 //____________________________________________________________________
339 void
340 AliFMDPattern::AddMarker(UShort_t det, Char_t rng, UShort_t sec, UShort_t str,
341                          TObject*, Float_t s, Float_t max)
342 {
343   // Add a marker to the display
344   //
345   //    det     Detector
346   //    rng     Ring
347   //    sec     Sector 
348   //    str     Strip
349   //    o       Object to refer to
350   //    s       Signal 
351   //    max     Maximum of signal 
352   //
353   Detector* d = 0;
354   switch (det) {
355   case 1: d = &fFMD1; break;
356   case 2: d = &fFMD2; break;
357   case 3: d = &fFMD3; break;
358   }
359   if (!d) return;
360   AliFMDGeometry*   geom = AliFMDGeometry::Instance();
361   Double_t x, y, z;
362   geom->Detector2XYZ(det, rng, sec, str, x, y, z);
363   if (true) {
364     AliFMDRing* r  = geom->GetRing(rng);
365     Double_t    t  = .9 * r->GetTheta() / 2;
366     Double_t    a  = gRandom->Uniform(-t,t) * TMath::Pi() / 180;
367     Double_t    x1 = x * TMath::Cos(a) - y * TMath::Sin(a);
368     Double_t    y1 = x * TMath::Sin(a) + y * TMath::Cos(a);
369     x = x1;
370     y = y1;
371   }
372   d->AddMarker(x, y, s, max);
373 }
374
375 //____________________________________________________________________
376 //
377 // EOF
378 //