CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTRunConditionVar Class Reference

#include <DTRunConditionVar.h>

Inheritance diagram for DTRunConditionVar:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 
 DTRunConditionVar (const edm::ParameterSet &pset)
 
 ~DTRunConditionVar () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Private Member Functions

void bookChamberHistos (DQMStore::IBooker &, const DTChamberId &dtCh, std::string histoType, int, float, float)
 

Private Attributes

std::map< uint32_t, std::map< std::string, MonitorElement * > > chamberHistos
 
bool debug
 
edm::EDGetTokenT< DTRecSegment4DCollectiondt4DSegmentsToken_
 
edm::ESHandle< DTGeometrydtGeom
 
double maxAnglePhiSegm
 
edm::ESHandle< DTMtimemTime
 
const DTMtimemTimeMap_
 
int nMinHitsPhi
 

Detailed Description

Description:

Author
: Paolo Bellan, Antonio Branca $date : 23/09/2011 15:42:04 CET $

Modification:

Definition at line 45 of file DTRunConditionVar.h.

Constructor & Destructor Documentation

DTRunConditionVar::DTRunConditionVar ( const edm::ParameterSet pset)

Definition at line 46 of file DTRunConditionVar.cc.

46  :
47  // Get the debug parameter for verbose output
48  debug(pSet.getUntrackedParameter<bool>("debug",false)),
49  nMinHitsPhi(pSet.getUntrackedParameter<int>("nMinHitsPhi")),
50  maxAnglePhiSegm(pSet.getUntrackedParameter<double>("maxAnglePhiSegm")),
51  dt4DSegmentsToken_(consumes<DTRecSegment4DCollection>(
52  pSet.getParameter<InputTag>("recoSegments")))
53 {
54 
55 }
edm::EDGetTokenT< DTRecSegment4DCollection > dt4DSegmentsToken_
DTRunConditionVar::~DTRunConditionVar ( )
override

Definition at line 57 of file DTRunConditionVar.cc.

References LogTrace.

58 {
59  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
60  << "DTRunConditionVar: destructor called";
61 
62  // free memory
63 }
#define LogTrace(id)

Member Function Documentation

void DTRunConditionVar::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
override

Definition at line 95 of file DTRunConditionVar.cc.

References chamberHistos, DTVelocityUnits::cm_per_ns, dt4DSegmentsToken_, dtGeom, Exception, HcalObjRepresent::Fill(), edm::EventSetup::get(), DTMtime::get(), LogTrace, maxAnglePhiSegm, mTime, mTimeMap_, nMinHitsPhi, Pi, DetId::rawId(), and xdir.

97 {
98 
99  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar") <<
100  "--- [DTRunConditionVar] Event analysed #Run: " <<
101  event.id().run() << " #Event: " << event.id().event() << endl;
102 
103  // Get the DT Geometry
105  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
106 
107  // Get the map of vdrift from the setup
108  eventSetup.get<DTMtimeRcd>().get(mTime);
109  mTimeMap_ = &*mTime;
110 
111  // Get the segment collection from the event
112  Handle<DTRecSegment4DCollection> all4DSegments;
113  event.getByToken(dt4DSegmentsToken_, all4DSegments);
114 
115  // Loop over the segments
116  for(DTRecSegment4DCollection::const_iterator segment = all4DSegments->begin();
117  segment != all4DSegments->end(); ++segment){
118 
119  // Get the chamber from the setup
120  DTChamberId DTid = (DTChamberId) segment->chamberId();
121  uint32_t indexCh = DTid.rawId();
122 
123  // Fill v-drift values
124  if( (*segment).hasPhi() ) {
125 
126  int nHitsPhi = (*segment).phiSegment()->degreesOfFreedom()+2;
127  double xdir = (*segment).phiSegment()->localDirection().x();
128  double zdir = (*segment).phiSegment()->localDirection().z();
129 
130  double anglePhiSegm = fabs(atan(xdir/zdir))*180./TMath::Pi();
131 
132  if( nHitsPhi >= nMinHitsPhi && anglePhiSegm <= maxAnglePhiSegm ) {
133 
134  double segmentVDrift = segment->phiSegment()->vDrift();
135 
136 
137  DTSuperLayerId indexSLPhi1(DTid,1);
138  DTSuperLayerId indexSLPhi2(DTid,3);
139 
140  float vDriftPhi1(0.), vDriftPhi2(0.);
141  float ResPhi1(0.), ResPhi2(0.);
142  int status1 = mTimeMap_->get(indexSLPhi1,vDriftPhi1,ResPhi1,DTVelocityUnits::cm_per_ns);
143  int status2 = mTimeMap_->get(indexSLPhi2,vDriftPhi2,ResPhi2,DTVelocityUnits::cm_per_ns);
144 
145  if(status1 != 0 || status2 != 0) {
146  DTSuperLayerId sl = (status1 != 0) ? indexSLPhi1 : indexSLPhi2;
147  throw cms::Exception("DTRunConditionVarClient") << "Could not find vDrift entry in DB for"
148  << sl << endl;
149  }
150 
151  float vDriftMed = (vDriftPhi1 + vDriftPhi2) / 2.;
152 
153  segmentVDrift = vDriftMed*(1. - segmentVDrift);
154 
155  double segmentT0 = segment->phiSegment()->t0();
156 
157  if( segment->phiSegment()->ist0Valid() ) (chamberHistos[indexCh])["T0_FromSegm"]->Fill(segmentT0);
158  if( segmentVDrift != vDriftMed ) (chamberHistos[indexCh])["VDrift_FromSegm"]->Fill(segmentVDrift);
159 
160  }
161  }
162 
163  } //end loop on segment
164 
165 } //end analyze
const double Pi
edm::ESHandle< DTMtime > mTime
edm::EDGetTokenT< DTRecSegment4DCollection > dt4DSegmentsToken_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
const DTMtime * mTimeMap_
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
#define LogTrace(id)
std::map< uint32_t, std::map< std::string, MonitorElement * > > chamberHistos
int get(int wheelId, int stationId, int sectorId, int slId, float &mTime, float &mTrms, DTTimeUnits::type unit) const
Definition: DTMtime.cc:82
T get() const
Definition: EventSetup.h:62
edm::ESHandle< DTGeometry > dtGeom
void DTRunConditionVar::bookChamberHistos ( DQMStore::IBooker ,
const DTChamberId dtCh,
std::string  histoType,
int  ,
float  ,
float   
)
private

