CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
HcalPedestalsAnalysis Class Reference

#include <HcalPedestalsAnalysis.h>

Inheritance diagram for HcalPedestalsAnalysis:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 
virtual void endJob ()
 
 HcalPedestalsAnalysis (const edm::ParameterSet &ps)
 
virtual ~HcalPedestalsAnalysis ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

std::vector< NewPedBunchBunches
 
bool dumpXML
 
bool firsttime
 
int firstTS
 
edm::InputTag hbheDigiCollectionTag_
 
TH1F * HBMeans
 
TH1F * HBWidths
 
TH1F * HEMeans
 
TH1F * HEWidths
 
edm::InputTag hfDigiCollectionTag_
 
TH1F * HFMeans
 
TH1F * HFWidths
 
bool hiSaveFlag
 
edm::InputTag hoDigiCollectionTag_
 
TH1F * HOMeans
 
TH1F * HOWidths
 
int ievt
 
int lastTS
 
std::string pedsADCfilename
 
std::string pedsfCfilename
 
HcalPedestalsrawPedsItem
 
HcalPedestalsrawPedsItemfc
 
HcalPedestalWidthsrawWidthsItem
 
HcalPedestalWidthsrawWidthsItemfc
 
std::string ROOTfilename
 
int runnum
 
TFile * theFile
 
HcalTopologytheTopology
 
bool verboseflag
 
std::string widthsADCfilename
 
std::string widthsfCfilename
 
std::string XMLfilename
 
std::string XMLtag
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 73 of file HcalPedestalsAnalysis.h.

Constructor & Destructor Documentation

HcalPedestalsAnalysis::HcalPedestalsAnalysis ( const edm::ParameterSet ps)

Definition at line 13 of file HcalPedestalsAnalysis.cc.

References gather_cfg::cout, dumpXML, firsttime, firstTS, edm::ParameterSet::getUntrackedParameter(), hiSaveFlag, ievt, lastTS, rawPedsItem, rawPedsItemfc, rawWidthsItem, rawWidthsItemfc, and verboseflag.

