CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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

virtual DTT0Data correction (const DTWireId &)
 
 DTT0FEBPathCorrection (const edm::ParameterSet &)
 
virtual void setES (const edm::EventSetup &setup)
 
float t0FEBPathCorrection (int wheel, int st, int sec, int sl, int l, int w)
 
virtual ~DTT0FEBPathCorrection ()
 
- 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 31 of file DTT0FEBPathCorrection.cc.

References calibChamber_, chosenChamberId_, and DTChamberId.

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

Definition at line 46 of file DTT0FEBPathCorrection.cc.

46  {
47 }

Member Function Documentation

DTT0Data DTT0FEBPathCorrection::correction ( const DTWireId wireId)
virtual

Implements dtCalibration::DTT0BaseCorrection.

Definition at line 57 of file DTT0FEBPathCorrection.cc.

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

57  {
58  // Compute for selected chamber (or All) correction using as reference chamber mean
59 
60  DTChamberId chamberId = wireId.layerId().superlayerId().chamberId();
61 
62  if( calibChamber_ == "" || calibChamber_ == "None" ) return defaultT0(wireId);
63  if( calibChamber_ != "All" && chamberId != chosenChamberId_ ) return defaultT0(wireId);
64 
65  // Access DB
66  float t0Mean,t0RMS;
67  int status = t0Map_->get(wireId,t0Mean,t0RMS,DTTimeUnits::counts);
68  if(status != 0)
69  throw cms::Exception("[DTT0FEBPathCorrection]") << "Could not find t0 entry in DB for"
70  << wireId << endl;
71  int wheel = chamberId.wheel();
72  int station = chamberId.station();
73  int sector = chamberId.sector();
74  int sl = wireId.layerId().superlayerId().superlayer();
75  int l = wireId.layerId().layer();
76  int wire = wireId.wire();
77  float t0MeanNew = t0Mean - t0FEBPathCorrection(wheel, station, sector, sl, l, wire);
78  float t0RMSNew = t0RMS;
79  return DTT0Data(t0MeanNew,t0RMSNew);
80 }
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
tuple status
Definition: ntuplemaker.py:245
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 82 of file DTT0FEBPathCorrection.cc.

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

Referenced by correction().

82  {
83  // Access default DB
84  float t0Mean,t0RMS;
85  int status = t0Map_->get(wireId,t0Mean,t0RMS,DTTimeUnits::counts);
86  if(!status){
87  return DTT0Data(t0Mean,t0RMS);
88  } else{
89  //...
90  throw cms::Exception("[DTT0FEBPathCorrection]") << "Could not find t0 entry in DB for"
91  << wireId << endl;
92  }
93 }
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
tuple status
Definition: ntuplemaker.py:245
Definition: DTT0.h:37
void DTT0FEBPathCorrection::setES ( const edm::EventSetup setup)
virtual

Implements dtCalibration::DTT0BaseCorrection.

Definition at line 49 of file DTT0FEBPathCorrection.cc.

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

49  {
50  // Get t0 record from DB
51  ESHandle<DTT0> t0H;
52  setup.get<DTT0Rcd>().get(t0H);
53  t0Map_ = &*t0H;
54  LogVerbatim("Calibration") << "[DTT0FEBPathCorrection] T0 version: " << t0H->version();
55 }
const T & get() const
Definition: EventSetup.h:56
Definition: DTT0Rcd.h:9
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 143 of file DTT0FEBPathCorrection.cc.

Referenced by correction().

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