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

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