source: trunk/MagicSoft/Mars/manalysis/MCT1PadSchweizer.cc@ 2174

Last change on this file since 2174 was 2167, checked in by wittek, 21 years ago
*** empty log message ***
File size: 25.1 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// MCT1PadSchweizer
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 "MCT1PadSchweizer.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(MCT1PadSchweizer);
91
92// --------------------------------------------------------------------------
93//
94// Default constructor.
95//
96MCT1PadSchweizer::MCT1PadSchweizer(const char *name, const char *title)
97{
98 fName = name ? name : "MCT1PadSchweizer";
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//
116MCT1PadSchweizer::~MCT1PadSchweizer()
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 MCT1PadSchweizer::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 MCT1PadSchweizer::SetPadFlag(Int_t padflag)
172{
173 fPadFlag = padflag;
174 *fLog << "MCT1PadSchweizer::SetPadFlag(); choose option " << fPadFlag << endl;
175}
176
177// --------------------------------------------------------------------------
178//
179// Set the pointers and prepare the histograms
180//
181Bool_t MCT1PadSchweizer::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 MCT1PadSchweizer::Process()
267{
268 //*fLog << "Entry MCT1PadSchweizer::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 << "MCT1PadSchweizer::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 << "MCT1PadSchweizer::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 << "MCT1PadSchweizer::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 if ( hblind->GetEntries() > 0.0 )
395 for (UInt_t i=0; i<numBlind; i++)
396 {
397 while(1)
398 {
399 idBlind = (Int_t) (hblind->GetRandom()+0.5);
400 equal = kFALSE;
401 for (UInt_t j=0; j<nlist; j++)
402 if (idBlind == listId[j])
403 {
404 equal = kTRUE;
405 break;
406 }
407 if (!equal) break;
408 }
409 listId[nlist] = idBlind;
410 nlist++;
411
412 fBlinds->SetPixelBlind(idBlind);
413 //*fLog << "idBlind = " << idBlind << endl;
414 }
415 fBlinds->SetReadyToSave();
416
417 delete hblind;
418
419
420 //-------------------------------------------
421 // for the current theta,
422 // generate a sigmabar according to the histogram fHSigmaTheta
423 //
424 Double_t sigmabar=0;
425 Int_t binNumber = fHSigmaTheta->GetXaxis()->FindBin(theta);
426
427 if (binPix != binNumber)
428 {
429 cout << "The binnings of the 2 histograms are not identical; aborting"
430 << endl;
431 return kERROR;
432 }
433
434 TH1D *hsigma;
435
436 hsigma = fHSigmaTheta->ProjectionY("", binNumber, binNumber, "");
437 if ( hsigma->GetEntries() == 0.0 )
438 {
439 *fLog << "MCT1PadSchweizer::Process(); projection for Theta bin "
440 << binNumber << " has no entries; Skip event " << endl;
441 // event cannot be padded; skip event
442 delete hsigma;
443
444 rc = 3;
445 fErrors[rc]++;
446 return kCONTINUE;
447 }
448 else
449 {
450 sigmabar = hsigma->GetRandom();
451 //*fLog << "Theta, bin number = " << theta << ", " << binNumber // << ", sigmabar = " << sigmabar << endl
452 }
453 delete hsigma;
454
455 const Double_t sigmabar2 = sigmabar*sigmabar;
456
457 //-------------------------------------------
458
459 //*fLog << "MCT1PadSchweizer::Process(); sigbarold, sigmabar = "
460 // << sigbarold << ", "<< sigmabar << endl;
461
462 // Skip event if target sigmabar is <= sigbarold
463 if (sigmabar <= sigbarold)
464 {
465 *fLog << "MCT1PadSchweizer::Process(); target sigmabar is less than sigbarold : "
466 << sigmabar << ", " << sigbarold << ", Skip event" << endl;
467
468 rc = 4;
469 fErrors[rc]++;
470 return kCONTINUE;
471 }
472
473
474 //-------------------------------------------
475 //
476 // Calculate average number of NSB photons to be added (lambdabar)
477 // from the value of sigmabar,
478 // - making assumptions about the average electronic noise (elNoise2) and
479 // - using a fixed value (F2excess) for the excess noise factor
480
481 Double_t elNoise2; // [photons]
482 Double_t F2excess = 1.3;
483 Double_t lambdabar; // [photons]
484
485
486
487 Int_t binTheta = fHDiffPixTheta->GetXaxis()->FindBin(theta);
488 if (binTheta != binNumber)
489 {
490 cout << "The binnings of the 2 histograms are not identical; aborting"
491 << endl;
492 return kERROR;
493 }
494
495 // Get RMS of (Sigma^2-sigmabar^2) in this Theta bin.
496 // The average electronic noise (to be added) has to be well above this RMS,
497 // otherwise the electronic noise of an individual pixel (elNoise2Pix)
498 // may become negative
499
500 TH1D *hnoise = fHDiffPixTheta->ProjectionZ("", binTheta, binTheta,
501 0, 9999, "");
502 Double_t RMS = hnoise->GetRMS(1);
503 delete hnoise;
504
505 elNoise2 = TMath::Min(RMS, sigmabar2 - sigbarold2);
506 //*fLog << "elNoise2 = " << elNoise2 << endl;
507
508 lambdabar = (sigmabar2 - sigbarold2 - elNoise2) / F2excess; // [photons]
509
510 // This value of lambdabar is the same for all pixels;
511 // note that lambdabar is normalized to the area of pixel 0
512
513 //---------- start loop over pixels ---------------------------------
514 // do the padding for each pixel
515 //
516 // pad only pixels - which are used (before image cleaning)
517 //
518 Double_t sig = 0;
519 Double_t sigma2 = 0;
520 Double_t diff = 0;
521 Double_t addSig2 = 0;
522 Double_t elNoise2Pix = 0;
523
524
525 for (UInt_t i=0; i<npix; i++)
526 {
527 MCerPhotPix &pix = (*fEvt)[i];
528 if ( !pix.IsPixelUsed() )
529 continue;
530
531 //if ( pix.GetNumPhotons() == 0.0)
532 //{
533 // *fLog << "MCT1PadSchweizer::Process(); no.of photons is 0 for used pixel"
534 // << endl;
535 // continue;
536 //}
537
538 Int_t j = pix.GetPixId();
539
540 Double_t ratioArea = fCam->GetPixRatio(j);
541
542 MPedestalPix &ppix = (*fPed)[j];
543 Double_t oldsigma = ppix.GetMeanRms();
544 Double_t oldsigma2 = oldsigma*oldsigma;
545
546 //---------------------------------
547 // throw the Sigma for this pixel
548 //
549 Int_t binPixel = fHDiffPixTheta->GetYaxis()->FindBin( (Double_t)j );
550
551 Int_t count;
552 Bool_t ok;
553
554 TH1D *hdiff;
555 TH1D *hsig;
556
557 switch (fPadFlag)
558 {
559 case 1 :
560 // throw the Sigma for this pixel from the distribution fHDiffPixTheta
561
562 hdiff = fHDiffPixTheta->ProjectionZ("", binTheta, binTheta,
563 binPixel, binPixel, "");
564 if ( hdiff->GetEntries() == 0 )
565 {
566 *fLog << "MCT1PadSchweizer::Process(); projection for Theta bin "
567 << binTheta << " and pixel bin " << binPixel
568 << " has no entries; aborting " << endl;
569 delete hdiff;
570
571 rc = 5;
572 fErrors[rc]++;
573 return kCONTINUE;
574 }
575
576 count = 0;
577 ok = kFALSE;
578 for (Int_t m=0; m<20; m++)
579 {
580 count += 1;
581 diff = hdiff->GetRandom();
582 // the following condition ensures that elNoise2Pix > 0.0
583
584 if ( (diff + sigmabar2 - oldsigma2/ratioArea
585 - lambdabar*F2excess) > 0.0 )
586 {
587 ok = kTRUE;
588 break;
589 }
590 }
591 if (!ok)
592 {
593
594 *fLog << "theta, j, count, sigmabar, diff = " << theta << ", "
595 << j << ", " << count << ", " << sigmabar << ", "
596 << diff << endl;
597 diff = lambdabar*F2excess + oldsigma2/ratioArea - sigmabar2;
598 }
599 delete hdiff;
600 sigma2 = diff + sigmabar2;
601 break;
602
603 case 2 :
604 // throw the Sigma for this pixel from the distribution fHSigmaPixTheta
605
606 hsig = fHSigmaPixTheta->ProjectionZ("", binTheta, binTheta,
607 binPixel, binPixel, "");
608 if ( hsig->GetEntries() == 0 )
609 {
610 *fLog << "MCT1PadSchweizer::Process(); projection for Theta bin "
611 << binTheta << " and pixel bin " << binPixel
612 << " has no entries; aborting " << endl;
613 delete hsig;
614
615 rc = 6;
616 fErrors[rc]++;
617 return kCONTINUE;
618 }
619
620 count = 0;
621 ok = kFALSE;
622 for (Int_t m=0; m<20; m++)
623 {
624 count += 1;
625
626 sig = hsig->GetRandom();
627 sigma2 = sig*sig/ratioArea;
628 // the following condition ensures that elNoise2Pix > 0.0
629
630 if ( (sigma2-oldsigma2/ratioArea-lambdabar*F2excess) > 0.0 )
631 {
632 ok = kTRUE;
633 break;
634 }
635 }
636 if (!ok)
637 {
638
639 *fLog << "theta, j, count, sigmabar, sig = " << theta << ", "
640 << j << ", " << count << ", " << sigmabar << ", "
641 << sig << endl;
642 sigma2 = lambdabar*F2excess + oldsigma2/ratioArea;
643 }
644 delete hsig;
645 break;
646 }
647
648 //---------------------------------
649 // get the additional sigma^2 for this pixel (due to the padding)
650
651 addSig2 = sigma2*ratioArea - oldsigma2;
652
653
654 //---------------------------------
655 // get the additional electronic noise for this pixel
656
657 elNoise2Pix = addSig2 - lambdabar*F2excess*ratioArea;
658
659
660 //---------------------------------
661 // throw actual number of additional NSB photons (NSB)
662 // and its RMS (sigmaNSB)
663
664 Double_t NSB0 = gRandom->Poisson(lambdabar*ratioArea);
665 Double_t arg = NSB0*(F2excess-1.0) + elNoise2Pix;
666 Double_t sigmaNSB0;
667
668 if (arg >= 0)
669 {
670 sigmaNSB0 = sqrt( arg );
671 }
672 else
673 {
674 *fLog << "MCT1PadSchweizer::Process(); argument of sqrt < 0 : "
675 << arg << endl;
676 sigmaNSB0 = 0.0000001;
677 }
678
679
680 //---------------------------------
681 // smear NSB0 according to sigmaNSB0
682 // and subtract lambdabar because of AC coupling
683
684 Double_t NSB = gRandom->Gaus(NSB0, sigmaNSB0) - lambdabar*ratioArea;
685
686 //---------------------------------
687
688 // add additional NSB to the number of photons
689 Double_t oldphotons = pix.GetNumPhotons();
690 Double_t newphotons = oldphotons + NSB;
691 pix.SetNumPhotons( newphotons );
692
693
694 fHNSB->Fill( NSB/sqrt(ratioArea) );
695 fHPhotons->Fill( oldphotons/sqrt(ratioArea), newphotons/sqrt(ratioArea) );
696
697
698 // error: add sigma of padded noise quadratically
699 Double_t olderror = pix.GetErrorPhot();
700 Double_t newerror = sqrt( olderror*olderror + addSig2 );
701 pix.SetErrorPhot( newerror );
702
703
704 Double_t newsigma = sqrt( oldsigma2 + addSig2 );
705 ppix.SetMeanRms( newsigma );
706
707 fHSigmaPedestal->Fill( oldsigma, newsigma );
708 }
709 //---------- end of loop over pixels ---------------------------------
710
711 // Calculate sigmabar again and crosscheck
712
713
714 //fSigmabar->Calc(*fCam, *fPed, *fEvt);
715 //*fLog << "after padding : " << endl;
716 //fSigmabar->Print("");
717
718
719 //*fLog << "Exit MCT1PadSchweizer::Process();" << endl;
720
721 rc = 0;
722 fErrors[rc]++;
723
724 return kTRUE;
725
726}
727
728// --------------------------------------------------------------------------
729//
730//
731Bool_t MCT1PadSchweizer::PostProcess()
732{
733 if (GetNumExecutions() != 0)
734 {
735
736 *fLog << inf << endl;
737 *fLog << GetDescriptor() << " execution statistics:" << endl;
738 *fLog << dec << setfill(' ');
739 *fLog << " " << setw(7) << fErrors[1] << " (" << setw(3)
740 << (int)(fErrors[1]*100/GetNumExecutions())
741 << "%) Evts skipped due to: Sigmabar_old > 0" << endl;
742
743 *fLog << " " << setw(7) << fErrors[2] << " (" << setw(3)
744 << (int)(fErrors[2]*100/GetNumExecutions())
745 << "%) Evts skipped due to: Zenith angle out of range" << endl;
746
747 *fLog << " " << setw(7) << fErrors[3] << " (" << setw(3)
748 << (int)(fErrors[3]*100/GetNumExecutions())
749 << "%) Evts skipped due to: No data for generating Sigmabar" << endl;
750
751 *fLog << " " << setw(7) << fErrors[4] << " (" << setw(3)
752 << (int)(fErrors[4]*100/GetNumExecutions())
753 << "%) Evts skipped due to: Target sigma <= Sigmabar_old" << endl;
754
755 *fLog << " " << setw(7) << fErrors[5] << " (" << setw(3)
756 << (int)(fErrors[5]*100/GetNumExecutions())
757 << "%) Evts skipped due to: No data for generating Sigma^2-Sigmabar^2" << endl;
758
759 *fLog << " " << setw(7) << fErrors[6] << " (" << setw(3)
760 << (int)(fErrors[6]*100/GetNumExecutions())
761 << "%) Evts skipped due to: No data for generating Sigma" << endl;
762
763 *fLog << " " << setw(7) << fErrors[7] << " (" << setw(3)
764 << (int)(fErrors[7]*100/GetNumExecutions())
765 << "%) Evts skipped due to: No data for generating Blind pixels" << endl;
766
767 *fLog << " " << fErrors[0] << " ("
768 << (int)(fErrors[0]*100/GetNumExecutions())
769 << "%) Evts survived the padding!" << endl;
770 *fLog << endl;
771
772 }
773
774 //---------------------------------------------------------------
775 TCanvas &c = *(MH::MakeDefCanvas("PadSchweizer", "", 900, 1200));
776 c.Divide(3, 4);
777 gROOT->SetSelectedPad(NULL);
778
779 c.cd(1);
780 fHNSB->SetDirectory(NULL);
781 fHNSB->DrawCopy();
782 fHNSB->SetBit(kCanDelete);
783
784 c.cd(2);
785 fHSigmaPedestal->SetDirectory(NULL);
786 fHSigmaPedestal->DrawCopy();
787 fHSigmaPedestal->SetBit(kCanDelete);
788
789 c.cd(3);
790 fHPhotons->SetDirectory(NULL);
791 fHPhotons->DrawCopy();
792 fHPhotons->SetBit(kCanDelete);
793
794 //--------------------------------------------------------------------
795
796
797 c.cd(4);
798 fHSigmaTheta->SetDirectory(NULL);
799 fHSigmaTheta->SetTitle("(Input) 2D : Sigmabar, \\Theta");
800 fHSigmaTheta->DrawCopy();
801 fHSigmaTheta->SetBit(kCanDelete);
802
803 //--------------------------------------------------------------------
804
805
806 c.cd(7);
807 fHBlindPixNTheta->SetDirectory(NULL);
808 fHBlindPixNTheta->SetTitle("(Input) 2D : no.of blind pixels, \\Theta");
809 fHBlindPixNTheta->DrawCopy();
810 fHBlindPixNTheta->SetBit(kCanDelete);
811
812 //--------------------------------------------------------------------
813
814
815 c.cd(10);
816 fHBlindPixIdTheta->SetDirectory(NULL);
817 fHBlindPixIdTheta->SetTitle("(Input) 2D : blind pixel Id, \\Theta");
818 fHBlindPixIdTheta->DrawCopy();
819 fHBlindPixIdTheta->SetBit(kCanDelete);
820
821
822 //--------------------------------------------------------------------
823 // draw the 3D histogram (input): Theta, pixel, Sigma^2-Sigmabar^2
824
825 c.cd(5);
826 TH2D *l1;
827 l1 = (TH2D*) ((TH3*)fHDiffPixTheta)->Project3D("zx");
828 l1->SetDirectory(NULL);
829 l1->SetTitle("(Input) Sigma^2-Sigmabar^2 vs. \\Theta (all pixels)");
830 l1->SetXTitle("\\Theta [\\circ]");
831 l1->SetYTitle("Sigma^2-Sigmabar^2");
832
833 l1->DrawCopy("box");
834 l1->SetBit(kCanDelete);;
835
836 c.cd(8);
837 TH2D *l2;
838 l2 = (TH2D*) ((TH3*)fHDiffPixTheta)->Project3D("zy");
839 l2->SetDirectory(NULL);
840 l2->SetTitle("(Input) Sigma^2-Sigmabar^2 vs. pixel number (all \\Theta)");
841 l2->SetXTitle("pixel");
842 l2->SetYTitle("Sigma^2-Sigmabar^2");
843
844 l2->DrawCopy("box");
845 l2->SetBit(kCanDelete);;
846
847 //--------------------------------------------------------------------
848 // draw the 3D histogram (input): Theta, pixel, Sigma
849
850 c.cd(6);
851 TH2D *k1;
852 k1 = (TH2D*) ((TH3*)fHSigmaPixTheta)->Project3D("zx");
853 k1->SetDirectory(NULL);
854 k1->SetTitle("(Input) Sigma vs. \\Theta (all pixels)");
855 k1->SetXTitle("\\Theta [\\circ]");
856 k1->SetYTitle("Sigma");
857
858 k1->DrawCopy("box");
859 k1->SetBit(kCanDelete);
860
861 c.cd(9);
862 TH2D *k2;
863 k2 = (TH2D*) ((TH3*)fHSigmaPixTheta)->Project3D("zy");
864 k2->SetDirectory(NULL);
865 k2->SetTitle("(Input) Sigma vs. pixel number (all \\Theta)");
866 k2->SetXTitle("pixel");
867 k2->SetYTitle("Sigma");
868
869 k2->DrawCopy("box");
870 k2->SetBit(kCanDelete);;
871
872
873 //--------------------------------------------------------------------
874
875
876 return kTRUE;
877}
878
879// --------------------------------------------------------------------------
880
881
882
883
884
885
886
887
888
889
890
891
892
893
Note: See TracBrowser for help on using the repository browser.