source: trunk/MagicSoft/Mars/mtemp/mifae/library/MHDisp.cc@ 5923

Last change on this file since 5923 was 5923, checked in by domingo, 20 years ago
*** empty log message ***
File size: 23.0 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Eva Domingo, 12/2004 <mailto:domingo@ifae.es>
19! Wolfgang Wittek, 12/2004 <mailto:wittek@mppmu.mpg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2005
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27// //
28// MHDisp //
29// //
30// container holding the histograms for Disp //
31// also computes the minimization parameter of the Disp optimization //
32// //
33/////////////////////////////////////////////////////////////////////////////
34#include "MHDisp.h"
35
36#include <math.h>
37
38#include <TH1.h>
39#include <TH2.h>
40#include <TProfile.h>
41#include <TArrayI.h>
42#include <TPad.h>
43#include <TCanvas.h>
44#include <TStyle.h>
45
46#include "MGeomCam.h"
47#include "MSrcPosCam.h"
48#include "MHillas.h"
49#include "MHillasExt.h"
50#include "MNewImagePar.h"
51#include "MMcEvt.hxx"
52#include "MPointingPos.h"
53#include "MImageParDisp.h"
54
55#include "MHMatrix.h"
56#include "MParameters.h"
57
58#include "MParList.h"
59
60#include "MLog.h"
61#include "MLogManip.h"
62
63ClassImp(MHDisp);
64
65using namespace std;
66
67enum dispVar_t {kX,kY,kMeanX,kMeanY,kDelta,kSize,kM3Long,kAsym,
68 kEnergy,kImpact,kLongitmax,kZd,kAz,kTotVar}; // enum variables for the
69 // matrix columns mapping
70
71// --------------------------------------------------------------------------
72//
73// Constructor
74//
75MHDisp::MHDisp(const char *imagepardispname,
76 const char *name, const char *title)
77 : fImageParDispName(imagepardispname)
78{
79 fName = name ? name : "MHDisp";
80 fTitle = title ? title : "Histograms for Disp";
81
82 fSelectedPos = 1; // default MHDisp flag (selects Correct Disp source position solution)
83
84 fMatrix = NULL;
85
86 // initialize mapping array dimension to the number of columns of fMatrix
87 fMap.Set(kTotVar);
88
89 //--------------------------------------------------
90 // creating the Disp related histograms
91
92 fHistEnergy = new TH1F("fHistEnergy",
93 "log10(Energy)", 50, 1., 3.);
94 fHistEnergy->SetDirectory(NULL);
95 fHistEnergy->UseCurrentStyle();
96 fHistEnergy->SetXTitle("log10(Energy) [GeV]");
97 fHistEnergy->SetYTitle("# events");
98
99 fHistSize = new TH1F("fHistSize",
100 "log10(Size)", 50, 2., 4.);
101 fHistSize->SetDirectory(NULL);
102 fHistSize->UseCurrentStyle();
103 fHistSize->SetXTitle("log10(Size) [#phot]");
104 fHistSize->SetYTitle("# events");
105
106 fHistcosZA = new TH1F("fHistcosZA",
107 "cos(Zenith Angle)", 10, 0., 1.);
108 fHistcosZA->SetDirectory(NULL);
109 fHistcosZA->UseCurrentStyle();
110 fHistcosZA->SetXTitle("cos(Theta)");
111 fHistcosZA->SetYTitle("# events");
112
113 fSkymapXY = new TH2F("fSkymapXY",
114 "Disp estimated source positions Skymap", 71, -2., 2., 71, -2., 2.);
115 fSkymapXY->SetDirectory(NULL);
116 fSkymapXY->UseCurrentStyle();
117 fSkymapXY->SetXTitle("X Disp source position [deg]");
118 fSkymapXY->SetYTitle("Y Disp source position [deg]");
119
120 fHistMinPar = new TH1F("fHistMinPar" ,
121 "Distance^2 between Disp and real srcpos", 100, 0., 2.);
122 fHistMinPar->SetDirectory(NULL);
123 fHistMinPar->UseCurrentStyle();
124 fHistMinPar->SetXTitle("Minimization parameter = d^2 Disp to real srcpos [deg^2]");
125 fHistMinPar->SetYTitle("# events");
126
127 fHistDuDv = new TH2F("fHistDuDv",
128 "Du vs. Dv (distances between Disp and real srcpos)",
129 100, -2., 2., 100, -2., 2.);
130 fHistDuDv->SetDirectory(NULL);
131 fHistDuDv->UseCurrentStyle();
132 fHistDuDv->SetXTitle("Dv = transveral distance [deg]");
133 fHistDuDv->SetYTitle("Du = longitudinal distance [deg]");
134
135 fHistMinParEnergy = new TH2F("fHistMinParEnergy",
136 "Minimization parameter vs. Energy", 50, 1., 3., 100, 0., 2.);
137 fHistMinParEnergy->SetDirectory(NULL);
138 fHistMinParEnergy->UseCurrentStyle();
139 fHistMinParEnergy->SetXTitle("log10(Energy) [GeV]");
140 fHistMinParEnergy->SetYTitle("Minimization parameter = d^2 Disp to real srcpos [deg^2]");
141
142 fHistMinParSize = new TH2F("fHistMinParSize",
143 "Minimization parameter vs. Size", 50, 2., 4., 100, 0., 2.);
144 fHistMinParSize->SetDirectory(NULL);
145 fHistMinParSize->UseCurrentStyle();
146 fHistMinParSize->SetXTitle("log10(Size) [#phot]");
147 fHistMinParSize->SetYTitle("Minimization parameter = d^2 Disp to real srcpos [deg^2]");
148
149 fHistDuEnergy = new TH2F("fHistDuEnergy",
150 "Du vs. Energy", 50, 1., 3., 100, -2., 2.);
151 fHistDuEnergy->SetDirectory(NULL);
152 fHistDuEnergy->UseCurrentStyle();
153 fHistDuEnergy->SetXTitle("log10(Energy) [GeV]");
154 fHistDuEnergy->SetYTitle("Du = longitudinal distance Disp to real srcpos [deg]");
155
156 fHistDuSize = new TH2F("fHistDuSize",
157 "Du vs. Size", 50, 2., 4., 100, -2., 2.);
158 fHistDuSize->SetDirectory(NULL);
159 fHistDuSize->UseCurrentStyle();
160 fHistDuSize->SetXTitle("log10(Size) [#phot]");
161 fHistDuSize->SetYTitle("Du = longitudinal distance Disp to real srcpos [deg]");
162
163 fHistDvEnergy = new TH2F("fHistDvEnergy",
164 "Dv vs. Energy", 50, 1., 3., 100, -2., 2.);
165 fHistDvEnergy->SetDirectory(NULL);
166 fHistDvEnergy->UseCurrentStyle();
167 fHistDvEnergy->SetXTitle("log10(Energy) [GeV]");
168 fHistDvEnergy->SetYTitle("Dv = transveral distance Disp to real srcpos [deg]");
169
170 fHistDvSize = new TH2F("fHistDvSize",
171 "Dv vs. Size", 50, 2., 4., 100, -2., 2.);
172 fHistDvSize->SetDirectory(NULL);
173 fHistDvSize->UseCurrentStyle();
174 fHistDvSize->SetXTitle("log10(Size) [#phot]");
175 fHistDvSize->SetYTitle("Dv = transveral distance Disp to real srcpos [deg]");
176
177 fEvCorrAssign = new TProfile("fEvCorrAssign",
178 "Fraction events srcpos well assign vs. Size", 50, 2., 4., 0., 1.);
179 fEvCorrAssign->SetDirectory(NULL);
180 fEvCorrAssign->SetStats(0);
181 fEvCorrAssign->SetXTitle("log10(Size) [#phot]");
182 fEvCorrAssign->SetYTitle("Fraction of events with srcpos well assign");
183}
184
185
186// --------------------------------------------------------------------------
187//
188// Destructor
189//
190MHDisp::~MHDisp()
191{
192 delete fHistEnergy;
193 delete fHistSize;
194 delete fHistcosZA;
195 delete fSkymapXY;
196 delete fHistMinPar;
197 delete fHistDuDv;
198 delete fHistMinParEnergy;
199 delete fHistMinParSize;
200 delete fHistDuEnergy;
201 delete fHistDuSize;
202 delete fHistDvEnergy;
203 delete fHistDvSize;
204 delete fEvCorrAssign;
205}
206
207// --------------------------------------------------------------------------
208//
209// Set the pointers to the containers
210//
211//
212Bool_t MHDisp::SetupFill(const MParList *pList)
213{
214 // reset all histograms and Minimization parameter computing variables
215 // before each new eventloop
216 fNumEv = 0;
217 fSumMinPar = 0.;
218 fMinPar = (MParameterD*)const_cast<MParList*>(pList)->FindCreateObj("MParameterD", "MinimizationParameter");
219 if (!fMinPar)
220 {
221 *fLog << err << "MParameterD (MinimizationParameter) not found and could not be created... aborting."
222 << endl;
223 return kFALSE;
224 }
225 fMinPar->SetVal(0);
226
227 fHistEnergy->Reset();
228 fHistSize->Reset();
229 fHistcosZA->Reset();
230 fSkymapXY->Reset();
231 fHistMinPar->Reset();
232 fHistDuDv->Reset();
233 fHistMinParEnergy->Reset();
234 fHistMinParSize->Reset();
235 fHistDuEnergy->Reset();
236 fHistDuSize->Reset();
237 fHistDvEnergy->Reset();
238 fHistDvSize->Reset();
239 fEvCorrAssign->Reset();
240
241
242 // look for the defined camera geometry to get mm to deg conversion factor
243 MGeomCam *cam = (MGeomCam*)pList->FindObject("MGeomCam");
244 if (!cam)
245 {
246 *fLog << err << "MGeomCam (Camera Geometry) not found... aborting."
247 << endl;
248 return kFALSE;
249 }
250 fMm2Deg = cam->GetConvMm2Deg();
251
252
253 // look for Disp related containers
254 fImageParDisp = (MImageParDisp*)pList->FindObject(fImageParDispName,
255 "MImageParDisp");
256 if (!fImageParDisp)
257 {
258 *fLog << err << fImageParDispName
259 << " [MImageParDisp] not found... aborting." << endl;
260 return kFALSE;
261 }
262
263
264 //-----------------------------------------------------------
265 // if fMatrix exists, skip preprocessing and just read events from matrix;
266 // if doesn't exist, read variables values from containers, so look for them
267 if (fMatrix)
268 return kTRUE;
269
270 fSrcPos = (MSrcPosCam*)pList->FindObject("MSrcPosCam");
271 if (!fSrcPos)
272 {
273 *fLog << err << "MSrcPosCam not found... aborting." << endl;
274 return kFALSE;
275 }
276
277 fHil = (MHillas*)pList->FindObject("MHillas");
278 if (!fHil)
279 {
280 *fLog << err << "MHillas not found... aborting." << endl;
281 return kFALSE;
282 }
283
284 fHilExt = (MHillasExt*)pList->FindObject("MHillasExt");
285 if (!fHilExt)
286 {
287 *fLog << err << "MHillasExt not found... aborting." << endl;
288 return kFALSE;
289 }
290
291 fNewPar = (MNewImagePar*)pList->FindObject("MNewImagePar");
292 if (!fNewPar)
293 {
294 *fLog << err << "MNewImagePar not found... aborting." << endl;
295 return kFALSE;
296 }
297
298 fMcEvt = (MMcEvt*)pList->FindObject("MMcEvt");
299 if (!fMcEvt)
300 {
301 *fLog << err << "MMcEvt not found... aborting." << endl;
302 return kFALSE;
303 }
304
305 fPointing = (MPointingPos*)pList->FindObject("MPointingPos");
306 if (!fPointing)
307 {
308 *fLog << err << "MPointingPos not found... aborting." << endl;
309 return kFALSE;
310 }
311
312 //------------------------------------------
313
314 return kTRUE;
315}
316
317
318// --------------------------------------------------------------------------
319//
320// Set which selection algorithm for the Disp estimated source position
321// solutions we want to follow when filling the histograms:
322// * iflag == 1 --> Correct source position, the closest to the real srcpos
323// (only applicable to MC or well known source real data)
324// * iflag == 2 --> Wrong source position, the furthest to the real srcpos
325// (same applicability than case 1)
326// * iflag == 3 --> Source position assigned as M3Long sign indicates
327// * iflag == 4 --> Source position assigned as Asym sign indicates
328//
329void MHDisp::SetSelectedPos(Int_t iflag)
330{
331 fSelectedPos = iflag;
332}
333
334
335// --------------------------------------------------------------------------
336//
337// Sets the Matrix and the array of mapped values for each Matrix column
338//
339void MHDisp::SetMatrixMap(MHMatrix *matrix, TArrayI &map)
340{
341 fMatrix = matrix;
342 fMap = map;
343}
344
345
346// --------------------------------------------------------------------------
347//
348// Returns the Matrix and the mapped value for each Matrix column
349//
350void MHDisp::GetMatrixMap(MHMatrix* &matrix, TArrayI &map)
351{
352 map = fMap;
353 matrix = fMatrix;
354}
355
356
357// --------------------------------------------------------------------------
358//
359// You can use this function if you want to use a MHMatrix instead of the
360// given containers for the Disp optimization. This function adds all
361// necessary columns to the given matrix. Afterwards, you should fill
362// the matrix with the corresponding data (eg from a file by using
363// MHMatrix::Fill). Then, if you loop through the matrix (eg using
364// MMatrixLoop), MHDisp::Fill will take the values from the matrix
365// instead of the containers.
366//
367void MHDisp::InitMapping(MHMatrix *mat)
368{
369 if (fMatrix)
370 return;
371
372 fMatrix = mat;
373
374 fMap[kX] = fMatrix->AddColumn("MSrcPosCam.fX");
375 fMap[kY] = fMatrix->AddColumn("MSrcPosCam.fY");
376 fMap[kMeanX] = fMatrix->AddColumn("MHillas.fMeanX");
377 fMap[kMeanY] = fMatrix->AddColumn("MHillas.fMeanY");
378 fMap[kDelta] = fMatrix->AddColumn("MHillas.fDelta");
379
380 fMap[kSize] = fMatrix->AddColumn("MHillas.fSize");
381
382 fMap[kM3Long] = fMatrix->AddColumn("MHillasExt.fM3Long");
383 fMap[kAsym] = fMatrix->AddColumn("MHillasExt.fAsym");
384
385 fMap[kEnergy] = fMatrix->AddColumn("MMcEvt.fEnergy");
386 fMap[kImpact] = fMatrix->AddColumn("MMcEvt.fImpact");
387 fMap[kLongitmax] = fMatrix->AddColumn("MMcEvt.fLongitmax");
388
389 fMap[kZd] = fMatrix->AddColumn("MPointingPos.fZd");
390 fMap[kAz] = fMatrix->AddColumn("MPointingPos.fAz");
391}
392
393
394// --------------------------------------------------------------------------
395//
396// Returns a mapped value from the Matrix
397//
398Double_t MHDisp::GetVal(Int_t i) const
399{
400 Double_t val = (*fMatrix)[fMap[i]];
401
402 //*fLog << "MHDisp::GetVal; i, fMatrix, fMap, val = "
403 // << i << ", " << fMatrix << ", " << fMap[i] << ", " << val << endl;
404
405 return val;
406}
407
408
409// --------------------------------------------------------------------------
410//
411// Fill the histograms and sum up the Minimization paramter outcome
412// of each processed event
413//
414Bool_t MHDisp::Fill(const MParContainer *par, const Stat_t w)
415{
416 Double_t energy = 0.;
417 Double_t impact = 0.;
418 Double_t xmax = 0.;
419
420 if ( fMatrix || (!fMatrix && fMcEvt) )
421 {
422 energy = fMatrix ? GetVal(kEnergy) : fMcEvt->GetEnergy();
423 impact = fMatrix ? GetVal(kImpact) : fMcEvt->GetImpact();
424 xmax = fMatrix ? GetVal(kLongitmax) : fMcEvt->GetLongitmax();
425 }
426 else
427 *fLog << "No MMcEvt container available (no Energy,ImpactPar or Xmax)" << endl;
428
429 Double_t theta = fMatrix ? GetVal(kZd) : fPointing->GetZd();
430 // Double_t phi = fMatrix ? GetVal(kAz) : fPointing->GetAz();
431
432 Double_t xsrcpos0 = fMatrix ? GetVal(kX) : fSrcPos->GetX();
433 Double_t ysrcpos0 = fMatrix ? GetVal(kY) : fSrcPos->GetY();
434 Double_t xmean0 = fMatrix ? GetVal(kMeanX) : fHil->GetMeanX();
435 Double_t ymean0 = fMatrix ? GetVal(kMeanY) : fHil->GetMeanY();
436 Double_t delta = fMatrix ? GetVal(kDelta) : fHil->GetDelta();
437
438 Double_t size = fMatrix ? GetVal(kSize) : fHil->GetSize();
439
440 Double_t m3long = fMatrix ? GetVal(kM3Long) : fHilExt->GetM3Long();
441 Double_t asym = fMatrix ? GetVal(kAsym) : fHilExt->GetAsym();
442
443 //------------------------------------------
444 // convert to deg
445 Double_t xsrcpos = xsrcpos0 * fMm2Deg;
446 Double_t ysrcpos = ysrcpos0 * fMm2Deg;
447 Double_t xmean = xmean0 * fMm2Deg;
448 Double_t ymean = ymean0 * fMm2Deg;
449
450 // calculate cosinus of the angle between d and a vectors
451 Double_t a = (xmean-xsrcpos)*cos(delta) + (ymean-ysrcpos)*sin(delta);
452 Double_t b = sqrt( (xmean-xsrcpos)*(xmean-xsrcpos) + (ymean-ysrcpos)*(ymean-ysrcpos) );
453 Double_t cosda = a/b;
454
455 // sign of cosda
456 Int_t s0 = cosda>0 ? 1 : -1;
457
458 // get the sign of M3Long and Asym variables
459 Int_t sm3 = m3long>0 ? 1 : -1;
460 Int_t sa = asym>0 ? 1 : -1;
461
462 // set the sign "s" that will select one Disp source position for each
463 // shower, according to each of the possible algorithms for solution selection:
464 // SelectedPos = 1 means choose the right source position
465 // 2 wrong
466 // 3 the position according to M3Long
467 // 4 the position according to Asym
468 Int_t s = s0;
469 if (fSelectedPos == 1)
470 s = s0;
471 else if (fSelectedPos == 2)
472 s = -s0;
473 else if (fSelectedPos == 3)
474 s = sm3;
475 else if (fSelectedPos == 4)
476 s = sa;
477 else
478 *fLog << "Illegal value for Disp srcpos selection algorithm: "
479 << " fSelectedPos = " << fSelectedPos << endl;
480
481 // count the number of events where source position has been correctly assigned
482 if (s == s0)
483 fEvCorrAssign->Fill(log10(size), 1.);
484 else
485 fEvCorrAssign->Fill(log10(size), 0.);
486
487 // get estimated Disp value
488 Double_t disp = fImageParDisp->GetDisp();
489
490 //------------------------------------------
491 // Disp estimated source position
492 Double_t xdisp = xmean - s*cos(delta)*disp;
493 Double_t ydisp = ymean - s*sin(delta)*disp;
494 fSkymapXY->Fill(xdisp,ydisp);
495
496 // Distance between estimated Disp and real source position
497 Double_t d2 = (xdisp-xsrcpos)*(xdisp-xsrcpos) + (ydisp-ysrcpos)*(ydisp-ysrcpos);
498 fHistMinPar->Fill(d2);
499
500 // Longitudinal and transversal distances between Disp and real source positon
501 Double_t Du = -s*( (xdisp-xsrcpos)*cos(delta) + (ydisp-ysrcpos)*sin(delta));
502 Double_t Dv = -s*(-(xdisp-xsrcpos)*sin(delta) + (ydisp-ysrcpos)*cos(delta));
503 fHistDuDv->Fill(Dv,Du);
504
505 // Fill Energy, Size and ZA distributions
506 fHistEnergy->Fill(log10(energy));
507 fHistSize->Fill(log10(size));
508 fHistcosZA->Fill(cos(theta/kRad2Deg));
509
510 // to check the size and energy dependence of the optimization
511 fHistMinParEnergy->Fill(log10(energy),d2);
512 fHistMinParSize->Fill(log10(size),d2);
513 fHistDuEnergy->Fill(log10(energy),Du);
514 fHistDuSize->Fill(log10(size),Du);
515 fHistDvEnergy->Fill(log10(energy),Dv);
516 fHistDvSize->Fill(log10(size),Dv);
517
518 // variables for the Minimization parameter calculation (= d^2)
519 fNumEv += 1;
520 fSumMinPar += d2;
521
522 return kTRUE;
523}
524
525
526// --------------------------------------------------------------------------
527//
528// Calculates the final Minimization parameter of the Disp optimization
529//
530Bool_t MHDisp::Finalize()
531{
532 fMinPar->SetVal(fSumMinPar/fNumEv);
533 *fLog << endl;
534 *fLog << "MHDisp::Finalize: SumMinPar, NumEv = " << fSumMinPar << ", " << fNumEv << endl;
535 *fLog << "MHDisp::Finalize: Final MinPar = " << fMinPar->GetVal() << endl;
536
537 fMinPar->SetVal(fHistMinPar->GetMean());
538 *fLog << "MHDisp::Finalize: Final MinPar = " << fMinPar->GetVal() << endl;
539
540 SetReadyToSave();
541
542 return kTRUE;
543}
544
545
546// --------------------------------------------------------------------------
547//
548// Creates a new canvas and draws the Disp related histograms into it.
549// Be careful: The histograms belongs to this object and won't get deleted
550// together with the canvas.
551//
552void MHDisp::Draw(Option_t *opt)
553{
554 gStyle->SetPalette(1);
555
556 TString title = GetName();
557 title += ": ";
558 title += GetTitle();
559
560 TCanvas *pad = new TCanvas(GetName(),title,0,0,900,1500);
561 pad->SetBorderMode(0);
562 pad->Divide(3,5);
563
564 pad->cd(1);
565 gPad->SetBorderMode(0);
566 fHistcosZA->SetTitleOffset(1.2,"Y");
567 fHistcosZA->DrawClone(opt);
568 fHistcosZA->SetBit(kCanDelete);
569 gPad->Modified();
570
571 pad->cd(2);
572 gPad->SetBorderMode(0);
573 fHistEnergy->SetTitleOffset(1.2,"Y");
574 fHistEnergy->DrawClone(opt);
575 fHistEnergy->SetBit(kCanDelete);
576 gPad->Modified();
577
578 pad->cd(3);
579 gPad->SetBorderMode(0);
580 fHistSize->SetTitleOffset(1.2,"Y");
581 fHistSize->DrawClone(opt);
582 fHistSize->SetBit(kCanDelete);
583 gPad->Modified();
584
585 pad->cd(4);
586 gPad->SetBorderMode(0);
587 fHistMinPar->SetTitleOffset(1.2,"Y");
588 fHistMinPar->DrawClone(opt);
589 fHistMinPar->SetBit(kCanDelete);
590 gPad->Modified();
591
592 TProfile *profMinParEnergy = fHistMinParEnergy->ProfileX("Minimization parameter vs. Energy",0,99999,"s");
593 profMinParEnergy->SetXTitle("log10(Energy) [GeV]");
594 profMinParEnergy->SetYTitle("Minimization parameter = d^2 Disp to real srcpos [deg^2]");
595 pad->cd(5);
596 gPad->SetBorderMode(0);
597 profMinParEnergy->SetTitleOffset(1.2,"Y");
598 profMinParEnergy->SetStats(0);
599 profMinParEnergy->DrawClone(opt);
600 profMinParEnergy->SetBit(kCanDelete);
601 gPad->Modified();
602
603 TProfile *profMinParSize = fHistMinParSize->ProfileX("Minimization parameter vs. Size",0,99999,"s");
604 profMinParSize->SetXTitle("log10(Size) [#phot]");
605 profMinParSize->SetYTitle("Minimization parameter = d^2 Disp to real srcpos [deg^2]");
606 pad->cd(6);
607 gPad->SetBorderMode(0);
608 profMinParSize->SetTitleOffset(1.2,"Y");
609 profMinParSize->SetStats(0);
610 profMinParSize->DrawClone(opt);
611 profMinParSize->SetBit(kCanDelete);
612 gPad->Modified();
613
614 pad->cd(7);
615 gPad->SetBorderMode(0);
616 fSkymapXY->SetTitleOffset(1.2,"Y");
617 fSkymapXY->DrawClone("COLZ");
618 fSkymapXY->SetBit(kCanDelete);
619 gPad->Modified();
620
621 pad->cd(8);
622 gPad->SetBorderMode(0);
623 fEvCorrAssign->SetTitleOffset(1.2,"Y");
624 fEvCorrAssign->DrawClone(opt);
625 fEvCorrAssign->SetBit(kCanDelete);
626 gPad->Modified();
627
628 pad->cd(9);
629 gPad->SetBorderMode(0);
630 fHistDuDv->SetTitleOffset(1.2,"Y");
631 fHistDuDv->DrawClone("COLZ");
632 fHistDuDv->SetBit(kCanDelete);
633 gPad->Modified();
634
635 TH1F *histDu = (TH1F*)fHistDuDv->ProjectionY("histDu");
636 histDu->SetTitle("Longitudinal distance Du");
637 histDu->SetXTitle("Du = longitudinal distance [deg]");
638 pad->cd(10);
639 gPad->SetBorderMode(0);
640 histDu->SetTitleOffset(1.2,"Y");
641 histDu->DrawClone(opt);
642 histDu->SetBit(kCanDelete);
643 gPad->Modified();
644
645 TProfile *profDuEnergy = fHistDuEnergy->ProfileX("Du vs. Energy",0,99999,"s");
646 profDuEnergy->SetXTitle("log10(Energy) [GeV]");
647 profDuEnergy->SetYTitle("Du = longitudinal distance [deg]");
648 pad->cd(11);
649 gPad->SetBorderMode(0);
650 profDuEnergy->SetTitleOffset(1.2,"Y");
651 profDuEnergy->SetStats(0);
652 profDuEnergy->DrawClone(opt);
653 profDuEnergy->SetBit(kCanDelete);
654 gPad->Modified();
655
656 TProfile *profDuSize = fHistDuSize->ProfileX("Du vs. Size",0,99999,"s");
657 profDuSize->SetXTitle("log10(Size) [#phot]");
658 profDuSize->SetYTitle("Du = longitudinal distance [deg]");
659 pad->cd(12);
660 gPad->SetBorderMode(0);
661 profDuSize->SetTitleOffset(1.2,"Y");
662 profDuSize->SetStats(0);
663 profDuSize->DrawClone(opt);
664 profDuSize->SetBit(kCanDelete);
665 gPad->Modified();
666
667 TH1F *histDv = (TH1F*)fHistDuDv->ProjectionX("histDv");
668 histDv->SetTitle("Transversal distance Dv");
669 histDv->SetXTitle("Dv = transversal distance [deg]");
670 pad->cd(13);
671 gPad->SetBorderMode(0);
672 histDv->SetTitleOffset(1.2,"Y");
673 histDv->DrawClone(opt);
674 histDv->SetBit(kCanDelete);
675 gPad->Modified();
676
677 TProfile *profDvEnergy = fHistDvEnergy->ProfileX("Dv vs. Energy",0,99999,"s");
678 profDvEnergy->SetXTitle("log10(Energy) [GeV]");
679 profDvEnergy->SetYTitle("Dv = transversal distance [deg]");
680 pad->cd(14);
681 gPad->SetBorderMode(0);
682 profDvEnergy->SetTitleOffset(1.2,"Y");
683 profDvEnergy->SetStats(0);
684 profDvEnergy->DrawClone(opt);
685 profDvEnergy->SetBit(kCanDelete);
686 gPad->Modified();
687
688 TProfile *profDvSize = fHistDvSize->ProfileX("Dv vs. Size",0,99999,"s");
689 profDvSize->SetXTitle("log10(Size) [#phot]");
690 profDvSize->SetYTitle("Dv = transversal distance [deg]");
691 pad->cd(15);
692 gPad->SetBorderMode(0);
693 profDvSize->SetTitleOffset(1.2,"Y");
694 profDvSize->SetStats(0);
695 profDvSize->DrawClone(opt);
696 profDvSize->SetBit(kCanDelete);
697 gPad->Modified();
698}
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
Note: See TracBrowser for help on using the repository browser.