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

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