Changeset 2365 for trunk/MagicSoft/Mars/mhistmc
- Timestamp:
- 09/29/03 17:39:49 (21 years ago)
- Location:
- trunk/MagicSoft/Mars/mhistmc
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/MagicSoft/Mars/mhistmc/MHMcTriggerLvl2.cc
r2173 r2365 33 33 #include "MHMcTriggerLvl2.h" 34 34 35 #include <math.h> 36 35 #include <TMath.h> 36 37 #include <TH2.h> 37 38 #include <TH1.h> 38 39 #include <TF1.h> … … 48 49 49 50 #include "MMcTriggerLvl2.h" 50 #include "MGeomCam.h"51 #include "MBinning.h"52 51 53 52 using namespace std; 54 55 /*56 Please, DON'T USE IFDEFS IN SUCH A CONTEXT, Thomas.57 --------------------------------------------------58 #ifndef COLOR_LINELPS59 #define COLOR_LINELPS Int_t colorlps = 160 COLOR_LINELPS;61 #endif62 63 #ifndef COLOR_LINESBC64 #define COLOR_LINESBC Int_t colorsbc = 165 COLOR_LINESBC;66 #endif67 68 #ifndef COLOR_LINEPS69 #define COLOR_LINEPS Int_t colorps = 170 COLOR_LINEPS;71 #endif72 */73 53 74 54 /* Use this insteadif you want to have some value which is the same for all … … 80 60 Int_t MHMcTriggerLvl2::fColorSbc = 1; 81 61 Int_t MHMcTriggerLvl2::fColorPs = 1; 62 Int_t MHMcTriggerLvl2::fColorPsE = 1; 82 63 83 64 ClassImp(MHMcTriggerLvl2); … … 95 76 fTitle = title ? title : "Trigger L2 image parameters"; 96 77 97 fHistLutPseudoSize = new TH1F("fHistLutPseudoSize", "number of compact pixels in one lut", 1 3, 0, 12);98 fHistPseudoSize = new TH1F("fHistPseudoSize", "Multiplicity of the cluster identified by the L2T", 41, 0, 40);99 fHistSizeBiggerCell = new TH1F("fHistSizeBiggerCell", "Number of fired pixel in bigger cell", 3 7, 0, 36);100 101 fHistLutPseudoSizeNorm = new TH1F("fHistLutPseudoSizeNorm", "Normalized Number of compact pixels in one lut", 1 3, 0, 12);102 fHistPseudoSizeNorm = new TH1F("fHistPseudoSizeNorm", "Normalized Multiplicity of the cluster identified by the L2T", 41, 0, 40);103 fHistSizeBiggerCellNorm = new TH1F("fHistSizeBiggerCellNorm", "Normalized Number of fired pixel in bigger cell", 3 7, 0, 36);78 fHistLutPseudoSize = new TH1F("fHistLutPseudoSize", "number of compact pixels in one lut", 12, 0, 12); 79 fHistPseudoSize = new TH1F("fHistPseudoSize", "Multiplicity of the cluster identified by the L2T", 397, 0, 397); 80 fHistSizeBiggerCell = new TH1F("fHistSizeBiggerCell", "Number of fired pixel in bigger cell", 36, 0, 36); 81 82 fHistLutPseudoSizeNorm = new TH1F("fHistLutPseudoSizeNorm", "Normalized Number of compact pixels in one lut", 12, 0, 12); 83 fHistPseudoSizeNorm = new TH1F("fHistPseudoSizeNorm", "Normalized Multiplicity of the cluster identified by the L2T", 397, 0, 397); 84 fHistSizeBiggerCellNorm = new TH1F("fHistSizeBiggerCellNorm", "Normalized Number of fired pixel in bigger cell", 36, 0, 36); 104 85 105 86 fHistLutPseudoSize->SetDirectory(NULL); … … 124 105 fHistSizeBiggerCellNorm->SetYTitle("Counts/Total Counts"); 125 106 126 fFNorm = new TF1("FNorm", "1", -1, 40); 107 fHistPseudoSizeEnergy = new TH2D("fHistPseudoSizeEnergy","Ps Size vs Energy", 40, 1, 5, 397, 0,397); 108 109 fHistPseudoSizeEnergy->SetName("fHistPseudoSizeEnergy"); 110 fHistPseudoSizeEnergy->SetTitle("PseudoSize vs. Energy"); 111 fHistPseudoSizeEnergy->SetXTitle("Log(E[GeV])"); 112 fHistPseudoSizeEnergy->SetYTitle("PseudoSize"); 113 114 fFNorm = new TF1("FNorm", "1", -1, 397); 127 115 } 128 116 … … 139 127 delete fHistPseudoSizeNorm; 140 128 delete fHistSizeBiggerCellNorm; 129 delete fHistPseudoSizeEnergy; 141 130 142 131 delete fFNorm; 143 132 } 133 144 134 145 135 // -------------------------------------------------------------------------- … … 155 145 fHistPseudoSize->Fill(h.GetPseudoSize()); 156 146 fHistSizeBiggerCell->Fill(h.GetSizeBiggerCell()); 147 //fHistPseudoSizeEnergy->Fill(TMath::Log10(h.GetEnergy()), h.GetPseudoSize()); 148 fHistPseudoSizeEnergy->Fill(TMath::Log10(h.GetEnergy()), h.GetLutPseudoSize()); 157 149 158 150 return kTRUE; 159 151 } 152 153 160 154 161 155 … … 186 180 histNorm.SetLineColor(col); 187 181 histNorm.DrawCopy(same?"same":""); 182 183 return c; 184 } 185 186 187 // -------------------------------------------------------------------------- 188 // 189 // This is the private function member which draw a clone of a 2D-histogram. 190 // This method is called by the DrawClone method. 191 // 192 TObject *MHMcTriggerLvl2::Draw2DHist(TH1 &hist, const TString &canvasname, Int_t &col) const 193 { 194 col++; 195 196 TCanvas *c = (TCanvas*)gROOT->FindObject(canvasname); 197 198 Bool_t same = kTRUE; 199 if (!c) 200 { 201 c = MakeDefCanvas(canvasname,canvasname, 800, 600); 202 same = kFALSE; 203 } 204 205 hist.SetLineColor(col); 206 hist.DrawCopy(same?"same":""); 188 207 189 208 return c; … … 213 232 if (!str.Contains("lps", TString::kIgnoreCase) && 214 233 !str.Contains("sbc", TString::kIgnoreCase) && 215 !str.Contains("ps", TString::kIgnoreCase)) 234 !str.Contains("ps", TString::kIgnoreCase) && 235 !str.Contains("energy", TString::kIgnoreCase)) 216 236 { 217 237 *fLog << "ARGH!@! Possible options are \"lps\", \"sbc\", \"ps\" or NULL!" <<endl; … … 219 239 } 220 240 221 TH1 *hist=NormalizeHist(fHistLutPseudoSizeNorm, fHistLutPseudoSize); 222 223 if (!hist) 224 return NULL; 225 226 if (str.Contains("lps",TString::kIgnoreCase)) 227 return DrawHist(*fHistLutPseudoSize, *hist, "CanvasLPS", fColorLps); 228 229 if (str.Contains("sbc",TString::kIgnoreCase)) 230 return DrawHist(*fHistSizeBiggerCell, *hist, "CanvasSBC", fColorSbc); 231 232 if (str.Contains("ps",TString::kIgnoreCase)) 233 return DrawHist(*fHistPseudoSize, *hist, "CanvasPS", fColorPs); 241 if (str.Contains("lps",TString::kIgnoreCase)){ 242 TH1 *hist=NormalizeHist(fHistLutPseudoSizeNorm, fHistLutPseudoSize); 243 return DrawHist(*fHistLutPseudoSize, *hist, "CanvasLPS", fColorLps); 244 } 245 246 if (str.Contains("sbc",TString::kIgnoreCase)){ 247 TH1 *hist=NormalizeHist(fHistSizeBiggerCellNorm, fHistSizeBiggerCell); 248 return DrawHist(*fHistSizeBiggerCell, *hist, "CanvasSBC", fColorSbc); 249 } 250 251 if (str.Contains("ps",TString::kIgnoreCase)){ 252 TH1 *hist=NormalizeHist(fHistPseudoSizeNorm, fHistPseudoSize); 253 return DrawHist(*fHistPseudoSize, *hist, "CanvasPS", fColorPs); 254 } 255 256 if (str.Contains("energy",TString::kIgnoreCase)) 257 return Draw2DHist(*fHistPseudoSizeEnergy, "CanvasPSE", fColorPsE); 234 258 235 259 return NULL; -
trunk/MagicSoft/Mars/mhistmc/MHMcTriggerLvl2.h
r2043 r2365 6 6 #endif 7 7 8 class TH2D; 8 9 class TH1F; 9 10 class TF1; … … 20 21 TH1F *fHistSizeBiggerCell; // Histogram of fSizeBiggerCell 21 22 TH1F *fHistSizeBiggerCellNorm; // Histogram of fSizeBiggerCell normalized on integral of distribution 22 23 TH2D *fHistPseudoSizeEnergy; // 2D-Histogram of fPseudoSize vs. Energy 23 24 TF1* fFNorm; // Function used to normalize histograms 24 25 … … 26 27 static Int_t fColorSbc; 27 28 static Int_t fColorPs; 29 static Int_t fColorPsE; 28 30 29 31 TObject *DrawHist(TH1 &hist, TH1 &histNorm, const TString &canvasname, Int_t &colore) const; 32 TObject *Draw2DHist(TH1 &hist, const TString &canvasname, Int_t &col) const; 30 33 31 34 public: … … 34 37 35 38 Bool_t Fill(const MParContainer *par, const Stat_t w=1); 36 39 // Bool_t Fill(const Double_t energy, const MParContainer *par, const Stat_t w=1); 37 40 TH1 *GetHistByName(const TString name); 38 41 … … 43 46 TH1F *GetHistSizeBiggerCell() const { return fHistSizeBiggerCell; } 44 47 TH1F *GetHistSizeBiggerCellNorm() const { return fHistSizeBiggerCellNorm; } 48 TH2D *GetHistPseudoSizeEnergy() const { return fHistPseudoSizeEnergy; } 45 49 46 50 void Draw(Option_t *opt=NULL);
Note:
See TracChangeset
for help on using the changeset viewer.