/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Josep Flix 04/2001 ! Author(s): Thomas Bretz 05/2001 ! Author(s): Sebastian Commichau 12/2003 ! Author(s): Javier Rico 01/2004 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MPedCalcPedRun // // This task takes a pedestal run file and fills MPedestalCam during // the Process() with the pedestal and rms computed in an event basis. // In the PostProcess() MPedestalCam is finally filled with the pedestal // mean and rms computed in a run basis. // More than one run (file) can be merged // // Input Containers: // MRawEvtData // // Output Containers: // MPedestalCam // // ///////////////////////////////////////////////////////////////////////////// #include "MPedCalcPedRun.h" #include "MParList.h" #include "MLog.h" #include "MLogManip.h" #include "MRawRunHeader.h" #include "MRawEvtPixelIter.h" #include "MRawEvtData.h" #include "MPedestalPix.h" #include "MPedestalCam.h" #include "MGeomCamMagic.h" ClassImp(MPedCalcPedRun); using namespace std; // -------------------------------------------------------------------------- // // default constructor // MPedCalcPedRun::MPedCalcPedRun(const char *name, const char *title) { fName = name ? name : "MPedCalcPedRun"; fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data"; AddToBranchList("fHiGainPixId"); AddToBranchList("fHiGainFadcSamples"); } // -------------------------------------------------------------------------- // // Look for the following input containers: // // - MRawEvtData // // The following output containers are also searched and created if // they were not found: // // - MPedestalCam // Int_t MPedCalcPedRun::PreProcess( MParList *pList ) { fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData"); if (!fRawEvt) { *fLog << err << "MRawEvtData not found... aborting." << endl; return kFALSE; } fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam"); if (!fPedestals) return kFALSE; fNumSamplesTot=0; return kTRUE; } // -------------------------------------------------------------------------- // // The ReInit searches for the following input containers: // - MRawRunHeader // // It also initializes the data arrays fSumx and fSumx2 // (only for the first read file) // Bool_t MPedCalcPedRun::ReInit(MParList *pList) { const MRawRunHeader *runheader = (MRawRunHeader*)pList->FindObject("MRawRunHeader"); if (!runheader) { *fLog << warn << dbginf; *fLog << "Warning - cannot check file type, MRawRunHeader not found." << endl; } else if (runheader->GetRunType() == kRTMonteCarlo) return kTRUE; fNumPixels = fPedestals->GetSize(); if(fSumx.GetSize()==0) { fSumx.Set(fNumPixels); fSumx2.Set(fNumPixels); fSumx.Reset(); fSumx2.Reset(); } // Calculate an even number for the hi gain samples to avoid // biases due to the fluctuation in pedestal from one slice to // the other one fNumHiGainSamples = runheader->GetNumSamplesHiGain() & ~1; return kTRUE; } // -------------------------------------------------------------------------- // // Fill the MPedestalCam container with the signal mean and rms for the event. // Store the measured signal in arrays fSumx and fSumx2 so that we can // calculate the overall mean and rms in the PostProcess() // Int_t MPedCalcPedRun::Process() { MRawEvtPixelIter pixel(fRawEvt); while (pixel.Next()) { Byte_t *ptr = pixel.GetHiGainSamples(); const Byte_t *end = ptr + fNumHiGainSamples; UInt_t sum = 0; UInt_t sqr = 0; do { sum += *ptr; sqr += *ptr * *ptr; } while (++ptr != end); const Float_t msum = (Float_t)sum; const Float_t msqr = (Float_t)sqr; const Float_t higainped = msum/fNumHiGainSamples; const Float_t higainrms = TMath::Sqrt((msqr-msum*msum/fNumHiGainSamples)/(fNumHiGainSamples-1.)); const UInt_t idx = pixel.GetPixelId(); (*fPedestals)[idx].Set(higainped, higainrms); fSumx[idx] += msum; // // The old version: // // fSumx2[idx] += msqr; // // The new version: // fSumx2[idx] += msum*msum; } fPedestals->SetReadyToSave(); fNumSamplesTot += fNumHiGainSamples; return kTRUE; } // -------------------------------------------------------------------------- // // Compute signal mean and rms in the whole run and store it in MPedestalCam // Int_t MPedCalcPedRun::PostProcess() { // Compute pedestals and rms from the whole run const ULong_t n = fNumSamplesTot; const ULong_t nevts = GetNumExecutions(); MRawEvtPixelIter pixel(fRawEvt); while (pixel.Next()) { const UInt_t pixid = pixel.GetPixelId(); const Float_t sum = fSumx.At(pixid); const Float_t sum2 = fSumx2.At(pixid); const Float_t higainped = sum/n; // // The old version: // // const Float_t higainrms = TMath::Sqrt((sum2-sum*sum/n)/(n-1.)); // // The new version: // const Float_t higainrms = TMath::Sqrt((sum2-sum*sum/nevts)/(nevts-1.)/(Float_t)fNumHiGainSamples); (*fPedestals)[pixid].Set(higainped, higainrms); } return kTRUE; }