13  :
14  hbheDigiCollectionTag_(ps.getParameter<edm::InputTag>("hbheDigiCollectionTag")),
15  hoDigiCollectionTag_(ps.getParameter<edm::InputTag>("hoDigiCollectionTag")),
16  hfDigiCollectionTag_(ps.getParameter<edm::InputTag>("hfDigiCollectionTag")) {
17 
18  std::cout << "Code version 10.6\n";
19  hiSaveFlag = ps.getUntrackedParameter<bool>("hiSaveFlag", false);
20  dumpXML = ps.getUntrackedParameter<bool>("dumpXML", true);
21  verboseflag = ps.getUntrackedParameter<bool>("verbose", false);
22  firstTS = ps.getUntrackedParameter<int>("firstTS", 0);
23  lastTS = ps.getUntrackedParameter<int>("lastTS", 9);
24  firsttime = true;
25  ievt = 0;
26 
27  rawPedsItem = 0;
28  rawWidthsItem = 0;
29  rawPedsItemfc = 0;
30  rawWidthsItemfc = 0;
31 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag hbheDigiCollectionTag_
HcalPedestalWidths * rawWidthsItem
tuple cout
Definition: gather_cfg.py:121
HcalPedestalWidths * rawWidthsItemfc
HcalPedestalsAnalysis::~HcalPedestalsAnalysis ( )
virtual

Definition at line 34 of file HcalPedestalsAnalysis.cc.

35 {}

Member Function Documentation

void HcalPedestalsAnalysis::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 206 of file HcalPedestalsAnalysis.cc.

References a, HcalQIESample::adc(), HcalElectronicsMap::allPrecisionId(), Bunches, NewPedBunch::cap, NewPedBunch::capfc, HcalQIESample::capid(), HcalQIECoder::charge(), NewPedBunch::detid, firsttime, firstTS, edm::EventSetup::get(), edm::Event::getByLabel(), hbheDigiCollectionTag_, HBMeans, HBWidths, HEMeans, HEWidths, hfDigiCollectionTag_, HFMeans, HFWidths, hoDigiCollectionTag_, HOMeans, HOWidths, i, HBHEDataFrame::id(), HFDataFrame::id(), HODataFrame::id(), edm::EventBase::id(), ievt, j, lastTS, NewPedBunch::num, pedsADCfilename, pedsfCfilename, NewPedBunch::prod, NewPedBunch::prodfc, edm::ESHandle< class >::product(), DetId::rawId(), rawPedsItem, rawPedsItemfc, rawWidthsItem, rawWidthsItemfc, ROOTfilename, edm::EventID::run(), runnum, HBHEDataFrame::sample(), HFDataFrame::sample(), HODataFrame::sample(), NewPedBunch::sig, NewPedBunch::sigfc, HBHEDataFrame::size(), HFDataFrame::size(), HODataFrame::size(), AlCaHLTBitMon_QueryRunRegistry::string, theFile, theTopology, NewPedBunch::usedflag, widthsADCfilename, widthsfCfilename, XMLfilename, and XMLtag.

207 {
208  using namespace edm;
209  using namespace std;
210 
214  edm::ESHandle<HcalDbService> conditions;
215  iSetup.get<HcalDbRecord>().get(conditions);
216 
217 
218  if(firsttime)
219  {
220 
222  iSetup.get<IdealGeometryRecord>().get( topology );
223  theTopology=new HcalTopology(*topology);
224 
229 
230  runnum = e.id().run();
231  std::string runnum_string;
232  std::stringstream tempstringout;
233  tempstringout << runnum;
234  runnum_string = tempstringout.str();
235  ROOTfilename = runnum_string + "-peds_ADC.root";
236  pedsADCfilename = runnum_string + "-peds_ADC_";
237  pedsfCfilename = runnum_string + "-peds_fC.txt";
238  widthsADCfilename = runnum_string + "-widths_ADC.txt";
239  widthsfCfilename = runnum_string + "-widths_fC.txt";
240  XMLfilename = runnum_string + "-peds_ADC_complete.xml";
241  XMLtag = "Hcal_pedestals_" + runnum_string;
242 
243  theFile = new TFile(ROOTfilename.c_str(), "RECREATE");
244  theFile->cd();
245  // Create sub-directories
246  theFile->mkdir("HB");
247  theFile->mkdir("HE");
248  theFile->mkdir("HF");
249  theFile->mkdir("HO");
250  theFile->cd();
251 
252  HBMeans = new TH1F("All Ped Means HB","All Ped Means HB", 100, 0, 9);
253  HBWidths = new TH1F("All Ped Widths HB","All Ped Widths HB", 100, 0, 3);
254  HEMeans = new TH1F("All Ped Means HE","All Ped Means HE", 100, 0, 9);
255  HEWidths = new TH1F("All Ped Widths HE","All Ped Widths HE", 100, 0, 3);
256  HFMeans = new TH1F("All Ped Means HF","All Ped Means HF", 100, 0, 9);
257  HFWidths = new TH1F("All Ped Widths HF","All Ped Widths HF", 100, 0, 3);
258  HOMeans = new TH1F("All Ped Means HO","All Ped Means HO", 100, 0, 9);
259  HOWidths = new TH1F("All Ped Widths HO","All Ped Widths HO", 100, 0, 3);
260 
262  iSetup.get<HcalElectronicsMapRcd>().get(refEMap);
263  const HcalElectronicsMap* myRefEMap = refEMap.product();
264  std::vector<HcalGenericDetId> listEMap = myRefEMap->allPrecisionId();
265  for (std::vector<HcalGenericDetId>::const_iterator it = listEMap.begin(); it != listEMap.end(); it++)
266  {
267  HcalGenericDetId mygenid(it->rawId());
268  if(mygenid.isHcalDetId())
269  {
270  NewPedBunch a;
271  HcalDetId chanid(mygenid.rawId());
272  a.detid = chanid;
273  a.usedflag = false;
274  for(int i = 0; i != 4; i++)
275  {
276  a.cap[i] = 0;
277  a.capfc[i] = 0;
278  for(int j = 0; j != 4; j++)
279  {
280  a.sig[i][j] = 0;
281  a.sigfc[i][j] = 0;
282  a.prod[i][j] = 0;
283  a.prodfc[i][j] = 0;
284  a.num[i][j] = 0;
285  }
286  }
287  Bunches.push_back(a);
288  }
289  }
290  firsttime = false;
291  }
292 
293  std::vector<NewPedBunch>::iterator bunch_it;
294 
295  for(HBHEDigiCollection::const_iterator j = hbhe->begin(); j != hbhe->end(); j++)
296  {
297  const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
298  for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
299  if(bunch_it->detid.rawId() == digi.id().rawId()) break;
300  bunch_it->usedflag = true;
301  for(int ts = firstTS; ts != lastTS+1; ts++)
302  {
303  const HcalQIECoder* coder = conditions->getHcalCoder(digi.id().rawId());
304  const HcalQIEShape* shape = conditions->getHcalShape(coder);
305  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts).capid()] += 1;
306  bunch_it->cap[digi.sample(ts).capid()] += digi.sample(ts).adc();
307  double charge1 = coder->charge(*shape, digi.sample(ts).adc(), digi.sample(ts).capid());
308  bunch_it->capfc[digi.sample(ts).capid()] += charge1;
309  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts).capid()] += (digi.sample(ts).adc() * digi.sample(ts).adc());
310  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts).capid()] += charge1 * charge1;
311  if((ts+1 < digi.size()) && (ts+1 < lastTS)){
312  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += digi.sample(ts).adc()*digi.sample(ts+1).adc();
313  double charge2 = coder->charge(*shape, digi.sample(ts+1).adc(), digi.sample(ts+1).capid());
314  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += charge1*charge2;
315  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += 1;
316  }
317  if((ts+2 < digi.size()) && (ts+2 < lastTS)){
318  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += digi.sample(ts).adc()*digi.sample(ts+2).adc();
319  double charge2 = coder->charge(*shape, digi.sample(ts+2).adc(), digi.sample(ts+2).capid());
320  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += charge1*charge2;
321  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += 1;
322  }
323  if((ts+3 < digi.size()) && (ts+3 < lastTS)){
324  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += digi.sample(ts).adc()*digi.sample(ts+3).adc();
325  double charge2 = coder->charge(*shape, digi.sample(ts+3).adc(), digi.sample(ts+3).capid());
326  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += charge1*charge2;
327  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += 1;
328  }
329  }
330  }
331 
332  for(HODigiCollection::const_iterator j = ho->begin(); j != ho->end(); j++)
333  {
334  const HODataFrame digi = (const HODataFrame)(*j);
335  for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
336  if(bunch_it->detid.rawId() == digi.id().rawId()) break;
337  bunch_it->usedflag = true;
338  for(int ts = firstTS; ts <= lastTS; ts++)
339  {
340  const HcalQIECoder* coder = conditions->getHcalCoder(digi.id().rawId());
341  const HcalQIEShape* shape = conditions->getHcalShape(coder);
342  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts).capid()] += 1;
343  bunch_it->cap[digi.sample(ts).capid()] += digi.sample(ts).adc();
344  double charge1 = coder->charge(*shape, digi.sample(ts).adc(), digi.sample(ts).capid());
345  bunch_it->capfc[digi.sample(ts).capid()] += charge1;
346  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts).capid()] += (digi.sample(ts).adc() * digi.sample(ts).adc());
347  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts).capid()] += charge1 * charge1;
348  if((ts+1 < digi.size()) && (ts+1 < lastTS)){
349  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += digi.sample(ts).adc()*digi.sample(ts+1).adc();
350  double charge2 = coder->charge(*shape, digi.sample(ts+1).adc(), digi.sample(ts+1).capid());
351  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += charge1*charge2;
352  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += 1;
353  }
354  if((ts+2 < digi.size()) && (ts+2 < lastTS)){
355  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += digi.sample(ts).adc()*digi.sample(ts+2).adc();
356  double charge2 = coder->charge(*shape, digi.sample(ts+2).adc(), digi.sample(ts+2).capid());
357  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += charge1*charge2;
358  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += 1;
359  }
360  if((ts+3 < digi.size()) && (ts+3 < lastTS)){
361  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += digi.sample(ts).adc()*digi.sample(ts+3).adc();
362  double charge2 = coder->charge(*shape, digi.sample(ts+3).adc(), digi.sample(ts+3).capid());
363  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += charge1*charge2;
364  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += 1;
365  }
366  }
367  }
368 
369  for(HFDigiCollection::const_iterator j = hf->begin(); j != hf->end(); j++)
370  {
371  const HFDataFrame digi = (const HFDataFrame)(*j);
372  for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
373  if(bunch_it->detid.rawId() == digi.id().rawId()) break;
374  bunch_it->usedflag = true;
375  for(int ts = firstTS; ts <= lastTS; ts++)
376  {
377  const HcalQIECoder* coder = conditions->getHcalCoder(digi.id().rawId());
378  const HcalQIEShape* shape = conditions->getHcalShape(coder);
379  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts).capid()] += 1;
380  bunch_it->cap[digi.sample(ts).capid()] += digi.sample(ts).adc();
381  double charge1 = coder->charge(*shape, digi.sample(ts).adc(), digi.sample(ts).capid());
382  bunch_it->capfc[digi.sample(ts).capid()] += charge1;
383  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts).capid()] += (digi.sample(ts).adc() * digi.sample(ts).adc());
384  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts).capid()] += charge1 * charge1;
385  if((ts+1 < digi.size()) && (ts+1 < lastTS)){
386  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += digi.sample(ts).adc()*digi.sample(ts+1).adc();
387  double charge2 = coder->charge(*shape, digi.sample(ts+1).adc(), digi.sample(ts+1).capid());
388  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += charge1*charge2;
389  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += 1;
390  }
391  if((ts+2 < digi.size()) && (ts+2 < lastTS)){
392  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += digi.sample(ts).adc()*digi.sample(ts+2).adc();
393  double charge2 = coder->charge(*shape, digi.sample(ts+2).adc(), digi.sample(ts+2).capid());
394  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += charge1*charge2;
395  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += 1;
396  }
397  if((ts+3 < digi.size()) && (ts+3 < lastTS)){
398  bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += digi.sample(ts).adc()*digi.sample(ts+3).adc();
399  double charge2 = coder->charge(*shape, digi.sample(ts+3).adc(), digi.sample(ts+3).capid());
400  bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += charge1*charge2;
401  bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += 1;
402  }
403  }
404  }
405 
406 //this is the last brace
407 ievt++;
408 }
int i
Definition: DBlmapReader.cc:9
edm::InputTag hbheDigiCollectionTag_
int adc() const
get the ADC sample
Definition: HcalQIESample.h:24
int size() const
total number of samples in the digi
Definition: HBHEDataFrame.h:26
std::vector< T >::const_iterator const_iterator
const HcalQIESample & sample(int i) const
access a sample
Definition: HODataFrame.h:40
HcalPedestalWidths * rawWidthsItem
const HcalDetId & id() const
Definition: HODataFrame.h:23
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
int size() const
total number of samples in the digi
Definition: HODataFrame.h:27
const HcalQIESample & sample(int i) const
access a sample
Definition: HFDataFrame.h:39
int j
Definition: DBlmapReader.cc:9
std::vector< HcalGenericDetId > allPrecisionId() const
HcalCastorDetId detid
int size() const
total number of samples in the digi
Definition: HFDataFrame.h:26
int capid() const
get the Capacitor id
Definition: HcalQIESample.h:28
const HcalQIESample & sample(int i) const
access a sample
Definition: HBHEDataFrame.h:39
T const * product() const
Definition: ESHandle.h:62
double a
Definition: hdecay.h:121
const HcalDetId & id() const
Definition: HBHEDataFrame.h:22
HcalPedestalWidths * rawWidthsItemfc
const HcalDetId & id() const
Definition: HFDataFrame.h:22
std::vector< NewPedBunch > Bunches
float charge(const HcalQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -&gt; fC conversion.
Definition: HcalQIECoder.cc:22
void HcalPedestalsAnalysis::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 37 of file HcalPedestalsAnalysis.cc.

References HcalCondObjectContainer< Item >::addValues(), Bunches, gather_cfg::cout, HcalDbASCIIIO::dumpObject(), HBMeans, HBWidths, HEMeans, HEWidths, HFMeans, HFWidths, hiSaveFlag, HOMeans, HOWidths, i, ievt, pedsADCfilename, rawPedsItem, rawPedsItemfc, rawWidthsItem, rawWidthsItemfc, HcalPedestalWidth::setSigma(), AlCaHLTBitMon_QueryRunRegistry::string, theFile, and verboseflag.

38 {
39 // std::cout << "NEvents " << ievt << std::endl;
40  if(ievt < 1000) return;
41  std::string ievt_string;
42  std::stringstream tempstringout;
43  tempstringout << ievt;
44  ievt_string = tempstringout.str();
45  pedsADCfilename += ievt_string;
46  pedsADCfilename += ".txt";
47  //Calculate pedestal constants
48  std::cout << "Calculating Pedestal constants...\n";
49  std::vector<NewPedBunch>::iterator bunch_it;
50  for(bunch_it=Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
51  {
52  if(bunch_it->usedflag){
53 
54  if(verboseflag) std::cout << "Analyzing channel " << bunch_it->detid << std::endl;
55  //pedestal constant is the mean
56  if(bunch_it->num[0][0]!=0) bunch_it->cap[0] /= bunch_it->num[0][0];
57  if(bunch_it->num[1][1]!=0) bunch_it->cap[1] /= bunch_it->num[1][1];
58  if(bunch_it->num[2][2]!=0) bunch_it->cap[2] /= bunch_it->num[2][2];
59  if(bunch_it->num[3][3]!=0) bunch_it->cap[3] /= bunch_it->num[3][3];
60  if(bunch_it->num[0][0]!=0) bunch_it->capfc[0] /= bunch_it->num[0][0];
61  if(bunch_it->num[1][1]!=0) bunch_it->capfc[1] /= bunch_it->num[1][1];
62  if(bunch_it->num[2][2]!=0) bunch_it->capfc[2] /= bunch_it->num[2][2];
63  if(bunch_it->num[3][3]!=0) bunch_it->capfc[3] /= bunch_it->num[3][3];
64  bunch_it->sig[0][0] = (bunch_it->prod[0][0]/bunch_it->num[0][0])-(bunch_it->cap[0]*bunch_it->cap[0]);
65  bunch_it->sig[0][1] = (bunch_it->prod[0][1]/bunch_it->num[0][1])-(bunch_it->cap[0]*bunch_it->cap[1]);
66  bunch_it->sig[0][2] = (bunch_it->prod[0][2]/bunch_it->num[0][2])-(bunch_it->cap[0]*bunch_it->cap[2]);
67  bunch_it->sig[0][3] = (bunch_it->prod[0][3]/bunch_it->num[0][3])-(bunch_it->cap[0]*bunch_it->cap[3]);
68  bunch_it->sig[1][0] = (bunch_it->prod[1][0]/bunch_it->num[1][0])-(bunch_it->cap[1]*bunch_it->cap[0]);
69  bunch_it->sig[1][1] = (bunch_it->prod[1][1]/bunch_it->num[1][1])-(bunch_it->cap[1]*bunch_it->cap[1]);
70  bunch_it->sig[1][2] = (bunch_it->prod[1][2]/bunch_it->num[1][2])-(bunch_it->cap[1]*bunch_it->cap[2]);
71  bunch_it->sig[1][3] = (bunch_it->prod[1][3]/bunch_it->num[1][3])-(bunch_it->cap[1]*bunch_it->cap[3]);
72  bunch_it->sig[2][0] = (bunch_it->prod[2][0]/bunch_it->num[2][0])-(bunch_it->cap[2]*bunch_it->cap[0]);
73  bunch_it->sig[2][1] = (bunch_it->prod[2][1]/bunch_it->num[2][1])-(bunch_it->cap[2]*bunch_it->cap[1]);
74  bunch_it->sig[2][2] = (bunch_it->prod[2][2]/bunch_it->num[2][2])-(bunch_it->cap[2]*bunch_it->cap[2]);
75  bunch_it->sig[2][3] = (bunch_it->prod[2][3]/bunch_it->num[2][3])-(bunch_it->cap[2]*bunch_it->cap[3]);
76  bunch_it->sig[3][0] = (bunch_it->prod[3][0]/bunch_it->num[3][0])-(bunch_it->cap[3]*bunch_it->cap[0]);
77  bunch_it->sig[3][1] = (bunch_it->prod[3][1]/bunch_it->num[3][1])-(bunch_it->cap[3]*bunch_it->cap[1]);
78  bunch_it->sig[3][2] = (bunch_it->prod[3][2]/bunch_it->num[3][2])-(bunch_it->cap[3]*bunch_it->cap[2]);
79  bunch_it->sig[3][3] = (bunch_it->prod[3][3]/bunch_it->num[3][3])-(bunch_it->cap[3]*bunch_it->cap[3]);
80 
81  bunch_it->sigfc[0][0] = (bunch_it->prodfc[0][0]/bunch_it->num[0][0])-(bunch_it->capfc[0]*bunch_it->capfc[0]);
82  bunch_it->sigfc[0][1] = (bunch_it->prodfc[0][1]/bunch_it->num[0][1])-(bunch_it->capfc[0]*bunch_it->capfc[1]);
83  bunch_it->sigfc[0][2] = (bunch_it->prodfc[0][2]/bunch_it->num[0][2])-(bunch_it->capfc[0]*bunch_it->capfc[2]);
84  bunch_it->sigfc[0][3] = (bunch_it->prodfc[0][3]/bunch_it->num[0][3])-(bunch_it->capfc[0]*bunch_it->capfc[3]);
85  bunch_it->sigfc[1][0] = (bunch_it->prodfc[1][0]/bunch_it->num[1][0])-(bunch_it->capfc[1]*bunch_it->capfc[0]);
86  bunch_it->sigfc[1][1] = (bunch_it->prodfc[1][1]/bunch_it->num[1][1])-(bunch_it->capfc[1]*bunch_it->capfc[1]);
87  bunch_it->sigfc[1][2] = (bunch_it->prodfc[1][2]/bunch_it->num[1][2])-(bunch_it->capfc[1]*bunch_it->capfc[2]);
88  bunch_it->sigfc[1][3] = (bunch_it->prodfc[1][3]/bunch_it->num[1][3])-(bunch_it->capfc[1]*bunch_it->capfc[3]);
89  bunch_it->sigfc[2][0] = (bunch_it->prodfc[2][0]/bunch_it->num[2][0])-(bunch_it->capfc[2]*bunch_it->capfc[0]);
90  bunch_it->sigfc[2][1] = (bunch_it->prodfc[2][1]/bunch_it->num[2][1])-(bunch_it->capfc[2]*bunch_it->capfc[1]);
91  bunch_it->sigfc[2][2] = (bunch_it->prodfc[2][2]/bunch_it->num[2][2])-(bunch_it->capfc[2]*bunch_it->capfc[2]);
92  bunch_it->sigfc[2][3] = (bunch_it->prodfc[2][3]/bunch_it->num[2][3])-(bunch_it->capfc[2]*bunch_it->capfc[3]);
93  bunch_it->sigfc[3][0] = (bunch_it->prodfc[3][0]/bunch_it->num[3][0])-(bunch_it->capfc[3]*bunch_it->capfc[0]);
94  bunch_it->sigfc[3][1] = (bunch_it->prodfc[3][1]/bunch_it->num[3][1])-(bunch_it->capfc[3]*bunch_it->capfc[1]);
95  bunch_it->sigfc[3][2] = (bunch_it->prodfc[3][2]/bunch_it->num[3][2])-(bunch_it->capfc[3]*bunch_it->capfc[2]);
96  bunch_it->sigfc[3][3] = (bunch_it->prodfc[3][3]/bunch_it->num[3][3])-(bunch_it->capfc[3]*bunch_it->capfc[3]);
97 
98  if(bunch_it->detid.subdet() == 1){
99  for(int i = 0; i != 4; i++){
100  HBMeans->Fill(bunch_it->cap[i]);
101  HBWidths->Fill(bunch_it->sig[i][i]);
102  }
103  }
104  if(bunch_it->detid.subdet() == 2){
105  for(int i = 0; i != 4; i++){
106  HEMeans->Fill(bunch_it->cap[i]);
107  HEWidths->Fill(bunch_it->sig[i][i]);
108  }
109  }
110  if(bunch_it->detid.subdet() == 3){
111  for(int i = 0; i != 4; i++){
112  HOMeans->Fill(bunch_it->cap[i]);
113  HOWidths->Fill(bunch_it->sig[i][i]);
114  }
115  }
116  if(bunch_it->detid.subdet() == 4){
117  for(int i = 0; i != 4; i++){
118  HFMeans->Fill(bunch_it->cap[i]);
119  HFWidths->Fill(bunch_it->sig[i][i]);
120  }
121  }
122 
123  const HcalPedestal item(bunch_it->detid, bunch_it->cap[0], bunch_it->cap[1], bunch_it->cap[2], bunch_it->cap[3],
124  bunch_it->sig[0][0], bunch_it->sig[1][1], bunch_it->sig[2][2], bunch_it->sig[3][3]);
125  rawPedsItem->addValues(item);
126  HcalPedestalWidth widthsp(bunch_it->detid);
127  widthsp.setSigma(0,0,bunch_it->sig[0][0]);
128  widthsp.setSigma(0,1,bunch_it->sig[0][1]);
129  widthsp.setSigma(0,2,bunch_it->sig[0][2]);
130  widthsp.setSigma(0,3,bunch_it->sig[0][3]);
131  widthsp.setSigma(1,1,bunch_it->sig[1][1]);
132  widthsp.setSigma(1,2,bunch_it->sig[1][2]);
133  widthsp.setSigma(1,3,bunch_it->sig[1][3]);
134  widthsp.setSigma(2,2,bunch_it->sig[2][2]);
135  widthsp.setSigma(2,3,bunch_it->sig[2][3]);
136  widthsp.setSigma(3,3,bunch_it->sig[3][3]);
137  rawWidthsItem->addValues(widthsp);
138 
139  const HcalPedestal itemfc(bunch_it->detid, bunch_it->capfc[0], bunch_it->capfc[1], bunch_it->capfc[2], bunch_it->capfc[3],
140  bunch_it->sigfc[0][0], bunch_it->sigfc[1][1], bunch_it->sigfc[2][2], bunch_it->sigfc[3][3]);
141  rawPedsItemfc->addValues(itemfc);
142  HcalPedestalWidth widthspfc(bunch_it->detid);
143  widthspfc.setSigma(0,0,bunch_it->sigfc[0][0]);
144  widthspfc.setSigma(0,1,bunch_it->sigfc[0][1]);
145  widthspfc.setSigma(0,2,bunch_it->sigfc[0][2]);
146  widthspfc.setSigma(0,3,bunch_it->sigfc[0][3]);
147  widthspfc.setSigma(1,1,bunch_it->sigfc[1][1]);
148  widthspfc.setSigma(1,2,bunch_it->sigfc[1][2]);
149  widthspfc.setSigma(1,3,bunch_it->sigfc[1][3]);
150  widthspfc.setSigma(2,2,bunch_it->sigfc[2][2]);
151  widthspfc.setSigma(2,3,bunch_it->sigfc[2][3]);
152  widthspfc.setSigma(3,3,bunch_it->sigfc[3][3]);
153  rawWidthsItemfc->addValues(widthspfc);
154  }
155  }
156 
157  // dump the resulting list of pedestals into a file
158  std::ofstream outStream1(pedsADCfilename.c_str());
159  HcalDbASCIIIO::dumpObject (outStream1, (*rawPedsItem) );
160 // std::ofstream outStream2(widthsADCfilename.c_str());
161 // HcalDbASCIIIO::dumpObject (outStream2, (*rawWidthsItem) );
162 //
163 // std::ofstream outStream3(pedsfCfilename.c_str());
164 // HcalDbASCIIIO::dumpObject (outStream3, (*rawPedsItemfc) );
165 // std::ofstream outStream4(widthsfCfilename.c_str());
166 // HcalDbASCIIIO::dumpObject (outStream4, (*rawWidthsItemfc) );
167 //
168 // if(dumpXML){
169 // std::ofstream outStream5(XMLfilename.c_str());
170 // HcalDbXml::dumpObject (outStream5, runnum, 0, 2147483647, XMLtag, 1, (*rawPedsItem), (*rawWidthsItem));
171 // }
172 //
173  if(hiSaveFlag){
174  theFile->Write();
175  }else{
176  theFile->cd();
177  theFile->cd("HB");
178  HBMeans->Write();
179  HBWidths->Write();
180  theFile->cd();
181  theFile->cd("HF");
182  HFMeans->Write();
183  HFWidths->Write();
184  theFile->cd();
185  theFile->cd("HE");
186  HEMeans->Write();
187  HEWidths->Write();
188  theFile->cd();
189  theFile->cd("HO");
190  HOMeans->Write();
191  HOWidths->Write();
192  }
193 
194  std::cout << "Writing ROOT file... ";
195  theFile->Close();
196  std::cout << "ROOT file closed.\n";
197 
198  delete rawPedsItem;
199  delete rawWidthsItem;
200  delete rawPedsItemfc;
201  delete rawWidthsItemfc;
202 }
int i
Definition: DBlmapReader.cc:9
HcalPedestalWidths * rawWidthsItem
bool dumpObject(std::ostream &fOutput, const HcalPedestals &fObject)
tuple cout
Definition: gather_cfg.py:121
HcalPedestalWidths * rawWidthsItemfc
bool addValues(const Item &myItem)
std::vector< NewPedBunch > Bunches

Member Data Documentation

std::vector<NewPedBunch> HcalPedestalsAnalysis::Bunches
private

Definition at line 86 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

bool HcalPedestalsAnalysis::dumpXML
private

Definition at line 89 of file HcalPedestalsAnalysis.h.

Referenced by HcalPedestalsAnalysis().

bool HcalPedestalsAnalysis::firsttime
private

Definition at line 113 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and HcalPedestalsAnalysis().

int HcalPedestalsAnalysis::firstTS
private

Definition at line 92 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and HcalPedestalsAnalysis().

edm::InputTag HcalPedestalsAnalysis::hbheDigiCollectionTag_
private

Definition at line 120 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

TH1F* HcalPedestalsAnalysis::HBMeans
private

Definition at line 103 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

TH1F* HcalPedestalsAnalysis::HBWidths
private

Definition at line 104 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

TH1F* HcalPedestalsAnalysis::HEMeans
private

Definition at line 105 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

TH1F* HcalPedestalsAnalysis::HEWidths
private

Definition at line 106 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

edm::InputTag HcalPedestalsAnalysis::hfDigiCollectionTag_
private

Definition at line 122 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

TH1F* HcalPedestalsAnalysis::HFMeans
private

Definition at line 107 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

TH1F* HcalPedestalsAnalysis::HFWidths
private

Definition at line 108 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

bool HcalPedestalsAnalysis::hiSaveFlag
private

Definition at line 88 of file HcalPedestalsAnalysis.h.

Referenced by endJob(), and HcalPedestalsAnalysis().

edm::InputTag HcalPedestalsAnalysis::hoDigiCollectionTag_
private

Definition at line 121 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

TH1F* HcalPedestalsAnalysis::HOMeans
private

Definition at line 109 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

TH1F* HcalPedestalsAnalysis::HOWidths
private

Definition at line 110 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

int HcalPedestalsAnalysis::ievt
private

Definition at line 94 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), endJob(), and HcalPedestalsAnalysis().

int HcalPedestalsAnalysis::lastTS
private

Definition at line 93 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and HcalPedestalsAnalysis().

std::string HcalPedestalsAnalysis::pedsADCfilename
private

Definition at line 96 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

std::string HcalPedestalsAnalysis::pedsfCfilename
private

Definition at line 97 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

HcalPedestals* HcalPedestalsAnalysis::rawPedsItem
private

Definition at line 114 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), endJob(), and HcalPedestalsAnalysis().

