CMS 3D CMS Logo

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

#include <DTLocalTriggerLutTest.h>

Inheritance diagram for DTLocalTriggerLutTest:
DTLocalTriggerBaseTest edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 DTLocalTriggerLutTest (const edm::ParameterSet &ps)
 Constructor. More...
 
virtual ~DTLocalTriggerLutTest ()
 Destructor. More...
 
- Public Member Functions inherited from DTLocalTriggerBaseTest
 DTLocalTriggerBaseTest ()
 Constructor. More...
 
virtual ~DTLocalTriggerBaseTest ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Protected Member Functions

void beginJob ()
 BeginJob. More...
 
void beginRun (const edm::Run &r, const edm::EventSetup &c)
 BeginRun. More...
 
void runClientDiagnostic ()
 Run client analysis. More...
 
- Protected Member Functions inherited from DTLocalTriggerBaseTest
void analyze (const edm::Event &e, const edm::EventSetup &c)
 Analyze. More...
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context)
 Perform begin lumiblock operations. More...
 
void bookCmsHistos (std::string hTag, std::string folder="", bool isGlb=false)
 Book the new MEs (CMS summary) More...
 
void bookSectorHistos (int wheel, int sector, std::string hTag, std::string folder="")
 Book the new MEs (for each sector) More...
 
void bookWheelHistos (int wheel, std::string hTag, std::string folder="")
 Book the new MEs (for each wheel) More...
 
std::string category ()
 Get message logger name. More...
 
void endJob ()
 EndJob. More...
 
void endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context)
 Perform client diagnostic in online. More...
 
void endRun (edm::Run const &run, edm::EventSetup const &context)
 Perform client diagnostic in offline. More...
 
std::string fullName (std::string htype)
 Create fullname from histo partial name. More...
 
template<class T >
TgetHisto (MonitorElement *me)
 Convert ME to Histogram fo type T. More...
 
std::string getMEName (std::string histoTag, std::string subfolder, const DTChamberId &chambid)
 Get the ME name (by chamber) More...
 
std::string getMEName (std::string histoTag, std::string subfolder, int wh)
 Get the ME name (by wheel) More...
 
std::pair< float, float > phiRange (const DTChamberId &id)
 Calculate phi range for histograms. More...
 
void setConfig (const edm::ParameterSet &ps, std::string name)
 Set configuration variables. More...
 
std::string & topFolder (bool isDCC)
 Get top folder name. More...
 
- 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)
 

Private Member Functions

void fillWhPlot (MonitorElement *plot, int sect, int stat, float value, bool lessIsBest=true)
 Fill summary plots managing double MB4 chambers. More...
 
int performLutTest (double mean, double RMS, double thresholdMean, double thresholdRMS)
 Perform Lut Test logical operations. More...
 

Private Attributes

bool doCorrStudy
 
double thresholdPhibMean
 
double thresholdPhibRMS
 
double thresholdPhiMean
 
double thresholdPhiRMS
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Attributes inherited from DTLocalTriggerBaseTest
std::string baseFolderDCC
 
std::string baseFolderDDU
 
std::map< std::string,
MonitorElement * > 
cmsME
 
DQMStoredbe
 
std::string hwSource
 
std::vector< std::string > hwSources
 
edm::ESHandle< DTGeometrymuonGeom
 
int nevents
 
unsigned int nLumiSegs
 
edm::ParameterSet parameters
 
int prescaleFactor
 
int run
 
bool runOnline
 
std::map< int, std::map
< std::string, MonitorElement * > > 
secME
 
std::string sourceFolder
 
std::string testName
 
std::string trigSource
 
std::vector< std::string > trigSources
 
std::map< int, std::map
< std::string, MonitorElement * > > 
whME
 

Detailed Description

Definition at line 18 of file DTLocalTriggerLutTest.h.

Constructor & Destructor Documentation

DTLocalTriggerLutTest::DTLocalTriggerLutTest ( const edm::ParameterSet ps)

Constructor.

Definition at line 34 of file DTLocalTriggerLutTest.cc.

References edm::ParameterSet::getUntrackedParameter().

