Changeset 4847 for trunk/MagicSoft/Mars
- Timestamp:
- 09/03/04 18:13:53 (20 years ago)
- Location:
- trunk/MagicSoft/Mars
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/MagicSoft/Mars/Changelog
r4846 r4847 25 25 - small modification in QE of new blind pixel 26 26 27 27 * mcalib/MCalibrationChargeCam.[h,cc] 28 - two new functions GetAveragedConvFADC2PhotPerArea and 29 GetAveragedConvFADC2PhotPerSector, to be used by the data check. 28 30 29 31 2004/09/03: Wolfgang Wittek -
trunk/MagicSoft/Mars/mcalib/MCalibrationChargeCam.cc
r4658 r4847 94 94 #include "MBadPixelsPix.h" 95 95 96 #include "MCalibrationQECam.h" 97 #include "MCalibrationQEPix.h" 98 96 99 #include "MCalibrationChargePix.h" 97 100 #include "MCalibrationChargeBlindPix.h" … … 562 565 563 566 567 // -------------------------------------------------------------------------- 568 // 569 // Calculates the average conversion factor FADC counts to photons for pixel sizes. 570 // The geometry container is used to get the necessary 571 // geometry information (area index). The MCalibrationQECam container is necessary for 572 // the quantum efficiency information. 573 // If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored 574 // in the calculation of the size average. 575 // 576 // Returns a TArrayF of dimension 2: 577 // arr[0]: averaged conversion factors (default: -1.) 578 // arr[1]: Error (rms) of averaged conversion factors (default: 0.) 579 // 580 // ATTENTION: THE USER HAS TO DELETE THE RETURNED TARRAYF ACCORDINGLY 581 // 582 TArrayF *MCalibrationChargeCam::GetAveragedConvFADC2PhotPerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam, 583 const UInt_t ai, MBadPixelsCam *bad) 584 { 585 586 const Int_t np = GetSize(); 587 588 Double_t mean = 0.; 589 Double_t mean2 = 0.; 590 Int_t nr = 0; 591 592 for (int i=0; i<np; i++) 593 { 594 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun)) 595 continue; 596 597 const UInt_t aidx = geom[i].GetAidx(); 598 599 if (ai != aidx) 600 continue; 601 602 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i]; 603 const Float_t conv = pix.GetMeanConvFADC2Phe(); 604 605 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i]; 606 const Float_t qe = qepix.GetQECascadesCombined(); 607 608 mean += conv/qe; 609 mean2 += conv*conv/qe/qe; 610 nr ++; 611 612 } 613 614 TArrayF *arr = new TArrayF(2); 615 arr->AddAt(nr ? mean/nr : -1.,0); 616 arr->AddAt(nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0. ,1); 617 618 return arr; 619 } 620 621 // -------------------------------------------------------------------------- 622 // 623 // Calculates the average conversion factor FADC counts to photons for camera sectors. 624 // The geometry container is used to get the necessary 625 // geometry information (area index). The MCalibrationQECam container is necessary for 626 // the quantum efficiency information. 627 // If the bad pixel container is given all pixels which have the flag 'kUnsuitableRun' are ignored 628 // in the calculation of the size average. 629 // 630 // Returns a TArrayF of dimension 2: 631 // arr[0]: averaged conversion factors (default: -1.) 632 // arr[1]: Error (rms) of averaged conversion factors (default: 0.) 633 // 634 // ATTENTION: THE USER HAS TO DELETE THE RETURNED TARRAYF ACCORDINGLY 635 // 636 TArrayF *MCalibrationChargeCam::GetAveragedConvFADC2PhotPerSector( const MGeomCam &geom, const MCalibrationQECam &qecam, 637 const UInt_t sec, MBadPixelsCam *bad) 638 { 639 const Int_t np = GetSize(); 640 641 Double_t mean = 0.; 642 Double_t mean2 = 0.; 643 Int_t nr = 0; 644 645 for (int i=0; i<np; i++) 646 { 647 if (bad && (*bad)[i].IsUnsuitable(MBadPixelsPix::kUnsuitableRun)) 648 continue; 649 650 const UInt_t sector = geom[i].GetSector(); 651 652 if (sector != sec) 653 continue; 654 655 const MCalibrationChargePix &pix = (MCalibrationChargePix&)(*this)[i]; 656 const Float_t conv = pix.GetMeanConvFADC2Phe(); 657 658 const MCalibrationQEPix &qepix = (MCalibrationQEPix&)qecam[i]; 659 const Float_t qe = qepix.GetQECascadesCombined(); 660 661 mean += conv/qe; 662 mean2 += conv*conv/qe/qe; 663 nr ++; 664 665 } 666 667 TArrayF *arr = new TArrayF(2); 668 arr->AddAt(nr ? mean/nr : -1.,0); 669 arr->AddAt(nr>1 ? TMath::Sqrt((mean2 - mean*mean/nr)/(nr-1)) : 0. ,1); 670 671 return arr; 672 } -
trunk/MagicSoft/Mars/mcalib/MCalibrationChargeCam.h
r4658 r4847 10 10 #endif 11 11 12 class MCalibrationQECam; 12 13 class MCalibrationChargeCam : public MCalibrationCam 13 14 { … … 46 47 Bool_t IsFFactorMethodValid () const; 47 48 49 TArrayF *GetAveragedConvFADC2PhotPerArea ( const MGeomCam &geom, const MCalibrationQECam &qecam, 50 const UInt_t ai=0, MBadPixelsCam *bad=NULL); 51 TArrayF *GetAveragedConvFADC2PhotPerSector( const MGeomCam &geom, const MCalibrationQECam &qecam, 52 const UInt_t sec=0, MBadPixelsCam *bad=NULL); 53 48 54 // Prints 49 55 void Print(Option_t *o="") const;
Note:
See TracChangeset
for help on using the changeset viewer.