00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023
00024
00025
00026 #include "FWCore/Framework/interface/EDAnalyzer.h"
00027
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030 #include "FWCore/Framework/interface/ESHandle.h"
00031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00032
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034 #include "DataFormats/HcalCalibObjects/interface/HOCalibVariables.h"
00035
00036 #include "DQMServices/Core/interface/DQMStore.h"
00037 #include "DQMServices/Core/interface/MonitorElement.h"
00038 #include "FWCore/ServiceRegistry/interface/Service.h"
00039 #include "DQMOffline/CalibCalo/src/DQMHOAlCaRecoStream.h"
00040
00041 #include <string>
00042
00043
00044
00045
00046
00047
00048
00049
00050 using namespace std;
00051 using namespace edm;
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066 DQMHOAlCaRecoStream::DQMHOAlCaRecoStream(const edm::ParameterSet& iConfig) :
00067 hoCalibVariableCollectionTag(iConfig.getParameter<edm::InputTag>("hoCalibVariableCollectionTag")) {
00068
00069
00070
00071 theRootFileName = iConfig.getUntrackedParameter<string>("RootFileName","tmp.root");
00072 folderName_ = iConfig.getUntrackedParameter<string>("folderName");
00073 m_sigmaValue = iConfig.getUntrackedParameter<double>("sigmaval",0.2);
00074 m_lowRadPosInMuch = iConfig.getUntrackedParameter<double>("lowradposinmuch",400.0);
00075 m_highRadPosInMuch = iConfig.getUntrackedParameter<double>("highradposinmuch",480.0);
00076 m_lowEdge = iConfig.getUntrackedParameter<double>("lowedge",-2.0);
00077 m_highEdge = iConfig.getUntrackedParameter<double>("highedge",6.0);
00078 m_nbins = iConfig.getUntrackedParameter<int>("nbins",40);
00079 saveToFile_ = iConfig.getUntrackedParameter<bool>("saveToFile",false);
00080 }
00081
00082
00083 DQMHOAlCaRecoStream::~DQMHOAlCaRecoStream()
00084 {
00085
00086
00087
00088
00089 }
00090
00091
00092
00093
00094
00095
00096
00097 void
00098 DQMHOAlCaRecoStream::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00099 {
00100 using namespace edm;
00101
00102 Nevents++;
00103
00104 edm::Handle<HOCalibVariableCollection>HOCalib;
00105 bool isCosMu = true;
00106
00107 iEvent.getByLabel(hoCalibVariableCollectionTag, HOCalib);
00108
00109 if(!HOCalib.isValid()){
00110 LogDebug("") << "DQMHOAlCaRecoStream:: Error! can't get HOCalib product!" << std::endl;
00111 return ;
00112 }
00113
00114
00115 if (isCosMu) {
00116 hMuonMultipl->Fill((*HOCalib).size(),1.);
00117 if ((*HOCalib).size() >0 ) {
00118 for (HOCalibVariableCollection::const_iterator hoC=(*HOCalib).begin(); hoC!=(*HOCalib).end(); hoC++){
00119
00120 float okt = 2.;
00121 double okx = std::pow((*hoC).trkvx,okt) + std::pow((*hoC).trkvy,okt);
00123 double dr=std::pow( okx, 0.5);
00124 if (dr <m_lowRadPosInMuch || dr >m_highRadPosInMuch) continue;
00125
00126 if ((*hoC).isect <0) continue;
00127 if (fabs((*hoC).trkth-acos(-1.)/2)<0.000001) continue;
00128 int ieta = int((std::abs((*hoC).isect)%10000)/100.)-30;
00129
00130 if (std::abs(ieta)>=16) continue;
00131
00132 Nmuons++;
00133
00134
00135 hMuonMom->Fill((*hoC).trkmm, 1.0);
00136 hMuonEta->Fill(-log(tan((*hoC).trkth/2.0)), 1.0);
00137 hMuonPhi->Fill((*hoC).trkph, 1.0);
00138 hDirCosine->Fill((*hoC).hoang, 1.0);
00139 hHOTime->Fill((*hoC).htime, 1.0);
00140
00141 double energy = (*hoC).hosig[4];
00142 double pedval = (*hoC).hocro;
00143 int iring = 0;
00144 if (ieta >=-15 && ieta <=-11) {iring = -2;}
00145 else if (ieta >=-10 && ieta <=-5) {iring = -1;}
00146 else if (ieta >= 5 && ieta <= 10) {iring = 1;}
00147 else if (ieta >= 11 && ieta <= 15) {iring = 2;}
00148
00149 hSigRing[iring+2]->Fill(energy,1.0);
00150 hPedRing[iring+2]->Fill(pedval,1.0);
00151
00152 for (int k=0; k<9; k++) {
00153 hSignal3x3[k]->Fill((*hoC).hosig[k]);
00154 }
00155 }
00156 }
00157 }
00158 }
00159
00160
00161
00162 void
00163 DQMHOAlCaRecoStream::beginJob()
00164 {
00165 dbe_ = edm::Service<DQMStore>().operator->();
00166 dbe_->setCurrentFolder(folderName_);
00167
00168 char title[200];
00169 char name[200];
00170
00171 hMuonMom = dbe_->book1D("hMuonMom", "Muon momentum (GeV)", 50, -100, 100);
00172 hMuonMom ->setAxisTitle("Muon momentum (GeV)",1);
00173
00174 hMuonEta = dbe_->book1D("hMuonEta", "Pseudo-rapidity of muon", 50, -1.5, 1.5);
00175 hMuonEta ->setAxisTitle("Pseudo-rapidity of muon",1);
00176
00177 hMuonPhi = dbe_->book1D("hMuonPhi", "Azimuthal angle of muon", 24, -acos(-1), acos(-1));
00178 hMuonPhi ->setAxisTitle("Azimuthal angle of muon",1);
00179
00180 hMuonMultipl = dbe_->book1D("hMuonMultipl", "Muon Multiplicity", 10, 0.5, 10.5);
00181 hMuonMultipl ->setAxisTitle("Muon Multiplicity",1);
00182
00183 hDirCosine = dbe_->book1D("hDirCosine", "Direction Cosine of muon at HO tower", 50, -1., 1.);
00184 hDirCosine ->setAxisTitle("Direction Cosine of muon at HO tower",1);
00185
00186 hHOTime = dbe_->book1D("hHOTime", "HO time distribution", 60, -20, 100.);
00187 hHOTime ->setAxisTitle("HO time distribution", 1);
00188
00189 for (int i=0; i<5; i++) {
00190 sprintf(name, "hSigRing_%i", i-2);
00191 sprintf(title, "HO signal in Ring_%i", i-2);
00192 hSigRing[i] = dbe_->book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
00193 hSigRing[i]->setAxisTitle(title,1);
00194
00195 sprintf(name, "hPedRing_%i", i-2);
00196 sprintf(title, "HO Pedestal in Ring_%i", i-2);
00197 hPedRing[i] = dbe_->book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
00198 hPedRing[i]->setAxisTitle(title,1);
00199 }
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219 for (int i=-1; i<=1; i++) {
00220 for (int j=-1; j<=1; j++) {
00221 int k = 3*(i+1)+j+1;
00222
00223 sprintf(title, "hSignal3x3_deta%i_dphi%i", i, j);
00224 hSignal3x3[k] = dbe_->book1D(title, title, m_nbins, m_lowEdge, m_highEdge);
00225 hSignal3x3[k]->setAxisTitle(title,1);
00226 }
00227 }
00228
00229 Nevents = 0;
00230 Nmuons = 0;
00231
00232 }
00233
00234
00235 void
00236 DQMHOAlCaRecoStream::endJob() {
00237 if (saveToFile_) {
00238
00239
00240 double scale = 1./max(1,Nmuons);
00241 hMuonMom->getTH1F()->Scale(scale);
00242 hMuonEta->getTH1F()->Scale(scale);
00243 hMuonPhi->getTH1F()->Scale(scale);
00244 hDirCosine->getTH1F()->Scale(scale);
00245 hHOTime->getTH1F()->Scale(scale);
00246
00247
00248 for (int k=0; k<5; k++) {
00249 hSigRing[k]->getTH1F()->Scale(scale);
00250 hPedRing[k]->getTH1F()->Scale(scale);
00251 }
00252
00253 for (int k=0; k<9; k++) {
00254 hSignal3x3[k]->getTH1F()->Scale(scale);
00255 }
00256
00257 dbe_->save(theRootFileName);
00258 }
00259
00260 }
00261