CMS 3D CMS Logo

DTTTrigWriter.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author S. Bolognesi
5  */
6 
10 
18 
22 
23 /* C++ Headers */
24 #include <vector>
25 #include <iostream>
26 #include <fstream>
27 #include <string>
28 #include <sstream>
29 #include "TFile.h"
30 #include "TH1.h"
31 
32 using namespace std;
33 using namespace edm;
34 
35 // Constructor
37  // get selected debug option
38  debug = pset.getUntrackedParameter<bool>("debug", false);
39 
40  // Open the root file which contains the histos
41  theRootInputFile = pset.getUntrackedParameter<string>("rootFileName");
42  theFile = new TFile(theRootInputFile.c_str(), "READ");
43  theFile->cd();
44  theFitter = new DTTimeBoxFitter();
45  if (debug)
46  theFitter->setVerbosity(1);
47 
48  double sigmaFit = pset.getUntrackedParameter<double>("sigmaTTrigFit", 10.);
49  theFitter->setFitSigma(sigmaFit);
50 
51  // the kfactor to be uploaded in the ttrig DB
52  kFactor = pset.getUntrackedParameter<double>("kFactor", -0.7);
53 
54  // Create the object to be written to DB
55  tTrig = new DTTtrig();
56 
57  if (debug)
58  cout << "[DTTTrigWriter]Constructor called!" << endl;
59 }
60 
61 // Destructor
63  if (debug)
64  cout << "[DTTTrigWriter]Destructor called!" << endl;
65  theFile->Close();
66  delete theFitter;
67 }
68 
69 // Do the job
70 void DTTTrigWriter::analyze(const Event& event, const EventSetup& eventSetup) {
71  if (debug)
72  cout << "[DTTTrigWriter]Analyzer called!" << endl;
73 
74  // Get the DT Geometry
75  ESHandle<DTGeometry> dtGeom;
76  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
77 
78  // Get all the sls from the setup
79  const vector<const DTSuperLayer*> superLayers = dtGeom->superLayers();
80 
81  // Loop over all SLs
82  for (auto sl = superLayers.begin(); sl != superLayers.end(); sl++) {
83  // Get the histo from file
84  DTSuperLayerId slId = (*sl)->id();
85  TH1F* histo = (TH1F*)theFile->Get((getTBoxName(slId)).c_str());
86  if (histo) { // Check that the histo exists
87  // Compute mean and sigma of the rising edge
88  pair<double, double> meanAndSigma = theFitter->fitTimeBox(histo);
89 
90  // Write them in DB object
91  tTrig->set(slId, meanAndSigma.first, meanAndSigma.second, kFactor, DTTimeUnits::ns);
92  if (debug) {
93  cout << " SL: " << slId << " mean = " << meanAndSigma.first << " sigma = " << meanAndSigma.second << endl;
94  }
95  }
96  }
97 }
98 
99 // Write objects to DB
101  if (debug)
102  cout << "[DTTTrigWriter]Writing ttrig object to DB!" << endl;
103 
104  // FIXME: to be read from cfg?
105  string tTrigRecord = "DTTtrigRcd";
106 
107  // Write the object to DB
108  DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
109 }
110 
111 // Compute the name of the time box histo
112 string DTTTrigWriter::getTBoxName(const DTSuperLayerId& slId) const {
113  string histoName = "Ch_" + std::to_string(slId.wheel()) + "_" + std::to_string(slId.station()) + "_" +
114  std::to_string(slId.sector()) + "_SL" + std::to_string(slId.superlayer()) + "_hTimeBox";
115  return histoName;
116 }
DTTTrigWriter::DTTTrigWriter
DTTTrigWriter(const edm::ParameterSet &pset)
Constructor.
Definition: DTTTrigWriter.cc:36
DTSuperLayerId
Definition: DTSuperLayerId.h:12
DTTimeBoxFitter.h
DTTtrig
Definition: DTTtrig.h:68
ESHandle.h
edm
HLT enums.
Definition: AlignableModifier.h:19
dttriganalyzer_cfi.tTrig
tTrig
Definition: dttriganalyzer_cfi.py:11
gather_cfg.cout
cout
Definition: gather_cfg.py:144
timingPdfMaker.histo
histo
Definition: timingPdfMaker.py:278
DTStatusFlagRcd.h
DTSuperLayerId::superlayer
int superlayer() const
Return the superlayer number (deprecated method name)
Definition: DTSuperLayerId.h:42
DTTTrigWriter::endJob
void endJob() override
Write ttrig in the DB.
Definition: DTTTrigWriter.cc:100
debug
#define debug
Definition: HDRShower.cc:19
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
interactiveExample.theFile
theFile
Definition: interactiveExample.py:17
DTTimeUnits::ns
Definition: DTTimeUnits.h:32
edm::ESHandle< DTGeometry >
DTStatusFlag.h
DTTTrigWriter::~DTTTrigWriter
~DTTTrigWriter() override
Destructor.
Definition: DTTTrigWriter.cc:62
DTGeometry.h
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
DTTTrigWriter.h
DTCalibDBUtils.h
DTGeometry::superLayers
const std::vector< const DTSuperLayer * > & superLayers() const
Return a vector of all SuperLayer.
Definition: DTGeometry.cc:86
edm::EventSetup
Definition: EventSetup.h:58
DTTtrig.h
get
#define get
DTChamberId::sector
int sector() const
Definition: DTChamberId.h:49
std
Definition: JetResolutionObject.h:76
dttriganalyzer_cfi.kFactor
kFactor
Definition: dttriganalyzer_cfi.py:7
DTTimeBoxFitter
Definition: DTTimeBoxFitter.h:17
DTTTrigWriter::getTBoxName
std::string getTBoxName(const DTSuperLayerId &slId) const
Definition: DTTTrigWriter.cc:112
HltBtagPostValidation_cff.histoName
histoName
Definition: HltBtagPostValidation_cff.py:17
EventSetup.h
ParameterSet.h
DTSuperLayerId.h
MuonGeometryRecord.h
event
Definition: event.py:1
edm::Event
Definition: Event.h:73
MuonGeometryRecord
Definition: MuonGeometryRecord.h:34
DTChamberId::wheel
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
DTCalibDBUtils::writeToDB
static void writeToDB(std::string record, T *payload)
Definition: DTCalibDBUtils.h:28
DTChamberId::station
int station() const
Return the station number.
Definition: DTChamberId.h:42
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
DTTTrigWriter::analyze
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
Compute the ttrig by fiting the TB rising edge.
Definition: DTTTrigWriter.cc:70