DTTTrigSyncTOFCorr

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208
/** \class DTTTrigSyncTOFCorr
 *  Concrete implementation of a DTTTrigBaseSync.
 *  This class define the offsets for RecHit building
 *  coherently to the digitization realized with the
 *  DTDigiSyncTOFCorr module.
 *  The offset is computes as:<br>
 *  offset = tTrig + wirePropCorr - tofCorr<br>
 *  where:<br>
 *      - tTrig is a fixed offset defined in tTrig parameter
 *        (default 500 ns)<br>
 *      - wirePropCorr is the correction for the signal propagation along the wire<br>
 *      - tofCorr is the correction for the TOF of the particle set according to
 *        tofCorrType parameter:<br>
 *        0: tofCorrType = TOF from IP to 3D Hit position (globPos)<br>
 *        1: tofCorrType = TOF correction for distance difference
 *                         between 3D center of the chamber and hit position<br>
 *        2: tofCorrType = TOF correction for distance difference
 *                         between 3D center of the wire and hit position
 *                         (This mode in available for backward compatibility)<br>
 *
 *  The emulatorOffset is computed as:
 *  <br>
 *  offset = int(ttrig/BXspace)*BXspace
 *  <br>
 *  where: <br>
 *     - ttrig from the fit of time box rising edge (taken from configuration, it is assumed to be in ns)
 *     - BXspace BX spacing (in ns). Taken from configuration (default 25ns).
 *   
 *  NOTE: this should approximate what is seen online by the BTI
 *
 *
 *
 *  \author G. Cerminara - INFN Torino
 */

#include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"
#include "FWCore/Framework/interface/ConsumesCollector.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/Exception.h"
#include "Geometry/DTGeometry/interface/DTLayer.h"
#include "Geometry/DTGeometry/interface/DTChamber.h"
#include "DataFormats/MuonDetId/interface/DTWireId.h"
#include <iostream>

class DTTTrigSyncTOFCorr : public DTTTrigBaseSync {
public:
  /// Constructor
  DTTTrigSyncTOFCorr(const edm::ParameterSet& config, edm::ConsumesCollector);

  /// Destructor
  ~DTTTrigSyncTOFCorr() override;

  // Operations

  /// Pass the Event Setup to the algo at each event
  void setES(const edm::EventSetup& setup) override {}

  /// Time (ns) to be subtracted to the digi time,
  /// Parameters are the layer and the wireId to which the
  /// digi is referred and the estimation of
  /// the 3D hit position (globPos)
  /// It also returns the different contributions separately:
  ///     - tTrig is the offset (t_trig)
  ///     - wirePropCorr is the delay for signal propagation along the wire
  ///     - tofCorr is the correction due to the particle TOF
  double offset(const DTLayer* layer,
                const DTWireId& wireId,
                const GlobalPoint& globPos,
                double& tTrig,
                double& wirePropCorr,
                double& tofCorr) const override;

  double offset(const DTWireId& wireId) const override;

  /// Time (ns) to be subtracted to the digi time for emulation purposes
  /// It does not take into account TOF and signal propagation along the wire
  /// It also returns the different contributions separately:
  ///     - tTrig is the offset (t_trig)
  ///     - t0cell is the t0 from pulses (always 0 in this case)
  double emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const override;

private:
  // The fixed t_trig to be subtracted to digi time (ns)
  const double theTTrig;
  // Velocity of signal propagation along the wire (cm/ns)
  // For the value
  // cfr. CMS-IN 2000-021:   (2.56+-0.17)x1e8 m/s
  //      CMS NOTE 2003-17:  (0.244)  m/ns = 24.4 cm/ns
  const double theVPropWire;

  // Select the mode for TOF correction:
  //     0: tofCorr = TOF from IP to 3D Hit position (globPos)
  //     1: tofCorr = TOF correction for distance difference
  //                  between 3D center of the chamber and hit position
  //     2: tofCorr = TOF correction for distance difference
  //                  between 3D center of the wire and hit position
  const int theTOFCorrType;

  // Set the verbosity level
  const bool debug;
  // spacing of BX in ns
  double theBXspace;
};

using namespace std;

