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