Definition at line 167 of file DTRunConditionVar.cc.

References DQMStore::IBooker::book1D(), chamberHistos, fftjetimagerecorder_cfi::histoLabel, LogTrace, DetId::rawId(), SimDataFormats::CaloAnalysis::sc, DTChamberId::sector(), DQMStore::IBooker::setCurrentFolder(), DTChamberId::station(), relativeConstraints::station, DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by bookHistograms().

167  {
168 
169  int wh = dtCh.wheel();
170  int sc = dtCh.sector();
171  int st = dtCh.station();
172  stringstream wheel; wheel << wh;
173  stringstream station; station << st;
174  stringstream sector; sector << sc;
175 
176  string bookingFolder = "DT/02-Segments/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str();
177  string histoTag = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
178 
179  ibooker.setCurrentFolder(bookingFolder);
180 
181  LogTrace ("DTDQM|DTMonitorModule|DTRunConditionVar")
182  << "[DTRunConditionVar]: booking histos in " << bookingFolder << endl;
183 
184  string histoName = histoType + histoTag;
185  string histoLabel = histoType;
186 
187  (chamberHistos[dtCh.rawId()])[histoType] =
188  ibooker.book1D(histoName,histoLabel,nbins,min,max);
189 
190 }
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
T min(T a, T b)
Definition: MathUtil.h:58
#define LogTrace(id)
std::map< uint32_t, std::map< std::string, MonitorElement * > > chamberHistos
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void DTRunConditionVar::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &   
)
override

Definition at line 65 of file DTRunConditionVar.cc.

References bookChamberHistos(), LogTrace, and makeMuonMisalignmentScenario::wheel.

67  {
68 
69  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
70  << "DTRunConditionVar: bookHistograms";
71 
72  for(int wheel=-2;wheel<=2;wheel++){
73  for(int sec=1; sec<=14; sec++) {
74  for(int stat=1; stat<=4; stat++) {
75 
76  bookChamberHistos(ibooker,DTChamberId(wheel,stat,sec),"VDrift_FromSegm",100,0.0043,0.0065);
77  bookChamberHistos(ibooker,DTChamberId(wheel,stat,sec),"T0_FromSegm",100,-25.,25.);
78 
79  }
80  }
81  }
82 
83 
84  return;
85 }
#define LogTrace(id)
void bookChamberHistos(DQMStore::IBooker &, const DTChamberId &dtCh, std::string histoType, int, float, float)
void DTRunConditionVar::dqmBeginRun ( const edm::Run run,
const edm::EventSetup setup 
)
override

Definition at line 87 of file DTRunConditionVar.cc.

References dtGeom, and edm::EventSetup::get().

88 {
89  // Get the DT Geometry
90  setup.get<MuonGeometryRecord>().get(dtGeom);
91 
92  return;
93 }
T get() const
Definition: EventSetup.h:62
edm::ESHandle< DTGeometry > dtGeom

Member Data Documentation

std::map<uint32_t, std::map<std::string, MonitorElement*> > DTRunConditionVar::chamberHistos
private

Definition at line 77 of file DTRunConditionVar.h.

Referenced by analyze(), and bookChamberHistos().

bool DTRunConditionVar::debug
private
edm::EDGetTokenT<DTRecSegment4DCollection> DTRunConditionVar::dt4DSegmentsToken_
private

Definition at line 70 of file DTRunConditionVar.h.

Referenced by analyze().

edm::ESHandle<DTGeometry> DTRunConditionVar::dtGeom
private

Definition at line 72 of file DTRunConditionVar.h.

Referenced by analyze(), and dqmBeginRun().

double DTRunConditionVar::maxAnglePhiSegm
private

Definition at line 68 of file DTRunConditionVar.h.

Referenced by analyze().

edm::ESHandle<DTMtime> DTRunConditionVar::mTime
private

Definition at line 74 of file DTRunConditionVar.h.

Referenced by analyze().

const DTMtime* DTRunConditionVar::mTimeMap_
private

Definition at line 75 of file DTRunConditionVar.h.

Referenced by analyze().

int DTRunConditionVar::nMinHitsPhi
private

Definition at line 67 of file DTRunConditionVar.h.

Referenced by analyze().