source: trunk/MagicSoft/Mars/manalysis/MPadSchweizer.cc@ 2037

Last change on this file since 2037 was 2025, checked in by wittek, 22 years ago
*** empty log message ***
File size: 21.2 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): Robert Wagner <mailto:magicsoft@rwagner.de> 10/2002
19! Author(s): Wolfgang Wittek <mailto:wittek@mppmu.mpg.de> 02/2003
20!
21! Copyright: MAGIC Software Development, 2000-2003
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MPadSchweizer
29//
30// This task applies padding such that for a given pixel and for a given
31// Theta bin the resulting distribution of the pedestal sigma is identical
32// to the distributions given by fHSigmaPixTheta and fHDiffPixTheta.
33//
34// The number of photons, its error and the pedestal sigmas are altered.
35// On average, the number of photons added is zero.
36//
37// The formulas used can be found in Thomas Schweizer's Thesis,
38// Section 2.2.1
39//
40// There are 2 options for the padding :
41//
42// 1) fPadFlag = 1 :
43// Generate first a Sigmabar using the 2D-histogram Sigmabar vs. Theta
44// (fHSigmaTheta). Then generate a pedestal sigma for each pixel using
45// the 3D-histogram Theta, pixel no., Sigma^2-Sigmabar^2
46// (fHDiffPixTheta).
47//
48// This is the preferred option as it takes into account the
49// correlations between the Sigma of a pixel and Sigmabar.
50//
51// 2) fPadFlag = 2 :
52// Generate a pedestal sigma for each pixel using the 3D-histogram
53// Theta, pixel no., Sigma (fHSigmaPixTheta).
54//
55//
56// The padding has to be done before the image cleaning because the
57// image cleaning depends on the pedestal sigmas.
58//
59// For random numbers gRandom is used.
60//
61// This implementation has been tested for CT1 data. For MAGIC some
62// modifications are necessary.
63//
64/////////////////////////////////////////////////////////////////////////////
65#include "MPadSchweizer.h"
66
67#include <stdio.h>
68
69#include <TH1.h>
70#include <TH2.h>
71#include <TH3.h>
72#include <TRandom.h>
73#include <TCanvas.h>
74
75#include "MBinning.h"
76#include "MSigmabar.h"
77#include "MMcEvt.hxx"
78#include "MLog.h"
79#include "MLogManip.h"
80#include "MParList.h"
81#include "MGeomCam.h"
82
83#include "MCerPhotPix.h"
84#include "MCerPhotEvt.h"
85
86#include "MPedestalCam.h"
87#include "MPedestalPix.h"
88
89ClassImp(MPadSchweizer);
90
91// --------------------------------------------------------------------------
92//
93// Default constructor.
94//
95MPadSchweizer::MPadSchweizer(const char *name, const char *title)
96{
97 fName = name ? name : "MPadSchweizer";
98 fTitle = title ? title : "Task for the padding (Schweizer)";
99
100 fHSigmaTheta = NULL;
101 fHSigmaPixTheta = NULL;
102 fHDiffPixTheta = NULL;
103}
104
105// --------------------------------------------------------------------------
106//
107// Destructor. STLL DOESN'T DESTRUCT EVERYTHING
108//
109MPadSchweizer::~MPadSchweizer()
110{
111 *fLog << all << "WARNING: ~MPadSchweizer called: Potential Memory Leak" << endl;
112}
113
114// --------------------------------------------------------------------------
115//
116// Set the references to the histograms to be used in the padding
117//
118// fHSigmaTheta 2D-histogram (Theta, sigmabar)
119// fHSigmaPixTheta 3D-hiostogram (Theta, pixel, sigma)
120// fHDiffPixTheta 3D-histogram (Theta, pixel, sigma^2-sigmabar^2)
121//
122//
123void MPadSchweizer::SetHistograms(TH2D *hist2, TH3D *hist3, TH3D *hist3Diff)
124{
125 fHSigmaTheta = hist2;
126 fHSigmaPixTheta = hist3;
127 fHDiffPixTheta = hist3Diff;
128
129 fHSigmaTheta->SetDirectory(NULL);
130 fHSigmaPixTheta->SetDirectory(NULL);
131 fHDiffPixTheta->SetDirectory(NULL);
132
133 Print();
134}
135
136// --------------------------------------------------------------------------
137//
138// Set the option for the padding
139//
140// There are 2 options for the padding :
141//
142// 1) fPadFlag = 1 :
143// Generate first a Sigmabar using the 2D-histogram Sigmabar vs. Theta
144// (fHSigmaTheta). Then generate a pedestal sigma for each pixel using
145// the 3D-histogram Theta, pixel no., Sigma^2-Sigmabar^2
146// (fHDiffPixTheta).
147//
148// This is the preferred option as it takes into account the
149// correlations between the Sigma of a pixel and Sigmabar.
150//
151// 2) fPadFlag = 2 :
152// Generate a pedestal sigma for each pixel using the 3D-histogram
153// Theta, pixel no., Sigma (fHSigmaPixTheta).
154//
155void MPadSchweizer::SetPadFlag(Int_t padflag)
156{
157 fPadFlag = padflag;
158 *fLog << "MPadSchweizer::SetPadFlag(); choose option " << fPadFlag << endl;
159}
160
161// --------------------------------------------------------------------------
162//
163// Set the pointers and prepare the histograms
164//
165Bool_t MPadSchweizer::PreProcess(MParList *pList)
166{
167 if ( !fHSigmaTheta || !fHSigmaPixTheta || !fHDiffPixTheta)
168 {
169 *fLog << err << "At least one of the histograms needed for the padding is not defined ... aborting." << endl;
170 return kFALSE;
171 }
172
173 fMcEvt = (MMcEvt*)pList->FindObject("MMcEvt");
174 if (!fMcEvt)
175 {
176 *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
177 return kFALSE;
178 }
179
180 fPed = (MPedestalCam*)pList->FindObject("MPedestalCam");
181 if (!fPed)
182 {
183 *fLog << err << "MPedestalCam not found... aborting." << endl;
184 return kFALSE;
185 }
186
187 fCam = (MGeomCam*)pList->FindObject("MGeomCam");
188 if (!fCam)
189 {
190 *fLog << err << "MGeomCam not found (no geometry information available)... aborting." << endl;
191 return kFALSE;
192 }
193
194 fEvt = (MCerPhotEvt*)pList->FindObject("MCerPhotEvt");
195 if (!fEvt)
196 {
197 *fLog << err << "MCerPhotEvt not found... aborting." << endl;
198 return kFALSE;
199 }
200
201 fSigmabar = (MSigmabar*)pList->FindCreateObj("MSigmabar");
202 if (!fSigmabar)
203 {
204 *fLog << err << "MSigmabar not found... aborting." << endl;
205 return kFALSE;
206 }
207
208
209 //--------------------------------------------------------------------
210 // histograms for checking the padding
211 //
212 // plot pedestal sigmas
213 fHSigmaPedestal = new TH2D("SigPed","Sigma: after vs. before padding",
214 100, 0.0, 5.0, 100, 0.0, 5.0);
215 fHSigmaPedestal->SetXTitle("Pedestal sigma before padding");
216 fHSigmaPedestal->SetYTitle("Pedestal sigma after padding");
217
218 // plot no.of photons (before vs. after padding)
219 fHPhotons = new TH2D("Photons","Photons: after vs.before padding",
220 100, -10.0, 90.0, 100, -10, 90);
221 fHPhotons->SetXTitle("no.of photons before padding");
222 fHPhotons->SetYTitle("no.of photons after padding");
223
224 // plot of added NSB
225 fHNSB = new TH1D("NSB","Distribution of added NSB", 100, -10.0, 10.0);
226 fHNSB->SetXTitle("no.of added NSB photons");
227 fHNSB->SetYTitle("no.of pixels");
228
229
230 //--------------------------------------------------------------------
231
232 memset(fErrors, 0, sizeof(fErrors));
233
234 return kTRUE;
235}
236
237// --------------------------------------------------------------------------
238//
239// Do the Padding
240// idealy the events to be padded should have been generated without NSB
241//
242Bool_t MPadSchweizer::Process()
243{
244 //*fLog << "Entry MPadSchweizer::Process();" << endl;
245
246 Int_t rc=0;
247
248 const UInt_t npix = fEvt->GetNumPixels();
249
250 //fSigmabar->Calc(*fCam, *fPed, *fEvt);
251 //*fLog << "before padding : " << endl;
252 //fSigmabar->Print("");
253
254
255 //$$$$$$$$$$$$$$$$$$$$$$$$$$
256 // to simulate the situation that before the padding the NSB and
257 // electronic noise are zero : set Sigma = 0 for all pixels
258 //for (UInt_t i=0; i<npix; i++)
259 //{
260 // MCerPhotPix &pix = fEvt->operator[](i);
261 // Int_t j = pix.GetPixId();
262
263 // MPedestalPix &ppix = fPed->operator[](j);
264 // ppix.SetMeanRms(0.0);
265 //}
266 //$$$$$$$$$$$$$$$$$$$$$$$$$$
267
268 //-------------------------------------------
269 // Calculate average sigma of the event
270 //
271 Double_t sigbarold = fSigmabar->Calc(*fCam, *fPed, *fEvt);
272 Double_t sigbarold2 = sigbarold*sigbarold;
273 //fSigmabar->Print("");
274
275 if (sigbarold > 0)
276 {
277 //*fLog << "MPadSchweizer::Process(); sigmabar of event to be padded is > 0 : "
278 // << sigbarold << ". Stop event loop " << endl;
279 // input data should have sigmabar = 0; stop event loop
280
281 rc = 1;
282 fErrors[rc]++;
283 return kCONTINUE;
284 }
285
286 const Double_t theta = kRad2Deg*fMcEvt->GetTelescopeTheta();
287 // *fLog << "theta = " << theta << endl;
288
289
290 //-------------------------------------------
291 // for the current theta,
292 // generate a sigmabar according to the histogram fHSigmaTheta
293 //
294 Double_t sigmabar=0;
295 Int_t binNumber = fHSigmaTheta->GetXaxis()->FindBin(theta);
296
297 TH1D *hsigma;
298
299 if ( binNumber < 1 || binNumber > fHSigmaTheta->GetNbinsX() )
300 {
301 //*fLog << "MPadSchweizer::Process(); binNumber out of range : theta, binNumber = "
302 // << theta << ", " << binNumber << "; Skip event " << endl;
303 // event cannot be padded; skip event
304
305 rc = 2;
306 fErrors[rc]++;
307 return kCONTINUE;
308 }
309 else
310 {
311 hsigma = fHSigmaTheta->ProjectionY("", binNumber, binNumber, "");
312 if ( hsigma->GetEntries() == 0.0 )
313 {
314 *fLog << "MPadSchweizer::Process(); projection for Theta bin "
315 << binNumber << " has no entries; Skip event " << endl;
316 // event cannot be padded; skip event
317 delete hsigma;
318
319 rc = 3;
320 fErrors[rc]++;
321 return kCONTINUE;
322 }
323 else
324 {
325 sigmabar = hsigma->GetRandom();
326
327 //*fLog << "Theta, bin number = " << theta << ", " << binNumber
328 // << ", sigmabar = " << sigmabar << endl;
329 }
330 delete hsigma;
331 }
332 const Double_t sigmabar2 = sigmabar*sigmabar;
333
334 //-------------------------------------------
335
336 //*fLog << "MPadSchweizer::Process(); sigbarold, sigmabar = "
337 // << sigbarold << ", "<< sigmabar << endl;
338
339 // Skip event if target sigmabar is <= sigbarold
340 if (sigmabar <= sigbarold)
341 {
342 *fLog << "MPadSchweizer::Process(); target sigmabar is less than sigbarold : "
343 << sigmabar << ", " << sigbarold << ", Skip event" << endl;
344
345 rc = 4;
346 fErrors[rc]++;
347 return kCONTINUE;
348 }
349
350
351 //-------------------------------------------
352 //
353 // Calculate average number of NSB photons to be added (lambdabar)
354 // from the value of sigmabar,
355 // - making assumptions about the average electronic noise (elNoise2) and
356 // - using a fixed value (F2excess) for the excess noise factor
357
358 Double_t elNoise2; // [photons]
359 Double_t F2excess = 1.3;
360 Double_t lambdabar; // [photons]
361
362
363
364 Int_t binTheta = fHDiffPixTheta->GetXaxis()->FindBin(theta);
365 if (binTheta != binNumber)
366 {
367 cout << "The binnings of the 2 histograms are not identical; aborting"
368 << endl;
369 return kERROR;
370 }
371
372 // Get RMS of (Sigma^2-sigmabar^2) in this Theta bin.
373 // The average electronic noise (to be added) has to be well above this RMS,
374 // otherwise the electronic noise of an individual pixel (elNoise2Pix)
375 // may become negative
376
377 TH1D *hnoise = fHDiffPixTheta->ProjectionZ("", binTheta, binTheta,
378 0, 9999, "");
379 Double_t RMS = hnoise->GetRMS(1);
380 delete hnoise;
381
382 elNoise2 = TMath::Min(RMS, sigmabar2 - sigbarold2);
383 //*fLog << "elNoise2 = " << elNoise2 << endl;
384
385 lambdabar = (sigmabar2 - sigbarold2 - elNoise2) / F2excess; // [photons]
386
387 // This value of lambdabar is the same for all pixels;
388 // note that lambdabar is normalized to the area of pixel 0
389
390 //---------- start loop over pixels ---------------------------------
391 // do the padding for each pixel
392 //
393 // pad only pixels - which are used (before image cleaning)
394 //
395 Double_t sig = 0;
396 Double_t sigma2 = 0;
397 Double_t diff = 0;
398 Double_t addSig2 = 0;
399 Double_t elNoise2Pix = 0;
400
401
402 for (UInt_t i=0; i<npix; i++)
403 {
404 MCerPhotPix &pix = (*fEvt)[i];
405 if ( !pix.IsPixelUsed() )
406 continue;
407
408 //if ( pix.GetNumPhotons() == 0.0)
409 //{
410 // *fLog << "MPadSchweizer::Process(); no.of photons is 0 for used pixel"
411 // << endl;
412 // continue;
413 //}
414
415 Int_t j = pix.GetPixId();
416
417 Double_t ratioArea = fCam->GetPixRatio(j);
418
419 MPedestalPix &ppix = (*fPed)[j];
420 Double_t oldsigma = ppix.GetMeanRms();
421 Double_t oldsigma2 = oldsigma*oldsigma;
422
423 //---------------------------------
424 // throw the Sigma for this pixel
425 //
426 Int_t binPixel = fHDiffPixTheta->GetYaxis()->FindBin( (Double_t)j );
427
428 Int_t count;
429 Bool_t ok;
430
431 TH1D *hdiff;
432 TH1D *hsig;
433
434 switch (fPadFlag)
435 {
436 case 1 :
437 // throw the Sigma for this pixel from the distribution fHDiffPixTheta
438
439 hdiff = fHDiffPixTheta->ProjectionZ("", binTheta, binTheta,
440 binPixel, binPixel, "");
441 if ( hdiff->GetEntries() == 0 )
442 {
443 *fLog << "MPadSchweizer::Process(); projection for Theta bin "
444 << binTheta << " and pixel bin " << binPixel
445 << " has no entries; aborting " << endl;
446 delete hdiff;
447
448 rc = 5;
449 fErrors[rc]++;
450 return kCONTINUE;
451 }
452
453 count = 0;
454 ok = kFALSE;
455 for (Int_t m=0; m<20; m++)
456 {
457 count += 1;
458 diff = hdiff->GetRandom();
459 // the following condition ensures that elNoise2Pix > 0.0
460
461 if ( (diff + sigmabar2 - oldsigma2/ratioArea
462 - lambdabar*F2excess) > 0.0 )
463 {
464 ok = kTRUE;
465 break;
466 }
467 }
468 if (!ok)
469 {
470
471 *fLog << "theta, j, count, sigmabar, diff = " << theta << ", "
472 << j << ", " << count << ", " << sigmabar << ", "
473 << diff << endl;
474 diff = lambdabar*F2excess + oldsigma2/ratioArea - sigmabar2;
475 }
476 delete hdiff;
477 sigma2 = diff + sigmabar2;
478 break;
479
480 case 2 :
481 // throw the Sigma for this pixel from the distribution fHSigmaPixTheta
482
483 hsig = fHSigmaPixTheta->ProjectionZ("", binTheta, binTheta,
484 binPixel, binPixel, "");
485 if ( hsig->GetEntries() == 0 )
486 {
487 *fLog << "MPadSchweizer::Process(); projection for Theta bin "
488 << binTheta << " and pixel bin " << binPixel
489 << " has no entries; aborting " << endl;
490 delete hsig;
491
492 rc = 6;
493 fErrors[rc]++;
494 return kCONTINUE;
495 }
496
497 count = 0;
498 ok = kFALSE;
499 for (Int_t m=0; m<20; m++)
500 {
501 count += 1;
502
503 sig = hsig->GetRandom();
504 sigma2 = sig*sig/ratioArea;
505 // the following condition ensures that elNoise2Pix > 0.0
506
507 if ( (sigma2-oldsigma2/ratioArea-lambdabar*F2excess) > 0.0 )
508 {
509 ok = kTRUE;
510 break;
511 }
512 }
513 if (!ok)
514 {
515
516 *fLog << "theta, j, count, sigmabar, sig = " << theta << ", "
517 << j << ", " << count << ", " << sigmabar << ", "
518 << sig << endl;
519 sigma2 = lambdabar*F2excess + oldsigma2/ratioArea;
520 }
521 delete hsig;
522 break;
523 }
524
525 //---------------------------------
526 // get the additional sigma^2 for this pixel (due to the padding)
527
528 addSig2 = sigma2*ratioArea - oldsigma2;
529
530
531 //---------------------------------
532 // get the additional electronic noise for this pixel
533
534 elNoise2Pix = addSig2 - lambdabar*F2excess*ratioArea;
535
536
537 //---------------------------------
538 // throw actual number of additional NSB photons (NSB)
539 // and its RMS (sigmaNSB)
540
541 Double_t NSB0 = gRandom->Poisson(lambdabar*ratioArea);
542 Double_t arg = NSB0*(F2excess-1.0) + elNoise2Pix;
543 Double_t sigmaNSB0;
544
545 if (arg >= 0)
546 {
547 sigmaNSB0 = sqrt( arg );
548 }
549 else
550 {
551 *fLog << "MPadSchweizer::Process(); argument of sqrt < 0 : "
552 << arg << endl;
553 sigmaNSB0 = 0.0000001;
554 }
555
556
557 //---------------------------------
558 // smear NSB0 according to sigmaNSB0
559 // and subtract lambdabar because of AC coupling
560
561 Double_t NSB = gRandom->Gaus(NSB0, sigmaNSB0) - lambdabar*ratioArea;
562
563 //---------------------------------
564
565 // add additional NSB to the number of photons
566 Double_t oldphotons = pix.GetNumPhotons();
567 Double_t newphotons = oldphotons + NSB;
568 pix.SetNumPhotons( newphotons );
569
570
571 fHNSB->Fill( NSB/sqrt(ratioArea) );
572 fHPhotons->Fill( oldphotons/sqrt(ratioArea), newphotons/sqrt(ratioArea) );
573
574
575 // error: add sigma of padded noise quadratically
576 Double_t olderror = pix.GetErrorPhot();
577 Double_t newerror = sqrt( olderror*olderror + addSig2 );
578 pix.SetErrorPhot( newerror );
579
580
581 Double_t newsigma = sqrt( oldsigma2 + addSig2 );
582 ppix.SetMeanRms( newsigma );
583
584 fHSigmaPedestal->Fill( oldsigma, newsigma );
585 }
586 //---------- end of loop over pixels ---------------------------------
587
588 // Calculate sigmabar again and crosscheck
589
590
591 //fSigmabar->Calc(*fCam, *fPed, *fEvt);
592 //*fLog << "after padding : " << endl;
593 //fSigmabar->Print("");
594
595
596 //*fLog << "Exit MPadSchweizer::Process();" << endl;
597
598 rc = 0;
599 fErrors[rc]++;
600
601 return kTRUE;
602
603}
604
605// --------------------------------------------------------------------------
606//
607//
608Bool_t MPadSchweizer::PostProcess()
609{
610 if (GetNumExecutions() != 0)
611 {
612
613 *fLog << inf << endl;
614 *fLog << GetDescriptor() << " execution statistics:" << endl;
615 *fLog << dec << setfill(' ');
616 *fLog << " " << setw(7) << fErrors[1] << " (" << setw(3)
617 << (int)(fErrors[1]*100/GetNumExecutions())
618 << "%) Evts skipped due to: Sigmabar_old > 0" << endl;
619
620 *fLog << " " << setw(7) << fErrors[2] << " (" << setw(3)
621 << (int)(fErrors[2]*100/GetNumExecutions())
622 << "%) Evts skipped due to: Zenith angle out of range" << endl;
623
624 *fLog << " " << setw(7) << fErrors[3] << " (" << setw(3)
625 << (int)(fErrors[3]*100/GetNumExecutions())
626 << "%) Evts skipped due to: No data for generating Sigmabar" << endl;
627
628 *fLog << " " << setw(7) << fErrors[4] << " (" << setw(3)
629 << (int)(fErrors[4]*100/GetNumExecutions())
630 << "%) Evts skipped due to: Target sigma <= Sigmabar_old" << endl;
631
632 *fLog << " " << setw(7) << fErrors[5] << " (" << setw(3)
633 << (int)(fErrors[5]*100/GetNumExecutions())
634 << "%) Evts skipped due to: No data for generating Sigma^2-Sigmabar^2" << endl;
635
636 *fLog << " " << setw(7) << fErrors[6] << " (" << setw(3)
637 << (int)(fErrors[6]*100/GetNumExecutions())
638 << "%) Evts skipped due to: No data for generating Sigma" << endl;
639
640 *fLog << " " << fErrors[0] << " ("
641 << (int)(fErrors[0]*100/GetNumExecutions())
642 << "%) Evts survived the padding!" << endl;
643 *fLog << endl;
644
645 }
646
647 //---------------------------------------------------------------
648 TCanvas &c = *(MH::MakeDefCanvas("PadSchweizer", "", 900, 900));
649 c.Divide(3, 3);
650 gROOT->SetSelectedPad(NULL);
651
652 c.cd(1);
653 fHNSB->SetDirectory(NULL);
654 fHNSB->DrawCopy();
655 fHNSB->SetBit(kCanDelete);
656
657 c.cd(2);
658 fHSigmaPedestal->SetDirectory(NULL);
659 fHSigmaPedestal->DrawCopy();
660 fHSigmaPedestal->SetBit(kCanDelete);
661
662 c.cd(3);
663 fHPhotons->SetDirectory(NULL);
664 fHPhotons->DrawCopy();
665 fHPhotons->SetBit(kCanDelete);
666
667 //--------------------------------------------------------------------
668
669
670 c.cd(4);
671 fHSigmaTheta->SetDirectory(NULL);
672 fHSigmaTheta->SetTitle("(Input) 2D : Sigmabar, \\Theta");
673 fHSigmaTheta->DrawCopy();
674 fHSigmaTheta->SetBit(kCanDelete);
675
676
677 //--------------------------------------------------------------------
678 // draw the 3D histogram (input): Theta, pixel, Sigma^2-Sigmabar^2
679
680 c.cd(5);
681 TH2D *l1;
682 l1 = (TH2D*) ((TH3*)fHDiffPixTheta)->Project3D("zx");
683 l1->SetDirectory(NULL);
684 l1->SetTitle("(Input) Sigma^2-Sigmabar^2 vs. \\Theta (all pixels)");
685 l1->SetXTitle("\\Theta [\\circ]");
686 l1->SetYTitle("Sigma^2-Sigmabar^2");
687
688 l1->DrawCopy("box");
689 l1->SetBit(kCanDelete);;
690
691 c.cd(6);
692 TH2D *l2;
693 l2 = (TH2D*) ((TH3*)fHDiffPixTheta)->Project3D("zy");
694 l2->SetDirectory(NULL);
695 l2->SetTitle("(Input) Sigma^2-Sigmabar^2 vs. pixel number (all \\Theta)");
696 l2->SetXTitle("pixel");
697 l2->SetYTitle("Sigma^2-Sigmabar^2");
698
699 l2->DrawCopy("box");
700 l2->SetBit(kCanDelete);;
701
702 //--------------------------------------------------------------------
703 // draw the 3D histogram (input): Theta, pixel, Sigma
704
705 c.cd(8);
706 TH2D *k1;
707 k1 = (TH2D*) ((TH3*)fHSigmaPixTheta)->Project3D("zx");
708 k1->SetDirectory(NULL);
709 k1->SetTitle("(Input) Sigma vs. \\Theta (all pixels)");
710 k1->SetXTitle("\\Theta [\\circ]");
711 k1->SetYTitle("Sigma");
712
713 k1->DrawCopy("box");
714 k1->SetBit(kCanDelete);
715
716 c.cd(9);
717 TH2D *k2;
718 k2 = (TH2D*) ((TH3*)fHSigmaPixTheta)->Project3D("zy");
719 k2->SetDirectory(NULL);
720 k2->SetTitle("(Input) Sigma vs. pixel number (all \\Theta)");
721 k2->SetXTitle("pixel");
722 k2->SetYTitle("Sigma");
723
724 k2->DrawCopy("box");
725 k2->SetBit(kCanDelete);;
726
727
728 //--------------------------------------------------------------------
729
730
731 return kTRUE;
732}
733
734// --------------------------------------------------------------------------
735
736
737
738
739
740
741
742
743
744
745
746
747
748
Note: See TracBrowser for help on using the repository browser.