DTTTrigSyncTOFCorr::DTTTrigSyncTOFCorr(const edm::ParameterSet& config, edm::ConsumesCollector)
    :                                                  // The fixed t0 (or t_trig) to be subtracted to digi time (ns)
      theTTrig(config.getParameter<double>("tTrig")),  // FIXME: Default was 500 ns
      // Velocity of signal propagation along the wire (cm/ns)
      theVPropWire(config.getParameter<double>("vPropWire")),  // FIXME: Default was 24.4 cm/ns
      // Select the mode for TOF correction:
      //     0: tofCorr = TOF from IP to 3D Hit position (globPos)
      //     1: tofCorr = TOF correction for distance difference
      //                  between 3D center of the chamber and hit position
      //                  (default)
      //     2: tofCorr = TOF correction for distance difference
      //                  between 3D center of the wire and hit position
      //                  (This mode in available for backward compatibility)
      theTOFCorrType(config.getParameter<int>("tofCorrType")),  // FIXME: Default was 1
      debug(config.getUntrackedParameter<bool>("debug")),
      theBXspace(config.getUntrackedParameter<double>("bxSpace", 25.))  // spacing of BX in ns
{}

DTTTrigSyncTOFCorr::~DTTTrigSyncTOFCorr() {}

double DTTTrigSyncTOFCorr::offset(const DTLayer* layer,
                                  const DTWireId& wireId,
                                  const GlobalPoint& globPos,
                                  double& tTrig,
                                  double& wirePropCorr,
                                  double& tofCorr) const {
  tTrig = offset(wireId);

  //Compute the time spent in signal propagation along wire.
  // NOTE: the FE is always at y>0
  float halfL = layer->specificTopology().cellLenght() / 2;
  float wireCoord = layer->toLocal(globPos).y();
  float propgL = halfL - wireCoord;
  wirePropCorr = propgL / theVPropWire;

  // Compute TOF correction treating it accordingly to
  // the tofCorrType card
  float flightToHit = globPos.mag();
  static const float cSpeed = 29.9792458;  // cm/ns
  tofCorr = 0.;
  switch (theTOFCorrType) {
    case 0: {
      // In this mode the subtraction of the TOF from IP to
      // estimate 3D hit digi position is done here
      // (No correction is needed anymore)
      tofCorr = -flightToHit / cSpeed;
      break;
    }
    case 1: {
      // Correction for TOF from the center of the chamber to hit position
      const DTChamber* chamber = layer->chamber();
      double flightToChamber = chamber->surface().position().mag();
      tofCorr = (flightToChamber - flightToHit) / cSpeed;
      break;
    }
    case 2: {
      // TOF from 3D center of the wire to hit position
      float flightToWire =
          layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
      tofCorr = (flightToWire - flightToHit) / cSpeed;
      break;
    }
    default: {
      throw cms::Exception("[DTTTrigSyncTOFCorr]")
          << " Invalid parameter: tofCorrType = " << theTOFCorrType << std::endl;
      break;
    }
  }

  if (debug) {
    cout << "[DTTTrigSyncTOFCorr] Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
         << "      various contributions are: " << endl
         << "      tTrig (ns):   " << tTrig << endl
         << "      Propagation along wire delay (ns): " << wirePropCorr << endl
         << "      TOF correction (ns): " << tofCorr << endl
         << endl;
  }
  //The global offset is the sum of various contributions
  return tTrig + wirePropCorr - tofCorr;
}

double DTTTrigSyncTOFCorr::offset(const DTWireId& wireId) const {
  // Correction for the float to int conversion
  // (half a bin on average) in DTDigi constructor
  static const float f2i_convCorr = (25. / 64.);  // ns
  //FIXME: This should be considered only if the Digi is constructed from a float....

  // The tTrig is taken from a parameter
  return theTTrig - f2i_convCorr;
}

double DTTTrigSyncTOFCorr::emulatorOffset(const DTWireId& wireId, double& tTrig, double& t0cell) const {
  tTrig = theTTrig;
  t0cell = 0.;

  return int(tTrig / theBXspace) * theBXspace;
}

#include "FWCore/PluginManager/interface/PluginFactory.h"
#include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"

DEFINE_EDM_PLUGIN(DTTTrigSyncFactory, DTTTrigSyncTOFCorr, "DTTTrigSyncTOFCorr");