#include <DTT0FEBPathCorrection.h>
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 DTT0 * | t0Map_ |
Definition at line 25 of file DTT0FEBPathCorrection.h.
DTT0FEBPathCorrection::DTT0FEBPathCorrection | ( | const edm::ParameterSet & | pset | ) |
Definition at line 27 of file DTT0FEBPathCorrection.cc.
References calibChamber_, and chosenChamberId_.
|
override |
Definition at line 42 of file DTT0FEBPathCorrection.cc.
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().
Definition at line 78 of file DTT0FEBPathCorrection.cc.
References DTTimeUnits::counts, Exception, DTT0::get(), mps_update::status, and t0Map_.
Referenced by correction().
|
overridevirtual |
Implements dtCalibration::DTT0BaseCorrection.
Definition at line 45 of file DTT0FEBPathCorrection.cc.
References edm::EventSetup::get(), t0Map_, and DTT0::version().
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)
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().
|
private |
Definition at line 40 of file DTT0FEBPathCorrection.h.
Referenced by correction(), and DTT0FEBPathCorrection().
|
private |
Definition at line 42 of file DTT0FEBPathCorrection.h.
Referenced by correction(), and DTT0FEBPathCorrection().
|
private |
Definition at line 43 of file DTT0FEBPathCorrection.h.
Referenced by correction(), defaultT0(), and setES().