HcalPedestals* HcalPedestalsAnalysis::rawPedsItemfc
private

Definition at line 116 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), endJob(), and HcalPedestalsAnalysis().

HcalPedestalWidths* HcalPedestalsAnalysis::rawWidthsItem
private

Definition at line 115 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), endJob(), and HcalPedestalsAnalysis().

HcalPedestalWidths* HcalPedestalsAnalysis::rawWidthsItemfc
private

Definition at line 117 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), endJob(), and HcalPedestalsAnalysis().

std::string HcalPedestalsAnalysis::ROOTfilename
private

Definition at line 95 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

int HcalPedestalsAnalysis::runnum
private

Definition at line 91 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

TFile* HcalPedestalsAnalysis::theFile
private

Definition at line 112 of file HcalPedestalsAnalysis.h.

Referenced by analyze(), and endJob().

HcalTopology* HcalPedestalsAnalysis::theTopology
private

Definition at line 118 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

bool HcalPedestalsAnalysis::verboseflag
private

Definition at line 90 of file HcalPedestalsAnalysis.h.

Referenced by endJob(), and HcalPedestalsAnalysis().

std::string HcalPedestalsAnalysis::widthsADCfilename
private

Definition at line 98 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

std::string HcalPedestalsAnalysis::widthsfCfilename
private

Definition at line 99 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

std::string HcalPedestalsAnalysis::XMLfilename
private

Definition at line 100 of file HcalPedestalsAnalysis.h.

Referenced by analyze().

std::string HcalPedestalsAnalysis::XMLtag
private

Definition at line 101 of file HcalPedestalsAnalysis.h.

Referenced by analyze().