CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/CalibMuon/DTDigiSync/src/DTTTrigSyncFromDB.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2009/12/07 17:22:18 $
00005  *  $Revision: 1.10 $
00006  *  \author G. Cerminara - INFN Torino
00007  */
00008 
00009 #include "DTTTrigSyncFromDB.h"
00010 
00011 
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "Geometry/DTGeometry/interface/DTLayer.h"
00016 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
00017 #include "DataFormats/MuonDetId/interface/DTWireId.h"
00018 #include "CondFormats/DTObjects/interface/DTT0.h"
00019 #include "CondFormats/DataRecord/interface/DTT0Rcd.h"
00020 #include "CondFormats/DTObjects/interface/DTTtrig.h"
00021 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
00022 
00023 #include <iostream>
00024 
00025 using namespace std;
00026 using namespace edm;
00027 
00028 
00029 DTTTrigSyncFromDB::DTTTrigSyncFromDB(const ParameterSet& config){
00030   debug = config.getUntrackedParameter<bool>("debug");
00031   // The velocity of signal propagation along the wire (cm/ns)
00032   theVPropWire = config.getParameter<double>("vPropWire");
00033   // Switch on/off the T0 correction from pulses
00034   doT0Correction = config.getParameter<bool>("doT0Correction");
00035   // Switch on/off the TOF correction for particles from IP
00036   doTOFCorrection = config.getParameter<bool>("doTOFCorrection");
00037   theTOFCorrType = config.getParameter<int>("tofCorrType");
00038   // Switch on/off the correction for the signal propagation along the wire
00039   doWirePropCorrection = config.getParameter<bool>("doWirePropCorrection");
00040   theWirePropCorrType = config.getParameter<int>("wirePropCorrType");
00041   // spacing of BX in ns
00042   theBXspace  = config.getUntrackedParameter<double>("bxSpace", 25.);
00043 
00044   thetTrigLabel = config.getParameter<string>("tTrigLabel");
00045 }
00046 
00047 
00048 
00049 DTTTrigSyncFromDB::~DTTTrigSyncFromDB(){}
00050 
00051 
00052 
00053 void DTTTrigSyncFromDB::setES(const EventSetup& setup) {
00054   if(doT0Correction)
00055     {
00056       // Get the map of t0 from pulses from the Setup
00057       ESHandle<DTT0> t0Handle;
00058       setup.get<DTT0Rcd>().get(t0Handle);
00059       tZeroMap = &*t0Handle;
00060       if(debug) {
00061         cout << "[DTTTrigSyncFromDB] t0 version: " << tZeroMap->version()<<endl;
00062           }
00063     }
00064 
00065   // Get the map of ttrig from the Setup
00066   ESHandle<DTTtrig> ttrigHandle;
00067   setup.get<DTTtrigRcd>().get(thetTrigLabel,ttrigHandle);
00068   tTrigMap = &*ttrigHandle;
00069     if(debug) {
00070       cout << "[DTTTrigSyncFromDB] ttrig version: " << tTrigMap->version() << endl;
00071   }
00072 }
00073 
00074 
00075 double DTTTrigSyncFromDB::offset(const DTLayer* layer,
00076                                   const DTWireId& wireId,
00077                                   const GlobalPoint& globPos,
00078                                   double& tTrig,
00079                                   double& wirePropCorr,
00080                                   double& tofCorr) {
00081   // Correction for the float to int conversion while writeing the ttrig in ns into an int variable
00082   // (half a bin on average)
00083   // FIXME: this should disappear as soon as the ttrig object will become a float
00084   //   static const float f2i_convCorr = (25./64.); // ns //FIXME: check how the conversion is performed
00085 
00086   tTrig = offset(wireId);
00087 
00088   // Compute the time spent in signal propagation along wire.
00089    // NOTE: the FE is always at y>0
00090   wirePropCorr = 0;
00091   if(doWirePropCorrection) {
00092     switch(theWirePropCorrType){
00093       // The ttrig computed from the timebox accounts on average for the signal propagation time
00094       // from the center of the wire to the frontend. Here we just have to correct for
00095       // the distance of the hit from the wire center.
00096     case 0: {
00097       float wireCoord = layer->toLocal(globPos).y();
00098       wirePropCorr = -wireCoord/theVPropWire;
00099       break;
00100       // FIXME: What if hits used for the time box are not distributed uniformly along the wire?
00101     }
00102     //On simulated data you need to subtract the total propagation time 
00103     case 1: {
00104       float halfL     = layer->specificTopology().cellLenght()/2;
00105       float wireCoord = layer->toLocal(globPos).y();
00106       float propgL    = halfL - wireCoord;
00107       wirePropCorr = propgL/theVPropWire;
00108       break;
00109     }
00110     default: {
00111     throw
00112       cms::Exception("[DTTTrigSyncFromDB]") << " Invalid parameter: wirePropCorrType = "
00113                                              << theWirePropCorrType 
00114                                              << std::endl;
00115     break;
00116     }
00117     }
00118   }
00119 
00120   // Compute TOF correction:
00121   tofCorr = 0.;
00122   // TOF Correction can be switched off with appropriate parameter
00123   if(doTOFCorrection) {
00124     float flightToHit = globPos.mag();
00125     static const float cSpeed = 29.9792458; // cm/ns
00126     switch(theTOFCorrType) {
00127     case 0: {
00128       // The ttrig computed from the real data accounts on average for the TOF correction 
00129       // Depending on the granularity used for the ttrig computation we just have to correct for the
00130       // TOF from the center of the chamber, SL, layer or wire to the hit position.
00131       // At the moment only SL granularity is considered
00132       // Correction for TOF from the center of the SL to hit position
00133       const DTSuperLayer *sl = layer->superLayer();
00134       double flightToSL = sl->surface().position().mag();
00135       tofCorr = (flightToSL-flightToHit)/cSpeed;
00136       break;
00137     }
00138     case 1: {
00139       // On simulated data you need to consider only the TOF from 3D center of the wire to hit position 
00140       // (because the TOF from the IP to the wire has been already subtracted in the digitization: 
00141       // SimMuon/DTDigitizer/DTDigiSyncTOFCorr.cc corrType=2)
00142       float flightToWire =
00143         layer->toGlobal(LocalPoint(layer->specificTopology().wirePosition(wireId.wire()), 0., 0.)).mag();
00144       tofCorr = (flightToWire-flightToHit)/cSpeed;
00145       break;
00146     }
00147     default: {
00148       throw
00149         cms::Exception("[DTTTrigSyncFromDB]") << " Invalid parameter: tofCorrType = "
00150                                               << theTOFCorrType 
00151                                               << std::endl;
00152       break;
00153     }
00154     }
00155   }
00156 
00157 
00158   if(debug) {
00159     cout << "[DTTTrigSyncFromDB] Channel: " << wireId << endl
00160          << "      Offset (ns): " << tTrig + wirePropCorr - tofCorr << endl
00161          << "      various contributions are: " << endl
00162          << "      tTrig + t0 (ns):   " << tTrig << endl
00163          //<< "      tZero (ns):   " << t0 << endl
00164          << "      Propagation along wire delay (ns): " <<  wirePropCorr << endl
00165          << "      TOF correction (ns): " << tofCorr << endl
00166          << endl;
00167   }
00168   //The global offset is the sum of various contributions
00169     return tTrig + wirePropCorr - tofCorr;
00170 }
00171 
00172 double DTTTrigSyncFromDB::offset(const DTWireId& wireId) {
00173   float t0 = 0;
00174   float t0rms = 0;
00175   if(doT0Correction)
00176     {
00177       // Read the t0 from pulses for this wire (ns)
00178        tZeroMap->get(wireId,
00179                      t0,
00180                      t0rms,
00181                      DTTimeUnits::ns);
00182 
00183     }
00184 
00185   // Read the ttrig for this wire
00186   float ttrigMean = 0;
00187   float ttrigSigma = 0;
00188   float kFactor = 0;
00189   // FIXME: should check the return value of the DTTtrigRcd::get(..) method
00190   if(tTrigMap->get(wireId.superlayerId(),
00191                    ttrigMean,
00192                    ttrigSigma,
00193                    kFactor,
00194                    DTTimeUnits::ns) != 0) {
00195     cout << "[DTTTrigSyncFromDB]*Error: ttrig not found for SL: " << wireId.superlayerId() << endl;
00196 //     FIXME: LogError.....
00197   }
00198   
00199   return t0 + ttrigMean + kFactor * ttrigSigma;
00200 
00201 }
00202 
00203 
00204 // Set the verbosity level
00205 bool DTTTrigSyncFromDB::debug = false;
00206 
00207 
00208 double DTTTrigSyncFromDB::emulatorOffset(const DTWireId& wireId,
00209                                          double &tTrig,
00210                                          double &t0cell) {
00211   float t0 = 0;
00212   float t0rms = 0;
00213   if(doT0Correction)
00214     {
00215       // Read the t0 from pulses for this wire (ns)
00216        tZeroMap->get(wireId,
00217                      t0,
00218                      t0rms,
00219                      DTTimeUnits::ns);
00220 
00221     }
00222 
00223   // Read the ttrig for this wire
00224   float ttrigMean = 0;
00225   float ttrigSigma = 0;
00226   float kFactor = 0;
00227   // FIXME: should check the return value of the DTTtrigRcd::get(..) method
00228   if(tTrigMap->get(wireId.superlayerId(),
00229                    ttrigMean,
00230                    ttrigSigma,
00231                    kFactor,
00232                    DTTimeUnits::ns) != 0) {
00233     cout << "[DTTTrigSyncFromDB]*Error: ttrig not found for SL: " << wireId.superlayerId() << endl;
00234 //     FIXME: LogError.....
00235   }
00236   
00237   tTrig = ttrigMean + kFactor * ttrigSigma;
00238   t0cell = t0;
00239 
00240   return int(tTrig/theBXspace)*theBXspace + t0cell;
00241 }