Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:27

0001 /*
0002  *  See header file for a description of this class.
0003  */
0004 
0005 #include "DTT0FEBPathCorrection.h"
0006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/ESHandle.h"
0009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0010 
0011 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0012 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0013 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
0014 
0015 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0016 #include "CondFormats/DTObjects/interface/DTT0.h"
0017 #include "CondFormats/DataRecord/interface/DTT0Rcd.h"
0018 #include "FWCore/Framework/interface/ConsumesCollector.h"
0019 
0020 #include <string>
0021 #include <sstream>
0022 
0023 using namespace std;
0024 using namespace edm;
0025 
0026 namespace dtCalibration {
0027 
0028   DTT0FEBPathCorrection::DTT0FEBPathCorrection(const ParameterSet& pset, edm::ConsumesCollector cc)
0029       : calibChamber_(pset.getParameter<string>("calibChamber")), t0Token_(cc.esConsumes<edm::Transition::BeginRun>()) {
0030     //DTChamberId chosenChamberId;
0031     if (!calibChamber_.empty() && calibChamber_ != "None" && calibChamber_ != "All") {
0032       stringstream linestr;
0033       int selWheel, selStation, selSector;
0034       linestr << calibChamber_;
0035       linestr >> selWheel >> selStation >> selSector;
0036       chosenChamberId_ = DTChamberId(selWheel, selStation, selSector);
0037       LogVerbatim("Calibration") << "[DTT0FEBPathCorrection] Chosen chamber: " << chosenChamberId_ << endl;
0038     }
0039     //FIXME: Check if chosen chamber is valid.
0040   }
0041 
0042   DTT0FEBPathCorrection::~DTT0FEBPathCorrection() {}
0043 
0044   void DTT0FEBPathCorrection::setES(const EventSetup& setup) {
0045     // Get t0 record from DB
0046     ESHandle<DTT0> t0H;
0047     t0H = setup.getHandle(t0Token_);
0048     t0Map_ = &setup.getData(t0Token_);
0049     LogVerbatim("Calibration") << "[DTT0FEBPathCorrection] T0 version: " << t0H->version();
0050   }
0051 
0052   DTT0Data DTT0FEBPathCorrection::correction(const DTWireId& wireId) {
0053     // Compute for selected chamber (or All) correction using as reference chamber mean
0054 
0055     DTChamberId chamberId = wireId.layerId().superlayerId().chamberId();
0056 
0057     if (calibChamber_.empty() || calibChamber_ == "None")
0058       return defaultT0(wireId);
0059     if (calibChamber_ != "All" && chamberId != chosenChamberId_)
0060       return defaultT0(wireId);
0061 
0062     // Access DB
0063     float t0Mean, t0RMS;
0064     int status = t0Map_->get(wireId, t0Mean, t0RMS, DTTimeUnits::counts);
0065     if (status != 0)
0066       throw cms::Exception("[DTT0FEBPathCorrection]") << "Could not find t0 entry in DB for" << wireId << endl;
0067     int wheel = chamberId.wheel();
0068     int station = chamberId.station();
0069     int sector = chamberId.sector();
0070     int sl = wireId.layerId().superlayerId().superlayer();
0071     int l = wireId.layerId().layer();
0072     int wire = wireId.wire();
0073     float t0MeanNew = t0Mean - t0FEBPathCorrection(wheel, station, sector, sl, l, wire);
0074     float t0RMSNew = t0RMS;
0075     return DTT0Data(t0MeanNew, t0RMSNew);
0076   }
0077 
0078   DTT0Data DTT0FEBPathCorrection::defaultT0(const DTWireId& wireId) {
0079     // Access default DB
0080     float t0Mean, t0RMS;
0081     int status = t0Map_->get(wireId, t0Mean, t0RMS, DTTimeUnits::counts);
0082     if (!status) {
0083       return DTT0Data(t0Mean, t0RMS);
0084     } else {
0085       //...
0086       throw cms::Exception("[DTT0FEBPathCorrection]") << "Could not find t0 entry in DB for" << wireId << endl;
0087     }
0088   }
0089 
0090   /**
0091 
0092   Return the value to be subtracted to t0s to correct for the difference of 
0093   path lenght for TP signals within the FEB, from the TP input connector 
0094   to the MAD.
0095 
0096   FEBs are alternately connected on the right (J9) and left (J8) TP 
0097   input connectors.
0098   Path lenghts (provided from Franco Gonella) for FEB-16 channels 
0099   for FEB-16s connected on the right TP connector (J9) are:
0100  
0101   CH0,1 = 32 mm 
0102   CH2,3 = 42 mm 
0103   CH4,5 = 52 mm 
0104   CH6,7 = 62 mm 
0105 
0106   CH8,9 = 111 mm 
0107   CH10,11 = 121 mm 
0108   CH12,13 = 131 mm 
0109   CH14,15 = 141 mm
0110 
0111   Given that ttrig calibration absorbs average offsets, 
0112   we assume thate only differences w.r.t. the average lenght (86.5 mm) remain.
0113 
0114   For FEBs connected on the right connector, values are swapped; so there is 
0115   a periodicity of 2 FEBS (8 cells)
0116 
0117   The mapping of FEB channels to wires, for the first two FEBs, is: 
0118 
0119          FEB 0 (J9)     FEB 1 (J8)     
0120   Wire | 1  2  3  4  |  5  6  7  8
0121   --------------------------------  
0122    L1  | 3  7 11 15  |  3  7 11 15 
0123    L2  | 1  5  9 13  |  1  5  9 13
0124    L3  | 2  6 10 14  |  2  6 10 14
0125    L4  | 0  4  8 12  |  0  4  8 12
0126 
0127 
0128   For FEB-20, distances from the left connector (J8) are:
0129   CH16,17 = 171 mm  
0130   CH18,19 = 181 mm
0131 
0132   We do not include a correction for this additional row of channels since
0133   they are at the edge of the SL so the effect cannot be seen on data (and moreover 
0134   the channel in L1 is usually not existing in the chamber)
0135 
0136 */
0137 
0138   float DTT0FEBPathCorrection::t0FEBPathCorrection(int wheel, int st, int sec, int sl, int l, int w) {
0139     // Skip correction for the last row of cells of FEB20 (see above)
0140     if ((st == 1 && ((sl != 2 && w == 49) || (sl == 2 && w == 57))) || ((st == 2 || st == 3) && (sl == 2 && w == 57)))
0141       return 0.;
0142 
0143     float dist[8] = {};
0144 
0145     // Path lenght differences for L1 and L3 (cm)
0146     if (l == 1 || l == 3) {
0147       dist[0] = +4.45;
0148       dist[1] = +2.45;
0149       dist[2] = -3.45;
0150       dist[3] = -5.45;
0151       dist[4] = -4.45;
0152       dist[5] = -2.45;
0153       dist[6] = +3.45;
0154       dist[7] = +5.45;
0155     }
0156 
0157     // Path lenght differences for L2 and L4 (cm)
0158     else {
0159       dist[0] = +5.45;
0160       dist[1] = +3.45;
0161       dist[2] = -2.45;
0162       dist[3] = -4.45;
0163       dist[4] = -5.45;
0164       dist[5] = -3.45;
0165       dist[6] = +2.45;
0166       dist[7] = +4.45;
0167     }
0168 
0169     // Wire position within the 8-cell period (2 FEBs).
0170     // Note that wire numbers start from 1.
0171     int pos = (w - 1) % 8;
0172 
0173     // Special case: in MB2 phi and MB4-10, the periodicity is broken at cell 49, as there is
0174     // an odd number of FEDs (15): the 13th FEB is connected on the left,
0175     // and not on the right; i.e. one FEB (4 colums of cells) is skipped from what would
0176     // be the regular structure.
0177     // The same happens in MB4-8/12 at cell 81.
0178     if ((st == 2 && sl != 2 && w >= 49) || (st == 4 && sec == 10 && w >= 49) ||
0179         (st == 4 && (sec == 8 || sec == 12) && w >= 81))
0180       pos = (w - 1 + 4) % 8;
0181 
0182     // Inverse of the signal propagation speed, determined from the
0183     // observed amplitude of the modulation. This matches what is found
0184     // with CAD simulation using reasonable assumptions on the PCB specs.
0185 
0186     return dist[pos] * 0.075;
0187   }
0188 
0189 }  // namespace dtCalibration