CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
dtCalibration::DTT0FEBPathCorrection Class Reference

#include <DTT0FEBPathCorrection.h>

Inheritance diagram for dtCalibration::DTT0FEBPathCorrection:
dtCalibration::DTT0BaseCorrection

Public Member Functions

DTT0Data correction (const DTWireId &) override
 
 DTT0FEBPathCorrection (const edm::ParameterSet &)
 
void setES (const edm::EventSetup &setup) override
 
float t0FEBPathCorrection (int wheel, int st, int sec, int sl, int l, int w)
 
 ~DTT0FEBPathCorrection () override
 
- Public Member Functions inherited from dtCalibration::DTT0BaseCorrection
 DTT0BaseCorrection ()
 
virtual ~DTT0BaseCorrection ()
 

Private Member Functions

DTT0Data defaultT0 (const DTWireId &)
 

Private Attributes

std::string calibChamber_
 
DTChamberId chosenChamberId_
 
const DTT0t0Map_
 

Detailed Description

Definition at line 25 of file DTT0FEBPathCorrection.h.

Constructor & Destructor Documentation

DTT0FEBPathCorrection::DTT0FEBPathCorrection ( const edm::ParameterSet pset)

Definition at line 27 of file DTT0FEBPathCorrection.cc.

References calibChamber_, and chosenChamberId_.

27  :
28  calibChamber_( pset.getParameter<string>("calibChamber") ) {
29 
30  //DTChamberId chosenChamberId;
31  if( calibChamber_ != "" && calibChamber_ != "None" && calibChamber_ != "All" ){
32  stringstream linestr;
33  int selWheel, selStation, selSector;
34  linestr << calibChamber_;
35  linestr >> selWheel >> selStation >> selSector;
36  chosenChamberId_ = DTChamberId(selWheel, selStation, selSector);
37  LogVerbatim("Calibration") << "[DTT0FEBPathCorrection] Chosen chamber: " << chosenChamberId_ << endl;
38  }
39  //FIXME: Check if chosen chamber is valid.
40 }
T getParameter(std::string const &) const
DTT0FEBPathCorrection::~DTT0FEBPathCorrection ( )
override

Definition at line 42 of file DTT0FEBPathCorrection.cc.

42  {
43 }

Member Function Documentation

DTT0Data DTT0FEBPathCorrection::correction ( const DTWireId wireId)
overridevirtual

Implements dtCalibration::DTT0BaseCorrection.

Definition at line 53 of file DTT0FEBPathCorrection.cc.

References calibChamber_, DTSuperLayerId::chamberId(), chosenChamberId_, DTTimeUnits::counts, defaultT0(), Exception, DTT0::get(), checklumidiff::l, DTLayerId::layer(), DTWireId::layerId(), DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, mps_update::status, DTSuperLayerId::superlayer(), DTLayerId::superlayerId(), t0FEBPathCorrection(), t0Map_, DTChamberId::wheel(), makeMuonMisalignmentScenario::wheel, and DTWireId::wire().

53  {
54  // Compute for selected chamber (or All) correction using as reference chamber mean
55 
56  DTChamberId chamberId = wireId.layerId().superlayerId().chamberId();
57 
58  if( calibChamber_ == "" || calibChamber_ == "None" ) return defaultT0(wireId);
59  if( calibChamber_ != "All" && chamberId != chosenChamberId_ ) return defaultT0(wireId);
60 
61  // Access DB
62  float t0Mean,t0RMS;
63  int status = t0Map_->get(wireId,t0Mean,t0RMS,DTTimeUnits::counts);
64  if(status != 0)
65  throw cms::Exception("[DTT0FEBPathCorrection]") << "Could not find t0 entry in DB for"
66  << wireId << endl;
67  int wheel = chamberId.wheel();
68  int station = chamberId.station();
69  int sector = chamberId.sector();
70  int sl = wireId.layerId().superlayerId().superlayer();
71  int l = wireId.layerId().layer();
72  int wire = wireId.wire();
73  float t0MeanNew = t0Mean - t0FEBPathCorrection(wheel, station, sector, sl, l, wire);
74  float t0RMSNew = t0RMS;
75  return DTT0Data(t0MeanNew,t0RMSNew);
76 }
DTChamberId chamberId() const
Return the corresponding ChamberId.
int layer() const
Return the layer number.
Definition: DTLayerId.h:53
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
int get(int wheelId, int stationId, int sectorId, int slId, int layerId, int cellId, float &t0mean, float &t0rms, DTTimeUnits::type unit) const
Definition: DTT0.cc:67
float t0FEBPathCorrection(int wheel, int st, int sec, int sl, int l, int w)
int wire() const
Return the wire number.
Definition: DTWireId.h:56
int superlayer() const
Return the superlayer number (deprecated method name)
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:62
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
Definition: DTT0.h:37
DTT0Data DTT0FEBPathCorrection::defaultT0 ( const DTWireId wireId)
private

Definition at line 78 of file DTT0FEBPathCorrection.cc.

References DTTimeUnits::counts, Exception, DTT0::get(), mps_update::status, and t0Map_.

Referenced by correction().

78  {
79  // Access default DB
80  float t0Mean,t0RMS;
81  int status = t0Map_->get(wireId,t0Mean,t0RMS,DTTimeUnits::counts);
82  if(!status){
83  return DTT0Data(t0Mean,t0RMS);
84  } else{
85  //...
86  throw cms::Exception("[DTT0FEBPathCorrection]") << "Could not find t0 entry in DB for"
87  << wireId << endl;
88  }
89 }
int get(int wheelId, int stationId, int sectorId, int slId, int layerId, int cellId, float &t0mean, float &t0rms, DTTimeUnits::type unit) const
Definition: DTT0.cc:67
Definition: DTT0.h:37
void DTT0FEBPathCorrection::setES ( const edm::EventSetup setup)
overridevirtual

