CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

SiPixelDaqInfo Class Reference

#include <SiPixelDaqInfo.h>

Inheritance diagram for SiPixelDaqInfo:
edm::EDAnalyzer

List of all members.

Public Member Functions

 SiPixelDaqInfo (const edm::ParameterSet &)
 ~SiPixelDaqInfo ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void beginLuminosityBlock (const edm::LuminosityBlock &, const edm::EventSetup &)
virtual void endJob ()
virtual void endLuminosityBlock (const edm::LuminosityBlock &, const edm::EventSetup &)
virtual void endRun (const edm::Run &, const edm::EventSetup &)

Private Attributes

std::string daqSource_
DQMStoredbe_
std::pair< int, int > FEDRange_
int FEDs_ [40]
MonitorElementFraction_
MonitorElementFractionBarrel_
MonitorElementFractionEndcap_
int NEvents_
int nFEDsBarrel_
int nFEDsEndcap_
int NumberOfFeds_

Detailed Description

Definition at line 25 of file SiPixelDaqInfo.h.


Constructor & Destructor Documentation

SiPixelDaqInfo::SiPixelDaqInfo ( const edm::ParameterSet ps) [explicit]

Definition at line 11 of file SiPixelDaqInfo.cc.

References edm::ParameterSet::getUntrackedParameter(), and i.

                                                        {
 
  FEDRange_.first  = ps.getUntrackedParameter<unsigned int>("MinimumPixelFEDId", 0);
  FEDRange_.second = ps.getUntrackedParameter<unsigned int>("MaximumPixelFEDId", 39);
  daqSource_       = ps.getUntrackedParameter<string>("daqSource",  "source");

  NumberOfFeds_ =FEDRange_.second -  FEDRange_.first +1;
  
  NEvents_ = 0;
  for(int i=0; i!=40; i++) FEDs_[i] = 0;

}
SiPixelDaqInfo::~SiPixelDaqInfo ( )

Definition at line 24 of file SiPixelDaqInfo.cc.

{}

Member Function Documentation

void SiPixelDaqInfo::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 165 of file SiPixelDaqInfo.cc.

References FEDRawDataCollection::FEDData(), edm::Event::getByLabel(), i, edm::HandleBase::isValid(), and FEDRawData::size().

                                                                               {
  NEvents_++;  
  //cout<<"in SiPixelDaqInfo::analyze now!"<<endl;
  if(NEvents_>=1 && NEvents_<=100){
    // check if any Pixel FED is in readout:
    edm::Handle<FEDRawDataCollection> rawDataHandle;
    iEvent.getByLabel(daqSource_, rawDataHandle);
    if(!rawDataHandle.isValid()){
      edm::LogInfo("SiPixelDaqInfo") << daqSource_ << " is empty!";
      return;
    }
    const FEDRawDataCollection& rawDataCollection = *rawDataHandle;
    nFEDsBarrel_ = 0; nFEDsEndcap_ = 0;
    for(int i = 0; i != 40; i++){
      if(rawDataCollection.FEDData(i).size() > 208 ) FEDs_[i]++;
    }
  }

}
void SiPixelDaqInfo::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 147 of file SiPixelDaqInfo.cc.

References DQMStore::bookFloat(), dbe_, cppFunctionSkipper::operator, and DQMStore::setCurrentFolder().

                             {

  dbe_ = 0;
  dbe_ = Service<DQMStore>().operator->();
  
 
  dbe_->setCurrentFolder("Pixel/EventInfo");
  Fraction_= dbe_->bookFloat("DAQSummary");  
  dbe_->setCurrentFolder("Pixel/EventInfo/DAQContents");
  FractionBarrel_= dbe_->bookFloat("PixelBarrelFraction");  
  FractionEndcap_= dbe_->bookFloat("PixelEndcapFraction");  
}
void SiPixelDaqInfo::beginLuminosityBlock ( const edm::LuminosityBlock lumiBlock,
const edm::EventSetup iSetup 
) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 26 of file SiPixelDaqInfo.cc.

{}
void SiPixelDaqInfo::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 161 of file SiPixelDaqInfo.cc.

{}
void SiPixelDaqInfo::endLuminosityBlock ( const edm::LuminosityBlock lumiBlock,
const edm::EventSetup iSetup 
) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 29 of file SiPixelDaqInfo.cc.

References edm::EventSetup::find(), edm::eventsetup::heterocontainer::HCTypeTag::findType(), edm::EventSetup::get(), and i.

                                                                                                         {
  edm::eventsetup::EventSetupRecordKey recordKey(edm::eventsetup::EventSetupRecordKey::TypeTag::findType("RunInfoRcd"));
  if(0 != iSetup.find( recordKey ) ) {
    // cout<<"record key found"<<endl;
    //get fed summary information
    ESHandle<RunInfo> sumFED;
    iSetup.get<RunInfoRcd>().get(sumFED);    
    vector<int> FedsInIds= sumFED->m_fed_in;   

    int FedCount=0;
    int FedCountBarrel=0;
    int FedCountEndcap=0;

    //loop on all active feds
    for(unsigned int fedItr=0;fedItr<FedsInIds.size(); ++fedItr) {
      int fedID=FedsInIds[fedItr];
      //make sure fed id is in allowed range  
      //cout<<fedID<<endl;   
      if(fedID>=FEDRange_.first && fedID<=FEDRange_.second){
        ++FedCount;
        if(fedID>=0 && fedID<=31) ++FedCountBarrel;
        else if(fedID>=32 && fedID<=39) ++FedCountEndcap;
      }
    }   
    
    //Fill active fed fraction ME
    if(FedCountBarrel<=32){
      FedCountBarrel = 0; FedCountEndcap = 0; FedCount = 0; NumberOfFeds_ = 40;
      for(int i=0; i!=40; i++){
        if(i<=31 && FEDs_[i]>0) FedCountBarrel++;
        if(i>=32 && FEDs_[i]>0) FedCountEndcap++;
        if(FEDs_[i]>0) FedCount++;
      }
    }
    if(NumberOfFeds_>0){
      //all Pixel:
      Fraction_->Fill( FedCount/NumberOfFeds_);
      //Barrel:
      FractionBarrel_->Fill( FedCountBarrel/32.);
      //Endcap:
      FractionEndcap_->Fill( FedCountEndcap/8.);
    }else{
      Fraction_->Fill(-1);
      FractionBarrel_->Fill(-1);
      FractionEndcap_->Fill(-1);
    } 
    
  }else{      
    Fraction_->Fill(-1);               
    FractionBarrel_->Fill(-1);
    FractionEndcap_->Fill(-1);
    return; 
  }
}
void SiPixelDaqInfo::endRun ( const edm::Run r,
const edm::EventSetup iSetup 
) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 84 of file SiPixelDaqInfo.cc.

References edm::EventSetup::find(), edm::eventsetup::heterocontainer::HCTypeTag::findType(), edm::EventSetup::get(), and i.

                                                                         {
  edm::eventsetup::EventSetupRecordKey recordKey(edm::eventsetup::EventSetupRecordKey::TypeTag::findType("RunInfoRcd"));
  if(0 != iSetup.find( recordKey ) ) {
    // cout<<"record key found"<<endl;
    //get fed summary information
    ESHandle<RunInfo> sumFED;
    iSetup.get<RunInfoRcd>().get(sumFED);    
    vector<int> FedsInIds= sumFED->m_fed_in;   

    int FedCount=0;
    int FedCountBarrel=0;
    int FedCountEndcap=0;

    //loop on all active feds
    for(unsigned int fedItr=0;fedItr<FedsInIds.size(); ++fedItr) {
      int fedID=FedsInIds[fedItr];
      //make sure fed id is in allowed range  
      //cout<<fedID<<endl;   
      if(fedID>=FEDRange_.first && fedID<=FEDRange_.second){
        ++FedCount;
        if(fedID>=0 && fedID<=31) ++FedCountBarrel;
        else if(fedID>=32 && fedID<=39) ++FedCountEndcap;
      }
    }   

    if(FedCountBarrel>32){
      FedCountBarrel = nFEDsBarrel_;
      FedCountEndcap = nFEDsEndcap_;
      FedCount = FedCountBarrel + FedCountEndcap;
      NumberOfFeds_ = 40;
    }

    //Fill active fed fraction ME
    if(FedCountBarrel<=32){
      FedCountBarrel = 0; FedCountEndcap = 0; FedCount = 0; NumberOfFeds_ = 40;
      for(int i=0; i!=40; i++){
        if(i<=31 && FEDs_[i]>0) FedCountBarrel++;
        if(i>=32 && FEDs_[i]>0) FedCountEndcap++;
        if(FEDs_[i]>0) FedCount++;
      }
    }
    if(NumberOfFeds_>0){
      //all Pixel:
      Fraction_->Fill( FedCount/NumberOfFeds_);
      //Barrel:
      FractionBarrel_->Fill( FedCountBarrel/32.);
      //Endcap:
      FractionEndcap_->Fill( FedCountEndcap/8.);
    }else{
      Fraction_->Fill(-1);
      FractionBarrel_->Fill(-1);
      FractionEndcap_->Fill(-1);
    } 
    
  }else{      
    Fraction_->Fill(-1);               
    FractionBarrel_->Fill(-1);
    FractionEndcap_->Fill(-1);
    return; 
  }
}

Member Data Documentation

std::string SiPixelDaqInfo::daqSource_ [private]

Definition at line 52 of file SiPixelDaqInfo.h.

Definition at line 39 of file SiPixelDaqInfo.h.

std::pair<int,int> SiPixelDaqInfo::FEDRange_ [private]

Definition at line 45 of file SiPixelDaqInfo.h.

int SiPixelDaqInfo::FEDs_[40] [private]

Definition at line 53 of file SiPixelDaqInfo.h.

Definition at line 41 of file SiPixelDaqInfo.h.

Definition at line 42 of file SiPixelDaqInfo.h.

Definition at line 43 of file SiPixelDaqInfo.h.

int SiPixelDaqInfo::NEvents_ [private]

Definition at line 49 of file SiPixelDaqInfo.h.

Definition at line 50 of file SiPixelDaqInfo.h.

Definition at line 51 of file SiPixelDaqInfo.h.

Definition at line 47 of file SiPixelDaqInfo.h.