Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /** \class DTTTrigSyncTOFCorr
0002  *  Concrete implementation of a DTTTrigBaseSync.
0003  *  This class define the offsets for RecHit building
0004  *  coherently to the digitization realized with the
0005  *  DTDigiSyncTOFCorr module.
0006  *  The offset is computes as:<br>
0007  *  offset = tTrig + wirePropCorr - tofCorr<br>
0008  *  where:<br>
0009  *      - tTrig is a fixed offset defined in tTrig parameter
0010  *        (default 500 ns)<br>
0011  *      - wirePropCorr is the correction for the signal propagation along the wire<br>
0012  *      - tofCorr is the correction for the TOF of the particle set according to
0013  *        tofCorrType parameter:<br>
0014  *        0: tofCorrType = TOF from IP to 3D Hit position (globPos)<br>
0015  *        1: tofCorrType = TOF correction for distance difference
0016  *                         between 3D center of the chamber and hit position<br>
0017  *        2: tofCorrType = TOF correction for distance difference
0018  *                         between 3D center of the wire and hit position
0019  *                         (This mode in available for backward compatibility)<br>
0020  *
0021  *  The emulatorOffset is computed as:
0022  *  <br>
0023  *  offset = int(ttrig/BXspace)*BXspace
0024  *  <br>
0025  *  where: <br>
0026  *     - ttrig from the fit of time box rising edge (taken from configuration, it is assumed to be in ns)
0027  *     - BXspace BX spacing (in ns). Taken from configuration (default 25ns).
0028  *   
0029  *  NOTE: this should approximate what is seen online by the BTI
0030  *
0031  *
0032  *
0033  *  \author G. Cerminara - INFN Torino
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   /// Constructor
0048   DTTTrigSyncTOFCorr(const edm::ParameterSet& config, edm::ConsumesCollector);
0049 
0050   /// Destructor
0051   ~DTTTrigSyncTOFCorr() override;
0052 
0053   // Operations
0054 
0055   /// Pass the Event Setup to the algo at each event
0056   void setES(const edm::EventSetup& setup) override {}
0057 
0058   /// Time (ns) to be subtracted to the digi time,
0059   /// Parameters are the layer and the wireId to which the
0060   /// digi is referred and the estimation of
0061   /// the 3D hit position (globPos)
0062   /// It also returns the different contributions separately:
0063   ///     - tTrig is the offset (t_trig)
0064   ///     - wirePropCorr is the delay for signal propagation along the wire
0065   ///     - tofCorr is the correction due to the particle TOF
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   /// Time (ns) to be subtracted to the digi time for emulation purposes
0076   /// It does not take into account TOF and signal propagation along the wire
0077   /// It also returns the different contributions separately:
0078   ///     - tTrig is the offset (t_trig)
0079   ///     - t0cell is the t0 from pulses (always 0 in this case)
0080   double emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const override;
0081 
0082 private:
0083   // The fixed t_trig to be subtracted to digi time (ns)
0084   const double theTTrig;
0085   // Velocity of signal propagation along the wire (cm/ns)
0086   // For the value
0087   // cfr. CMS-IN 2000-021:   (2.56+-0.17)x1e8 m/s
0088   //      CMS NOTE 2003-17:  (0.244)  m/ns = 24.4 cm/ns
0089   const double theVPropWire;
0090 
0091   // Select the mode for TOF correction:
0092   //     0: tofCorr = TOF from IP to 3D Hit position (globPos)
0093   //     1: tofCorr = TOF correction for distance difference
0094   //                  between 3D center of the chamber and hit position
0095   //     2: tofCorr = TOF correction for distance difference
0096   //                  between 3D center of the wire and hit position
0097   const int theTOFCorrType;
0098 
0099   // Set the verbosity level
0100   const bool debug;
0101   // spacing of BX in ns
0102   double theBXspace;
0103 };
0104 
0105 using namespace std;
0106 
0107 DTTTrigSyncTOFCorr::DTTTrigSyncTOFCorr(const edm::ParameterSet& config, edm::ConsumesCollector)
0108     :                                                  // The fixed t0 (or t_trig) to be subtracted to digi time (ns)
0109       theTTrig(config.getParameter<double>("tTrig")),  // FIXME: Default was 500 ns
0110       // Velocity of signal propagation along the wire (cm/ns)
0111       theVPropWire(config.getParameter<double>("vPropWire")),  // FIXME: Default was 24.4 cm/ns
0112       // Select the mode for TOF correction:
0113       //     0: tofCorr = TOF from IP to 3D Hit position (globPos)
0114       //     1: tofCorr = TOF correction for distance difference
0115       //                  between 3D center of the chamber and hit position
0116       //                  (default)
0117       //     2: tofCorr = TOF correction for distance difference
0118       //                  between 3D center of the wire and hit position
0119       //                  (This mode in available for backward compatibility)
0120       theTOFCorrType(config.getParameter<int>("tofCorrType")),  // FIXME: Default was 1
0121       debug(config.getUntrackedParameter<bool>("debug")),
0122       theBXspace(config.getUntrackedParameter<double>("bxSpace", 25.))  // spacing of BX in ns
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   //Compute the time spent in signal propagation along wire.
0136   // NOTE: the FE is always at y>0
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   // Compute TOF correction treating it accordingly to
0143   // the tofCorrType card
0144   float flightToHit = globPos.mag();
0145   static const float cSpeed = 29.9792458;  // cm/ns
0146   tofCorr = 0.;
0147   switch (theTOFCorrType) {
0148     case 0: {
0149       // In this mode the subtraction of the TOF from IP to
0150       // estimate 3D hit digi position is done here
0151       // (No correction is needed anymore)
0152       tofCorr = -flightToHit / cSpeed;
0153       break;
0154     }
0155     case 1: {
0156       // Correction for TOF from the center of the chamber to hit position
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       // TOF from 3D center of the wire to hit position
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   //The global offset is the sum of various contributions
0185   return tTrig + wirePropCorr - tofCorr;
0186 }
0187 
0188 double DTTTrigSyncTOFCorr::offset(const DTWireId& wireId) const {
0189   // Correction for the float to int conversion
0190   // (half a bin on average) in DTDigi constructor
0191   static const float f2i_convCorr = (25. / 64.);  // ns
0192   //FIXME: This should be considered only if the Digi is constructed from a float....
0193 
0194   // The tTrig is taken from a parameter
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");