#include <DQMOffline/Muon/src/MuonTrackResidualsTest.h>
Public Member Functions | |
MuonTrackResidualsTest (const edm::ParameterSet &ps) | |
Constructor. | |
virtual | ~MuonTrackResidualsTest () |
Destructor. | |
Protected Member Functions | |
void | analyze (const edm::Event &e, const edm::EventSetup &c) |
Analyze. | |
void | beginJob (const edm::EventSetup &c) |
BeginJob. | |
void | beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) |
void | endJob () |
Endjob. | |
void | endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c) |
DQM Client Diagnostic. | |
Private Attributes | |
std::map< std::string, std::vector< std::string > > | histoNames |
std::map< std::string, MonitorElement * > | MeanHistos |
std::string | metname |
int | nevents |
unsigned int | nLumiSegs |
edm::ParameterSet | parameters |
int | prescaleFactor |
int | run |
std::map< std::string, MonitorElement * > | SigmaHistos |
DQMStore * | theDbe |
Definition at line 39 of file MuonTrackResidualsTest.h.
MuonTrackResidualsTest::MuonTrackResidualsTest | ( | const edm::ParameterSet & | ps | ) |
Constructor.
Definition at line 36 of file MuonTrackResidualsTest.cc.
References edm::ParameterSet::getUntrackedParameter(), parameters, prescaleFactor, and theDbe.
00036 { 00037 parameters = ps; 00038 00039 theDbe = edm::Service<DQMStore>().operator->(); 00040 00041 prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1); 00042 00043 }
MuonTrackResidualsTest::~MuonTrackResidualsTest | ( | ) | [virtual] |
void MuonTrackResidualsTest::analyze | ( | const edm::Event & | e, | |
const edm::EventSetup & | c | |||
) | [protected, virtual] |
Analyze.
Implements edm::EDAnalyzer.
Definition at line 112 of file MuonTrackResidualsTest.cc.
References LogTrace, metname, and nevents.
00112 { 00113 00114 nevents++; 00115 LogTrace(metname)<< "[MuonTrackResidualsTest]: "<<nevents<<" events"; 00116 00117 }
void MuonTrackResidualsTest::beginJob | ( | const edm::EventSetup & | c | ) | [protected, virtual] |
BeginJob.
Reimplemented from edm::EDAnalyzer.
Definition at line 53 of file MuonTrackResidualsTest.cc.
References DQMStore::book1D(), c, lat::endl(), histoNames, LogTrace, MeanHistos, metname, nevents, DQMStore::setCurrentFolder(), SigmaHistos, and theDbe.
00053 { 00054 00055 metname = "trackResidualsTest"; 00056 theDbe->setCurrentFolder("Muons/Tests/trackResidualsTest"); 00057 00058 LogTrace(metname) << "[MuonTrackResidualsTest] Parameters initialization"<<endl; 00059 00060 00061 string histName, MeanHistoName, SigmaHistoName, MeanHistoTitle, SigmaHistoTitle; 00062 vector<string> type; 00063 type.push_back("eta"); 00064 type.push_back("theta"); 00065 type.push_back("phi"); 00066 00067 00068 for(unsigned int c=0; c<type.size(); c++){ 00069 00070 MeanHistoName = "MeanTest_" + type[c]; 00071 SigmaHistoName = "SigmaTest_" + type[c]; 00072 00073 MeanHistoTitle = "Mean of the #" + type[c] + " residuals distribution"; 00074 SigmaHistoTitle = "Sigma of the #" + type[c] + " residuals distribution"; 00075 00076 histName = "Res_GlbSta_"+type[c]; 00077 histoNames[type[c]].push_back(histName); 00078 histName = "Res_TkGlb_"+type[c]; 00079 histoNames[type[c]].push_back(histName); 00080 histName = "Res_TkSta_"+type[c]; 00081 histoNames[type[c]].push_back(histName); 00082 00083 00084 MeanHistos[type[c]] = theDbe->book1D(MeanHistoName.c_str(),MeanHistoTitle.c_str(),3,0.5,3.5); 00085 (MeanHistos[type[c]])->setBinLabel(1,"Res_StaGlb",1); 00086 (MeanHistos[type[c]])->setBinLabel(2,"Res_TkGlb",1); 00087 (MeanHistos[type[c]])->setBinLabel(3,"Res_TkSta",1); 00088 00089 00090 SigmaHistos[type[c]] = theDbe->book1D(SigmaHistoName.c_str(),SigmaHistoTitle.c_str(),3,0.5,3.5); 00091 (SigmaHistos[type[c]])->setBinLabel(1,"Res_StaGlb",1); 00092 (SigmaHistos[type[c]])->setBinLabel(2,"Res_TkGlb",1); 00093 (SigmaHistos[type[c]])->setBinLabel(3,"Res_TkSta",1); 00094 00095 } 00096 00097 nevents = 0; 00098 00099 }
void MuonTrackResidualsTest::beginLuminosityBlock | ( | edm::LuminosityBlock const & | lumiSeg, | |
edm::EventSetup const & | context | |||
) | [protected, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 102 of file MuonTrackResidualsTest.cc.
References LogTrace, metname, run, and edm::LuminosityBlock::run().
00102 { 00103 00104 LogTrace(metname)<<"[MuonTrackResidualsTest]: Begin of LS transition"; 00105 00106 // Get the run number 00107 run = lumiSeg.run(); 00108 00109 }
Endjob.
Reimplemented from edm::EDAnalyzer.
Definition at line 196 of file MuonTrackResidualsTest.cc.
References LogTrace, metname, DQMStore::rmdir(), and theDbe.
00196 { 00197 00198 LogTrace(metname)<< "[MuonTrackResidualsTest] endjob called!"; 00199 theDbe->rmdir("Muons/Tests/trackResidualsTest"); 00200 00201 }
void MuonTrackResidualsTest::endLuminosityBlock | ( | edm::LuminosityBlock const & | lumiSeg, | |
edm::EventSetup const & | c | |||
) | [protected, virtual] |
DQM Client Diagnostic.
Reimplemented from edm::EDAnalyzer.
Definition at line 121 of file MuonTrackResidualsTest.cc.
References lat::endl(), DQMStore::get(), QReport::getBadChannels(), QReport::getMessage(), MonitorElement::getQReport(), QReport::getStatus(), edm::ParameterSet::getUntrackedParameter(), histo, histoNames, edm::LuminosityBlock::id(), LogTrace, edm::LuminosityBlockID::luminosityBlock(), mean(), MeanHistos, metname, nLumiSegs, parameters, path(), prescaleFactor, SigmaHistos, and theDbe.
00121 { 00122 00123 LogTrace(metname)<<"[MuonTrackResidualsTest]: End of LS transition, performing the DQM client operation"<<endl; 00124 00125 // counts number of lumiSegs 00126 nLumiSegs = lumiSeg.id().luminosityBlock(); 00127 00128 // prescale factor 00129 if ( nLumiSegs%prescaleFactor != 0 ) return; 00130 00131 00132 for(map<string, vector<string> > ::const_iterator histo = histoNames.begin(); 00133 histo != histoNames.end(); 00134 histo++) { 00135 00136 for (unsigned int type=0; type< (*histo).second.size(); type++){ 00137 00138 string path = "Muons/MuonRecoAnalyzer/" + (*histo).second[type]; 00139 MonitorElement * res_histo = theDbe->get(path); 00140 if (res_histo) { 00141 00142 // gaussian test 00143 string GaussianCriterionName = 00144 parameters.getUntrackedParameter<string>("resDistributionTestName", 00145 "ResidualsDistributionGaussianTest"); 00146 const QReport * GaussianReport = res_histo->getQReport(GaussianCriterionName); 00147 if(GaussianReport){ 00148 LogTrace(metname) << "-------- histo : "<<(*histo).second[type]<<" "<<GaussianReport->getMessage()<<" ------- "<<GaussianReport->getStatus(); 00149 } 00150 int BinNumber = type+1; 00151 float mean = (*res_histo).getMean(1); 00152 float sigma = (*res_histo).getRMS(1); 00153 MeanHistos.find((*histo).first)->second->setBinContent(BinNumber, mean); 00154 SigmaHistos.find((*histo).first)->second->setBinContent(BinNumber, sigma); 00155 } 00156 } 00157 } 00158 00159 00160 // Mean test 00161 string MeanCriterionName = parameters.getUntrackedParameter<string>("meanTestName","ResidualsMeanInRange"); 00162 for(map<string, MonitorElement*>::const_iterator hMean = MeanHistos.begin(); 00163 hMean != MeanHistos.end(); 00164 hMean++) { 00165 const QReport * theMeanQReport = (*hMean).second->getQReport(MeanCriterionName); 00166 if(theMeanQReport) { 00167 vector<dqm::me_util::Channel> badChannels = theMeanQReport->getBadChannels(); 00168 for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); 00169 channel != badChannels.end(); channel++) { 00170 LogTrace(metname)<< "type:"<<(*hMean).first<<" Bad mean channels: "<<(*channel).getBin()<<" Contents : "<<(*channel).getContents()<<endl; 00171 } 00172 LogTrace(metname)<< "-------- type: "<<(*hMean).first<<" "<<theMeanQReport->getMessage()<<" ------- "<<theMeanQReport->getStatus()<<endl; 00173 } 00174 } 00175 00176 // Sigma test 00177 string SigmaCriterionName = parameters.getUntrackedParameter<string>("sigmaTestName","ResidualsSigmaInRange"); 00178 for(map<string, MonitorElement*>::const_iterator hSigma = SigmaHistos.begin(); 00179 hSigma != SigmaHistos.end(); 00180 hSigma++) { 00181 const QReport * theSigmaQReport = (*hSigma).second->getQReport(SigmaCriterionName); 00182 if(theSigmaQReport) { 00183 vector<dqm::me_util::Channel> badChannels = theSigmaQReport->getBadChannels(); 00184 for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); 00185 channel != badChannels.end(); channel++) { 00186 LogTrace(metname)<< "type:"<<(*hSigma).first<<" Bad sigma channels: "<<(*channel).getBin()<<" Contents : "<<(*channel).getContents()<<endl; 00187 } 00188 LogTrace(metname) << "-------- type: "<<(*hSigma).first<<" "<<theSigmaQReport->getMessage()<<" ------- "<<theSigmaQReport->getStatus()<<endl; 00189 } 00190 } 00191 }
std::map< std::string, std::vector<std::string> > MuonTrackResidualsTest::histoNames [private] |
Definition at line 80 of file MuonTrackResidualsTest.h.
Referenced by beginJob(), and endLuminosityBlock().
std::map< std::string, MonitorElement* > MuonTrackResidualsTest::MeanHistos [private] |
Definition at line 83 of file MuonTrackResidualsTest.h.
Referenced by beginJob(), and endLuminosityBlock().
std::string MuonTrackResidualsTest::metname [private] |
Definition at line 74 of file MuonTrackResidualsTest.h.
Referenced by analyze(), beginJob(), beginLuminosityBlock(), endJob(), endLuminosityBlock(), and ~MuonTrackResidualsTest().
int MuonTrackResidualsTest::nevents [private] |
Definition at line 69 of file MuonTrackResidualsTest.h.
Referenced by analyze(), beginJob(), and ~MuonTrackResidualsTest().
unsigned int MuonTrackResidualsTest::nLumiSegs [private] |
Definition at line 77 of file MuonTrackResidualsTest.h.
Referenced by endLuminosityBlock(), and MuonTrackResidualsTest().
int MuonTrackResidualsTest::prescaleFactor [private] |
Definition at line 71 of file MuonTrackResidualsTest.h.
Referenced by endLuminosityBlock(), and MuonTrackResidualsTest().
int MuonTrackResidualsTest::run [private] |
std::map< std::string ,MonitorElement* > MuonTrackResidualsTest::SigmaHistos [private] |
Definition at line 84 of file MuonTrackResidualsTest.h.
Referenced by beginJob(), and endLuminosityBlock().
DQMStore* MuonTrackResidualsTest::theDbe [private] |
Definition at line 76 of file MuonTrackResidualsTest.h.
Referenced by beginJob(), endJob(), endLuminosityBlock(), and MuonTrackResidualsTest().