34  {
35 
36  setConfig(ps,"DTLocalTriggerLut");
37  baseFolderDCC = "DT/03-LocalTrigger-DCC/";
38  baseFolderDDU = "DT/04-LocalTrigger-DDU/";
39  thresholdPhiMean = ps.getUntrackedParameter<double>("thresholdPhiMean",1.5);
40  thresholdPhiRMS = ps.getUntrackedParameter<double>("thresholdPhiRMS",.5);
41  thresholdPhibMean = ps.getUntrackedParameter<double>("thresholdPhibMean",1.5);
42  thresholdPhibRMS = ps.getUntrackedParameter<double>("thresholdPhibRMS",.8);
43  doCorrStudy = ps.getUntrackedParameter<bool>("doCorrelationStudy",false);
44 
45 
46 }
T getUntrackedParameter(std::string const &, T const &) const
void setConfig(const edm::ParameterSet &ps, std::string name)
Set configuration variables.
DTLocalTriggerLutTest::~DTLocalTriggerLutTest ( )
virtual

Destructor.

Definition at line 49 of file DTLocalTriggerLutTest.cc.

49  {
50 
51 }

Member Function Documentation

void DTLocalTriggerLutTest::beginJob ( void  )
protectedvirtual

BeginJob.

Reimplemented from DTLocalTriggerBaseTest.

Definition at line 54 of file DTLocalTriggerLutTest.cc.

References DTLocalTriggerBaseTest::beginJob(), and Parameters::parameters.

54  {
55 
57 
58  vector<string>::const_iterator iTr = trigSources.begin();
59  vector<string>::const_iterator trEnd = trigSources.end();
60  vector<string>::const_iterator iHw = hwSources.begin();
61  vector<string>::const_iterator hwEnd = hwSources.end();
62 
63  //Booking
64  if(parameters.getUntrackedParameter<bool>("staticBooking", true)){
65  for (; iTr != trEnd; ++iTr){
66  trigSource = (*iTr);
67  for (; iHw != hwEnd; ++iHw){
68  hwSource = (*iHw);
69  // Loop over the TriggerUnits
70  for (int wh=-2; wh<=2; ++wh){
71  bookWheelHistos(wh,"PhiResidualMean");
72  bookWheelHistos(wh,"PhiResidualRMS");
73  bookWheelHistos(wh,"PhibResidualMean");
74  bookWheelHistos(wh,"PhibResidualRMS");
75  if (doCorrStudy) {
76  bookWheelHistos(wh,"PhiTkvsTrigSlope");
77  bookWheelHistos(wh,"PhiTkvsTrigIntercept");
78  bookWheelHistos(wh,"PhiTkvsTrigCorr");
79  bookWheelHistos(wh,"PhibTkvsTrigSlope");
80  bookWheelHistos(wh,"PhibTkvsTrigIntercept");
81  bookWheelHistos(wh,"PhibTkvsTrigCorr");
82  }
83  }
84  }
85  }
86  }
87 
88  // Summary test histo booking (only static)
89  for (iTr = trigSources.begin(); iTr != trEnd; ++iTr){
90  trigSource = (*iTr);
91  for (iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw){
92  hwSource = (*iHw);
93  // Loop over the TriggerUnits
94  for (int wh=-2; wh<=2; ++wh){
95  bookWheelHistos(wh,"PhiLutSummary");
96  bookWheelHistos(wh,"PhibLutSummary");
97  }
98  bookCmsHistos("PhiLutSummary");
99  bookCmsHistos("PhibLutSummary");
100  }
101  }
102 
103 }
T getUntrackedParameter(std::string const &, T const &) const
std::vector< std::string > trigSources
void bookCmsHistos(std::string hTag, std::string folder="", bool isGlb=false)
Book the new MEs (CMS summary)
void bookWheelHistos(int wheel, std::string hTag, std::string folder="")
Book the new MEs (for each wheel)
std::vector< std::string > hwSources
void DTLocalTriggerLutTest::beginRun ( const edm::Run r,
const edm::EventSetup c 
)
protectedvirtual

BeginRun.

Reimplemented from DTLocalTriggerBaseTest.

Definition at line 106 of file DTLocalTriggerLutTest.cc.

References DTLocalTriggerBaseTest::beginRun().

106  {
107 
109 
110 }
void beginRun(edm::Run const &run, edm::EventSetup const &context)
BeginRun.
void DTLocalTriggerLutTest::fillWhPlot ( MonitorElement plot,
int  sect,
int  stat,
float  value,
bool  lessIsBest = true 
)
private

Fill summary plots managing double MB4 chambers.

Definition at line 331 of file DTLocalTriggerLutTest.cc.

References MonitorElement::getBinContent(), and MonitorElement::setBinContent().

