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

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