Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class DTTTrigSyncFromDB
0002  *  Concrete implementation of a DTTTrigBaseSync.
0003  *  This class define the offset for RecHit building
0004  *  of data and simulation.
0005  *  The offset is computes as: 
0006  *  <br>
0007  *  offset = t0 + tTrig + wirePropCorr - tofCorr 
0008  *  <br>
0009  *  where: <br>
0010  *     - t0 from test pulses (taken from DB, it is assumed to be in ns; can be switched off)
0011  *     - ttrig from the fit of time boxrising edge (taken from DB, it is assumed to be in ns)
0012  *       (At the moment a single value is read for ttrig offset 
0013  *       but this may change in the future)
0014  *     - signal propagation along the wire (can be switched off):
0015  *       it is assumed the ttrig accounts on average for
0016  *       correction from the center of the wire to the frontend.
0017  *       Here we just have to correct for the distance of the hit from the wire center.
0018  *     - TOF correction (can be switched off for cosmics):
0019  *       the ttrig already accounts for average TOF correction, 
0020  *       depending on the granularity used for the ttrig computation we just have to correct for the
0021  *       TOF from the center of the chamber, SL, layer or wire to the hit position.
0022  *       NOTE: particles are assumed as coming from the IP.
0023  *
0024  *  The emulatorOffset is computed as:
0025  *  <br>
0026  *  offset = int(ttrig/BXspace)*BXspace + t0
0027  *  <br>
0028  *  where: <br>
0029  *     - t0 from test pulses (taken from DB, it is assumed to be in ns; can be switched off)
0030  *     - ttrig from the fit of time box rising edge (taken from DB, it is assumed to be in ns)
0031  *     - BXspace BX spacing (in ns). Can be configured.
0032  *   
0033  *  NOTE: this should approximate what is seen online by the BTI
0034  *
0035  *  \author G. Cerminara - INFN Torino
0036  */
0037 
0038 #include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"
0039 #include "FWCore/Framework/interface/ConsumesCollector.h"
0040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0041 #include "FWCore/Framework/interface/EventSetup.h"
0042 #include "FWCore/Framework/interface/ESHandle.h"
0043 #include "Geometry/DTGeometry/interface/DTLayer.h"
0044 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
0045 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0046 #include "CondFormats/DTObjects/interface/DTT0.h"
0047 #include "CondFormats/DataRecord/interface/DTT0Rcd.h"
0048 #include "CondFormats/DTObjects/interface/DTTtrig.h"
0049 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
0050 
0051 #include <iostream>
0052 
0053 namespace edm {
0054   class ParameterSet;
0055 }
0056 
0057 class DTTTrigSyncFromDB : public DTTTrigBaseSync {
0058 public:
0059   /// Constructor
0060   DTTTrigSyncFromDB(const edm::ParameterSet& config, edm::ConsumesCollector);
0061 
0062   /// Destructor
0063   ~DTTTrigSyncFromDB() override;
0064 
0065   // Operations
0066 
0067   /// Pass the Event Setup to the algo at each event
0068   void setES(const edm::EventSetup& setup) override;
0069 
0070   /// Time (ns) to be subtracted to the digi time,
0071   /// Parameters are the layer and the wireId to which the
0072   /// digi is referred and the estimation of
0073   /// the 3D hit position (globPos)
0074   double offset(const DTLayer* layer,
0075                 const DTWireId& wireId,
0076                 const GlobalPoint& globPos,
0077                 double& tTrig,
0078                 double& wirePropCorr,
0079                 double& tofCorr) const override;
0080 
0081   /// Time (ns) to be subtracted to the digi time.
0082   /// It does not take into account TOF and signal propagation along the wire
0083   double offset(const DTWireId& wireId) const override;
0084 
0085   /// Time (ns) to be subtracted to the digi time for emulation purposes
0086   /// It does not take into account TOF and signal propagation along the wire
0087   /// It also returns the different contributions separately:
0088   ///     - tTrig is the offset (t_trig)
0089   ///     - t0cell is the t0 from pulses
0090   double emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const override;
0091 
0092 private:
0093   edm::ESGetToken<DTT0, DTT0Rcd> t0Token_;
0094   const edm::ESGetToken<DTTtrig, DTTtrigRcd> ttrigToken_;
0095   const DTT0* tZeroMap;
0096   const DTTtrig* tTrigMap;
0097   // Set the verbosity level
0098   const bool debug;
0099   // The velocity of signal propagation along the wire (cm/ns)
0100   double theVPropWire;
0101   // Switch on/off the T0 correction from pulses
0102   bool doT0Correction;
0103   // Switch on/off the TOF correction for particles from IP
0104   bool doTOFCorrection;
0105   int theTOFCorrType;
0106   // Switch on/off the correction for the signal propagation along the wire
0107   bool doWirePropCorrection;
0108   int theWirePropCorrType;
0109   // spacing of BX in ns
0110   double theBXspace;
0111 };
0112 
0113 using namespace std;
0114 using namespace edm;
0115 
0116 DTTTrigSyncFromDB::DTTTrigSyncFromDB(const ParameterSet& config, edm::ConsumesCollector cc)
0117     : ttrigToken_(cc.esConsumes(edm::ESInputTag("", config.getParameter<string>("tTrigLabel")))),
0118       debug(config.getUntrackedParameter<bool>("debug")),
0119       // The velocity of signal propagation along the wire (cm/ns)
0120       theVPropWire(config.getParameter<double>("vPropWire")),
0121       // Switch on/off the T0 correction from pulses
0122       doT0Correction(config.getParameter<bool>("doT0Correction")),
0123       // Switch on/off the TOF correction for particles from IP
0124       doTOFCorrection(config.getParameter<bool>("doTOFCorrection")),
0125       theTOFCorrType(config.getParameter<int>("tofCorrType")),
0126       // Switch on/off the correction for the signal propagation along the wire
0127       doWirePropCorrection(config.getParameter<bool>("doWirePropCorrection")),
0128       theWirePropCorrType(config.getParameter<int>("wirePropCorrType")),
0129       // spacing of BX in ns
0130       theBXspace(config.getUntrackedParameter<double>("bxSpace", 25.)) {
0131   if (doT0Correction) {
0132     t0Token_ = cc.esConsumes(edm::ESInputTag("", config.getParameter<string>("t0Label")));
0133   }
0134 }
0135 
0136 DTTTrigSyncFromDB::~DTTTrigSyncFromDB() {}
0137 
0138 void DTTTrigSyncFromDB::setES(const EventSetup& setup) {
0139   if (doT0Correction) {
0140     // Get the map of t0 from pulses from the Setup
0141     tZeroMap = &setup.getData(t0Token_);
0142     if (debug) {
0143       cout << "[DTTTrigSyncFromDB] t0 version: " << tZeroMap->version() << endl;
0144     }
0145   }
0146 
0147   // Get the map of ttrig from the Setup
0148   tTrigMap = &setup.getData(ttrigToken_);
0149   if (debug) {
0150     cout << "[DTTTrigSyncFromDB] ttrig version: " << tTrigMap->version() << endl;
0151   }
0152 }
0153 
0154 double DTTTrigSyncFromDB::offset(const DTLayer* layer,
0155                                  const DTWireId& wireId,
0156                                  const GlobalPoint& globPos,
0157                                  double& tTrig,
0158                                  double& wirePropCorr,
0159                                  double& tofCorr) const {
0160   // Correction for the float to int conversion while writeing the ttrig in ns into an int variable
0161   // (half a bin on average)
0162   // FIXME: this should disappear as soon as the ttrig object will become a float
0163   //   static const float f2i_convCorr = (25./64.); // ns //FIXME: check how the conversion is performed
0164 
0165   tTrig = offset(wireId);
0166 
0167   // Compute the time spent in signal propagation along wire.
0168   // NOTE: the FE is always at y>0
0169   wirePropCorr = 0;
0170   if (doWirePropCorrection) {
0171     switch (theWirePropCorrType) {
0172         // The ttrig computed from the timebox accounts on average for the signal propagation time
0173         // from the center of the wire to the frontend. Here we just have to correct for
0174         // the distance of the hit from the wire center.
0175       case 0: {
0176         float wireCoord = layer->toLocal(globPos).y();
0177         wirePropCorr = -wireCoord / theVPropWire;
0178         break;
0179         // FIXME: What if hits used for the time box are not distributed uniformly along the wire?
0180       }
0181       //On simulated data you need to subtract the total propagation time
0182       case 1: {
0183         float halfL = layer->specificTopology().cellLenght() / 2;
0184         float wireCoord = layer->toLocal(globPos).y();
0185         float propgL = halfL - wireCoord;
0186         wirePropCorr = propgL / theVPropWire;
0187         break;
0188       }
0189       default: {
0190         throw cms::Exception("[DTTTrigSyncFromDB]")
0191             << " Invalid parameter: wirePropCorrType = " << theWirePropCorrType << std::endl;
0192         break;
0193       }
0194     }
0195   }
0196 
0197   // Compute TOF correction:
0198   tofCorr = 0.;
0199   // TOF Correction can be switched off with appropriate parameter
0200   if (doTOFCorrection) {
0201     float flightToHit = globPos.mag();
0202     static const float cSpeed = 29.9792458;  // cm/ns
0203     switch (theTOFCorrType) {
0204       case 0: {
0205         // The ttrig computed from the real data accounts on average for the TOF correction
0206         // Depending on the granularity used for the ttrig computation we just have to correct for the
0207         // TOF from the center of the chamber, SL, layer or wire to the hit position.
0208         // At the moment only SL granularity is considered
0209         // Correction for TOF from the center of the SL to hit position
0210         const DTSuperLayer* sl = layer->superLayer();
0211         double flightToSL = sl->surface().position().mag();
0212         tofCorr = (flightToSL - flightToHit) / cSpeed;
0213         break;
0214       }
0215       case 1: {
0216         // On simulated data you need to consider only the TOF from 3D center of the wire to hit position
0217         // (because the TOF from the IP to the wire has been already subtracted in the digitization:
0218         // SimMuon/DTDigitizer/DTDigiSyncTOFCorr.cc corrType=2)
0219         float flightToWire =
0220             layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
0221         tofCorr = (flightToWire - flightToHit) / cSpeed;
0222         break;
0223       }
0224       default: {
0225         throw cms::Exception("[DTTTrigSyncFromDB]")
0226             << " Invalid parameter: tofCorrType = " << theTOFCorrType << std::endl;
0227         break;
0228       }
0229     }
0230   }
0231 
0232   if (debug) {
0233     cout << "[DTTTrigSyncFromDB] Channel: " << wireId << endl
0234          << "      Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
0235          << "      various contributions are: " << endl
0236          << "      tTrig + t0 (ns):   " << tTrig
0237          << endl
0238          //<< "      tZero (ns):   " << t0 << endl
0239          << "      Propagation along wire delay (ns): " << wirePropCorr << endl
0240          << "      TOF correction (ns): " << tofCorr << endl
0241          << endl;
0242   }
0243   //The global offset is the sum of various contributions
0244   return tTrig + wirePropCorr - tofCorr;
0245 }
0246 
0247 double DTTTrigSyncFromDB::offset(const DTWireId& wireId) const {
0248   float t0 = 0;
0249   float t0rms = 0;
0250   if (doT0Correction) {
0251     // Read the t0 from pulses for this wire (ns)
0252     tZeroMap->get(wireId, t0, t0rms, DTTimeUnits::ns);
0253   }
0254 
0255   // Read the ttrig for this wire
0256   float ttrigMean = 0;
0257   float ttrigSigma = 0;
0258   float kFactor = 0;
0259   // FIXME: should check the return value of the DTTtrigRcd::get(..) method
0260   if (tTrigMap->get(wireId.superlayerId(), ttrigMean, ttrigSigma, kFactor, DTTimeUnits::ns) != 0) {
0261     cout << "[DTTTrigSyncFromDB]*Error: ttrig not found for SL: " << wireId.superlayerId() << endl;
0262     //     FIXME: LogError.....
0263   }
0264 
0265   return t0 + ttrigMean + kFactor * ttrigSigma;
0266 }
0267 
0268 double DTTTrigSyncFromDB::emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const {
0269   float t0 = 0;
0270   float t0rms = 0;
0271   if (doT0Correction) {
0272     // Read the t0 from pulses for this wire (ns)
0273     tZeroMap->get(wireId, t0, t0rms, DTTimeUnits::ns);
0274   }
0275 
0276   // Read the ttrig for this wire
0277   float ttrigMean = 0;
0278   float ttrigSigma = 0;
0279   float kFactor = 0;
0280   // FIXME: should check the return value of the DTTtrigRcd::get(..) method
0281   if (tTrigMap->get(wireId.superlayerId(), ttrigMean, ttrigSigma, kFactor, DTTimeUnits::ns) != 0) {
0282     cout << "[DTTTrigSyncFromDB]*Error: ttrig not found for SL: " << wireId.superlayerId() << endl;
0283     //     FIXME: LogError.....
0284   }
0285 
0286   tTrig = ttrigMean + kFactor * ttrigSigma;
0287   t0cell = t0;
0288 
0289   return int(tTrig / theBXspace) * theBXspace + t0cell;
0290 }
0291 
0292 #include "FWCore/PluginManager/interface/PluginFactory.h"
0293 #include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"
0294 
0295 DEFINE_EDM_PLUGIN(DTTTrigSyncFactory, DTTTrigSyncFromDB, "DTTTrigSyncFromDB");