CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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)
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  dtGeom = eventSetup.getHandle(dtGeomToken_);
76 
77  // Get all the sls from the setup
78  const vector<const DTSuperLayer*> superLayers = dtGeom->superLayers();
79 
80  // Loop over all SLs
81  for (auto sl = superLayers.begin(); sl != superLayers.end(); sl++) {
82  // Get the histo from file
83  DTSuperLayerId slId = (*sl)->id();
84  TH1F* histo = (TH1F*)theFile->Get((getTBoxName(slId)).c_str());
85  if (histo) { // Check that the histo exists
86  // Compute mean and sigma of the rising edge
87  pair<double, double> meanAndSigma = theFitter->fitTimeBox(histo);
88 
89  // Write them in DB object
90  tTrig->set(slId, meanAndSigma.first, meanAndSigma.second, kFactor, DTTimeUnits::ns);
91  if (debug) {
92  cout << " SL: " << slId << " mean = " << meanAndSigma.first << " sigma = " << meanAndSigma.second << endl;
93  }
94  }
95  }
96 }
97 
98 // Write objects to DB
100  if (debug)
101  cout << "[DTTTrigWriter]Writing ttrig object to DB!" << endl;
102 
103  // FIXME: to be read from cfg?
104  string tTrigRecord = "DTTtrigRcd";
105 
106  // Write the object to DB
107  DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
108 }
109 
110 // Compute the name of the time box histo
111 string DTTTrigWriter::getTBoxName(const DTSuperLayerId& slId) const {
112  string histoName = "Ch_" + std::to_string(slId.wheel()) + "_" + std::to_string(slId.station()) + "_" +
113  std::to_string(slId.sector()) + "_SL" + std::to_string(slId.superlayer()) + "_hTimeBox";
114  return histoName;
115 }
DTTimeBoxFitter * theFitter
Definition: DTTTrigWriter.h:62
T getUntrackedParameter(std::string const &, T const &) const
std::string theRootInputFile
Definition: DTTTrigWriter.h:59
int set(int wheelId, int stationId, int sectorId, int slId, float tTrig, float tTrms, float kFact, DTTimeUnits::type unit)
Definition: DTTtrig.cc:172
void setVerbosity(unsigned int lvl)
Set the verbosity of the output: 0 = silent, 1 = info, 2 = debug.
TFile * theFile
Definition: DTTTrigWriter.h:56
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
Compute the ttrig by fiting the TB rising edge.
DTTTrigWriter(const edm::ParameterSet &pset)
Constructor.
std::string getTBoxName(const DTSuperLayerId &slId) const
edm::ESHandle< DTGeometry > dtGeom
Definition: DTTTrigWriter.h:68
void setFitSigma(double sigma)
int superlayer() const
Return the superlayer number (deprecated method name)
DTTtrig * tTrig
Definition: DTTTrigWriter.h:65
int sector() const
Definition: DTChamberId.h:49
std::pair< double, double > fitTimeBox(TH1F *hTimeBox)
Fit the rising edge of the time box returning mean value and sigma (first and second respectively) ...
tuple cout
Definition: gather_cfg.py:144
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
static void writeToDB(std::string record, T *payload)
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
Definition: DTTTrigWriter.h:69
void endJob() override
Write ttrig in the DB.
~DTTTrigWriter() override
Destructor.