Implements dtCalibration::DTT0BaseCorrection.

Definition at line 45 of file DTT0FEBPathCorrection.cc.

References edm::EventSetup::get(), t0Map_, and DTT0::version().

45  {
46  // Get t0 record from DB
47  ESHandle<DTT0> t0H;
48  setup.get<DTT0Rcd>().get(t0H);
49  t0Map_ = &*t0H;
50  LogVerbatim("Calibration") << "[DTT0FEBPathCorrection] T0 version: " << t0H->version();
51 }
Definition: DTT0Rcd.h:9
T get() const
Definition: EventSetup.h:71
const std::string & version() const
access version
Definition: DTT0.cc:118
float DTT0FEBPathCorrection::t0FEBPathCorrection ( int  wheel,
int  st,
int  sec,
int  sl,
int  l,
int  w 
)

Return the value to be subtracted to t0s to correct for the difference of path lenght for TP signals within the FEB, from the TP input connector to the MAD.

FEBs are alternately connected on the right (J9) and left (J8) TP input connectors. Path lenghts (provided from Franco Gonella) for FEB-16 channels for FEB-16s connected on the right TP connector (J9) are:

CH0,1 = 32 mm CH2,3 = 42 mm CH4,5 = 52 mm CH6,7 = 62 mm

CH8,9 = 111 mm CH10,11 = 121 mm CH12,13 = 131 mm CH14,15 = 141 mm

Given that ttrig calibration absorbs average offsets, we assume thate only differences w.r.t. the average lenght (86.5 mm) remain.

For FEBs connected on the right connector, values are swapped; so there is a periodicity of 2 FEBS (8 cells)

The mapping of FEB channels to wires, for the first two FEBs, is:

   FEB 0 (J9)     FEB 1 (J8)     

Wire | 1 2 3 4 | 5 6 7 8

L1 | 3 7 11 15 | 3 7 11 15 L2 | 1 5 9 13 | 1 5 9 13 L3 | 2 6 10 14 | 2 6 10 14 L4 | 0 4 8 12 | 0 4 8 12

For FEB-20, distances from the left connector (J8) are: CH16,17 = 171 mm CH18,19 = 181 mm

We do not include a correction for this additional row of channels since they are at the edge of the SL so the effect cannot be seen on data (and moreover the channel in L1 is usually not existing in the chamber)

Definition at line 139 of file DTT0FEBPathCorrection.cc.

Referenced by correction().

139  {
140 
141  // Skip correction for the last row of cells of FEB20 (see above)
142  if( (st==1 && ((sl!=2 && w ==49) || (sl==2 && w ==57))) ||
143  ((st==2||st==3)&& (sl==2 && w ==57)) ) return 0.;
144 
145 
146  float dist[8]={};
147 
148  // Path lenght differences for L1 and L3 (cm)
149  if(l==1 || l ==3){
150 
151  dist[0] = +4.45;
152  dist[1] = +2.45;
153  dist[2] = -3.45;
154  dist[3] = -5.45;
155  dist[4] = -4.45;
156  dist[5] = -2.45;
157  dist[6] = +3.45;
158  dist[7] = +5.45;
159  }
160 
161  // Path lenght differences for L2 and L4 (cm)
162  else {
163 
164  dist[0] = +5.45;
165  dist[1] = +3.45;
166  dist[2] = -2.45;
167  dist[3] = -4.45;
168  dist[4] = -5.45;
169  dist[5] = -3.45;
170  dist[6] = +2.45;
171  dist[7] = +4.45;
172  }
173 
174 
175  // Wire position within the 8-cell period (2 FEBs).
176  // Note that wire numbers start from 1.
177  int pos = (w-1)%8;
178 
179  // Special case: in MB2 phi and MB4-10, the periodicity is broken at cell 49, as there is
180  // an odd number of FEDs (15): the 13th FEB is connected on the left,
181  // and not on the right; i.e. one FEB (4 colums of cells) is skipped from what would
182  // be the regular structure.
183  // The same happens in MB4-8/12 at cell 81.
184  if ((st==2 && sl!=2 && w>=49) ||
185  (st==4 && sec==10 && w>=49) ||
186  (st==4 && (sec==8||sec==12) && w>=81) ) pos =(w-1+4)%8;
187 
188  // Inverse of the signal propagation speed, determined from the
189  // observed amplitude of the modulation. This matches what is found
190  // with CAD simulation using reasonable assumptions on the PCB specs.
191 
192  return dist[pos]*0.075;
193 
194 }
const double w
Definition: UKUtility.cc:23

Member Data Documentation

std::string dtCalibration::DTT0FEBPathCorrection::calibChamber_
private

Definition at line 40 of file DTT0FEBPathCorrection.h.

Referenced by correction(), and DTT0FEBPathCorrection().

DTChamberId dtCalibration::DTT0FEBPathCorrection::chosenChamberId_
private

Definition at line 42 of file DTT0FEBPathCorrection.h.

Referenced by correction(), and DTT0FEBPathCorrection().

const DTT0* dtCalibration::DTT0FEBPathCorrection::t0Map_
private

Definition at line 43 of file DTT0FEBPathCorrection.h.

Referenced by correction(), defaultT0(), and setES().