00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #include <DQMOffline/Muon/src/MuonTrackResidualsTest.h>
00012
00013
00014 #include <FWCore/Framework/interface/Event.h>
00015 #include "DataFormats/Common/interface/Handle.h"
00016 #include <FWCore/Framework/interface/ESHandle.h>
00017 #include <FWCore/Framework/interface/MakerMacros.h>
00018 #include <FWCore/Framework/interface/EventSetup.h>
00019 #include <FWCore/ParameterSet/interface/ParameterSet.h>
00020
00021
00022 #include "DQMServices/Core/interface/DQMStore.h"
00023 #include "DQMServices/Core/interface/MonitorElement.h"
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025
00026 #include <iostream>
00027 #include <stdio.h>
00028 #include <string>
00029 #include <math.h>
00030 #include "TF1.h"
00031
00032 using namespace edm;
00033 using namespace std;
00034
00035
00036 MuonTrackResidualsTest::MuonTrackResidualsTest(const edm::ParameterSet& ps){
00037 parameters = ps;
00038
00039 theDbe = edm::Service<DQMStore>().operator->();
00040
00041 prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);
00042
00043 }
00044
00045
00046 MuonTrackResidualsTest::~MuonTrackResidualsTest(){
00047
00048 LogTrace(metname) << "DTResolutionTest: analyzed " << nevents << " events";
00049
00050 }
00051
00052
00053 void MuonTrackResidualsTest::beginJob(const edm::EventSetup& context){
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 }
00100
00101
00102 void MuonTrackResidualsTest::beginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
00103
00104 LogTrace(metname)<<"[MuonTrackResidualsTest]: Begin of LS transition";
00105
00106
00107 run = lumiSeg.run();
00108
00109 }
00110
00111
00112 void MuonTrackResidualsTest::analyze(const edm::Event& e, const edm::EventSetup& context){
00113
00114 nevents++;
00115 LogTrace(metname)<< "[MuonTrackResidualsTest]: "<<nevents<<" events";
00116
00117 }
00118
00119
00120
00121 void MuonTrackResidualsTest::endLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
00122
00123 LogTrace(metname)<<"[MuonTrackResidualsTest]: End of LS transition, performing the DQM client operation"<<endl;
00124
00125
00126 nLumiSegs = lumiSeg.id().luminosityBlock();
00127
00128
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
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
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
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 }
00192
00193
00194
00195
00196 void MuonTrackResidualsTest::endJob(){
00197
00198 LogTrace(metname)<< "[MuonTrackResidualsTest] endjob called!";
00199 theDbe->rmdir("Muons/Tests/trackResidualsTest");
00200
00201 }
00202
00203
00204