00001
00002
00003
00004 #include <sstream>
00005 #include <string>
00006 #include <stdio.h>
00007
00008
00009
00010 #include "TGraph.h"
00011 #include "TH2F.h"
00012 #include "TH1F.h"
00013 #include "TH1I.h"
00014 #include <TProfile.h>
00015
00016
00017
00018 #include "DaqDataHistogram.h"
00019
00020
00021
00022 #include "Event/ReadoutHeader.h"
00023 #include "Event/Readout.h"
00024 #include "Event/DaqCrate.h"
00025
00026 #include "EventReadoutFormat/EventHeader.h"
00027 #include "EventReadoutFormat/EventReadout.h"
00028 #include "EventReadoutFormat/EventTraits.h"
00029 #include "EventReadoutFormat/RomFragment.h"
00030 #include "EventReadoutFormat/RomHeader.h"
00031 #include "FeeReadoutFormat/FeeHit.h"
00032 #include "FeeReadoutFormat/FeeReadout.h"
00033 #include "FileReadoutFormat/DataSeparatorRecord.h"
00034 #include "LtbReadoutFormat/LtbFrame.h"
00035 #include "LtbReadoutFormat/LtbReadout.h"
00036 #include "LtbReadoutFormat/LtbTime.h"
00037 #include "FadcReadoutFormat/FadcBuffer.h"
00038 #include "FadcReadoutFormat/FadcReadout.h"
00039 #include "FadcReadoutFormat/FadcTraits.h"
00040 #include "FadcReadoutFormat/FadcData.h"
00041 #include "StatisticsSvc/IStatisticsSvc.h"
00042
00043 #include "FecReadoutFormat/FecReadout.h"
00044 #include "FecReadoutFormat/FecData.h"
00045 #include "RtmReadoutFormat/RtmReadout.h"
00046 #include "RtmReadoutFormat/RtmData.h"
00047
00048 #include "DaqReadoutFormat/RomData.h"
00049
00050 using namespace std;
00051 using namespace DayaBay;
00052 using namespace DybDaq;
00053
00054 DaqDataHistogram::DaqDataHistogram(const string& name,
00055 ISvcLocator* pSvcLocator) :
00056 GaudiAlgorithm(name,
00057 pSvcLocator),
00058 m_statsSvc(0),
00059 m_daqCrateOrNot(0),
00060 m_firstTriggerTime(0),
00061 m_eventCount(0)
00062 {
00063 declareProperty("PrintFreq", m_printFreq=0);
00064 declareProperty("AdcSumMax", m_adcSumMax=100000);
00065 declareProperty("HighGainFactor", m_highGainFactor=10);
00066 }
00067
00068 DaqDataHistogram::~DaqDataHistogram()
00069 {
00070 }
00071
00072 StatusCode DaqDataHistogram::initialize()
00073 {
00074 info() << "Initializing " << endreq;
00075
00076 StatusCode sc = service("StatisticsSvc", m_statsSvc, true);
00077 if (sc.isFailure()) {
00078 error() << "Failed to get StatisticsSvc" << endreq;
00079 return sc;
00080 }
00081
00082 verbose()<< "done init" << endreq;
00083 return sc;
00084 }
00085
00086 StatusCode DaqDataHistogram::finalize()
00087 {
00088 info() << "Finalizing " << endreq;
00089
00090 if(m_daqCrateOrNot){
00091 double scale = 1/(m_lastTriggerTime - m_firstTriggerTime);
00092 m_scale->Scale(scale);
00093 }
00094
00095 if( m_statsSvc ) {
00096 m_statsSvc->release();
00097 }
00098
00099 return StatusCode::SUCCESS;
00100 }
00101
00102 StatusCode DaqDataHistogram::execute()
00103 {
00104 verbose() << "executing: " << m_eventCount << endreq;
00105 const ReadoutHeader* readoutHeader = get<ReadoutHeader>("/Event/Readout/ReadoutHeader");
00106 if (0 == readoutHeader) {
00107 error() << "Failed to get readout header." << endreq;
00108 return StatusCode::FAILURE;
00109 }
00110
00111 const DayaBay::DaqCrate* daqCrate = readoutHeader->daqCrate();
00112
00113 const DybDaq::EventReadout* event = 0;
00114
00115 if ( daqCrate == 0 ) {
00116 const Readout* readout = readoutHeader->readout();
00117 if ( readout == 0 ) {
00118 error() << "Failed to get DAQ readout or readout from header" << endreq;
00119 return StatusCode::FAILURE;
00120 }
00121 event = readout->eventReadout();
00122 if ( event == 0 ) {
00123 error() << "Failed to get event readout from readout" << endreq;
00124 return StatusCode::FAILURE;
00125 }
00126 } else {
00127 event = &(daqCrate->eventReadout());
00128 }
00129
00130 if(!event){
00131 error() << "Failed to get event both from DaqCrate and Readout" << endreq;
00132 return StatusCode::SUCCESS;
00133 } else {
00134 const DybDaq::EventHeader* eventHeader = &(event->header());
00135 int run = eventHeader->run();
00136 int detectorInt = eventHeader->detector();
00137 if(detectorInt == 7){
00138 TH1F* h_fecId = dynamic_cast<TH1F*>(rpcGetOrMakeHist(1, detectorInt, 0, 0, kFecId));
00139
00140 static const EventTraits& eventTraits = dynamic_cast<const EventTraits&>(event->header().daqTraits());
00141 static const unsigned int kRpcRomModuleType = eventTraits.moduleType(EventTraits::kRpcRomModule);
00142 static const unsigned int kRpcRtmModuleType = eventTraits.moduleType(EventTraits::kRpcRtmModule);
00143
00144 const EventReadout::RomFragmentPtrList& fragments = event->romFragments();
00145 for(EventReadout::RomFragmentPtrList::const_iterator fragment = fragments.begin();
00146 fragment != fragments.end(); ++fragment) {
00147
00148 const RomHeader& header = (*fragment)->header();
00149
00150 if ( kRpcRomModuleType == header.moduleType() ) {
00151 handleFecReadout( dynamic_cast<const FecReadout&>((*fragment)->unwrappedData()) );
00152
00153 } else if ( kRpcRtmModuleType == header.moduleType() ) {
00154 handleRtmReadout( dynamic_cast<const RtmReadout&>((*fragment)->unwrappedData()) );
00155
00156 } else {
00157 std::cout << "Current fragment data is not RPC!" << std::endl;
00158 }
00159 }
00160 }
00161 }
00162
00163 if(daqCrate){
00164 m_daqCrateOrNot = 1;
00165
00166 const DayaBay::DaqPmtCrate* pmtCrate = daqCrate->asPmtCrate();
00167 if(!pmtCrate){
00168
00169 return StatusCode::SUCCESS;
00170 }
00171
00172 const DayaBay::Detector& detector = pmtCrate->detector();
00173 info() << "Detector = " << detector << endreq;
00174
00175 if(!m_firstTriggerTime){
00176 m_firstTriggerTime = pmtCrate->triggerTime().GetSeconds();
00177 }
00178
00179 int runNumber = pmtCrate->runNumber();
00180 TH1F* h_triggerType = dynamic_cast<TH1F*>
00181 (getOrMakeHist(1, detector, 0, 0, kTriggerType));
00182 TH1F* h_rateVsTriggerType = dynamic_cast<TH1F*>
00183 (getOrMakeHist(1, detector, 0, 0, kRateVsTriggerType));
00184 TH1F* h_adcSum = dynamic_cast<TH1F*>
00185 (getOrMakeHist(1, detector, 0, 0, kAdcSum));
00186 TH1F* h_nChannel = dynamic_cast<TH1F*>
00187 (getOrMakeHist(1, detector, 0, 0, kNchannel));
00188 TH1F* h_hitNumber = dynamic_cast<TH1F*>
00189 (getOrMakeHist(1, detector, 0, 0, kHitNumber));
00190 TH1F* h_hitRate = dynamic_cast<TH1F*>
00191 (getOrMakeHist(1, detector, 0, 0, kHitRate));
00192 TH1F* h_adcMean = dynamic_cast<TH1F*>
00193 (getOrMakeHist(1, detector, 0, 0, kAdcMean));
00194 TH1F* h_adcRMS = dynamic_cast<TH1F*>
00195 (getOrMakeHist(1, detector, 0, 0, kAdcRMS));
00196 TH1F* h_tdcMean = dynamic_cast<TH1F*>
00197 (getOrMakeHist(1, detector, 0, 0, kTdcMean));
00198 TH1F* h_tdcRMS = dynamic_cast<TH1F*>
00199 (getOrMakeHist(1, detector, 0, 0, kTdcRMS));
00200 TH1F* h_preAdcMean = dynamic_cast<TH1F*>
00201 (getOrMakeHist(1, detector, 0, 0, kPreAdcMean));
00202 TH1F* h_preAdcRMS = dynamic_cast<TH1F*>
00203 (getOrMakeHist(1, detector, 0, 0, kPreAdcRMS));
00204 TH2F* h_adcVsChannel = dynamic_cast<TH2F*>
00205 (getOrMakeHist(1, detector, 0, 0, kAdcVsChannel));
00206 TH2F* h_tdcVsChannel = dynamic_cast<TH2F*>
00207 (getOrMakeHist(1, detector, 0, 0, kTdcVsChannel));
00208 TH2F* h_preAdcVsChannel = dynamic_cast<TH2F*>
00209 (getOrMakeHist(1, detector, 0, 0, kPreAdcVsChannel));
00210 TH1F* h_darkRate = dynamic_cast<TH1F*>
00211 (getOrMakeHist(1, detector, 0, 0, kDarkRate));
00212 TH2F* h_adcSumVsTriggerType = dynamic_cast<TH2F*>
00213 (getOrMakeHist(1, detector, 0, 0, kAdcSumVsTriggerType));
00214 TProfile* h_darkNoiseVsChannel = dynamic_cast<TProfile*>
00215 (getOrMakeHist(1, detector, 0, 0, kDarkNoiseVsChannel));
00216
00217 h_hitNumber->Fill(-1);
00218 h_darkRate->Fill(-1);
00219
00220 const DayaBay::DaqPmtCrate::PmtChannelPtrList& channels
00221 = pmtCrate->channelReadouts();
00222
00223 h_nChannel->Fill( channels.size() );
00224
00225 DayaBay::DaqPmtCrate::PmtChannelPtrList::const_iterator channelIter,
00226 channelEnd = channels.end();
00227
00228 float adcSum = 0;
00229 for(channelIter = channels.begin(); channelIter!=channelEnd; channelIter++) {
00230 const DayaBay::DaqPmtChannel& channel = *(*channelIter);
00231
00232 int board = channel.channelId().board();
00233 int connector = channel.channelId().connector();
00234 int localId = board*16+connector;
00235 info() << "(board,connector) = (" << board << "," << connector << ")" << endreq;
00236
00237 h_hitNumber->Fill(localId);
00238
00239 TH1I* h_tdc = dynamic_cast<TH1I*>
00240 (getOrMakeHist(1, detector, board, connector, kTdc));
00241 TH1I* h_peakCycle = dynamic_cast<TH1I*>
00242 (getOrMakeHist(1, detector, board, connector, kPeakCycle));
00243 TH1I* h_range = dynamic_cast<TH1I*>
00244 (getOrMakeHist(1, detector, board, connector, kRange));
00245 TH1I* h_hitCount = dynamic_cast<TH1I*>
00246 (getOrMakeHist(1, detector, board, connector, kHitCount));
00247 TH1I* h_adcFine = dynamic_cast<TH1I*>
00248 (getOrMakeHist(1, detector, board, connector, kAdcFine));
00249 TH1I* h_adcCoarse = dynamic_cast<TH1I*>
00250 (getOrMakeHist(1, detector, board, connector, kAdcCoarse));
00251 TH1I* h_preAdcFine = dynamic_cast<TH1I*>
00252 (getOrMakeHist(1, detector, board, connector, kPreAdcFine));
00253 TH1I* h_preAdcCoarse = dynamic_cast<TH1I*>
00254 (getOrMakeHist(1, detector, board, connector, kPreAdcCoarse));
00255 TH1I* h_darkNoise = dynamic_cast<TH1I*>
00256 (getOrMakeHist(1, detector, board, connector, kDarkNoise));
00257
00258 for(unsigned int i=0; i<channel.hitCount(); i++) {
00259 int tdc = channel.tdc(i);
00260 int adc = channel.adc(i);
00261 float preAdc = channel.preAdcAvg(i);
00262 float deltaAdc = adc - preAdc;
00263 debug() << "tdc=" << tdc << ", adc=" << adc << endreq;
00264
00265 h_tdc->Fill(tdc);
00266 h_peakCycle->Fill(channel.peakCycle(i));
00267 h_hitCount->Fill(i+1);
00268 if(channel.isHighGainAdc(i)) {
00269 h_range->Fill(1);
00270 h_adcFine->Fill(adc);
00271 h_preAdcFine->Fill(preAdc);
00272 if(tdc>=1070&&tdc<1170) {
00273 h_darkNoise->Fill(deltaAdc);
00274 h_darkNoiseVsChannel->Fill(localId,deltaAdc);
00275 }
00276 adcSum += deltaAdc;
00277 h_adcVsChannel->Fill(localId,deltaAdc);
00278 h_preAdcVsChannel->Fill(localId,preAdc);
00279 } else {
00280 h_range->Fill(2);
00281 h_adcCoarse->Fill(adc);
00282 h_preAdcCoarse->Fill(preAdc);
00283 adcSum += (deltaAdc)*m_highGainFactor;
00284 h_adcVsChannel->Fill(localId,(deltaAdc)*m_highGainFactor);
00285 h_adcVsChannel->Fill(localId,(preAdc)*m_highGainFactor);
00286 }
00287 h_tdcVsChannel->Fill(localId,tdc);
00288 }
00289 }
00290
00291 h_adcSum->Fill(adcSum);
00292
00293 int trigType = daqTrigType(pmtCrate->triggerType());
00294 for(unsigned i=0; i<16; i++) {
00295 if(trigType & 1<<i) {
00296 h_triggerType->Fill(i);
00297 h_rateVsTriggerType->Fill(i);
00298 h_adcSumVsTriggerType->Fill(i,adcSum);
00299 }
00300 }
00301
00302 m_lastTriggerTime = pmtCrate->triggerTime().GetSeconds();
00303 ++m_eventCount;
00304
00305 if (0 == (m_eventCount % 1000)) {
00306 info() << " Processed " << m_eventCount << " events" << endreq;
00307 }
00308 }
00309 return StatusCode::SUCCESS;
00310 }
00311
00312 void DaqDataHistogram::handleFecReadout(const FecReadout& fecReadout){
00313 const FecReadout::FecDataPtrList& fecList = fecReadout.fecDataList();
00314
00315 FecReadout::FecDataPtrList::const_iterator it = fecList.begin();
00316 while ( it != fecList.end() ) {
00317 handleFecData( *it );
00318
00319 ++it;
00320 }
00321 }
00322
00323 void DaqDataHistogram::handleFecData(const FecData* fecData){
00324 TH1F* h_fecId = m_statsSvc->getTH1F(m_rpcHistPath[kFecId].c_str());
00325 h_fecId->Fill(fecData->rpcFecId());
00326 }
00327
00328 void DaqDataHistogram::handleRtmReadout(const RtmReadout& rtmReadout){
00329 const RtmReadout::RtmDataPtrList& rtmList = rtmReadout.rtmDataList();
00330
00331 RtmReadout::RtmDataPtrList::const_iterator it = rtmList.begin();
00332 while ( it != rtmList.end() ) {
00333 handleRtmData( *it );
00334
00335 const RtmData* rtmData = *it;
00336 ++it;
00337 }
00338 }
00339
00340 void DaqDataHistogram::handleRtmData(const RtmData* rtmData){
00341 cout << " ** RtmData: " << std::endl;
00342
00343
00344
00345
00346 }
00347
00348 TH1* DaqDataHistogram::rpcGetOrMakeHist(int run, int detector, int row, int column, int histIndex)
00349 {
00350 debug() << "Detector=" << detector << endreq;
00351
00352 std::map<int, TH1**>::iterator histIter = m_rpcHist.find(detector);
00353 TH1** histograms = 0;
00354 if( histIter == m_rpcHist.end() ){
00355 info() << "Create " << MAXRPCHISTS << " histograms for detector " << detector << endreq;
00356
00357 histograms = new TH1*[MAXRPCHISTS];
00358 for(int i = 0; i < MAXRPCHISTS; i++) histograms[i] = 0;
00359 m_rpcHist[detector] = histograms;
00360 } else {
00361
00362 histograms = histIter->second;
00363 }
00364
00365 int iHist = row * NCOLUMN * NRPCHISTOGRAMS + column * NRPCHISTOGRAMS + histIndex;
00366 info() << "Histogram index = " << iHist << endreq;
00367 TH1* hist = histograms[iHist];
00368 if(!hist){
00369
00370 std::string histName;
00371 if(row){
00372 switch(histIndex){
00373 case kHitMap:
00374 hist = new TH1F("HitMap", "HitMap", 32, -0.5, 31.5);
00375 hist->GetXaxis()->SetTitle("ChannelId");
00376 hist->GetYaxis()->SetTitle("Entries");
00377 break;
00378 default:
00379 error() << "Unknown histogram for single channel: " << histIndex << endreq;
00380 return 0;
00381 }
00382 } else {
00383 switch(histIndex){
00384 case kFecId:
00385 hist = new TH1F("h_sum_FecId", "h_sum_FecId", 20, 0, 20);
00386 hist->GetXaxis()->SetTitle("FecId");
00387 hist->GetYaxis()->SetTitle("Entries");
00388 break;
00389 default:
00390 error() << "Unknown histogram for single channel: " << histIndex << endreq;
00391 return 0;
00392 }
00393 }
00394
00395 info() << "Making histogram: " << hist->GetName() << endreq;
00396 std::string rpcHistPath = rpcGetPath(run, detector, row, column, hist->GetName());
00397 m_rpcHistPath[histIndex] = rpcHistPath;
00398 m_statsSvc->put( rpcHistPath, hist );
00399 histograms[iHist] = hist;
00400 }
00401
00402 return hist;
00403 }
00404
00405 std::string DaqDataHistogram::rpcGetPath(int run, int detector,
00406 int row, int column, const char* histName)
00407 {
00408
00409 std::stringstream path;
00410
00411
00412 path << "/file1/EH1/RPC/";
00413 if(row){
00414 path << "row" << std::setfill('0') << std::setw(2)
00415 << row << "_column" << std::setfill('0') << std::setw(2)
00416 << column;
00417 } else {
00418 path << "Raw";
00419 }
00420 path << "/" << histName;
00421 debug() << "Histogram path = " << path << endreq;
00422 return path.str();
00423 }
00424
00425 TH1* DaqDataHistogram::getOrMakeHist(int run, const DayaBay::Detector& detector,
00426 int board, int connector, int histIndex)
00427 {
00428 debug() << "Detector=" << detector << "," << detector.fullPackedData() << endreq;
00429
00430 int detectorInt = detector.fullPackedData();
00431 std::map<int,TH1**>::iterator histIter = m_hist.find(detectorInt);
00432 TH1** histograms = 0;
00433 if( histIter == m_hist.end() ){
00434 info() << "Create " << MAXHISTS << " histograms for " << detector << endreq;
00435
00436 histograms = new TH1*[MAXHISTS];
00437 for(int i=0; i<MAXHISTS; i++) histograms[i] = 0;
00438 m_hist[detectorInt] = histograms;
00439 }else{
00440
00441 histograms = histIter->second;
00442 }
00443
00444 int iHist = board * NCONNECTORS * NHISTOGRAMS + connector * NHISTOGRAMS + histIndex;
00445 info() << "Histogram index = " << iHist << endreq;
00446 TH1* hist = histograms[board * NCONNECTORS * NHISTOGRAMS + connector * NHISTOGRAMS + histIndex];
00447 if(!hist){
00448
00449 std::string histName;
00450 if(detector.detectorId()){
00451 if(board){
00452
00453 switch(histIndex){
00454 case kTdc:
00455 hist = new TH1I("TDC","TDC",1000,300,1300);
00456 hist->GetXaxis()->SetTitle("TDC");
00457 hist->GetYaxis()->SetTitle("Entries");
00458 break;
00459 case kAdcFine:
00460
00461 hist = new TH1I("ADCFine","ADC (Fine Range)",4096,0,4096);
00462 hist->GetXaxis()->SetTitle("ADC Fine");
00463 hist->GetYaxis()->SetTitle("Entries");
00464 break;
00465 case kAdcCoarse:
00466
00467 hist = new TH1I("ADCCoarse","ADC (Coarse Range)",4096,0,4096);
00468 hist->GetXaxis()->SetTitle("ADC Coarse");
00469 hist->GetYaxis()->SetTitle("Entries");
00470 break;
00471 case kPeakCycle:
00472 hist = new TH1I("PeakCycle","PeakCycle",14,0,14);
00473 hist->GetXaxis()->SetTitle("PeakCycle");
00474 hist->GetYaxis()->SetTitle("Entries");
00475 break;
00476 case kRange:
00477 hist = new TH1I("Range","Range",3,0,3);
00478 hist->GetXaxis()->SetTitle("Range");
00479 hist->GetYaxis()->SetTitle("Entries");
00480 break;
00481 case kHitCount:
00482 hist = new TH1I("HitCount","HitCount",20,0,20);
00483 hist->GetXaxis()->SetTitle("HitCount");
00484 hist->GetYaxis()->SetTitle("Entries");
00485 break;
00486 case kPreAdcFine:
00487
00488 hist = new TH1I("PreADCFine","PreADC (Fine Range)",200,0,200);
00489 hist->GetXaxis()->SetTitle("PreADC Fine");
00490 hist->GetYaxis()->SetTitle("Entries");
00491 break;
00492 case kPreAdcCoarse:
00493
00494 hist = new TH1I("PreADCCoarse","PreADC (Coarse Range)",200,0,200);
00495 hist->GetXaxis()->SetTitle("PreADC Coarse");
00496 hist->GetYaxis()->SetTitle("Entries");
00497 break;
00498 case kDarkNoise:
00499 hist = new TH1I("DarkNoise","DarkNoise",250,-50,200);
00500 hist->GetXaxis()->SetTitle("Dark Noise");
00501 hist->GetYaxis()->SetTitle("Entries");
00502 break;
00503 default:
00504 error() << "Unknown histogram for single channel: " << histIndex << endreq;
00505 return 0;
00506 }
00507 } else {
00508
00509 switch(histIndex){
00510 case kTriggerType:
00511 hist = new TH1F("h_sum_triggerType", "h_triggerType;triggerType", 16, 0, 16);
00512 break;
00513 case kRateVsTriggerType:
00514 hist = new TH1F("h_sum_rateVsTriggerType", "h_rateVsTriggerType;triggerType", 16, 0, 16);
00515 hist->GetYaxis()->SetTitle("Rate");
00516 m_scale = hist;
00517 break;
00518 case kAdcSum:
00519 hist = new TH1F("h_sum_FEEADCSUM", "h_FEEADCSUM", 1000 , 0, m_adcSumMax);
00520 break;
00521 case kNchannel:
00522 hist = new TH1F("h_sum_nChannel","Number of Channels in Event",200,0,200);
00523 hist->GetXaxis()->SetTitle("Number of Channels");
00524 break;
00525 case kHitRate:
00526 hist = new TH1F("h_sum_hitRate", "h_hitRate;channelIndex;hit_rate", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00527 break;
00528 case kHitNumber:
00529 hist = new TH1F("h_tmp_hitNumber", "h_hitNumber;channelIndex;hit_number", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00530 break;
00531 case kAdcMean:
00532 hist = new TH1F("h_sum_ADCMean", "h_ADCMean;channelIndex", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00533 break;
00534 case kAdcRMS:
00535 hist = new TH1F("h_sum_ADCRMS", "h_ADCRMS;channelIndex", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00536 break;
00537 case kTdcMean:
00538 hist = new TH1F("h_sum_TDCMean", "h_TDCMean;channelIndex", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00539 break;
00540 case kTdcRMS:
00541 hist = new TH1F("h_sum_TDCRMS", "h_TDCRMS;channelIndex", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00542 break;
00543 case kPreAdcMean:
00544 hist = new TH1F("h_sum_preAdcMean", "h_preAdcMean;channelIndex", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00545 break;
00546 case kPreAdcRMS:
00547 hist = new TH1F("h_sum_preAdcRMS", "h_preAdcRMS;channelIndex", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00548 break;
00549 case kAdcVsChannel:
00550 hist = new TH2F("h_sum_ADC_Channel", "h_sum_ADC_Channel", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5, (int)(4096*m_highGainFactor), 0, 4096*m_highGainFactor);
00551 hist->GetXaxis()->SetTitle("Board#times16+Channel");
00552 hist->GetYaxis()->SetTitle("ADC");
00553 break;
00554 case kTdcVsChannel:
00555 hist = new TH2F("h_sum_TDC_Channel", "h_sum_TDC_Channel", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5, 1000, 300, 1300);
00556 hist->GetXaxis()->SetTitle("Board#times16+Channel");
00557 hist->GetYaxis()->SetTitle("TDC");
00558 break;
00559 case kPreAdcVsChannel:
00560 hist = new TH2F("h_sum_preADC_Channel", "h_sum_preADC_Channel", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5, (int)(4096*m_highGainFactor), 0, 4096*m_highGainFactor);
00561 hist->GetXaxis()->SetTitle("Board#times16+Channel");
00562 hist->GetYaxis()->SetTitle("preADC");
00563 break;
00564 case kDarkRate:
00565 hist = new TH1F("h_sum_DarkRate", "h_DarkRate", MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5);
00566 hist->GetXaxis()->SetTitle("Board#times16+Channel");
00567 hist->GetYaxis()->SetTitle("Dark Rate (kHz)");
00568 break;
00569 case kDarkNoiseVsChannel:
00570 hist = new TProfile("h_sum_DarkNoise","h_sum_DarkNoise",MAXMODNUM*MAXCHNNUM , 0.5, MAXMODNUM*MAXCHNNUM+0.5,-50,250);
00571 hist->GetXaxis()->SetTitle("Board#times16+Channel");
00572 hist->GetYaxis()->SetTitle("Dark Noise (ADC)");
00573 break;
00574 case kAdcSumVsTriggerType:
00575 hist = new TH2F("h_sum_AdcSum_TriggerType","h_sum_AdcSum_TriggerType", 16, 0, 16, 1000 , 0, m_adcSumMax);
00576 hist->GetXaxis()->SetTitle("Trigger Type");
00577 hist->GetYaxis()->SetTitle("ADCSum");
00578 break;
00579 default:
00580 error() << "Unknown histogram for event: " << histIndex << endreq;
00581 return 0;
00582 }
00583 }
00584 }
00585 debug() << "Making histogram: " << hist->GetName() << endreq;
00586 m_statsSvc->put( getPath(run, detector, board, connector, hist->GetName()), hist );
00587 histograms[ board * NCONNECTORS * NHISTOGRAMS + connector * NHISTOGRAMS + histIndex] = hist;
00588 }
00589
00590 return hist;
00591 }
00592
00593 std::string DaqDataHistogram::getPath(int run, const Detector &detector,
00594 int board, int connector, const char* histName)
00595 {
00596
00597 std::stringstream path;
00598
00599
00600 path << "/file1/";
00601 path << Site::AsString(detector.site()) << "/" << DetectorId::AsString(detector.detectorId());
00602 if(board){
00603 path << "/board" << std::setfill('0') << std::setw(2)
00604 << board << "_connector" << std::setfill('0') << std::setw(2)
00605 << connector;
00606 } else {
00607 path << "/Raw";
00608 }
00609 path << "/" << histName;
00610 debug() << "Histogram path = " << path << endreq;
00611 return path.str();
00612 }
00613
00614 int DaqDataHistogram::daqTrigType(int offTrigType)
00615 {
00616 int trigType = 0;
00617 if(offTrigType & DayaBay::Trigger::kMult) {
00618 trigType = trigType | (1<<8);
00619 }
00620 if(offTrigType & DayaBay::Trigger::kExternal) {
00621 trigType = trigType | (1<<1);
00622 }
00623 if(offTrigType & DayaBay::Trigger::kManual) {
00624 trigType = trigType | (1<<0);
00625 }
00626 if(offTrigType & DayaBay::Trigger::kPeriodic) {
00627 trigType = trigType | (1<<2);
00628 }
00629 if(offTrigType & DayaBay::Trigger::kESumAll) {
00630 trigType = trigType | (1<<12);
00631 }
00632 if(offTrigType & DayaBay::Trigger::kESumHigh) {
00633 trigType = trigType | (1<<10);
00634 }
00635 if(offTrigType & DayaBay::Trigger::kESumLow) {
00636 trigType = trigType | (1<<11);
00637 }
00638 if(offTrigType & DayaBay::Trigger::kESumADC) {
00639 trigType = trigType | (1<<9);
00640 }
00641 if(offTrigType & DayaBay::Trigger::kCross) {
00642 trigType = trigType | (1<<1);
00643 }
00644 if(offTrigType & DayaBay::Trigger::kESum) {
00645 trigType = trigType | (1<<12);
00646 }
00647 return trigType;
00648 }