| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

In This Package:

PmtCalibFullModel Class Reference

#include <PmtCalibFullModel.h>

Inheritance diagram for PmtCalibFullModel:

[legend]
Collaboration diagram for PmtCalibFullModel:
[legend]
List of all members.

Public Types

 SUCCESS
 NO_INTERFACE
 VERSMISMATCH
 LAST_ERROR
 SUCCESS
 NO_INTERFACE
 VERSMISMATCH
 LAST_ERROR
enum  Status

Public Member Functions

 PmtCalibFullModel (const std::string &type, const std::string &name, const IInterface *parent)
virtual ~PmtCalibFullModel ()
virtual StatusCode initialize ()
virtual StatusCode finalize ()
virtual StatusCode process (const DayaBay::ReadoutHeader &readout)
 process(const DayaBay::ReadoutPmtCrate&) process a single readout, extracting the information needed for calibration
virtual StatusCode calibrate ()
 calibrate() is called after processing many readouts.
INTupleSvcntupleSvc () const
INTupleSvcevtColSvc () const
IDataProviderSvcdetSvc () const
IDataProviderSvcevtSvc () const
IIncidentSvcincSvc () const
IChronoStatSvcchronoSvc () const
IHistogramSvchistoSvc () const
IAlgContextSvccontextSvc () const
DataObjectput (IDataProviderSvc *svc, DataObject *object, const std::string &address, const bool useRootInTES=true) const
DataObjectput (DataObject *object, const std::string &address, const bool useRootInTES=true) const
Gaudi::Utils::GetData< TYPE
>::return_type 
get (IDataProviderSvc *svc, const std::string &location, const bool useRootInTES=true) const
Gaudi::Utils::GetData< TYPE
>::return_type 
get (const std::string &location, const bool useRootInTES=true) const
TYPE * getDet (IDataProviderSvc *svc, const std::string &location) const
TYPE * getDet (const std::string &location) const
bool exist (IDataProviderSvc *svc, const std::string &location, const bool useRootInTES=true) const
bool exist (const std::string &location, const bool useRootInTES=true) const
bool existDet (IDataProviderSvc *svc, const std::string &location) const
bool existDet (const std::string &location) const
TYPE * getOrCreate (IDataProviderSvc *svc, const std::string &location, const bool useRootInTES=true) const
TYPE * getOrCreate (const std::string &location, const bool useRootInTES=true) const
TOOL * tool (const std::string &type, const std::string &name, const IInterface *parent=0, bool create=true) const
TOOL * tool (const std::string &type, const IInterface *parent=0, bool create=true) const
SERVICE * svc (const std::string &name, const bool create=true) const
IUpdateManagerSvcupdMgrSvc () const
IDataProviderSvcfastContainersSvc () const
StatusCode Error (const std::string &msg, const StatusCode st=StatusCode::FAILURE, const size_t mx=10) const
StatusCode Warning (const std::string &msg, const StatusCode st=StatusCode::FAILURE, const size_t mx=10) const
StatusCode Print (const std::string &msg, const StatusCode st=StatusCode::SUCCESS, const MSG::Level lev=MSG::INFO) const
StatusCode Assert (const bool ok, const std::string &message="", const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Assert (const bool ok, const char *message, const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Exception (const std::string &msg, const GaudiException &exc, const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Exception (const std::string &msg, const std::exception &exc, const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Exception (const std::string &msg="no message", const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
MsgStreammsgStream (const MSG::Level level) const
MsgStreamalways () const
MsgStreamfatal () const
MsgStreamerr () const
MsgStreamerror () const
MsgStreamwarning () const
MsgStreaminfo () const
MsgStreamdebug () const
MsgStreamverbose () const
MsgStreammsg () const
const Statisticscounters () const
StatEntitycounter (const std::string &tag) const
MSG::Level msgLevel () const
bool msgLevel (const MSG::Level level) const
void resetMsgStream () const
bool typePrint () const
bool propsPrint () const
bool statPrint () const
bool errorsPrint () const
long printStat (const MSG::Level level=MSG::ALWAYS) const
long printErrors (const MSG::Level level=MSG::ALWAYS) const
long printProps (const MSG::Level level=MSG::ALWAYS) const
void registerCondition (const std::string &condition, StatusCode(CallerClass::*mf)()=NULL)
void registerCondition (const std::string &condition, CondType *&condPtrDest, StatusCode(CallerClass::*mf)()=NULL)
void registerCondition (char *condition, StatusCode(CallerClass::*mf)()=NULL)
void registerCondition (TargetClass *condition, StatusCode(CallerClass::*mf)()=NULL)
StatusCode runUpdate ()
TransientFastContainer< T > * getFastContainer (const std::string &location, typename TransientFastContainer< T >::size_type initial=0)
StatusCode release (const IInterface *interface) const
virtual unsigned long release ()
const std::string & context () const
const std::string & rootInTES () const
double globalTimeOffset () const
virtual StatusCode queryInterface (const InterfaceID &riid, void **ppvUnknown)
virtual unsigned long addRef ()
virtual const std::string & name () const
virtual const std::string & type () const
virtual const IInterfaceparent () const
virtual StatusCode configure ()
virtual StatusCode start ()
virtual StatusCode stop ()
virtual StatusCode terminate ()
virtual StatusCode reinitialize ()
virtual StatusCode restart ()
virtual Gaudi::StateMachine::State FSMState () const
virtual Gaudi::StateMachine::State targetFSMState () const
virtual StatusCode sysInitialize ()
virtual StatusCode sysStart ()
virtual StatusCode sysStop ()
virtual StatusCode sysFinalize ()
virtual StatusCode sysReinitialize ()
virtual StatusCode sysRestart ()
virtual StatusCode setProperty (const Property &p)
virtual StatusCode setProperty (const std::string &s)
virtual StatusCode setProperty (const std::string &n, const std::string &v)
StatusCode setProperty (const std::string &name, const TYPE &value)
virtual StatusCode getProperty (Property *p) const
virtual const PropertygetProperty (const std::string &name) const
virtual StatusCode getProperty (const std::string &n, std::string &v) const
virtual const std::vector<
Property * > & 
getProperties () const
PropertyMgrgetPropertyMgr ()
ISvcLocatorserviceLocator () const
ISvcLocatorsvcLoc () const
IMessageSvcmsgSvc () const
IToolSvctoolSvc () const
StatusCode setProperties ()
StatusCode service (const std::string &name, T *&svc, bool createIf=true) const
StatusCode service (const std::string &type, const std::string &name, T *&svc) const
void declInterface (const InterfaceID &, void *)
PropertydeclareProperty (const std::string &name, T &property, const std::string &doc="none") const
PropertydeclareRemoteProperty (const std::string &name, IProperty *rsvc, const std::string &rname="") const
IAuditorSvcauditorSvc () const
IMonitorSvcmonitorSvc () const
void declareInfo (const std::string &name, const T &var, const std::string &desc) const
void declareInfo (const std::string &name, const std::string &format, const void *var, int size, const std::string &desc) const
virtual const std::string & type () const =0
virtual const IInterfaceparent () const =0
virtual StatusCode configure ()=0
virtual StatusCode start ()=0
virtual StatusCode stop ()=0
virtual StatusCode terminate ()=0
virtual StatusCode reinitialize ()=0
virtual StatusCode restart ()=0
virtual Gaudi::StateMachine::State FSMState () const =0
virtual StatusCode sysInitialize ()=0
virtual StatusCode sysStart ()=0
virtual StatusCode sysStop ()=0
virtual StatusCode sysFinalize ()=0
virtual StatusCode sysReinitialize ()=0
virtual StatusCode sysRestart ()=0
virtual unsigned long refCount () const =0
virtual const std::string & name () const =0
virtual StatusCode queryInterface (const InterfaceID &riid, void **ppvInterface)=0
virtual unsigned long addRef ()=0
virtual unsigned long release ()=0

Static Public Member Functions

static const InterfaceIDinterfaceID ()
static const InterfaceIDinterfaceID ()
static const InterfaceIDinterfaceID ()
static const InterfaceIDinterfaceID ()
 Retrieve interface ID.

Public Attributes

 SUCCESS
 NO_INTERFACE
 VERSMISMATCH
 LAST_ERROR

Protected Types

typedef std::map< std::string,
StatEntity
Statistics
typedef std::map< std::string,
unsigned int > 
Counter
typedef std::vector< IAlgTool * > AlgTools
typedef std::pair< IInterface *,
std::string > 
ServiceEntry
typedef std::vector< ServiceEntryServices

Protected Member Functions

StatusCode releaseTool (const IAlgTool *tool) const
StatusCode releaseSvc (const IInterface *svc) const
int outputLevel () const
virtual unsigned long refCount () const
IntegerPropertyoutputLevelProperty ()
void initOutputLevel (Property &prop)

Static Protected Attributes

static const bool IgnoreRootInTES
static const bool UseRootInTES

Private Member Functions

bool hasStats (const DayaBay::Detector &detector)
 Check if the statistics for this detector have been initialized.
StatusCode prepareStats (const Context &context)
 Get the current statistics for a detector, or create if needed.
std::string getPath (const DayaBay::Detector &detector)
std::string getPath (const DayaBay::Detector &detector, int ring)
std::string getPath (const DayaBay::FeeChannelId &channelId)

Private Attributes

std::string m_cableSvcName
std::string m_calibSvcName
std::string m_floatFeePedesSvcName
std::string m_filepath
ICableSvcm_cableSvc
ICalibDataSvcm_calibSvc
IStatisticsSvcm_statsSvc
IFloatingFeePedestalSvcm_floatFeePedesSvc
bool m_useFloatFeePedes
std::vector< DayaBay::Detectorm_processedDetectors

Detailed Description

Definition at line 72 of file PmtCalibFullModel.h.


Constructor & Destructor Documentation

PmtCalibFullModel::PmtCalibFullModel ( const std::string &  type,
const std::string &  name,
const IInterface parent 
)

Definition at line 17 of file PmtCalibFullModel.cc.

00020   : GaudiTool(type,name,parent)
00021   , m_cableSvc(0)
00022   , m_statsSvc(0)
00023 {
00024   declareInterface< IPmtCalibParamTool >(this) ;   
00025   declareProperty("CableSvcName",m_cableSvcName="StaticCableSvc",
00026                   "Name of service to map between detector, hardware, and electronic IDs");
00027   declareProperty("CalibSvcName",m_calibSvcName="StaticCalibDataSvc",
00028                   "Name of service to access FEE channel baselines");
00029   declareProperty("FloatingFeePedestalSvcName", m_floatFeePedesSvcName="FloatingFeePedestalSvc",
00030                   "Name of service to access floating FEE channel baselines");
00031   declareProperty("UseFloatFeePedes", m_useFloatFeePedes=true,
00032                   "Use FloatFeePedes or not? In case of false, it will use CalibSvc.");
00033   declareProperty("FilePath",m_filepath="/file0/PmtCalibFullModel",
00034                   "File path of registered histograms.");
00035 }

PmtCalibFullModel::~PmtCalibFullModel (  )  [virtual]

Definition at line 37 of file PmtCalibFullModel.cc.

00037 {}


Member Function Documentation

StatusCode PmtCalibFullModel::initialize (  )  [virtual]

Reimplemented from GaudiTool.

Definition at line 39 of file PmtCalibFullModel.cc.

00040 {
00041   info() << "initialize()" << endreq;
00042   StatusCode sc = this->GaudiTool::initialize();
00043   if( sc != StatusCode::SUCCESS ) return sc;
00044 
00045   // Get Cable Service
00046   m_cableSvc = svc<ICableSvc>(m_cableSvcName,true);
00047   if(!m_cableSvc){
00048     error() << "Failed to access cable svc: " << m_cableSvcName << endreq;
00049     return StatusCode::FAILURE;
00050   }
00051 
00052   // Get Calibration Service
00053   m_calibSvc = svc<ICalibDataSvc>(m_calibSvcName,true);
00054   if(!m_calibSvc){
00055     error() << "Failed to access calibration svc: " << m_calibSvcName << endreq;
00056     return StatusCode::FAILURE;
00057   }
00058 
00059   // Get floating fee pedestal service
00060   m_floatFeePedesSvc = svc<IFloatingFeePedestalSvc>(m_floatFeePedesSvcName,true);
00061   if(!m_floatFeePedesSvc){
00062     error() << "Failed to access FloatingFeePedestalSvc: " << m_floatFeePedesSvcName << endreq;
00063     return StatusCode::FAILURE;
00064   }
00065 
00066   // Get the histogram service
00067   m_statsSvc = svc<IStatisticsSvc>("StatisticsSvc",true);
00068   if(!m_statsSvc) {
00069     error() << " No StatisticsSvc available." << endreq;
00070     return StatusCode::FAILURE;
00071   }
00072 
00073   // Initialize data from input file.
00074   // Processed detectors are identified by directory name.
00075   std::vector<std::string> detDirs = m_statsSvc->getSubFolders(m_filepath);
00076   std::vector<std::string>::iterator dirIter, dirEnd = detDirs.end();
00077   for(dirIter=detDirs.begin(); dirIter != dirEnd; dirIter++){
00078     info() << "Processing input directory: " << m_filepath << "/" 
00079            << *dirIter << endreq;
00080     DayaBay::Detector detector( 
00081                  DayaBay::Detector::siteDetPackedFromString(*dirIter) );
00082     if( detector.site() == Site::kUnknown 
00083         || detector.detectorId() == DetectorId::kUnknown ){
00084       warning() << "Input directory: " << m_filepath << "/" 
00085                 << *dirIter << " not processed." << endreq;
00086     }
00087     m_processedDetectors.push_back(detector);
00088   }
00089 
00090   return StatusCode::SUCCESS;
00091 }

StatusCode PmtCalibFullModel::finalize (  )  [virtual]

Reimplemented from GaudiTool.

Definition at line 93 of file PmtCalibFullModel.cc.

00094 {
00095   StatusCode sc = this->GaudiTool::finalize();
00096   if( sc != StatusCode::SUCCESS ) return sc;
00097 
00098   return StatusCode::SUCCESS;
00099 }

StatusCode PmtCalibFullModel::process ( const DayaBay::ReadoutHeader readout  )  [virtual]

process(const DayaBay::ReadoutPmtCrate&) process a single readout, extracting the information needed for calibration

Implements IPmtCalibParamTool.

Definition at line 103 of file PmtCalibFullModel.cc.

00103                                                                               {
00104   
00105   const DayaBay::DaqCrate* daqCrate = readoutHeader.daqCrate();
00106   if(!daqCrate){
00107     error() << "Failed to get DAQ crate from header" << endreq;
00108     return StatusCode::FAILURE;
00109   }
00110   
00111   const DayaBay::Detector& detector = daqCrate->detector();
00112   
00113   // Convert to PMT crate readout
00114   const DayaBay::DaqPmtCrate* pmtCrate
00115     = daqCrate->asPmtCrate();
00116   if(!pmtCrate){
00117     // Not an AD, continue
00118     debug() << "process(): Do not process detector: " << detector.detName() 
00119             << endreq;
00120     return StatusCode::SUCCESS;
00121   }
00122   
00123   Context context = readoutHeader.context();
00124   if( !hasStats(detector) ){
00125     // New detector; initialize all the detector statistics
00126     StatusCode sc = prepareStats(context);
00127     if( sc != StatusCode::SUCCESS ) return sc;
00128   }
00129   ServiceMode svcMode(context, 0);
00130   
00131   // Retrieve the parameters and histograms
00132   TObject* nReadoutsObj = 
00133     m_statsSvc->get(this->getPath(detector) + "/nReadouts");
00134   TParameter<int>* nReadoutsPar= dynamic_cast<TParameter<int>*>(nReadoutsObj);
00135   if(!nReadoutsPar) return StatusCode::FAILURE;
00136   TH1F* adcSumH = m_statsSvc->getTH1F( this->getPath(detector) + "/adcSum");
00137   if(!adcSumH) return StatusCode::FAILURE;
00138   TH1F* tdcMedianH = m_statsSvc->getTH1F( this->getPath(detector)+"/tdcMedian");
00139   if(!tdcMedianH) return StatusCode::FAILURE;
00140   
00141   // Add the current readout to the calibration statistics
00142   std::vector<int> readoutTdc;
00143   double adcSum = 0;
00144   // Loop over channels
00145   const DayaBay::DaqPmtCrate::PmtChannelPtrList& channels
00146     = pmtCrate->channelReadouts();
00147   
00148   DayaBay::DaqPmtCrate::PmtChannelPtrList::const_iterator channelIter, 
00149     channelEnd = channels.end();
00150   for(channelIter = channels.begin(); channelIter!=channelEnd; channelIter++) { 
00151     const DayaBay::DaqPmtChannel& channel = *(*channelIter);
00152     const DayaBay::FeeChannelId& chanId = channel.channelId();
00153     //const DayaBay::FeeCalibData* feeCalib = m_calibSvc->feeCalibData(chanId,svcMode);
00154     
00155     
00156     // Get Parameters and Histograms
00157     TObject* nHitsObj = m_statsSvc->get(this->getPath(chanId) + "/nHits");
00158     TParameter<int>* nHitsPar = dynamic_cast<TParameter<int>*>(nHitsObj);
00159     if(!nHitsPar) return StatusCode::FAILURE;
00160     TH1F* tdcRawH = m_statsSvc->getTH1F( this->getPath(chanId) + "/tdcRaw");
00161     if(!tdcRawH) return StatusCode::FAILURE;
00162     TH1F* adcRawH = m_statsSvc->getTH1F( this->getPath(chanId) + "/adcRaw");
00163     if(!adcRawH) return StatusCode::FAILURE;
00164     TH1F* adcH = m_statsSvc->getTH1F( this->getPath(chanId) + "/adc");
00165     if(!adcH) return StatusCode::FAILURE;
00166     TH1F* adcByClockH = m_statsSvc->getTH1F( this->getPath(chanId) 
00167                                              + "/adcByClock");
00168     if(!adcByClockH) return StatusCode::FAILURE;
00169     
00170     //std::vector<int>::const_iterator tdcIter, tdcEnd = pmtChan.tdc().end();
00171     // Loop over tdc values
00172     for(unsigned int i=0; i<channel.hitCount(); i++) {
00173       int tdc=channel.tdc(i);
00174       readoutTdc.push_back(tdc);
00175       tdcRawH->Fill(tdc);
00176     }
00177     
00178     // Loop over adc values
00179     //double adcBaseline = feeCalib->m_adcBaselineLow;
00180     for(unsigned int i=0; i<channel.hitCount(); i++) {
00181       int adcClock = channel.peakCycle(i);
00182       int adc = channel.adc(i);
00183 
00184       double adcBaseline;
00185       if( m_useFloatFeePedes ) {
00186         if(channel.isHighGainAdc(i)){
00187           adcBaseline = m_floatFeePedesSvc->pedestal( chanId, DayaBay::FeeGain::kHigh);
00188         }else{
00189           adcBaseline = m_floatFeePedesSvc->pedestal( chanId, DayaBay::FeeGain::kLow);
00190         }
00191         
00192         // In case floatFeePedesSvc didn't get enough statistics yet, use preAdc instead.
00193         if( adcBaseline<0 ) adcBaseline = channel.preAdcAvg(i);
00194       } else {
00195         // adcBaseline = feeCalib->m_adcBaselineLow;  // The pedestal run result is not good enoughreating CalibParamConfDbMerge.lock and Starting CalibParamConfDbMerge
00196         adcBaseline = channel.preAdcAvg(i);  // Use preAdc is much better
00197       }
00198   
00199       // Restrict to fine range gain measurement
00200       if( adc>0 && channel.isHighGainAdc(i) ) {
00201         adcRawH->Fill(adc);
00202         adcH->Fill(adc-adcBaseline);
00203         adcByClockH->Fill(adcClock,adc-adcBaseline);
00204         adcSum += adc-adcBaseline;
00205       }
00206     }
00207     // Increment number of hits on this channel
00208     nHitsPar->SetVal( nHitsPar->GetVal() + 1 );
00209 
00210   }
00211   if(readoutTdc.size() > 0){
00212     // Find the median TDC
00213     std::sort(readoutTdc.begin(), readoutTdc.end());
00214     int medianIdx = readoutTdc.size()/2;
00215     int medianTdc = readoutTdc[medianIdx];
00216     for(channelIter = channels.begin(); channelIter!=channelEnd; channelIter++) { 
00217       const DayaBay::DaqPmtChannel& channel = *(*channelIter);
00218       const DayaBay::FeeChannelId& chanId = channel.channelId();
00219       TH1F* tdcByMedianH = 
00220         m_statsSvc->getTH1F( this->getPath(chanId) + "/tdcByMedian");
00221       if(!tdcByMedianH) return StatusCode::FAILURE;
00222       //std::vector<int>::const_iterator tdcIter, tdcEnd = pmtChan.tdc().end();
00223       // Loop over tdc values
00224       for(unsigned int i=0; i<channel.hitCount(); i++) {
00225         int tdc = channel.tdc(i);
00226         tdcByMedianH->Fill(tdc-medianTdc);
00227       }
00228     }
00229     tdcMedianH->Fill(medianTdc);
00230   }
00231 
00232   adcSumH->Fill(adcSum);
00233   // Increment number of readouts from this detector
00234   nReadoutsPar->SetVal( nReadoutsPar->GetVal() + 1 );
00235   return StatusCode::SUCCESS;  
00236 }

StatusCode PmtCalibFullModel::calibrate (  )  [virtual]

calibrate() is called after processing many readouts.

This method is responsible for calculating the calibration parameters.

Implements IPmtCalibParamTool.

Definition at line 238 of file PmtCalibFullModel.cc.

00238                                        {
00239   // Loop over detectors in summary data, and calculate parameters
00240 
00241   std::vector<DayaBay::Detector>::iterator detIter, 
00242     detEnd = m_processedDetectors.end();
00243   for(detIter = m_processedDetectors.begin(); detIter != detEnd; detIter++){
00244     DayaBay::Detector detector = *detIter;
00245 
00246     // Retrieve the data description and histograms
00247     TObject* nReadoutsObj = 
00248       m_statsSvc->get(this->getPath(detector) + "/nReadouts");
00249     TParameter<int>* nReadoutsPar= dynamic_cast<TParameter<int>*>(nReadoutsObj);
00250     if(!nReadoutsPar) return StatusCode::FAILURE;
00251 
00252     TH1F* meanOccupancyH = m_statsSvc->getTH1F( this->getPath(detector)
00253                                                 +"/meanOccupancy");
00254     if(!meanOccupancyH) return StatusCode::FAILURE;
00255 
00256     TH1F* meanTdcOffsetH = m_statsSvc->getTH1F( this->getPath(detector)
00257                                                 +"/meanTdcOffset");
00258     if(!meanTdcOffsetH) return StatusCode::FAILURE;
00259 
00260     TH1F* gainChannelH = m_statsSvc->getTH1F( this->getPath(detector)
00261                                                 +"/gainChannel");
00262     if(!gainChannelH) return StatusCode::FAILURE;
00263 
00264     TH1F* tdcOffsetChannelH = m_statsSvc->getTH1F( this->getPath(detector)
00265                                                 +"/tdcOffsetChannel");
00266     if(!tdcOffsetChannelH) return StatusCode::FAILURE;
00267 
00268     TH1F* occupancyChannelH = m_statsSvc->getTH1F( this->getPath(detector)
00269                                                 +"/occupancyChannel");
00270     if(!occupancyChannelH) return StatusCode::FAILURE;
00271 
00272     TH2F* occupancyH2D = m_statsSvc->getTH2F( this->getPath(detector) 
00273                                               + "/occupancy");
00274     if(!occupancyH2D) return StatusCode::FAILURE;
00275     TH2F* gainH2D = m_statsSvc->getTH2F( this->getPath(detector) 
00276                                               + "/gain2D");
00277     if(!gainH2D) return StatusCode::FAILURE;
00278     TH1F* gainH1D = m_statsSvc->getTH1F( this->getPath(detector) 
00279                                               + "/gain1D");
00280     if(!gainH1D) return StatusCode::FAILURE;
00281 
00282     TObject* simFlagObj = m_statsSvc->get(this->getPath(detector) +"/simFlag");
00283     TParameter<int>* simFlagPar = dynamic_cast<TParameter<int>*>(simFlagObj);
00284     if(!simFlagPar) return StatusCode::FAILURE;
00285 
00286     TObject* timeSecObj = m_statsSvc->get(this->getPath(detector) +"/timeSec");
00287     TParameter<int>* timeSecPar = dynamic_cast<TParameter<int>*>(timeSecObj);
00288     if(!timeSecPar) return StatusCode::FAILURE;
00289     
00290     TObject* timeNanoSecObj = 
00291       m_statsSvc->get(this->getPath(detector) +"/timeNanoSec");
00292     TParameter<int>* timeNanoSecPar = 
00293       dynamic_cast<TParameter<int>*>(timeNanoSecObj);
00294     if(!timeNanoSecPar) return StatusCode::FAILURE;
00295 
00296     Site::Site_t site = detector.site();
00297     DetectorId::DetectorId_t detId = detector.detectorId();
00298     SimFlag::SimFlag_t simFlag = (SimFlag::SimFlag_t)(simFlagPar->GetVal());
00299     time_t timeSec = (time_t)(timeSecPar->GetVal());
00300     int timeNanoSec = timeNanoSecPar->GetVal();
00301     int nReadouts = nReadoutsPar->GetVal();
00302 
00303     // Context, ServiceMode
00304     Context context(site,simFlag,TimeStamp(timeSec,timeNanoSec),detId);
00305     ServiceMode svcMode(context, 0);
00306 
00307     // Prepare data for ring-to-ring comparison
00308     int nRings = 8;
00309     std::map<int,double> meanOccRing;
00310     std::map<int,double> meanOccWeight;
00311     std::map<int,double> meanTdcRing;
00312     std::map<int,double> meanTdcWeight;
00313     for(int ring = 1; ring <= nRings; ring++){
00314       meanOccRing[ring] = 0.0;
00315       meanOccWeight[ring] = 0.0;
00316       meanTdcRing[ring] = 0.0;
00317       meanTdcWeight[ring] = 0.0;
00318     }
00319     double minValidOcc = 0.01;  // Cut on PMTs with bad/low occupancy
00320     std::vector<DayaBay::FeeChannelId> badChannels;
00321 
00322     // Get list of all FEE channels
00323     const std::vector<DayaBay::FeeChannelId>& channelList 
00324       = m_cableSvc->feeChannelIds( svcMode );
00325     std::vector<DayaBay::FeeChannelId>::const_iterator chanIter, 
00326       chanEnd = channelList.end();
00327     // Initialize statistics for each channel
00328     for(chanIter = channelList.begin(); chanIter != chanEnd; chanIter++){
00329       DayaBay::FeeChannelId chanId = *chanIter;
00330 
00331       // Analyze the channel data
00332       // Determine ring and column
00333       DayaBay::AdPmtSensor pmtId = 
00334         m_cableSvc->adPmtSensor(chanId, svcMode);
00335       int ring = pmtId.ring();
00336       int column = pmtId.column();
00337 
00338       // Channel status
00339       bool occupancyOk = true;
00340       bool gainOk = true;
00341 
00342       // Get Parameters and Histograms
00343       TObject* nHitsObj = m_statsSvc->get(this->getPath(chanId) + "/nHits");
00344       TParameter<int>* nHitsPar = dynamic_cast<TParameter<int>*>(nHitsObj);
00345       if(!nHitsPar) return StatusCode::FAILURE;
00346       TH1F* tdcByMedianH= m_statsSvc->getTH1F( this->getPath(chanId)
00347                                                +"/tdcByMedian");
00348       if(!tdcByMedianH) return StatusCode::FAILURE;
00349       TH1F* adcH = m_statsSvc->getTH1F( this->getPath(chanId) + "/adc");
00350       if(!adcH) return StatusCode::FAILURE;
00351       TH1F* occupancyH = m_statsSvc->getTH1F( this->getPath(detector,ring)
00352                                               +"/occupancy");
00353       if(!occupancyH) return StatusCode::FAILURE;
00354       TH1F* tdcOffsetH = m_statsSvc->getTH1F( this->getPath(detector,ring)
00355                                               +"/tdcOffset");
00356       if(!tdcOffsetH) return StatusCode::FAILURE;
00357       TH1F* adcMedianH = m_statsSvc->getTH1F( this->getPath(detector,ring)
00358                                               +"/adcMedian");
00359       if(!adcMedianH) return StatusCode::FAILURE;
00360       TH1F* adcSigmaH = m_statsSvc->getTH1F( this->getPath(detector,ring)
00361                                              +"/adcSigma");
00362       if(!adcSigmaH) return StatusCode::FAILURE;
00363 
00364 
00365       // ADC median and sigma
00366       double adcMean = 0;
00367       double adcMeanUncert = 0;
00368       double adcSigma = 0;
00369       double adcSigmaUncert = 0;
00370 
00371       double occupancytosave=0;
00372       double occupancyerrtosave=0;
00373       {
00374         {
00375           int nHits = nHitsPar->GetVal();
00376           if(nReadouts > 0){
00377             occupancytosave = nHits / ((double) nReadouts);
00378             occupancyerrtosave = sqrt(occupancytosave*(1-occupancytosave)*1./nReadouts);
00379           }
00380         }
00381 
00382         //Fit code modified by jpochoa on 03/18/11
00383         TF1 *nimmodelfunc = new TF1("nimmodelfunc",NIMmodel,5,300,8);
00384         nimmodelfunc->SetParNames(  "N","Q_{0}","#sigma_{0}","Q_{1}","#sigma_{1}","w", "#alpha","#mu");
00385         double AreaGuess, Q1MeanGuess, Q1SigmaGuess;
00386         float q0=0;//for a pedestal substracted distribution
00387         float sigma0=0;//for a pedestal substracted distribution
00388         AreaGuess = adcH->GetEntries();
00389         Q1MeanGuess = adcH->GetMean() - q0;
00390         Q1SigmaGuess = adcH->GetRMS();
00391         nimmodelfunc->SetParameters(AreaGuess, 90,1.9,Q1MeanGuess, Q1SigmaGuess, 0.01, 0.03,    0.1);
00392         //fix four parameters
00393         nimmodelfunc->FixParameter(1,q0);
00394         nimmodelfunc->FixParameter(2,sigma0);
00395         nimmodelfunc->FixParameter(5,0.);
00396         nimmodelfunc->FixParameter(6,0.);
00397        
00398         adcH->Fit(nimmodelfunc,"R");
00399         if(nimmodelfunc!=NULL){
00400           
00401           adcMean=nimmodelfunc->GetParameter(3);
00402           adcMeanUncert =nimmodelfunc->GetParError(3);
00403           adcSigma = nimmodelfunc->GetParameter(4);
00404           adcSigmaUncert = nimmodelfunc->GetParError(3);
00405                   
00406           if(fabs(adcH->GetMean()-nimmodelfunc->GetParameter(3)) > adcH->GetRMS()){
00407             warning() << "Issues with fit on ring, column: " << ring << " / " 
00408                       << column << " mean/peak: " << adcMean << " / " << adcH->GetMean() 
00409                       << endreq;
00410             gainOk = false;
00411           }//print warning if issues
00412           
00413           
00414         } else {
00415 
00416           adcMean=0;
00417           adcMeanUncert=0;
00418           adcSigma = 0;
00419           adcSigmaUncert = 0;
00420 
00421           warning() << "PMT full model fit error on ring/column: " << ring << "," << column << endreq;
00422           gainOk = false;
00423         }
00424                
00425         adcMedianH->SetBinContent(adcMedianH->FindBin(column),adcMean);
00426         adcMedianH->SetBinError(adcMedianH->FindBin(column),adcMeanUncert);
00427         adcSigmaH->SetBinContent(adcSigmaH->FindBin(column),adcSigma);
00428         adcSigmaH->SetBinError(adcSigmaH->FindBin(column),adcSigmaUncert);
00429 
00430         gainChannelH->SetBinContent(ring*24+column,adcMean);
00431         gainChannelH->SetBinError(ring*24+column,adcMeanUncert);
00432 
00433         gainH2D->SetBinContent(gainH2D->GetXaxis()->FindBin(column),
00434                                gainH2D->GetYaxis()->FindBin(ring),
00435                                adcMean);
00436         gainH2D->SetBinError(gainH2D->GetXaxis()->FindBin(column),
00437                                gainH2D->GetYaxis()->FindBin(ring),
00438                                adcMeanUncert);
00439         gainH1D->Fill(adcMean);
00440 
00441       }
00442 
00443 
00444       // Channel Occupancy
00445       double occupancy = occupancytosave;
00446       double occupancyUncert = occupancyerrtosave;
00447       {
00448         occupancyH->SetBinContent(occupancyH->FindBin(column),occupancy);
00449         occupancyH->SetBinError(occupancyH->FindBin(column),occupancyUncert);
00450         if(occupancy < minValidOcc) occupancyOk = false;
00451 
00452         occupancyChannelH->SetBinContent(ring*24+column,occupancy);
00453         occupancyChannelH->SetBinError(ring*24+column,occupancyUncert);
00454         occupancyH2D->SetBinContent(occupancyH2D->GetXaxis()->FindBin(column),
00455                                     occupancyH2D->GetYaxis()->FindBin(ring),
00456                                     occupancy);
00457         occupancyH2D->SetBinError(occupancyH2D->GetXaxis()->FindBin(column),
00458                                     occupancyH2D->GetYaxis()->FindBin(ring),
00459                                     occupancyUncert);
00460 
00461       }
00462       // TDC offset
00463       double tdcOffset = 0;
00464       double tdcOffsetUncert = 0;
00465       {
00466         // Use a fit to the leading edge to determine the time offset
00467         int maxBin = tdcByMedianH->GetMaximumBin();
00468         double fitMin = tdcByMedianH->GetBinCenter(maxBin - 2);
00469         double fitMax = tdcByMedianH->GetBinCenter(maxBin + 10);
00470         tdcByMedianH->Fit("gaus","","",fitMin, fitMax);
00471         TF1* fitFunc = tdcByMedianH->GetFunction("gaus");
00472         if( fitFunc ){
00473           tdcOffset = fitFunc->GetParameter(1);
00474           tdcOffsetUncert = fitFunc->GetParError(1);
00475           
00476           if( fitFunc->GetNDF() < 3 ) 
00477             // Catch bad fits
00478             tdcOffsetUncert = 0;
00479           else
00480             // Scale uncertainty by the quality of the fit
00481             tdcOffsetUncert *= sqrt(fitFunc->GetChisquare() / fitFunc->GetNDF());
00482           
00483           tdcOffsetH->SetBinContent(tdcOffsetH->FindBin(column),tdcOffset);
00484           tdcOffsetH->SetBinError(tdcOffsetH->FindBin(column),tdcOffsetUncert);
00485 
00486           tdcOffsetChannelH->SetBinContent(ring*24+column,tdcOffset);
00487           tdcOffsetChannelH->SetBinError(ring*24+column,tdcOffsetUncert);
00488 
00489         }else{
00490           warning() << "PMT  " << pmtId << " has no TDC data for fitting." 
00491                     << endreq;
00492         }
00493       }
00494 
00495 
00496 
00497       // Record ring averages, ignoring bad channels
00498       if( occupancyOk && gainOk ){
00499         if( occupancyUncert > 0 ){
00500           meanOccRing[ring] += occupancy / (occupancyUncert*occupancyUncert); 
00501           meanOccWeight[ring] += 1.0 / (occupancyUncert*occupancyUncert); 
00502         }
00503         if( tdcOffsetUncert > 0 ){
00504           meanTdcRing[ring] += tdcOffset / (tdcOffsetUncert*tdcOffsetUncert); 
00505           meanTdcWeight[ring] += 1.0 / (tdcOffsetUncert*tdcOffsetUncert); 
00506         }
00507       }
00508 
00509     }
00510 
00511     // Record mean time offset and mean occupancy by ring
00512     for(int ring = 1; ring <= nRings; ring++){
00513       if( meanOccWeight[ring] > 0 ){
00514         meanOccupancyH->SetBinContent(ring,
00515                                       meanOccRing[ring] / meanOccWeight[ring]);
00516         meanOccupancyH->SetBinError(ring, 
00517                                     sqrt( 1.0 / meanOccWeight[ring] ));
00518       }
00519       if( meanTdcWeight[ring] > 0 ){
00520         meanTdcOffsetH->SetBinContent(ring,
00521                                       meanTdcRing[ring] / meanTdcWeight[ring]);
00522         meanTdcOffsetH->SetBinError(ring, 
00523                                     sqrt( 1.0 / meanTdcWeight[ring] ));
00524       }
00525     } // Loop over rings
00526  
00527   } // Loop over detectors
00528  
00529   return StatusCode::SUCCESS;
00530 }

bool PmtCalibFullModel::hasStats ( const DayaBay::Detector detector  )  [private]

Check if the statistics for this detector have been initialized.

Definition at line 533 of file PmtCalibFullModel.cc.

00533                                                                {
00534   // Check if statistics have been initialized for this detector
00535   return (std::find(m_processedDetectors.begin(), 
00536                     m_processedDetectors.end(), 
00537                     detector) 
00538           != m_processedDetectors.end());
00539 }

StatusCode PmtCalibFullModel::prepareStats ( const Context context  )  [private]

Get the current statistics for a detector, or create if needed.

Definition at line 541 of file PmtCalibFullModel.cc.

00541                                                                 {
00542   // Create the histograms and parameters for this detector's statistics
00543   DayaBay::Detector detector(context.GetSite(), context.GetDetId());
00544 
00545   // Site
00546   {
00547     std::string name = "site";
00548     std::ostringstream path;
00549     path << this->getPath(detector) << "/" << name;
00550     TParameter<int>* par = new TParameter<int>(name.c_str(), 
00551                                                context.GetSite());
00552     if( m_statsSvc->put(path.str(),par).isFailure() ) {
00553       error() << "prepareStats(): Could not register " << name 
00554               << " at " << path << endreq;
00555       delete par;
00556       par = 0;
00557       return StatusCode::FAILURE;
00558     }
00559   }
00560 
00561   // Detector ID
00562   {
00563     std::string name = "detectorId";
00564     std::ostringstream path;
00565     path << this->getPath(detector) << "/" << name;
00566     TParameter<int>* par = new TParameter<int>(name.c_str(), 
00567                                                context.GetDetId());
00568     if( m_statsSvc->put(path.str(),par).isFailure() ) {
00569       error() << "prepareStats(): Could not register " << name 
00570               << " at " << path << endreq;
00571       delete par;
00572       par = 0;
00573       return StatusCode::FAILURE;
00574     }
00575   }
00576 
00577   // SimFlag
00578   {
00579     std::string name = "simFlag";
00580     std::ostringstream path;
00581     path << this->getPath(detector) << "/" << name;
00582     TParameter<int>* par = new TParameter<int>(name.c_str(), 
00583                                                context.GetSimFlag());
00584     if( m_statsSvc->put(path.str(),par).isFailure() ) {
00585       error() << "prepareStats(): Could not register " << name 
00586               << " at " << path << endreq;
00587       delete par;
00588       par = 0;
00589       return StatusCode::FAILURE;
00590     }
00591   }
00592 
00593   // Timestamp: seconds
00594   {
00595     std::string name = "timeSec";
00596     std::ostringstream path;
00597     path << this->getPath(detector) << "/" << name;
00598     TParameter<int>* par = new TParameter<int>(name.c_str(), 
00599                                                context.GetTimeStamp().GetSec());
00600     if( m_statsSvc->put(path.str(),par).isFailure() ) {
00601       error() << "prepareStats(): Could not register " << name 
00602               << " at " << path << endreq;
00603       delete par;
00604       par = 0;
00605       return StatusCode::FAILURE;
00606     }
00607   }
00608 
00609   // Timestamp: nanoseconds
00610   {
00611     std::string name = "timeNanoSec";
00612     std::ostringstream path;
00613     path << this->getPath(detector) << "/" << name;
00614     TParameter<int>* par = 
00615       new TParameter<int>(name.c_str(), context.GetTimeStamp().GetNanoSec());
00616     if( m_statsSvc->put(path.str(),par).isFailure() ) {
00617       error() << "prepareStats(): Could not register " << name 
00618               << " at " << path << endreq;
00619       delete par;
00620       par = 0;
00621       return StatusCode::FAILURE;
00622     }
00623   }
00624 
00625   // Number of Readouts
00626   {
00627     std::string name = "nReadouts";
00628     std::ostringstream path;
00629     path << this->getPath(detector) << "/" << name;
00630     TParameter<int>* par = new TParameter<int>(name.c_str(), 0);
00631     if( m_statsSvc->put(path.str(),par).isFailure() ) {
00632       error() << "prepareStats(): Could not register " << name 
00633               << " at " << path << endreq;
00634       delete par;
00635       par = 0;
00636       return StatusCode::FAILURE;
00637     }
00638   }
00639 
00640   // ADC Sum spectrum
00641   {
00642     std::ostringstream title, path;
00643     std::string name = "adcSum";
00644     title << "ADC Sum for readouts in " << detector.detName(); 
00645     path << this->getPath(detector) << "/" << name;
00646     TH1F* adcSum = new TH1F(name.c_str(),title.str().c_str(),
00647                             2000,0,20000);
00648     adcSum->GetXaxis()->SetTitle("ADC Sum value");
00649     adcSum->GetYaxis()->SetTitle("Entries");
00650     // DEBUG
00651     info() << "name= " << adcSum->GetName() 
00652            << " at path = \"" << path.str() << "\"" << endreq;
00653     if( m_statsSvc->put(path.str(),adcSum).isFailure() ) {
00654       error() << "prepareStats(): Could not register " << path << endreq;
00655       delete adcSum;
00656       return StatusCode::FAILURE;
00657     }
00658   }
00659 
00660   // channel Gain
00661   {
00662     std::ostringstream title, path;
00663     std::string name = "gainChannel";
00664     title << "gainChannel " << detector.detName(); 
00665     path << this->getPath(detector) << "/" << name;
00666     TH1F* gainChannel = new TH1F(name.c_str(),title.str().c_str(),
00667                             300,0,300);
00668     gainChannel->GetXaxis()->SetTitle("channel ID");
00669     gainChannel->GetYaxis()->SetTitle("gain in ADC bin");
00670     // DEBUG
00671     info() << "name= " << gainChannel->GetName() 
00672            << " at path = \"" << path.str() << "\"" << endreq;
00673     if( m_statsSvc->put(path.str(),gainChannel).isFailure() ) {
00674       error() << "prepareStats(): Could not register " << path << endreq;
00675       delete gainChannel;
00676       return StatusCode::FAILURE;
00677     }
00678   }
00679 
00680   // channel tdc offset
00681   {
00682     std::ostringstream title, path;
00683     std::string name = "tdcOffsetChannel";
00684     title << "tdcOffsetChannel " << detector.detName(); 
00685     path << this->getPath(detector) << "/" << name;
00686     TH1F* tdcOffsetChannel = new TH1F(name.c_str(),title.str().c_str(),
00687                             300,0,300);
00688     tdcOffsetChannel->GetXaxis()->SetTitle("channel ID");
00689     tdcOffsetChannel->GetYaxis()->SetTitle("tdc offset in TDC bin");
00690     // DEBUG
00691     info() << "name= " << tdcOffsetChannel->GetName() 
00692            << " at path = \"" << path.str() << "\"" << endreq;
00693     if( m_statsSvc->put(path.str(),tdcOffsetChannel).isFailure() ) {
00694       error() << "prepareStats(): Could not register " << path << endreq;
00695       delete tdcOffsetChannel;
00696       return StatusCode::FAILURE;
00697     }
00698   }
00699 
00700   // channel occupancy
00701   {
00702     std::ostringstream title, path;
00703     std::string name = "occupancyChannel";
00704     title << "occupancyChannel " << detector.detName(); 
00705     path << this->getPath(detector) << "/" << name;
00706     TH1F* occupancyChannel = new TH1F(name.c_str(),title.str().c_str(),
00707                             300,0,300);
00708     occupancyChannel->GetXaxis()->SetTitle("channel ID");
00709     occupancyChannel->GetYaxis()->SetTitle("occupancy in mean pe");
00710     // DEBUG
00711     info() << "name= " << occupancyChannel->GetName() 
00712            << " at path = \"" << path.str() << "\"" << endreq;
00713     if( m_statsSvc->put(path.str(),occupancyChannel).isFailure() ) {
00714       error() << "prepareStats(): Could not register " << path << endreq;
00715       delete occupancyChannel;
00716       return StatusCode::FAILURE;
00717     }
00718   }
00719 
00720 
00721   // TDC Median spectrum
00722   {
00723     std::ostringstream title, path;
00724     std::string name = "tdcMedian";
00725     title << "Median TDC for readouts in " << detector.detName(); 
00726     path << this->getPath(detector) << "/" << name;
00727     TH1F* tdcMedian = new TH1F(name.c_str(),title.str().c_str(),
00728                                300,0,300);
00729     tdcMedian->GetXaxis()->SetTitle("TDC value");
00730     tdcMedian->GetYaxis()->SetTitle("Entries");
00731     if( m_statsSvc->put(path.str(),tdcMedian).isFailure() ) {
00732       error() << "prepareStats(): Could not register " << path << endreq;
00733       delete tdcMedian;
00734       return StatusCode::FAILURE;
00735     }
00736   }
00737   // Occupancy for each channel (jpochoa)
00738     {
00739       std::ostringstream title, path;
00740       std::string name = "occupancy";
00741       title << "Occupancy per channel " << detector.detName();
00742       path << this->getPath(detector) << "/" << name;
00743       TH2F* occupancy = new TH2F(name.c_str(),title.str().c_str(),49,0.25,24.75,19,-0.75,8.75);
00744       occupancy->GetXaxis()->SetTitle("Column");
00745       occupancy->GetYaxis()->SetTitle("Ring");
00746       if( m_statsSvc->put(path.str(),occupancy).isFailure() ) {
00747         error() << "prepareStats(): Could not register " << path << endreq;
00748         delete occupancy;
00749         return StatusCode::FAILURE;
00750       }
00751     }
00752   // Gain for each channel (jpochoa)
00753   {
00754       std::ostringstream title, path;
00755       std::string name = "gain2D";
00756       title << "Gain per channel " << detector.detName();
00757       path << this->getPath(detector) << "/" << name;
00758       TH2F* gain2D = new TH2F(name.c_str(),title.str().c_str(),49,0.25,24.75,19,-0.75,8.75);
00759       gain2D->GetXaxis()->SetTitle("Column");
00760       gain2D->GetYaxis()->SetTitle("Ring");
00761       if( m_statsSvc->put(path.str(),gain2D).isFailure() ) {
00762         error() << "prepareStats(): Could not register " << path << endreq;
00763         delete gain2D;
00764         return StatusCode::FAILURE;
00765       }
00766     }
00767     // Gain for all channels (jpochoa)
00768   {
00769       std::ostringstream title, path;
00770       std::string name = "gain1D";
00771       title << "Gain for all channels " << detector.detName();
00772       path << this->getPath(detector) << "/" << name;
00773       TH1F* gain1D = new TH1F(name.c_str(),title.str().c_str(),50,20,60);
00774       gain1D->GetXaxis()->SetTitle("ADC/PE");
00775       gain1D->GetYaxis()->SetTitle("Number of entries");
00776       if( m_statsSvc->put(path.str(),gain1D).isFailure() ) {
00777         error() << "prepareStats(): Could not register " << path << endreq;
00778         delete gain1D;
00779         return StatusCode::FAILURE;
00780       }
00781     }
00782     
00783 
00784   // Mean TDC Offset by AD ring
00785   {
00786     std::ostringstream title, path;
00787     std::string name = "meanTdcOffset";
00788     title << "Mean TDC offset by ring in " << detector.detName(); 
00789     path << this->getPath(detector) << "/" << name;
00790     TH1F* meanTdcOffset = new TH1F(name.c_str(),title.str().c_str(),
00791                                    8,1,9);
00792     meanTdcOffset->GetXaxis()->SetTitle("AD Ring");
00793     meanTdcOffset->GetYaxis()->SetTitle("Mean TDC Offset");
00794     if( m_statsSvc->put(path.str(),meanTdcOffset).isFailure() ) {
00795       error() << "prepareStats(): Could not register " << path << endreq;
00796       delete meanTdcOffset;
00797       return StatusCode::FAILURE;
00798     }
00799   }
00800   // Mean Occupancy by AD ring
00801   {
00802     std::ostringstream title, path;
00803     std::string name = "meanOccupancy";
00804     title << "Mean occupancy by ring in " << detector.detName(); 
00805     path << this->getPath(detector) << "/" << name;
00806     TH1F* meanOccupancy = new TH1F(name.c_str(),title.str().c_str(),
00807                                    8,1,9);
00808     meanOccupancy->GetXaxis()->SetTitle("AD Ring");
00809     meanOccupancy->GetYaxis()->SetTitle("Mean Occupancy");
00810     if( m_statsSvc->put(path.str(),meanOccupancy).isFailure() ) {
00811       error() << "prepareStats(): Could not register " << path << endreq;
00812       delete meanOccupancy;
00813       return StatusCode::FAILURE;
00814     }
00815   }
00816 
00817   // Make sure summary histograms for each ring exist in detector statistics
00818   for(int ring = 0; ring <= 8; ring++){
00819     // Occupancy by ring
00820     {
00821       std::ostringstream title, path;
00822       std::string name = "occupancy";
00823       title << "Occupancy for PMTs in " << detector.detName() 
00824             << " ring " << ring; 
00825       path << this->getPath(detector, ring) << "/" << name;
00826       TH1F* occupancy = new TH1F(name.c_str(),title.str().c_str(),
00827                                  24,1,25);
00828       occupancy->GetXaxis()->SetTitle("Column");
00829       occupancy->GetYaxis()->SetTitle("Occupancy");
00830       if( m_statsSvc->put(path.str(),occupancy).isFailure() ) {
00831         error() << "prepareStats(): Could not register " << path << endreq;
00832         delete occupancy;
00833         return StatusCode::FAILURE;
00834       }
00835     }
00836     // TDC offset by ring
00837     {
00838       std::ostringstream title, path;
00839       std::string name = "tdcOffset";
00840       title << "Time Offset for PMTs in " << detector.detName() 
00841             << " ring " << ring; 
00842       path << this->getPath(detector, ring) << "/" << name;
00843       TH1F* tdcOffset = new TH1F(name.c_str(),title.str().c_str(),
00844                            24,1,25);
00845       tdcOffset->GetXaxis()->SetTitle("Column");
00846       tdcOffset->GetYaxis()->SetTitle("Time Offset [tdc]");
00847       if( m_statsSvc->put(path.str(),tdcOffset).isFailure() ) {
00848         error() << "prepareStats(): Could not register " << path << endreq;
00849         delete tdcOffset;
00850         return StatusCode::FAILURE;
00851       }
00852     }
00853     // ADC Median by ring
00854     {
00855       std::ostringstream title, path;
00856       std::string name = "adcMedian";
00857       title << "ADC Median for PMTs in " << detector.detName() 
00858             << " ring " << ring; 
00859       path << this->getPath(detector, ring) << "/" << name;
00860       TH1F* adcMedian = new TH1F(name.c_str(),title.str().c_str(),
00861                                  24,1,25);
00862       adcMedian->GetXaxis()->SetTitle("Column");
00863       adcMedian->GetYaxis()->SetTitle("Median ADC");
00864       if( m_statsSvc->put(path.str(),adcMedian).isFailure() ) {
00865         error() << "prepareStats(): Could not register " << path << endreq;
00866         delete adcMedian;
00867         return StatusCode::FAILURE;
00868       }
00869     }
00870     // ADC Sigma by ring
00871     {
00872       std::ostringstream title, path;
00873       std::string name = "adcSigma";
00874       title << "ADC Sigma for PMTs in " << detector.detName() 
00875             << " ring " << ring;
00876       path << this->getPath(detector, ring) << "/" << name;
00877       TH1F* adcSigma = new TH1F(name.c_str(),title.str().c_str(),
00878                                 24,1,25);
00879       adcSigma->GetXaxis()->SetTitle("Column");
00880       adcSigma->GetYaxis()->SetTitle("ADC Sigma");
00881       if( m_statsSvc->put(path.str(),adcSigma).isFailure() ) {
00882         error() << "prepareStats(): Could not register " << path << endreq;
00883         delete adcSigma;
00884         return StatusCode::FAILURE;
00885       }
00886     }
00887   }
00888 
00889   ServiceMode svcMode(context, 0);
00890 
00891   // Get list of all FEE channels
00892   const std::vector<DayaBay::FeeChannelId>& channelList 
00893     = m_cableSvc->feeChannelIds( svcMode );
00894   std::vector<DayaBay::FeeChannelId>::const_iterator chanIter, 
00895     chanEnd = channelList.end();
00896   // Initialize statistics for each channel
00897   for(chanIter = channelList.begin(); chanIter != chanEnd; chanIter++){
00898     DayaBay::FeeChannelId chanId = *chanIter;
00899 
00900     // Board Number
00901     {
00902       std::string name = "board";
00903       std::ostringstream path;
00904       path << this->getPath(chanId) << "/" << name;
00905       TParameter<int>* par = new TParameter<int>(name.c_str(), chanId.board());
00906       if( m_statsSvc->put(path.str(),par).isFailure() ) {
00907         error() << "prepareStats(): Could not register " << name 
00908                 << " at " << path << endreq;
00909         delete par;
00910         par = 0;
00911         return StatusCode::FAILURE;
00912       }
00913     }
00914 
00915     // Connector Number
00916     {
00917       std::string name = "connector";
00918       std::ostringstream path;
00919       path << this->getPath(chanId) << "/" << name;
00920       TParameter<int>* par = new TParameter<int>(name.c_str(), 
00921                                                  chanId.connector());
00922       if( m_statsSvc->put(path.str(),par).isFailure() ) {
00923         error() << "prepareStats(): Could not register " << name 
00924                 << " at " << path << endreq;
00925         delete par;
00926         par = 0;
00927         return StatusCode::FAILURE;
00928       }
00929     }
00930 
00931     // Number of Hits
00932     {
00933       std::string name = "nHits";
00934       std::ostringstream path;
00935       path << this->getPath(chanId) << "/" << name;
00936       TParameter<int>* par = new TParameter<int>(name.c_str(), 0);
00937       if( m_statsSvc->put(path.str(),par).isFailure() ) {
00938         error() << "prepareStats(): Could not register " << name 
00939                 << " at " << path << endreq;
00940         delete par;
00941         par = 0;
00942         return StatusCode::FAILURE;
00943       }
00944     }
00945 
00946     // Raw TDC spectrum
00947     {
00948       std::string name = "tdcRaw";
00949       std::ostringstream title, path;
00950       title << "Raw TDC values for " << detector.detName() 
00951             << " channel " << chanId.board() << "_"
00952             << chanId.connector();
00953       path << this->getPath(chanId) << "/" << name;
00954       TH1F* tdcRaw = new TH1F(name.c_str(),title.str().c_str(),
00955                               300,0,300);
00956       tdcRaw->GetXaxis()->SetTitle("TDC value");
00957       tdcRaw->GetYaxis()->SetTitle("Entries");
00958       if( m_statsSvc->put(path.str(),tdcRaw).isFailure() ) {
00959         error() << "prepareStats(): Could not register " << path << endreq;
00960         delete tdcRaw;
00961         return StatusCode::FAILURE;
00962       }
00963     }
00964 
00965     // TDC spectrum, median corrected
00966     {
00967       std::string name = "tdcByMedian";
00968       std::ostringstream title, path;
00969       title << "Median-corrected TDC values for " << detector.detName() 
00970             << " channel " << chanId.board() << "_"
00971             << chanId.connector();
00972       path << this->getPath(chanId) << "/" << name;
00973       TH1F* tdcByMedian = new TH1F(name.c_str(),title.str().c_str(),
00974                                    600,-300,300);
00975       tdcByMedian->GetXaxis()->SetTitle("TDC value");
00976       tdcByMedian->GetYaxis()->SetTitle("Entries");
00977       if( m_statsSvc->put(path.str(),tdcByMedian).isFailure() ) {
00978         error() << "prepareStats(): Could not register " << path << endreq;
00979         delete tdcByMedian;
00980         return StatusCode::FAILURE;
00981       }
00982     }
00983 
00984     // Raw ADC spectrum
00985     {
00986       std::ostringstream title, path;
00987       std::string name = "adcRaw";
00988       title << "ADC values for " << detector.detName() 
00989             << " channel " << chanId.board() << "_"
00990             << chanId.connector();
00991       path << this->getPath(chanId) << "/" << name;
00992       TH1F* adcRaw = new TH1F(name.c_str(),title.str().c_str(),
00993                               4096,0,4096);
00994       adcRaw->GetXaxis()->SetTitle("ADC value");
00995       adcRaw->GetYaxis()->SetTitle("Entries");
00996       if( m_statsSvc->put(path.str(),adcRaw).isFailure() ) {
00997         error() << "prepareStats(): Could not register " << path << endreq;
00998         delete adcRaw;
00999         return StatusCode::FAILURE;
01000       }
01001     }
01002     
01003     // ADC spectrum, baseline subtracted
01004     {
01005       std::ostringstream title, path;
01006       std::string name = "adc";
01007       title << "Baseline-subtracted ADC values for " << detector.detName() 
01008             << " channel " << chanId.board() << "_"
01009             << chanId.connector();
01010       path << this->getPath(chanId) << "/" << name;
01011       TH1F* adc = new TH1F(name.c_str(),title.str().c_str(),
01012                            5096,-1000,4096);
01013       adc->GetXaxis()->SetTitle("ADC value");
01014       adc->GetYaxis()->SetTitle("Entries");
01015       if( m_statsSvc->put(path.str(),adc).isFailure() ) {
01016         error() << "prepareStats(): Could not register " << path << endreq;
01017         delete adc;
01018         return StatusCode::FAILURE;
01019       }
01020     }
01021 
01022     // ADC spectrum by Clock Cycle
01023     {
01024       std::ostringstream title, path;
01025       std::string name = "adcByClock";
01026       title << "ADC values by clock cycle for " << detector.detName() 
01027             << " channel " << chanId.board() << "_"
01028             << chanId.connector();
01029       path << this->getPath(chanId) << "/" << name;
01030       TH1F* adcByClock = new TH1F(name.c_str(),title.str().c_str(),
01031                                   20,0,20);
01032       adcByClock->GetXaxis()->SetTitle("ADC Clock Cycle");
01033       adcByClock->GetYaxis()->SetTitle("ADC");
01034       if( m_statsSvc->put(path.str(),adcByClock).isFailure() ) {
01035         error() << "prepareStats(): Could not register " << path << endreq;
01036         delete adcByClock;
01037         return StatusCode::FAILURE;
01038       }
01039     }
01040   }  
01041 
01042   m_processedDetectors.push_back(detector);
01043 
01044   return StatusCode::SUCCESS;
01045 }

std::string PmtCalibFullModel::getPath ( const DayaBay::Detector detector  )  [private]

Definition at line 1047 of file PmtCalibFullModel.cc.

01047                                                                      {
01048   return m_filepath + "/" + detector.detName();
01049 }

std::string PmtCalibFullModel::getPath ( const DayaBay::Detector detector,
int  ring 
) [private]

Definition at line 1051 of file PmtCalibFullModel.cc.

01052                                                   {
01053   std::ostringstream path;
01054   path << m_filepath << "/" << detector.detName() << "/ring_" << ring; 
01055   return path.str();
01056 }

std::string PmtCalibFullModel::getPath ( const DayaBay::FeeChannelId channelId  )  [private]

Definition at line 1058 of file PmtCalibFullModel.cc.

01058                                                                      {
01059   std::ostringstream path;
01060   path << m_filepath << "/" << chanId.detName() 
01061        << "/chan_" << chanId.board() << "_" << chanId.connector(); 
01062   return path.str();
01063 }

const InterfaceID & IPmtCalibParamTool::interfaceID (  )  [static, inherited]

Retrieve interface ID.

Reimplemented from IAlgTool.

Definition at line 8 of file IPmtCalibParamTool.cc.

00009 { 
00010     return IID_IPmtCalibParamTool; 
00011 }


Member Data Documentation

std::string PmtCalibFullModel::m_cableSvcName [private]

Definition at line 111 of file PmtCalibFullModel.h.

std::string PmtCalibFullModel::m_calibSvcName [private]

Definition at line 114 of file PmtCalibFullModel.h.

std::string PmtCalibFullModel::m_floatFeePedesSvcName [private]

Definition at line 117 of file PmtCalibFullModel.h.

std::string PmtCalibFullModel::m_filepath [private]

Definition at line 120 of file PmtCalibFullModel.h.

ICableSvc* PmtCalibFullModel::m_cableSvc [private]

Definition at line 123 of file PmtCalibFullModel.h.

ICalibDataSvc* PmtCalibFullModel::m_calibSvc [private]

Definition at line 126 of file PmtCalibFullModel.h.

IStatisticsSvc* PmtCalibFullModel::m_statsSvc [private]

Definition at line 129 of file PmtCalibFullModel.h.

IFloatingFeePedestalSvc* PmtCalibFullModel::m_floatFeePedesSvc [private]

Definition at line 132 of file PmtCalibFullModel.h.

bool PmtCalibFullModel::m_useFloatFeePedes [private]

Definition at line 135 of file PmtCalibFullModel.h.

std::vector<DayaBay::Detector> PmtCalibFullModel::m_processedDetectors [private]

Definition at line 138 of file PmtCalibFullModel.h.


The documentation for this class was generated from the following files:
| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

Generated on Mon Apr 11 20:29:41 2011 for CalibParam by doxygen 1.4.7