Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:10

0001 
0002 /******* \class DTRunConditionVar *******
0003  *
0004  * Description:
0005  *
0006  *  detailed description
0007  *
0008  * \author : Paolo Bellan, Antonio Branca
0009  * $date   : 23/09/2011 15:42:04 CET $
0010  *
0011  * Modification:
0012  *
0013  *********************************/
0014 
0015 #include "DTRunConditionVar.h"
0016 
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/EventSetup.h"
0019 #include "FWCore/ServiceRegistry/interface/Service.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 
0022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0023 
0024 #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
0025 
0026 #include "DQMServices/Core/interface/DQMStore.h"
0027 
0028 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0029 
0030 #include "RecoMuon/Navigation/interface/DirectMuonNavigation.h"
0031 
0032 #include "DataFormats/DetId/interface/DetId.h"
0033 #include "TrackingTools/DetLayers/interface/DetLayer.h"
0034 #include "DataFormats/Common/interface/RefToBase.h"
0035 
0036 #include "TMath.h"
0037 #include <cmath>
0038 
0039 using namespace std;
0040 using namespace edm;
0041 
0042 DTRunConditionVar::DTRunConditionVar(const ParameterSet& pSet)
0043     :  // Get the debug parameter for verbose output
0044       debug(pSet.getUntrackedParameter<bool>("debug", false)),
0045       nMinHitsPhi(pSet.getUntrackedParameter<int>("nMinHitsPhi")),
0046       maxAnglePhiSegm(pSet.getUntrackedParameter<double>("maxAnglePhiSegm")),
0047       dt4DSegmentsToken_(consumes<DTRecSegment4DCollection>(pSet.getUntrackedParameter<InputTag>("recoSegments"))),
0048       muonGeomToken_(esConsumes<edm::Transition::BeginRun>()),
0049       readLegacyVDriftDB(pSet.getParameter<bool>("readLegacyVDriftDB")) {
0050   if (readLegacyVDriftDB) {
0051     mTimeToken_ = esConsumes();
0052   } else {
0053     vDriftToken_ = esConsumes();
0054   }
0055 }
0056 
0057 DTRunConditionVar::~DTRunConditionVar() {
0058   LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar") << "DTRunConditionVar: destructor called";
0059 
0060   // free memory
0061 }
0062 
0063 void DTRunConditionVar::bookHistograms(DQMStore::IBooker& ibooker,
0064                                        edm::Run const& iRun,
0065                                        edm::EventSetup const& /* iSetup */) {
0066   LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar") << "DTRunConditionVar: bookHistograms";
0067 
0068   for (int wheel = -2; wheel <= 2; wheel++) {
0069     for (int sec = 1; sec <= 14; sec++) {
0070       for (int stat = 1; stat <= 4; stat++) {
0071         bookChamberHistos(ibooker, DTChamberId(wheel, stat, sec), "VDrift_FromSegm", 100, 0.0043, 0.0065);
0072         bookChamberHistos(ibooker, DTChamberId(wheel, stat, sec), "T0_FromSegm", 100, -25., 25.);
0073       }
0074     }
0075   }
0076 
0077   return;
0078 }
0079 
0080 void DTRunConditionVar::dqmBeginRun(const Run& run, const EventSetup& setup) {
0081   // Get the DT Geometry
0082   dtGeom = &setup.getData(muonGeomToken_);
0083   return;
0084 }
0085 
0086 void DTRunConditionVar::analyze(const Event& event, const EventSetup& eventSetup) {
0087   LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
0088       << "--- [DTRunConditionVar] Event analysed #Run: " << event.id().run() << " #Event: " << event.id().event()
0089       << endl;
0090 
0091   // Get the map of vdrift from the setup
0092   if (readLegacyVDriftDB) {
0093     mTimeMap_ = &eventSetup.getData(mTimeToken_);
0094     vDriftMap_ = nullptr;
0095   } else {
0096     vDriftMap_ = &eventSetup.getData(vDriftToken_);
0097     mTimeMap_ = nullptr;
0098     // Consistency check: no parametrization is implemented for the time being
0099     int version = vDriftMap_->version();
0100     if (version != 1) {
0101       throw cms::Exception("Configuration") << "only version 1 is presently supported for VDriftDB";
0102     }
0103   }
0104 
0105   // Get the segment collection from the event
0106   Handle<DTRecSegment4DCollection> all4DSegments;
0107   event.getByToken(dt4DSegmentsToken_, all4DSegments);
0108 
0109   // Loop over the segments
0110   for (DTRecSegment4DCollection::const_iterator segment = all4DSegments->begin(); segment != all4DSegments->end();
0111        ++segment) {
0112     // Get the chamber from the setup
0113     DTChamberId DTid = (DTChamberId)segment->chamberId();
0114     uint32_t indexCh = DTid.rawId();
0115 
0116     // Fill v-drift values
0117     if ((*segment).hasPhi()) {
0118       int nHitsPhi = (*segment).phiSegment()->degreesOfFreedom() + 2;
0119       double xdir = (*segment).phiSegment()->localDirection().x();
0120       double zdir = (*segment).phiSegment()->localDirection().z();
0121 
0122       double anglePhiSegm = fabs(atan(xdir / zdir)) * 180. / TMath::Pi();
0123 
0124       if (nHitsPhi >= nMinHitsPhi && anglePhiSegm <= maxAnglePhiSegm) {
0125         double segmentVDrift = segment->phiSegment()->vDrift();
0126 
0127         DTSuperLayerId indexSLPhi1(DTid, 1);
0128         DTSuperLayerId indexSLPhi2(DTid, 3);
0129 
0130         float vDriftPhi1(0.), vDriftPhi2(0.);
0131         float ResPhi1(0.), ResPhi2(0.);
0132         if (readLegacyVDriftDB) {  // Legacy format
0133           int status1 = mTimeMap_->get(indexSLPhi1, vDriftPhi1, ResPhi1, DTVelocityUnits::cm_per_ns);
0134           int status2 = mTimeMap_->get(indexSLPhi2, vDriftPhi2, ResPhi2, DTVelocityUnits::cm_per_ns);
0135 
0136           if (status1 != 0 || status2 != 0) {
0137             DTSuperLayerId sl = (status1 != 0) ? indexSLPhi1 : indexSLPhi2;
0138             throw cms::Exception("DTRunConditionVarClient") << "Could not find vDrift entry in DB for" << sl << endl;
0139           }
0140         } else {
0141           vDriftPhi1 = vDriftMap_->get(DTWireId(indexSLPhi1.rawId()));
0142           vDriftPhi2 = vDriftMap_->get(DTWireId(indexSLPhi2.rawId()));
0143         }
0144 
0145         float vDriftMed = (vDriftPhi1 + vDriftPhi2) / 2.;
0146 
0147         segmentVDrift = vDriftMed * (1. - segmentVDrift);
0148 
0149         double segmentT0 = segment->phiSegment()->t0();
0150 
0151         if (segment->phiSegment()->ist0Valid())
0152           (chamberHistos[indexCh])["T0_FromSegm"]->Fill(segmentT0);
0153         if (segmentVDrift != vDriftMed)
0154           (chamberHistos[indexCh])["VDrift_FromSegm"]->Fill(segmentVDrift);
0155       }
0156     }
0157 
0158   }  //end loop on segment
0159 
0160 }  //end analyze
0161 
0162 void DTRunConditionVar::bookChamberHistos(
0163     DQMStore::IBooker& ibooker, const DTChamberId& dtCh, string histoType, int nbins, float min, float max) {
0164   int wh = dtCh.wheel();
0165   int sc = dtCh.sector();
0166   int st = dtCh.station();
0167   stringstream wheel;
0168   wheel << wh;
0169   stringstream station;
0170   station << st;
0171   stringstream sector;
0172   sector << sc;
0173 
0174   string bookingFolder = "DT/02-Segments/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str();
0175   string histoTag = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
0176 
0177   ibooker.setCurrentFolder(bookingFolder);
0178 
0179   LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
0180       << "[DTRunConditionVar]: booking histos in " << bookingFolder << endl;
0181 
0182   string histoName = histoType + histoTag;
0183   const string& histoLabel = histoType;
0184 
0185   (chamberHistos[dtCh.rawId()])[histoType] = ibooker.book1D(histoName, histoLabel, nbins, min, max);
0186 }
0187 
0188 // Local Variables:
0189 // show-trailing-whitespace: t
0190 // truncate-lines: t
0191 // End: