CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

APVCyclePhaseMonitor Class Reference

#include <DPGAnalysis/SiStripTools/plugins/APVCyclePhaseMonitor.cc>

Inheritance diagram for APVCyclePhaseMonitor:
edm::EDAnalyzer

List of all members.

Public Member Functions

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

Private Member Functions

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

Private Attributes

const edm::InputTag _apvphasecollection
std::map< std::string, TH1F * > _hphases
std::map< std::string, TProfile * > _hphasevsorbit
std::map< std::string, TH1F * > _hselectedphases
std::map< std::string, TH1F * > _hselectedphasessize
std::map< std::string, TH1F * > _hselectedphasesvector
std::map< std::string, TProfile * > _hselectedphasevectorvsorbit
std::map< std::string, TProfile * > _hselectedphasevsorbit
unsigned int _nevents
std::vector< std::string > _selectedparts
std::vector< std::string > _selectedvectorparts

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 49 of file APVCyclePhaseMonitor.cc.


Constructor & Destructor Documentation

APVCyclePhaseMonitor::APVCyclePhaseMonitor ( const edm::ParameterSet iConfig) [explicit]

Definition at line 88 of file APVCyclePhaseMonitor.cc.

References _apvphasecollection.

                                                                        :
  _apvphasecollection(iConfig.getParameter<edm::InputTag>("apvCyclePhaseCollection")),
  _selectedparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedPartitions",std::vector<std::string>())),
  _selectedvectorparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedVectorPartitions",std::vector<std::string>())),
  _hphases(),_hselectedphases(),_hselectedphasesvector(),_hselectedphasessize(),
  _hphasevsorbit(),_hselectedphasevsorbit(),_hselectedphasevectorvsorbit(),
  _nevents(0)
{
   //now do what ever initialization is needed

  edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecollection << " used";

}
APVCyclePhaseMonitor::~APVCyclePhaseMonitor ( )

Definition at line 103 of file APVCyclePhaseMonitor.cc.

{
 
   // do anything here that needs to be done at desctruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 118 of file APVCyclePhaseMonitor.cc.

References _apvphasecollection, _hphases, _hphasevsorbit, _hselectedphases, _hselectedphasessize, _hselectedphasesvector, _hselectedphasevectorvsorbit, _hselectedphasevsorbit, _nevents, _selectedparts, _selectedvectorparts, combineCards::dirname, edm::Event::getByLabel(), TFileDirectory::mkdir(), edm::EventBase::orbitNumber(), and edm::Event::run().

{
   using namespace edm;

   _nevents++;

   edm::Handle<APVCyclePhaseCollection> apvphases;
   iEvent.getByLabel(_apvphasecollection,apvphases);

   // improve the matchin between default and actual partitions
   
   edm::Service<TFileService> tfserv;

   for(std::map<std::string,int>::const_iterator phase = apvphases->get().begin(); phase != apvphases->get().end(); ++phase) {

     if(_hphases.find(phase->first)==_hphases.end()) {
       char dirname[300];
       sprintf(dirname,"run_%d",iEvent.run());
       TFileDirectory subrun = tfserv->mkdir(dirname);

       char hname[300];

       sprintf(hname,"phase_%s",phase->first.c_str());
       edm::LogInfo("TH1FBeingBooked") << "TH1F " << hname << " being booked" ;
       _hphases[phase->first] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
       _hphases[phase->first]->GetXaxis()->SetTitle("BX mod 70"); _hphases[phase->first]->GetYaxis()->SetTitle("Events");

       sprintf(hname,"phasevsorbit_%s",phase->first.c_str());
       edm::LogInfo("TProfileBeingBooked") << "TProfile " << hname << " being booked" ;
       _hphasevsorbit[phase->first] = subrun.make<TProfile>(hname,hname,1800,0.,1800*11223);
       _hphasevsorbit[phase->first]->SetBit(TH1::kCanRebin);
       _hphasevsorbit[phase->first]->GetXaxis()->SetTitle("time [orbit#]"); _hphasevsorbit[phase->first]->GetYaxis()->SetTitle("Phase");
       
     }
     _hphases[phase->first]->Fill(phase->second);
     _hphasevsorbit[phase->first]->Fill(iEvent.orbitNumber(),phase->second);
   }

   // selected partitions
   
   for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
     if(_hselectedphases.find(*part)!=_hselectedphases.end()) {
       _hselectedphases[*part]->Fill(apvphases->getPhase(*part));
     }
     if(_hselectedphasevsorbit.find(*part)!=_hselectedphasevsorbit.end()) {
       _hselectedphasevsorbit[*part]->Fill(iEvent.orbitNumber(),apvphases->getPhase(*part));
     }
   }

   for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
       part!=_selectedvectorparts.end();++part) {

     const std::vector<int> phases = apvphases->getPhases(*part);

      if(_hselectedphasessize.find(*part)!=_hselectedphasessize.end()) {
        _hselectedphasessize[*part]->Fill(phases.size());
      }

     for(std::vector<int>::const_iterator phase=phases.begin();phase!=phases.end();++phase) {
       if(_hselectedphasesvector.find(*part)!=_hselectedphasesvector.end()) {
         _hselectedphasesvector[*part]->Fill(*phase);
       }
       if(_hselectedphasevectorvsorbit.find(*part)!=_hselectedphasevectorvsorbit.end()) {
         _hselectedphasevectorvsorbit[*part]->Fill(iEvent.orbitNumber(),*phase);
       }
     }

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

Reimplemented from edm::EDAnalyzer.

Definition at line 257 of file APVCyclePhaseMonitor.cc.

{

}
void APVCyclePhaseMonitor::beginRun ( const edm::Run iRun,
const edm::EventSetup  
) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 189 of file APVCyclePhaseMonitor.cc.

References _hphases, _hphasevsorbit, _hselectedphases, _hselectedphasessize, _hselectedphasesvector, _hselectedphasevectorvsorbit, _hselectedphasevsorbit, _selectedparts, _selectedvectorparts, combineCards::dirname, and edm::RunBase::run().

{
  
  _hphases.clear();
  _hselectedphases.clear();
  _hselectedphasesvector.clear();
  _hselectedphasessize.clear();
  _hphasevsorbit.clear();
  _hselectedphasevsorbit.clear();
  _hselectedphasevectorvsorbit.clear();

  edm::Service<TFileService> tfserv;

  char dirname[300];
  sprintf(dirname,"run_%d",iRun.run());
  TFileDirectory subrun = tfserv->mkdir(dirname);
  
  for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {

    char hname[300];
    
    sprintf(hname,"selected_phase_%s",part->c_str());
    edm::LogInfo("SelectedTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
    _hselectedphases[*part] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
    _hselectedphases[*part]->GetXaxis()->SetTitle("BX mod 70"); _hselectedphases[*part]->GetYaxis()->SetTitle("Events");
    
    sprintf(hname,"selected_phasevsorbit_%s",part->c_str());
    edm::LogInfo("SelectedTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
    _hselectedphasevsorbit[*part] = subrun.make<TProfile>(hname,hname,1800,0.,1800*11223);
    _hselectedphasevsorbit[*part]->SetBit(TH1::kCanRebin);
    _hselectedphasevsorbit[*part]->GetXaxis()->SetTitle("time [orbit#]"); 
    _hselectedphasevsorbit[*part]->GetYaxis()->SetTitle("Phase");
  }

  for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
      part!=_selectedvectorparts.end();++part) {

    char hname[300];
    
    sprintf(hname,"selected_phase_vector_%s",part->c_str());
    edm::LogInfo("SelectedVectTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
    _hselectedphasesvector[*part] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
    _hselectedphasesvector[*part]->GetXaxis()->SetTitle("BX mod 70"); 
    _hselectedphasesvector[*part]->GetYaxis()->SetTitle("Events");
    
    sprintf(hname,"selected_phase_size_%s",part->c_str());
    edm::LogInfo("SelectedVectSizeTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
    _hselectedphasessize[*part] = subrun.make<TH1F>(hname,hname,10,-0.5,9.5);
    _hselectedphasessize[*part]->GetXaxis()->SetTitle("Number of Phases"); 
    _hselectedphasessize[*part]->GetYaxis()->SetTitle("Events");
    
    sprintf(hname,"selected_phasevectorvsorbit_%s",part->c_str());
    edm::LogInfo("SelectedVectTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
    _hselectedphasevectorvsorbit[*part] = subrun.make<TProfile>(hname,hname,1800,0.,1800*11223);
    _hselectedphasevectorvsorbit[*part]->SetBit(TH1::kCanRebin);
    _hselectedphasevectorvsorbit[*part]->GetXaxis()->SetTitle("time [orbit#]"); 
    _hselectedphasevectorvsorbit[*part]->GetYaxis()->SetTitle("Phase");
  }
}
void APVCyclePhaseMonitor::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 264 of file APVCyclePhaseMonitor.cc.

References _nevents.

                             {

  edm::LogInfo("EndOfJob") << _nevents << " analyzed events";

}
void APVCyclePhaseMonitor::endRun ( const edm::Run iRun,
const edm::EventSetup  
) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 250 of file APVCyclePhaseMonitor.cc.

{
}

Member Data Documentation

Definition at line 64 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and APVCyclePhaseMonitor().

std::map<std::string,TH1F*> APVCyclePhaseMonitor::_hphases [private]

Definition at line 67 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TProfile*> APVCyclePhaseMonitor::_hphasevsorbit [private]

Definition at line 71 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TH1F*> APVCyclePhaseMonitor::_hselectedphases [private]

Definition at line 68 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TH1F*> APVCyclePhaseMonitor::_hselectedphasessize [private]

Definition at line 70 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TH1F*> APVCyclePhaseMonitor::_hselectedphasesvector [private]

Definition at line 69 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TProfile*> APVCyclePhaseMonitor::_hselectedphasevectorvsorbit [private]

Definition at line 73 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TProfile*> APVCyclePhaseMonitor::_hselectedphasevsorbit [private]

Definition at line 72 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

unsigned int APVCyclePhaseMonitor::_nevents [private]

Definition at line 74 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and endJob().

std::vector<std::string> APVCyclePhaseMonitor::_selectedparts [private]

Definition at line 65 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::vector<std::string> APVCyclePhaseMonitor::_selectedvectorparts [private]

Definition at line 66 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().