331  {
332 
333  if (sect>12) {
334  int scsect = sect==13 ? 4 : 10;
335  if ( (fabs(value)>fabs(plot->getBinContent(scsect,stat)))==lessIsBest) {
336  plot->setBinContent(scsect,stat,value);
337  }
338  }
339  else {
340  plot->setBinContent(sect,stat,value);
341  }
342 
343  return;
344 
345 }
void setBinContent(int binx, double content)
set content of bin (1-D)
double getBinContent(int binx) const
get content of bin (1-D)
int DTLocalTriggerLutTest::performLutTest ( double  mean,
double  RMS,
double  thresholdMean,
double  thresholdRMS 
)
private

Perform Lut Test logical operations.

Definition at line 322 of file DTLocalTriggerLutTest.cc.

322  {
323 
324  bool meanErr = fabs(mean)>thresholdMean;
325  bool rmsErr = RMS>thresholdRMS;
326 
327  return (meanErr || rmsErr) ? 2+(meanErr!=rmsErr) : 0 ;
328 
329 }
void DTLocalTriggerLutTest::runClientDiagnostic ( )
protectedvirtual

Run client analysis.

Implements DTLocalTriggerBaseTest.

Definition at line 113 of file DTLocalTriggerLutTest.cc.

References python.rootplot.argparse::category, spr::find(), newFWLiteAna::fullName, edm::second(), and DTChamberId::wheel().

