CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DQMHOAlCaRecoStream.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: DQMHOAlCaRecoStream
4 // Class: DQMHOAlCaRecoStream
5 //
13 //
14 // Original Author: Gobinda Majumder
15 // Created: Mon Mar 2 12:33:08 CET 2009
16 // $Id: DQMHOAlCaRecoStream.cc,v 1.11 2012/09/26 21:08:48 wdd Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
25 
27 
32 
35 
40 
41 #include <string>
42 
43 
44 
45 
46 
47 //
48 // class decleration
49 //
50 using namespace std;
51 using namespace edm;
52 
53 
54 
55 //
56 // constants, enums and typedefs
57 //
58 
59 //
60 // static data member definitions
61 //
62 
63 //
64 // constructors and destructor
65 //
67  hoCalibVariableCollectionTag(iConfig.getParameter<edm::InputTag>("hoCalibVariableCollectionTag")) {
68 
69  //now do what ever initialization is needed
70 
71  theRootFileName = iConfig.getUntrackedParameter<string>("RootFileName","tmp.root");
72  folderName_ = iConfig.getUntrackedParameter<string>("folderName");
73  m_sigmaValue = iConfig.getUntrackedParameter<double>("sigmaval",0.2);
74  m_lowRadPosInMuch = iConfig.getUntrackedParameter<double>("lowradposinmuch",400.0);
75  m_highRadPosInMuch = iConfig.getUntrackedParameter<double>("highradposinmuch",480.0);
76  m_lowEdge = iConfig.getUntrackedParameter<double>("lowedge",-2.0);
77  m_highEdge = iConfig.getUntrackedParameter<double>("highedge",6.0);
78  m_nbins = iConfig.getUntrackedParameter<int>("nbins",40);
79  saveToFile_ = iConfig.getUntrackedParameter<bool>("saveToFile",false);
80 }
81 
82 
84 {
85 
86  // do anything here that needs to be done at desctruction time
87  // (e.g. close files, deallocate resources etc.)
88 
89 }
90 
91 
92 //
93 // member functions
94 //
95 
96 // ------------ method called to for each event ------------
97 void
99 {
100  using namespace edm;
101 
102  Nevents++;
103 
105  bool isCosMu = true;
106 
107  iEvent.getByLabel(hoCalibVariableCollectionTag, HOCalib);
108 
109  if(!HOCalib.isValid()){
110  LogDebug("") << "DQMHOAlCaRecoStream:: Error! can't get HOCalib product!" << std::endl;
111  return ;
112  }
113 
114 
115  if (isCosMu) {
116  hMuonMultipl->Fill((*HOCalib).size(),1.);
117  if ((*HOCalib).size() >0 ) {
118  for (HOCalibVariableCollection::const_iterator hoC=(*HOCalib).begin(); hoC!=(*HOCalib).end(); hoC++){
119 // OK!!!!
120  float okt = 2.;
121  double okx = std::pow((*hoC).trkvx,okt) + std::pow((*hoC).trkvy,okt);
123  double dr=std::pow( okx, 0.5);
124  if (dr <m_lowRadPosInMuch || dr >m_highRadPosInMuch) continue;
125 
126  if ((*hoC).isect <0) continue;
127  if (fabs((*hoC).trkth-acos(-1.)/2)<0.000001) continue;
128  int ieta = int((std::abs((*hoC).isect)%10000)/100.)-30;
129 
130  if (std::abs(ieta)>=16) continue;
131 
132  Nmuons++;
133 
134 
135  hMuonMom->Fill((*hoC).trkmm, 1.0);
136  hMuonEta->Fill(-log(tan((*hoC).trkth/2.0)), 1.0);
137  hMuonPhi->Fill((*hoC).trkph, 1.0);
138  hDirCosine->Fill((*hoC).hoang, 1.0);
139  hHOTime->Fill((*hoC).htime, 1.0);
140 
141  double energy = (*hoC).hosig[4];
142  double pedval = (*hoC).hocro;
143  int iring = 0;
144  if (ieta >=-15 && ieta <=-11) {iring = -2;}
145  else if (ieta >=-10 && ieta <=-5) {iring = -1;}
146  else if (ieta >= 5 && ieta <= 10) {iring = 1;}
147  else if (ieta >= 11 && ieta <= 15) {iring = 2;}
148 
149  hSigRing[iring+2]->Fill(energy,1.0);
150  hPedRing[iring+2]->Fill(pedval,1.0);
151 
152  for (int k=0; k<9; k++) {
153  hSignal3x3[k]->Fill((*hoC).hosig[k]);
154  }
155  } //for (HOCalibVariableCollection::const_iterator hoC=(*HOCalib).begin()
156  } // if ((*HOCalib).size() >0 ) {
157  } // if (isCosMu) {
158 }
159 
160 
161 // ------------ method called once each job just before starting event loop ------------
162 void
164 {
167 
168  char title[200];
169  char name[200];
170 
171  hMuonMom = dbe_->book1D("hMuonMom", "Muon momentum (GeV)", 50, -100, 100);
172  hMuonMom ->setAxisTitle("Muon momentum (GeV)",1);
173 
174  hMuonEta = dbe_->book1D("hMuonEta", "Pseudo-rapidity of muon", 50, -1.5, 1.5);
175  hMuonEta ->setAxisTitle("Pseudo-rapidity of muon",1);
176 
177  hMuonPhi = dbe_->book1D("hMuonPhi", "Azimuthal angle of muon", 24, -acos(-1), acos(-1));
178  hMuonPhi ->setAxisTitle("Azimuthal angle of muon",1);
179 
180  hMuonMultipl = dbe_->book1D("hMuonMultipl", "Muon Multiplicity", 10, 0.5, 10.5);
181  hMuonMultipl ->setAxisTitle("Muon Multiplicity",1);
182 
183  hDirCosine = dbe_->book1D("hDirCosine", "Direction Cosine of muon at HO tower", 50, -1., 1.);
184  hDirCosine ->setAxisTitle("Direction Cosine of muon at HO tower",1);
185 
186  hHOTime = dbe_->book1D("hHOTime", "HO time distribution", 60, -20, 100.);
187  hHOTime ->setAxisTitle("HO time distribution", 1);
188 
189  for (int i=0; i<5; i++) {
190  sprintf(name, "hSigRing_%i", i-2);
191  sprintf(title, "HO signal in Ring_%i", i-2);
192  hSigRing[i] = dbe_->book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
193  hSigRing[i]->setAxisTitle(title,1);
194 
195  sprintf(name, "hPedRing_%i", i-2);
196  sprintf(title, "HO Pedestal in Ring_%i", i-2);
197  hPedRing[i] = dbe_->book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
198  hPedRing[i]->setAxisTitle(title,1);
199  }
200 
201  // hSigRingm1 = dbe_->book1D("hSigRingm1", "HO signal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
202  // hSigRingm1->setAxisTitle("HO signal in Ring-1",1);
203 
204  // hSigRing00 = dbe_->book1D("hSigRing00", "HO signal in Ring_0", m_nbins, m_lowEdge, m_highEdge);
205  // hSigRing00->setAxisTitle("HO signal in Ring_0",1);
206 
207  // hSigRingp1 = dbe_->book1D("hSigRingp1", "HO signal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
208  // hSigRingp1->setAxisTitle("HO signal in Ring+1",1);
209 
210  // hSigRingp2 = dbe_->book1D("hSigRingp2", "HO signal in Ring-2", m_nbins, m_lowEdge, m_highEdge);
211  // hSigRingp2->setAxisTitle("HO signal in Ring+2",1);
212 
213  // hPedRingm2 = dbe_->book1D("hPedRingm2", "HO pedestal in Ring-2", m_nbins, m_lowEdge, m_highEdge);
214  // hPedRingm1 = dbe_->book1D("hPedRingm1", "HO pedestal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
215  // hPedRing00 = dbe_->book1D("hPedRing00", "HO pedestal in Ring_0", m_nbins, m_lowEdge, m_highEdge);
216  // hPedRingp1 = dbe_->book1D("hPedRingp1", "HO pedestal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
217  // hPedRingp2 = dbe_->book1D("hPedRingp2", "HO pedestal in Ring-2", m_nbins, m_lowEdge, m_highEdge);
218 
219  for (int i=-1; i<=1; i++) {
220  for (int j=-1; j<=1; j++) {
221  int k = 3*(i+1)+j+1;
222 
223  sprintf(title, "hSignal3x3_deta%i_dphi%i", i, j);
224  hSignal3x3[k] = dbe_->book1D(title, title, m_nbins, m_lowEdge, m_highEdge);
225  hSignal3x3[k]->setAxisTitle(title,1);
226  }
227  }
228 
229  Nevents = 0;
230  Nmuons = 0;
231 
232 }
233 
234 // ------------ method called once each job just after ending the event loop ------------
235 void
237  if (saveToFile_) {
238 
239  // double scale = 1./max(1,Nevents);
240  double scale = 1./max(1,Nmuons);
241  hMuonMom->getTH1F()->Scale(scale);
242  hMuonEta->getTH1F()->Scale(scale);
243  hMuonPhi->getTH1F()->Scale(scale);
244  hDirCosine->getTH1F()->Scale(scale);
245  hHOTime->getTH1F()->Scale(scale);
246 
247  // scale = 1./max(1,Nmuons);
248  for (int k=0; k<5; k++) {
249  hSigRing[k]->getTH1F()->Scale(scale);
250  hPedRing[k]->getTH1F()->Scale(scale);
251  }
252 
253  for (int k=0; k<9; k++) {
254  hSignal3x3[k]->getTH1F()->Scale(scale);
255  }
256 
258  }
259 
260 }
261 
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
MonitorElement * hMuonMultipl
MonitorElement * hPedRing[5]
MonitorElement * hMuonEta
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:722
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2118
DQMHOAlCaRecoStream(const edm::ParameterSet &)
#define abs(x)
Definition: mlp_lapack.h:159
MonitorElement * hSigRing[5]
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:243
const T & max(const T &a, const T &b)
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
int j
Definition: DBlmapReader.cc:9
edm::InputTag hoCalibVariableCollectionTag
MonitorElement * hHOTime
MonitorElement * hSignal3x3[9]
MonitorElement * hMuonMom
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
MonitorElement * hMuonPhi
int k[5][pyjets_maxn]
TH1F * getTH1F(void) const
virtual void analyze(const edm::Event &, const edm::EventSetup &)
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
MonitorElement * hDirCosine
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434