]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDFancy.h
DP:Misalignment of CPV added
[u/mrichter/AliRoot.git] / FMD / AliFMDFancy.h
1 #ifndef AliFMDFANCY_H
2 #define AliFMDFANCY_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights
4  * reserved. 
5  *
6  * See cxx source for full Copyright notice                               
7  */
8 /** @file    AliFMDFancy.h
9     @author  Christian Holm Christensen <cholm@nbi.dk>
10     @date    Mon Mar 27 12:39:09 2006
11     @brief   FMD Event display (as fancys)
12 */
13 //___________________________________________________________________
14 //
15 // The classes defined here, are utility classes for reading in data
16 // for the FMD.  They are  put in a seperate library to not polute the
17 // normal libraries.  The classes are intended to be used as base
18 // classes for customized class that do some sort of analysis on the
19 // various types of data produced by the FMD. 
20 //
21 #include "AliFMDDisplay.h"
22 #include <TObjArray.h>
23 #include <TGraph2D.h>
24 #include <TLatex.h>
25 #include <TLine.h>
26 class TCanvas;
27 class TPad;
28 class TH1;
29 class TH2;
30 class TH3;
31
32
33 //___________________________________________________________________
34 /** @class AliFMDFancy 
35     @brief Utility class to visualize FMD data in 2D. 
36     @ingroup FMD_util
37  */
38 class AliFMDFancy : public AliFMDDisplay
39 {
40 public:
41   struct Detector 
42   {
43     Detector(UShort_t id);
44     ~Detector();
45     
46     void Init();
47     void Begin(Int_t event=0);
48     void Clear(Int_t event=0);
49     void End();
50     void AddMarker(Char_t rng, UShort_t sec, UShort_t str, 
51                    Float_t v, Float_t max);
52     TH1*      fFrame;
53     Int_t     fId;
54     TObjArray fShapes;
55     Int_t     fNInnerHits;
56     TGraph2D  fInnerHits;
57     Int_t     fNOuterHits;
58     TGraph2D  fOuterHits;
59     Double_t  fMaxR;
60     Double_t  fMinZ;
61     Double_t  fMaxZ;
62   private:
63     void AddHistogram(TGraph2D& g, const char* opt="");
64     Detector(const Detector& );
65     Detector& operator=(const Detector& ) { return *this; }
66   };
67     
68   /** Constructor
69       @param gAliceFile galice file*/
70   AliFMDFancy(const char* gAliceFile="galice.root");
71   /** DTOR */
72   virtual ~AliFMDFancy();
73
74   /** Initialize
75       @return  @c false on error */
76   virtual Bool_t Init();
77   /** Called at beginning of an event 
78       @param event Event number
79       @return @c false on error  */
80   virtual Bool_t Begin(Int_t event);
81   /** Called at end of an event 
82       @return @c false on error  */
83   virtual Bool_t End();
84  protected:
85   AliFMDFancy(const AliFMDFancy& );
86   AliFMDFancy& operator=(const AliFMDFancy& ) { return *this; }
87   /** Add a marker to the display
88       @param det Detector
89       @param rng Ring
90       @param sec Sector 
91       @param str Strip
92       @param o   Object to refer to
93       @param s   Signal 
94       @param max Maximum of signal */
95   virtual void AddMarker(UShort_t det, Char_t rng, UShort_t sec, UShort_t str, 
96                          TObject* o, Float_t s, Float_t max);
97
98   virtual Bool_t ProcessHit(AliFMDHit* hit, TParticle*);
99   /** FMD1 Pad */
100   TPad*  fFMD1Pad;
101   /** FMD1 Frame */ 
102   Detector fFMD1;
103   /** FMD2 Pad  */
104   TPad*  fFMD2Pad;
105   /** FMD2 Frame */ 
106   Detector fFMD2;
107   /** FMD3 Pad */
108   TPad*  fFMD3Pad;
109   /** FMD3 Frame */ 
110   Detector fFMD3;
111   /** Summary pad */
112   TPad*    fSummary;
113   /** Text fields */
114   TLatex fEvent;
115   TLatex fFMD1IHits;
116   TLatex fFMD2IHits;
117   TLatex fFMD2OHits;
118   TLatex fFMD3IHits;
119   TLatex fFMD3OHits;
120   TLine  fLine;
121   TLatex fTotal;
122   
123   ClassDef(AliFMDFancy,0)
124 };
125
126
127 #endif
128 //____________________________________________________________________
129 //
130 // Local Variables:
131 //   mode: C++
132 // End:
133 //
134 // EOF
135 //