Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "DTTTrigSyncTOFCorr.h"
00010
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012 #include "FWCore/Utilities/interface/Exception.h"
00013 #include "Geometry/DTGeometry/interface/DTLayer.h"
00014 #include "Geometry/DTGeometry/interface/DTChamber.h"
00015 #include "DataFormats/MuonDetId/interface/DTWireId.h"
00016 #include <iostream>
00017
00018 using namespace std;
00019
00020
00021
00022 DTTTrigSyncTOFCorr::DTTTrigSyncTOFCorr(const edm::ParameterSet& config){
00023 theTTrig = config.getParameter<double>("tTrig");
00024 theVPropWire = config.getParameter<double>("vPropWire");
00025 theTOFCorrType = config.getParameter<int>("tofCorrType");
00026 debug = config.getUntrackedParameter<bool>("debug");
00027
00028 theBXspace = config.getUntrackedParameter<double>("bxSpace", 25.);
00029
00030 }
00031
00032
00033
00034 DTTTrigSyncTOFCorr::~DTTTrigSyncTOFCorr(){}
00035
00036
00037
00038 double DTTTrigSyncTOFCorr::offset(const DTLayer* layer,
00039 const DTWireId& wireId,
00040 const GlobalPoint& globPos,
00041 double& tTrig,
00042 double& wirePropCorr,
00043 double& tofCorr) {
00044 tTrig = offset(wireId);
00045
00046
00047
00048 float halfL = layer->specificTopology().cellLenght()/2;
00049 float wireCoord = layer->toLocal(globPos).y();
00050 float propgL = halfL - wireCoord;
00051 wirePropCorr = propgL/theVPropWire;
00052
00053
00054
00055
00056 float flightToHit = globPos.mag();
00057 static const float cSpeed = 29.9792458;
00058 tofCorr = 0.;
00059 switch(theTOFCorrType) {
00060 case 0: {
00061
00062
00063
00064 tofCorr = -flightToHit/cSpeed;
00065 break;
00066 }
00067 case 1: {
00068
00069 const DTChamber * chamber = layer->chamber();
00070 double flightToChamber = chamber->surface().position().mag();
00071 tofCorr = (flightToChamber-flightToHit)/cSpeed;
00072 break;
00073 }
00074 case 2: {
00075
00076 float flightToWire =
00077 layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
00078 tofCorr = (flightToWire-flightToHit)/cSpeed;
00079 break;
00080 }
00081 default: {
00082 throw
00083 cms::Exception("[DTTTrigSyncTOFCorr]") << " Invalid parameter: tofCorrType = "
00084 << theTOFCorrType
00085 << std::endl;
00086 break;
00087 }
00088 }
00089
00090 if(debug) {
00091 cout << "[DTTTrigSyncTOFCorr] Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
00092 << " various contributions are: " << endl
00093 << " tTrig (ns): " << tTrig << endl
00094 << " Propagation along wire delay (ns): " << wirePropCorr << endl
00095 << " TOF correction (ns): " << tofCorr << endl
00096 << endl;
00097 }
00098
00099 return tTrig + wirePropCorr - tofCorr;
00100 }
00101
00102 double DTTTrigSyncTOFCorr::offset(const DTWireId& wireId) {
00103
00104
00105 static const float f2i_convCorr = (25./64.);
00106
00107
00108
00109
00110 return theTTrig - f2i_convCorr;
00111 }
00112
00113
00114
00115 double DTTTrigSyncTOFCorr::theTTrig;
00116
00117
00118
00119
00120 double DTTTrigSyncTOFCorr::theVPropWire;
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132 int DTTTrigSyncTOFCorr::theTOFCorrType;
00133
00134
00135
00136
00137 bool DTTTrigSyncTOFCorr::debug;
00138
00139
00140
00141 double DTTTrigSyncTOFCorr::emulatorOffset(const DTWireId& wireId,
00142 double &tTrig,
00143 double &t0cell) {
00144 tTrig = theTTrig;
00145 t0cell = 0.;
00146
00147 return int(tTrig/theBXspace)*theBXspace;
00148 }
00149