File indexing completed on 2024-04-06 11:58:31
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036 #include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"
0037 #include "FWCore/Framework/interface/ConsumesCollector.h"
0038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0039 #include "FWCore/Utilities/interface/Exception.h"
0040 #include "Geometry/DTGeometry/interface/DTLayer.h"
0041 #include "Geometry/DTGeometry/interface/DTChamber.h"
0042 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0043 #include <iostream>
0044
0045 class DTTTrigSyncTOFCorr : public DTTTrigBaseSync {
0046 public:
0047
0048 DTTTrigSyncTOFCorr(const edm::ParameterSet& config, edm::ConsumesCollector);
0049
0050
0051 ~DTTTrigSyncTOFCorr() override;
0052
0053
0054
0055
0056 void setES(const edm::EventSetup& setup) override {}
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066 double offset(const DTLayer* layer,
0067 const DTWireId& wireId,
0068 const GlobalPoint& globPos,
0069 double& tTrig,
0070 double& wirePropCorr,
0071 double& tofCorr) const override;
0072
0073 double offset(const DTWireId& wireId) const override;
0074
0075
0076
0077
0078
0079
0080 double emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const override;
0081
0082 private:
0083
0084 const double theTTrig;
0085
0086
0087
0088
0089 const double theVPropWire;
0090
0091
0092
0093
0094
0095
0096
0097 const int theTOFCorrType;
0098
0099
0100 const bool debug;
0101
0102 double theBXspace;
0103 };
0104
0105 using namespace std;
0106
0107 DTTTrigSyncTOFCorr::DTTTrigSyncTOFCorr(const edm::ParameterSet& config, edm::ConsumesCollector)
0108 :
0109 theTTrig(config.getParameter<double>("tTrig")),
0110
0111 theVPropWire(config.getParameter<double>("vPropWire")),
0112
0113
0114
0115
0116
0117
0118
0119
0120 theTOFCorrType(config.getParameter<int>("tofCorrType")),
0121 debug(config.getUntrackedParameter<bool>("debug")),
0122 theBXspace(config.getUntrackedParameter<double>("bxSpace", 25.))
0123 {}
0124
0125 DTTTrigSyncTOFCorr::~DTTTrigSyncTOFCorr() {}
0126
0127 double DTTTrigSyncTOFCorr::offset(const DTLayer* layer,
0128 const DTWireId& wireId,
0129 const GlobalPoint& globPos,
0130 double& tTrig,
0131 double& wirePropCorr,
0132 double& tofCorr) const {
0133 tTrig = offset(wireId);
0134
0135
0136
0137 float halfL = layer->specificTopology().cellLenght() / 2;
0138 float wireCoord = layer->toLocal(globPos).y();
0139 float propgL = halfL - wireCoord;
0140 wirePropCorr = propgL / theVPropWire;
0141
0142
0143
0144 float flightToHit = globPos.mag();
0145 static const float cSpeed = 29.9792458;
0146 tofCorr = 0.;
0147 switch (theTOFCorrType) {
0148 case 0: {
0149
0150
0151
0152 tofCorr = -flightToHit / cSpeed;
0153 break;
0154 }
0155 case 1: {
0156
0157 const DTChamber* chamber = layer->chamber();
0158 double flightToChamber = chamber->surface().position().mag();
0159 tofCorr = (flightToChamber - flightToHit) / cSpeed;
0160 break;
0161 }
0162 case 2: {
0163
0164 float flightToWire =
0165 layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
0166 tofCorr = (flightToWire - flightToHit) / cSpeed;
0167 break;
0168 }
0169 default: {
0170 throw cms::Exception("[DTTTrigSyncTOFCorr]")
0171 << " Invalid parameter: tofCorrType = " << theTOFCorrType << std::endl;
0172 break;
0173 }
0174 }
0175
0176 if (debug) {
0177 cout << "[DTTTrigSyncTOFCorr] Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
0178 << " various contributions are: " << endl
0179 << " tTrig (ns): " << tTrig << endl
0180 << " Propagation along wire delay (ns): " << wirePropCorr << endl
0181 << " TOF correction (ns): " << tofCorr << endl
0182 << endl;
0183 }
0184
0185 return tTrig + wirePropCorr - tofCorr;
0186 }
0187
0188 double DTTTrigSyncTOFCorr::offset(const DTWireId& wireId) const {
0189
0190
0191 static const float f2i_convCorr = (25. / 64.);
0192
0193
0194
0195 return theTTrig - f2i_convCorr;
0196 }
0197
0198 double DTTTrigSyncTOFCorr::emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const {
0199 tTrig = theTTrig;
0200 t0cell = 0.;
0201
0202 return int(tTrig / theBXspace) * theBXspace;
0203 }
0204
0205 #include "FWCore/PluginManager/interface/PluginFactory.h"
0206 #include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"
0207
0208 DEFINE_EDM_PLUGIN(DTTTrigSyncFactory, DTTTrigSyncTOFCorr, "DTTTrigSyncTOFCorr");