source: trunk/MagicSoft/Mars/mhflux/MHAlpha.h@ 7144

Last change on this file since 7144 was 7122, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 4.0 KB
Line 
1#ifndef MARS_MHAlpha
2#define MARS_MHAlpha
3
4#ifndef MARS_MH
5#include "MH.h"
6#endif
7
8#ifndef MARS_MTime
9#include "MTime.h"
10#endif
11
12#ifndef MARS_MAlphaFitter
13#include "MAlphaFitter.h"
14#endif
15
16#ifndef ROOT_TH3
17#include <TH3.h>
18#endif
19
20class MParList;
21class MParameterD;
22class MHillas;
23class MHMatrix;
24class MPointingPos;
25
26class MHAlpha : public MH
27{
28protected:
29 TH3D fHist; // Alpha vs. theta and energy
30 TH1D fHistTime; //! temporary histogram to get alpha vs. time
31
32 TString fNameParameter;
33
34 MParContainer *fParameter; //!
35
36 const TH3D *fOffData;
37
38 MAlphaFitter fFit; // SEEMS THAT STREAMER HAS SOME PROBLEMS... MAYBE IF FUNC IS USED AT THE SAME TIME FOR FITS (PAINT)
39
40private:
41 TH1D fHEnergy; // excess events vs energy
42 TH1D fHTheta; // excess events vs theta
43 TH1D fHTime; // excess events vs time;
44
45 MParameterD *fResult; //!
46 MParameterD *fEnergy; //!
47 MHillas *fHillas; //!
48 MPointingPos *fPointPos; //!
49
50 MTime *fTimeEffOn; //! Time to steer filling of fHTime
51 MTime *fTime; //! Event-Time to finalize fHTime
52 MTime fLastTime; //! Last fTimeEffOn
53
54 UShort_t fNumTimeBins; // Number of time bins to fill together
55 UShort_t fNumRebin; //!
56
57 //const TString fNameProjAlpha; //! This should make sure, that gROOT doen't confuse the projection with something else
58protected:
59 MHMatrix *fMatrix; //!
60 Int_t fMap[5]; //!
61
62 Bool_t fSkipHistTime;
63 Bool_t fSkipHistTheta;
64 Bool_t fSkipHistEnergy;
65 Bool_t fForceUsingSize;
66
67private:
68 Float_t FitEnergyBins(Bool_t paint=kFALSE);
69 void FitThetaBins(Bool_t paint=kFALSE);
70
71 void UpdateAlphaTime(Bool_t final=kFALSE);
72 void InitAlphaTime(const MTime &t);
73 void FinalAlphaTime(MBinning &bins);
74
75 void PaintText(Double_t val, Double_t error) const;
76
77 Int_t DistancetoPrimitive(Int_t px, Int_t py);
78
79 virtual Bool_t GetParameter(const MParList &pl);
80 virtual Double_t GetVal() const;
81 virtual const char *GetParameterRule() const
82 {
83 return "MHillasSrc.fAlpha";
84 }
85
86public:
87 MHAlpha(const char *name=NULL, const char *title=NULL);
88
89 // MH
90 Bool_t SetupFill(const MParList *pl);
91 Bool_t Fill(const MParContainer *par, const Stat_t w=1);
92 Bool_t Finalize();
93
94 // MHAlpha
95 const TH3D &GetHist() const { return fHist; }
96 const MAlphaFitter &GetAlphaFitter() const { return fFit; }
97
98 const TH1D &GetHEnergy() const { return fHEnergy; }
99
100 // Setter
101 void SetNameParameter(const char *name) { fNameParameter=name; }
102 void SetOffData(const MHAlpha &h)
103 {
104 fOffData = &h.fHist;
105 fForceUsingSize = h.fForceUsingSize;
106 fNumTimeBins = h.fNumTimeBins;
107 }
108 void SetNumTimeBins(UShort_t n) { fNumTimeBins = n; }
109/*
110 void SetSizeCuts(Float_t min, Float_t max) { fSizeMin=min; fSizeMax=max; }
111 void SetSizeMin(Float_t min) { fSizeMin=min; }
112 void SetSizeMax(Float_t max) { fSizeMax=max; }
113 void SetEnergyCuts(Float_t min, Float_t max) { fEnergyMin=min; fEnergyMax=max; }
114 void SetEnergyMin(Float_t min) { fEnergyMin=min; }
115 void SetEnergyMax(Float_t max) { fEnergyMax=max; }
116
117 void SetCuts(const MHAlpha &h) {
118 fSizeMin = h.fSizeMin; fEnergyMin = h.fEnergyMin;
119 fSizeMax = h.fSizeMax; fEnergyMax = h.fEnergyMax;
120 }
121 */
122
123 void SkipHistTime(Bool_t b=kTRUE) { fSkipHistTime=b; }
124 void SkipHistTheta(Bool_t b=kTRUE) { fSkipHistTheta=b; }
125 void SkipHistEnergy(Bool_t b=kTRUE) { fSkipHistEnergy=b; }
126 void ForceUsingSize(Bool_t b=kTRUE) { fForceUsingSize=b; }
127
128 void DrawAll(); //*MENU*
129
130 virtual void InitMapping(MHMatrix *mat, Int_t type=0);
131 void StopMapping();
132
133 // TObject
134 void Paint(Option_t *opt="");
135 void Draw(Option_t *option="");
136
137 // MParContainer
138 Int_t ReadEnv(const TEnv &env, TString prefix, Bool_t print=kFALSE);
139
140 ClassDef(MHAlpha, 1) // Alpha-Plot which is fitted online
141};
142
143#endif
Note: See TracBrowser for help on using the repository browser.