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 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
24 
26 
31 
33 
38 
39 #include <string>
40 
41 
42 
43 
44 
45 //
46 // class decleration
47 //
48 using namespace std;
49 using namespace edm;
50 
51 
52 
53 //
54 // constants, enums and typedefs
55 //
56 
57 //
58 // static data member definitions
59 //
60 
61 //
62 // constructors and destructor
63 //
65  hoCalibVariableCollectionTag(consumes<HOCalibVariableCollection>(iConfig.getParameter<edm::InputTag>("hoCalibVariableCollectionTag"))) {
66 
67  //now do what ever initialization is needed
68 
69  theRootFileName = iConfig.getUntrackedParameter<string>("RootFileName","tmp.root");
70  folderName_ = iConfig.getUntrackedParameter<string>("folderName");
71  m_sigmaValue = iConfig.getUntrackedParameter<double>("sigmaval",0.2);
72  m_lowRadPosInMuch = iConfig.getUntrackedParameter<double>("lowradposinmuch",400.0);
73  m_highRadPosInMuch = iConfig.getUntrackedParameter<double>("highradposinmuch",480.0);
74  m_lowEdge = iConfig.getUntrackedParameter<double>("lowedge",-2.0);
75  m_highEdge = iConfig.getUntrackedParameter<double>("highedge",6.0);
76  m_nbins = iConfig.getUntrackedParameter<int>("nbins",40);
77  saveToFile_ = iConfig.getUntrackedParameter<bool>("saveToFile",false);
78 }
79 
80 
82 {
83 
84  // do anything here that needs to be done at desctruction time
85  // (e.g. close files, deallocate resources etc.)
86 
87 }
88 
89 
90 //
91 // member functions
92 //
93 
94 // ------------ method called to for each event ------------
95 void
97 {
98  using namespace edm;
99 
100  Nevents++;
101 
103  bool isCosMu = true;
104 
105  iEvent.getByToken(hoCalibVariableCollectionTag, HOCalib);
106 
107  if(!HOCalib.isValid()){
108  LogDebug("") << "DQMHOAlCaRecoStream:: Error! can't get HOCalib product!" << std::endl;
109  return ;
110  }
111 
112 
113  if (isCosMu) {
114  hMuonMultipl->Fill((*HOCalib).size(),1.);
115  if ((*HOCalib).size() >0 ) {
116  for (HOCalibVariableCollection::const_iterator hoC=(*HOCalib).begin(); hoC!=(*HOCalib).end(); hoC++){
117 // OK!!!!
118  float okt = 2.;
119  double okx = std::pow((*hoC).trkvx,okt) + std::pow((*hoC).trkvy,okt);
121  double dr=std::pow( okx, 0.5);
122  if (dr <m_lowRadPosInMuch || dr >m_highRadPosInMuch) continue;
123 
124  if ((*hoC).isect <0) continue;
125  if (fabs((*hoC).trkth-acos(-1.)/2)<0.000001) continue;
126  int ieta = int((std::abs((*hoC).isect)%10000)/100.)-30;
127 
128  if (std::abs(ieta)>=16) continue;
129 
130  Nmuons++;
131 
132 
133  hMuonMom->Fill((*hoC).trkmm, 1.0);
134  hMuonEta->Fill(-log(tan((*hoC).trkth/2.0)), 1.0);
135  hMuonPhi->Fill((*hoC).trkph, 1.0);
136  hDirCosine->Fill((*hoC).hoang, 1.0);
137  hHOTime->Fill((*hoC).htime, 1.0);
138 
139  double energy = (*hoC).hosig[4];
140  double pedval = (*hoC).hocro;
141  int iring = 0;
142  if (ieta >=-15 && ieta <=-11) {iring = -2;}
143  else if (ieta >=-10 && ieta <=-5) {iring = -1;}
144  else if (ieta >= 5 && ieta <= 10) {iring = 1;}
145  else if (ieta >= 11 && ieta <= 15) {iring = 2;}
146 
147  hSigRing[iring+2]->Fill(energy,1.0);
148  hPedRing[iring+2]->Fill(pedval,1.0);
149 
150  for (int k=0; k<9; k++) {
151  hSignal3x3[k]->Fill((*hoC).hosig[k]);
152  }
153  } //for (HOCalibVariableCollection::const_iterator hoC=(*HOCalib).begin()
154  } // if ((*HOCalib).size() >0 ) {
155  } // if (isCosMu) {
156 }
157 
158 
159 // ------------ method called once each job just before starting event loop ------------
160 void
162 {
163  ibooker.setCurrentFolder(folderName_);
164 
165  char title[200];
166  char name[200];
167 
168  hMuonMom = ibooker.book1D("hMuonMom", "Muon momentum (GeV)", 50, -100, 100);
169  hMuonMom ->setAxisTitle("Muon momentum (GeV)",1);
170 
171  hMuonEta = ibooker.book1D("hMuonEta", "Pseudo-rapidity of muon", 50, -1.5, 1.5);
172  hMuonEta ->setAxisTitle("Pseudo-rapidity of muon",1);
173 
174  hMuonPhi = ibooker.book1D("hMuonPhi", "Azimuthal angle of muon", 24, -acos(-1), acos(-1));
175  hMuonPhi ->setAxisTitle("Azimuthal angle of muon",1);
176 
177  hMuonMultipl = ibooker.book1D("hMuonMultipl", "Muon Multiplicity", 10, 0.5, 10.5);
178  hMuonMultipl ->setAxisTitle("Muon Multiplicity",1);
179 
180  hDirCosine = ibooker.book1D("hDirCosine", "Direction Cosine of muon at HO tower", 50, -1., 1.);
181  hDirCosine ->setAxisTitle("Direction Cosine of muon at HO tower",1);
182 
183  hHOTime = ibooker.book1D("hHOTime", "HO time distribution", 60, -20, 100.);
184  hHOTime ->setAxisTitle("HO time distribution", 1);
185 
186  for (int i=0; i<5; i++) {
187  sprintf(name, "hSigRing_%i", i-2);
188  sprintf(title, "HO signal in Ring_%i", i-2);
189  hSigRing[i] = ibooker.book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
190  hSigRing[i]->setAxisTitle(title,1);
191 
192  sprintf(name, "hPedRing_%i", i-2);
193  sprintf(title, "HO Pedestal in Ring_%i", i-2);
194  hPedRing[i] = ibooker.book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
195  hPedRing[i]->setAxisTitle(title,1);
196  }
197 
198  // hSigRingm1 = ibooker.book1D("hSigRingm1", "HO signal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
199  // hSigRingm1->setAxisTitle("HO signal in Ring-1",1);
200 
201  // hSigRing00 = ibooker.book1D("hSigRing00", "HO signal in Ring_0", m_nbins, m_lowEdge, m_highEdge);
202  // hSigRing00->setAxisTitle("HO signal in Ring_0",1);
203 
204  // hSigRingp1 = ibooker.book1D("hSigRingp1", "HO signal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
205  // hSigRingp1->setAxisTitle("HO signal in Ring+1",1);
206 
207  // hSigRingp2 = ibooker.book1D("hSigRingp2", "HO signal in Ring-2", m_nbins, m_lowEdge, m_highEdge);
208  // hSigRingp2->setAxisTitle("HO signal in Ring+2",1);
209 
210  // hPedRingm2 = ibooker.book1D("hPedRingm2", "HO pedestal in Ring-2", m_nbins, m_lowEdge, m_highEdge);
211  // hPedRingm1 = ibooker.book1D("hPedRingm1", "HO pedestal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
212  // hPedRing00 = ibooker.book1D("hPedRing00", "HO pedestal in Ring_0", m_nbins, m_lowEdge, m_highEdge);
213  // hPedRingp1 = ibooker.book1D("hPedRingp1", "HO pedestal in Ring-1", m_nbins, m_lowEdge, m_highEdge);
214  // hPedRingp2 = ibooker.book1D("hPedRingp2", "HO pedestal in Ring-2", m_nbins, m_lowEdge, m_highEdge);
215 
216  for (int i=-1; i<=1; i++) {
217  for (int j=-1; j<=1; j++) {
218  int k = 3*(i+1)+j+1;
219 
220  sprintf(title, "hSignal3x3_deta%i_dphi%i", i, j);
221  hSignal3x3[k] = ibooker.book1D(title, title, m_nbins, m_lowEdge, m_highEdge);
222  hSignal3x3[k]->setAxisTitle(title,1);
223  }
224  }
225 
226  Nevents = 0;
227  Nmuons = 0;
228 
229 }
230 
231 
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< HOCalibVariableCollection > hoCalibVariableCollectionTag
MonitorElement * hMuonMultipl
MonitorElement * hPedRing[5]
MonitorElement * hMuonEta
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
DQMHOAlCaRecoStream(const edm::ParameterSet &)
MonitorElement * hSigRing[5]
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
MonitorElement * hHOTime
MonitorElement * hSignal3x3[9]
MonitorElement * hMuonMom
MonitorElement * hMuonPhi
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
std::vector< HOCalibVariables > HOCalibVariableCollection
collection of HOcalibration variabale
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
return(e1-e2)*(e1-e2)+dp *dp
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
Definition: Run.h:42