113  {
114 
115  // Loop over Trig & Hw sources
116  for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr){
117  trigSource = (*iTr);
118  for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw){
119  hwSource = (*iHw);
120  vector<DTChamber*>::const_iterator chIt = muonGeom->chambers().begin();
121  vector<DTChamber*>::const_iterator chEnd = muonGeom->chambers().end();
122  for (; chIt != chEnd; ++chIt) {
123  DTChamberId chId((*chIt)->id());
124  int wh = chId.wheel();
125  int sect = chId.sector();
126  int stat = chId.station();
127 
128 
129  if (doCorrStudy) {
130  // Perform Correlation Plots analysis (DCC + segment Phi)
131  TH2F * TrackPhitkvsPhitrig = getHisto<TH2F>(dbe->get(getMEName("PhitkvsPhitrig","Segment", chId)));
132 
133  if (TrackPhitkvsPhitrig && TrackPhitkvsPhitrig->GetEntries()>10) {
134 
135  // Fill client histos
136  if( whME[wh].find(fullName("PhiTkvsTrigCorr")) == whME[wh].end() ){
137  bookWheelHistos(wh,"PhiTkvsTrigSlope");
138  bookWheelHistos(wh,"PhiTkvsTrigIntercept");
139  bookWheelHistos(wh,"PhiTkvsTrigCorr");
140  }
141 
142  TProfile* PhitkvsPhitrigProf = TrackPhitkvsPhitrig->ProfileX();
143  double phiInt = 0;
144  double phiSlope = 0;
145  double phiCorr = 0;
146  try {
147  PhitkvsPhitrigProf->Fit("pol1","CQO");
148  TF1 *ffPhi= PhitkvsPhitrigProf->GetFunction("pol1");
149  if (ffPhi) {
150  phiInt = ffPhi->GetParameter(0);
151  phiSlope = ffPhi->GetParameter(1);
152  phiCorr = TrackPhitkvsPhitrig->GetCorrelationFactor();
153  }
154  } catch (cms::Exception& iException) {
155  edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhitkvsPhitrig for Wheel " << wh
156  <<" Sector " << sect << " Station " << stat;
157  }
158 
159  std::map<std::string,MonitorElement*> &innerME = whME[wh];
160  fillWhPlot(innerME.find(fullName("PhiTkvsTrigSlope"))->second,sect,stat,phiSlope-1);
161  fillWhPlot(innerME.find(fullName("PhiTkvsTrigIntercept"))->second,sect,stat,phiInt);
162  fillWhPlot(innerME.find(fullName("PhiTkvsTrigCorr"))->second,sect,stat,phiCorr,false);
163 
164  }
165 
166  // Perform Correlation Plots analysis (DCC + segment Phib)
167  TH2F * TrackPhibtkvsPhibtrig = getHisto<TH2F>(dbe->get(getMEName("PhibtkvsPhibtrig","Segment", chId)));
168 
169  if (stat != 3 && TrackPhibtkvsPhibtrig && TrackPhibtkvsPhibtrig->GetEntries()>10) {// station 3 has no meaningful MB3 phi bending information
170 
171  // Fill client histos
172  if( whME[wh].find(fullName("PhibTkvsTrigCorr")) == whME[wh].end() ){
173  bookWheelHistos(wh,"PhibTkvsTrigSlope");
174  bookWheelHistos(wh,"PhibTkvsTrigIntercept");
175  bookWheelHistos(wh,"PhibTkvsTrigCorr");
176  }
177 
178  TProfile* PhibtkvsPhibtrigProf = TrackPhibtkvsPhibtrig->ProfileX();
179  double phibInt = 0;
180  double phibSlope = 0;
181  double phibCorr = 0;
182  try {
183  PhibtkvsPhibtrigProf->Fit("pol1","CQO");
184  TF1 *ffPhib= PhibtkvsPhibtrigProf->GetFunction("pol1");
185  if (ffPhib) {
186  phibInt = ffPhib->GetParameter(0);
187  phibSlope = ffPhib->GetParameter(1);
188  phibCorr = TrackPhibtkvsPhibtrig->GetCorrelationFactor();
189  }
190  } catch (cms::Exception& iException) {
191  edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhibtkvsPhibtrig for Wheel " << wh
192  <<" Sector " << sect << " Station " << stat;
193  }
194 
195  std::map<std::string,MonitorElement*> &innerME = whME[wh];
196  fillWhPlot(innerME.find(fullName("PhibTkvsTrigSlope"))->second,sect,stat,phibSlope-1);
197  fillWhPlot(innerME.find(fullName("PhibTkvsTrigIntercept"))->second,sect,stat,phibInt);
198  fillWhPlot(innerME.find(fullName("PhibTkvsTrigCorr"))->second,sect,stat,phibCorr,false);
199 
200  }
201 
202  }
203 
204  // Make Phi Residual Summary
205  TH1F * PhiResidual = getHisto<TH1F>(dbe->get(getMEName("PhiResidual","Segment", chId)));
206  int phiSummary = 1;
207 
208  if (PhiResidual && PhiResidual->GetEffectiveEntries()>10) {
209 
210  // Fill client histos
211  if( whME[wh].find(fullName("PhiResidualMean")) == whME[wh].end() ){
212  bookWheelHistos(wh,"PhiResidualMean");
213  bookWheelHistos(wh,"PhiResidualRMS");
214  }
215 
216  double peak = PhiResidual->GetBinCenter(PhiResidual->GetMaximumBin());
217  double phiMean = 0;
218  double phiRMS = 0;
219  try {
220  PhiResidual->Fit("gaus","CQO","",peak-5,peak+5);
221  TF1 *ffPhi = PhiResidual->GetFunction("gaus");
222  if ( ffPhi ) {
223  phiMean = ffPhi->GetParameter(1);
224  phiRMS = ffPhi->GetParameter(2);
225  }
226  } catch (cms::Exception& iException) {
227  edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhiResidual for Wheel " << wh
228  <<" Sector " << sect << " Station " << stat;
229  }
230 
231  std::map<std::string,MonitorElement*> &innerME = whME[wh];
232  fillWhPlot(innerME.find(fullName("PhiResidualMean"))->second,sect,stat,phiMean);
233  fillWhPlot(innerME.find(fullName("PhiResidualRMS"))->second,sect,stat,phiRMS);
234 
235  phiSummary = performLutTest(phiMean,phiRMS,thresholdPhiMean,thresholdPhiRMS);
236 
237  }
238  fillWhPlot(whME[wh].find(fullName("PhiLutSummary"))->second,sect,stat,phiSummary);
239 
240  // Make Phib Residual Summary
241  TH1F * PhibResidual = getHisto<TH1F>(dbe->get(getMEName("PhibResidual","Segment", chId)));
242  int phibSummary = stat==3 ? 0 : 1; // station 3 has no meaningful MB3 phi bending information
243 
244  if (stat != 3 && PhibResidual && PhibResidual->GetEffectiveEntries()>10) {// station 3 has no meaningful MB3 phi bending information
245 
246  // Fill client histos
247  if( whME[wh].find(fullName("PhibResidualMean")) == whME[wh].end() ){
248  bookWheelHistos(wh,"PhibResidualMean");
249  bookWheelHistos(wh,"PhibResidualRMS");
250  }
251 
252  double peak = PhibResidual->GetBinCenter(PhibResidual->GetMaximumBin());
253  double phibMean = 0;
254  double phibRMS = 0;
255  try {
256  PhibResidual->Fit("gaus","CQO","",peak-5,peak+5);
257  TF1 *ffPhib = PhibResidual->GetFunction("gaus");
258  if ( ffPhib ) {
259  phibMean = ffPhib->GetParameter(1);
260  phibRMS = ffPhib->GetParameter(2);
261  }
262  } catch (cms::Exception& iException) {
263  edm::LogError(category()) << "[" << testName << "Test]: Error fitting PhibResidual for Wheel " << wh
264  <<" Sector " << sect << " Station " << stat;
265  }
266 
267  std::map<std::string,MonitorElement*> &innerME = whME[wh];
268  fillWhPlot(innerME.find(fullName("PhibResidualMean"))->second,sect,stat,phibMean);
269  fillWhPlot(innerME.find(fullName("PhibResidualRMS"))->second,sect,stat,phibRMS);
270 
271  phibSummary = performLutTest(phibMean,phibRMS,thresholdPhibMean,thresholdPhibRMS);
272 
273  }
274  fillWhPlot(whME[wh].find(fullName("PhibLutSummary"))->second,sect,stat,phibSummary);
275 
276  }
277  }
278  }
279 
280  // Barrel Summary Plots
281  for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr){
282  trigSource = (*iTr);
283  for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw){
284  hwSource = (*iHw);
285  for (int wh=-2; wh<=2; ++wh){
286  std::map<std::string,MonitorElement*> *innerME = &(whME[wh]);
287 
288  TH2F* phiWhSummary = getHisto<TH2F>(innerME->find(fullName("PhiLutSummary"))->second);
289  TH2F* phibWhSummary = getHisto<TH2F>(innerME->find(fullName("PhibLutSummary"))->second);
290  for (int sect=1; sect<=12; ++sect){
291  int phiErr = 0;
292  int phibErr = 0;
293  int phiNoData = 0;
294  int phibNoData = 0;
295  for (int stat=1; stat<=4; ++stat){
296  switch (static_cast<int>(phiWhSummary->GetBinContent(sect,stat))) {
297  case 1:
298  phiNoData++;
299  case 2:
300  case 3:
301  phiErr++;
302  }
303  switch (static_cast<int>(phibWhSummary->GetBinContent(sect,stat))) {
304  case 1:
305  phibNoData++;
306  case 2:
307  case 3:
308  phibErr++;
309  }
310  }
311  if (phiNoData == 4) phiErr = 5;
312  if (phibNoData == 3) phibErr = 5; // MB3 has no phib information
313  cmsME.find(fullName("PhiLutSummary"))->second->setBinContent(sect,wh+3,phiErr);
314  cmsME.find(fullName("PhibLutSummary"))->second->setBinContent(sect,wh+3,phibErr);
315  }
316  }
317  }
318  }
319 
320 }
std::vector< std::string > trigSources
std::map< int, std::map< std::string, MonitorElement * > > whME
int performLutTest(double mean, double RMS, double thresholdMean, double thresholdRMS)
Perform Lut Test logical operations.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
edm::ESHandle< DTGeometry > muonGeom
void fillWhPlot(MonitorElement *plot, int sect, int stat, float value, bool lessIsBest=true)
Fill summary plots managing double MB4 chambers.
U second(std::pair< T, U > const &p)
void bookWheelHistos(int wheel, std::string hTag, std::string folder="")
Book the new MEs (for each wheel)
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1623
std::string category()
Get message logger name.
std::string getMEName(std::string histoTag, std::string subfolder, const DTChamberId &chambid)
Get the ME name (by chamber)
std::map< std::string, MonitorElement * > cmsME
std::vector< std::string > hwSources
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
std::string fullName(std::string htype)
Create fullname from histo partial name.

Member Data Documentation

bool DTLocalTriggerLutTest::doCorrStudy
private

Definition at line 50 of file DTLocalTriggerLutTest.h.

double DTLocalTriggerLutTest::thresholdPhibMean
private

Definition at line 48 of file DTLocalTriggerLutTest.h.

double DTLocalTriggerLutTest::thresholdPhibRMS
private

Definition at line 49 of file DTLocalTriggerLutTest.h.

double DTLocalTriggerLutTest::thresholdPhiMean
private

Definition at line 48 of file DTLocalTriggerLutTest.h.

double DTLocalTriggerLutTest::thresholdPhiRMS
private

Definition at line 49 of file DTLocalTriggerLutTest.h.