Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-03 00:12:44

0001 // -*- C++ -*-
0002 //
0003 // Package:    Validation/MtdValidation
0004 // Class:      BtlLocalRecoValidation
0005 //
0006 /**\class BtlLocalRecoValidation BtlLocalRecoValidation.cc Validation/MtdValidation/plugins/BtlLocalRecoValidation.cc
0007 
0008  Description: BTL RECO hits and clusters validation
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 
0014 #include <string>
0015 
0016 #include "FWCore/Framework/interface/Frameworkfwd.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 
0021 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0022 #include "DQMServices/Core/interface/DQMStore.h"
0023 
0024 #include "DataFormats/Common/interface/ValidHandle.h"
0025 #include "DataFormats/Math/interface/GeantUnits.h"
0026 #include "DataFormats/ForwardDetId/interface/BTLDetId.h"
0027 #include "DataFormats/FTLRecHit/interface/FTLRecHitCollections.h"
0028 #include "DataFormats/FTLRecHit/interface/FTLClusterCollections.h"
0029 #include "DataFormats/TrackerRecHit2D/interface/MTDTrackingRecHit.h"
0030 
0031 #include "SimDataFormats/CaloAnalysis/interface/MtdSimLayerCluster.h"
0032 #include "SimDataFormats/Associations/interface/MtdRecoClusterToSimLayerClusterAssociationMap.h"
0033 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
0034 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
0035 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0036 
0037 #include "Geometry/Records/interface/MTDDigiGeometryRecord.h"
0038 #include "Geometry/Records/interface/MTDTopologyRcd.h"
0039 #include "Geometry/MTDGeometryBuilder/interface/MTDGeometry.h"
0040 #include "Geometry/MTDGeometryBuilder/interface/MTDTopology.h"
0041 
0042 #include "Geometry/MTDGeometryBuilder/interface/ProxyMTDTopology.h"
0043 #include "Geometry/MTDGeometryBuilder/interface/RectangularMTDTopology.h"
0044 
0045 #include "Geometry/MTDCommonData/interface/MTDTopologyMode.h"
0046 
0047 #include "RecoLocalFastTime/Records/interface/MTDCPERecord.h"
0048 #include "RecoLocalFastTime/FTLClusterizer/interface/MTDClusterParameterEstimator.h"
0049 
0050 #include "MTDHit.h"
0051 
0052 class BtlLocalRecoValidation : public DQMEDAnalyzer {
0053 public:
0054   explicit BtlLocalRecoValidation(const edm::ParameterSet&);
0055   ~BtlLocalRecoValidation() override;
0056 
0057   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0058 
0059 private:
0060   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0061 
0062   void analyze(const edm::Event&, const edm::EventSetup&) override;
0063 
0064   bool isSameCluster(const FTLCluster&, const FTLCluster&);
0065 
0066   // ------------ member data ------------
0067 
0068   const std::string folder_;
0069   const double hitMinEnergy_;
0070   const bool optionalPlots_;
0071   const bool uncalibRecHitsPlots_;
0072   const double hitMinAmplitude_;
0073 
0074   edm::EDGetTokenT<FTLRecHitCollection> btlRecHitsToken_;
0075   edm::EDGetTokenT<FTLUncalibratedRecHitCollection> btlUncalibRecHitsToken_;
0076   edm::EDGetTokenT<CrossingFrame<PSimHit>> btlSimHitsToken_;
0077   edm::EDGetTokenT<FTLClusterCollection> btlRecCluToken_;
0078   edm::EDGetTokenT<MTDTrackingDetSetVector> mtdTrackingHitToken_;
0079   edm::EDGetTokenT<MtdRecoClusterToSimLayerClusterAssociationMap> r2sAssociationMapToken_;
0080 
0081   const edm::ESGetToken<MTDGeometry, MTDDigiGeometryRecord> mtdgeoToken_;
0082   const edm::ESGetToken<MTDTopology, MTDTopologyRcd> mtdtopoToken_;
0083   const edm::ESGetToken<MTDClusterParameterEstimator, MTDCPERecord> cpeToken_;
0084 
0085   // --- histograms declaration
0086 
0087   MonitorElement* meNevents_;
0088 
0089   MonitorElement* meNhits_;
0090 
0091   static constexpr int nRU_ = BTLDetId::kRUPerRod;
0092 
0093   MonitorElement* meHitEnergy_;
0094   MonitorElement* meHitEnergyRUSlice_[nRU_];
0095   MonitorElement* meHitLogEnergy_;
0096   MonitorElement* meHitTime_;
0097   MonitorElement* meHitTimeError_;
0098 
0099   MonitorElement* meOccupancy_;
0100 
0101   //local position monitoring
0102   MonitorElement* meLocalOccupancy_;
0103   MonitorElement* meHitXlocal_;
0104   MonitorElement* meHitYlocal_;
0105   MonitorElement* meHitZlocal_;
0106 
0107   MonitorElement* meHitZ_;
0108   MonitorElement* meHitPhi_;
0109   MonitorElement* meHitEta_;
0110 
0111   MonitorElement* meHitTvsE_;
0112   MonitorElement* meHitEvsPhi_;
0113   MonitorElement* meHitEvsEta_;
0114   MonitorElement* meHitEvsZ_;
0115   MonitorElement* meHitTvsPhi_;
0116   MonitorElement* meHitTvsEta_;
0117   MonitorElement* meHitTvsZ_;
0118   MonitorElement* meHitLongPos_;
0119 
0120   MonitorElement* meTimeRes_;
0121   MonitorElement* meTimeResVsE_;
0122   MonitorElement* meEnergyRes_;
0123   MonitorElement* meEnergyRelResVsE_;
0124 
0125   MonitorElement* meLongPosPull_;
0126   MonitorElement* meLongPosPullvsE_;
0127   MonitorElement* meLongPosPullvsEta_;
0128 
0129   MonitorElement* meTPullvsE_;
0130   MonitorElement* meTPullvsEta_;
0131   MonitorElement* meUnmatchedRecHit_;
0132 
0133   MonitorElement* meNclusters_;
0134 
0135   MonitorElement* meCluTime_;
0136   MonitorElement* meCluTimeError_;
0137   MonitorElement* meCluEnergy_;
0138   MonitorElement* meCluPhi_;
0139   MonitorElement* meCluEta_;
0140   MonitorElement* meCluHits_;
0141   MonitorElement* meCluZvsPhi_;
0142   MonitorElement* meCluEnergyvsEta_;
0143   MonitorElement* meCluHitsvsEta_;
0144 
0145   MonitorElement* meCluTimeRes_;
0146   MonitorElement* meCluEnergyRes_;
0147   MonitorElement* meCluTResvsE_;
0148   MonitorElement* meCluTResvsEta_;
0149   MonitorElement* meCluTPullvsE_;
0150   MonitorElement* meCluTPullvsEta_;
0151   MonitorElement* meCluRhoRes_;
0152   MonitorElement* meCluPhiRes_;
0153   MonitorElement* meCluLocalXRes_;
0154 
0155   MonitorElement* meCluLocalYResZGlobPlus_;
0156   MonitorElement* meCluLocalYResZGlobMinus_;
0157 
0158   MonitorElement* meCluZRes_;
0159   MonitorElement* meCluLocalXPull_;
0160 
0161   MonitorElement* meCluLocalYPullZGlobPlus_;
0162   MonitorElement* meCluLocalYPullZGlobMinus_;
0163 
0164   MonitorElement* meCluSingCrystalLocalYRes_;
0165   MonitorElement* meCluSingCrystalLocalYResZGlobPlus_;
0166   MonitorElement* meCluSingCrystalLocalYResZGlobMinus_;
0167 
0168   MonitorElement* meCluMultiCrystalLocalYRes_;
0169   MonitorElement* meCluMultiCrystalLocalYResZGlobPlus_;
0170   MonitorElement* meCluMultiCrystalLocalYResZGlobMinus_;
0171 
0172   MonitorElement* meCluCentralLocalYRes_;
0173   MonitorElement* meCluCentralLocalYResZGlobPlus_;
0174   MonitorElement* meCluCentralLocalYResZGlobMinus_;
0175 
0176   MonitorElement* meCluForwardLocalYRes_;
0177   MonitorElement* meCluForwardPlusLocalYRes_;
0178   MonitorElement* meCluForwardMinusLocalYRes_;
0179 
0180   MonitorElement* meCluZPull_;
0181   MonitorElement* meCluYXLocal_;
0182   MonitorElement* meCluYXLocalSim_;
0183   MonitorElement* meCluXLocalErr_;
0184   MonitorElement* meCluYLocalErr_;
0185 
0186   // resolution wrt to MtdSimLayerClusters
0187   MonitorElement* meCluTrackIdOffset_;
0188 
0189   MonitorElement* meCluTimeRes_simLC_;
0190   MonitorElement* meCluEnergyRes_simLC_;
0191   MonitorElement* meCluTResvsE_simLC_;
0192   MonitorElement* meCluTResvsEta_simLC_;
0193   MonitorElement* meCluTPullvsE_simLC_;
0194   MonitorElement* meCluTPullvsEta_simLC_;
0195   MonitorElement* meCluRhoRes_simLC_;
0196   MonitorElement* meCluPhiRes_simLC_;
0197   MonitorElement* meCluLocalXRes_simLC_;
0198 
0199   MonitorElement* meCluLocalYResZGlobPlus_simLC_;
0200   MonitorElement* meCluLocalYResZGlobMinus_simLC_;
0201 
0202   MonitorElement* meCluZRes_simLC_;
0203   MonitorElement* meCluLocalXPull_simLC_;
0204 
0205   MonitorElement* meCluLocalYPullZGlobPlus_simLC_;
0206   MonitorElement* meCluLocalYPullZGlobMinus_simLC_;
0207 
0208   MonitorElement* meCluSingCrystalLocalYRes_simLC_;
0209   MonitorElement* meCluSingCrystalLocalYResZGlobPlus_simLC_;
0210   MonitorElement* meCluSingCrystalLocalYResZGlobMinus_simLC_;
0211 
0212   MonitorElement* meCluMultiCrystalLocalYRes_simLC_;
0213   MonitorElement* meCluMultiCrystalLocalYResZGlobPlus_simLC_;
0214   MonitorElement* meCluMultiCrystalLocalYResZGlobMinus_simLC_;
0215 
0216   MonitorElement* meCluCentralLocalYRes_simLC_;
0217   MonitorElement* meCluCentralLocalYResZGlobPlus_simLC_;
0218   MonitorElement* meCluCentralLocalYResZGlobMinus_simLC_;
0219 
0220   MonitorElement* meCluForwardLocalYRes_simLC_;
0221   MonitorElement* meCluForwardPlusLocalYRes_simLC_;
0222   MonitorElement* meCluForwardMinusLocalYRes_simLC_;
0223 
0224   MonitorElement* meCluZPull_simLC_;
0225   MonitorElement* meCluYXLocalSim_simLC_;
0226 
0227   MonitorElement* meCluTimeRes_simLC_fromIndirectHits_;
0228   MonitorElement* meCluEnergyRes_simLC_fromIndirectHits_;
0229   MonitorElement* meCluTResvsE_simLC_fromIndirectHits_;
0230   MonitorElement* meCluTResvsEta_simLC_fromIndirectHits_;
0231   MonitorElement* meCluTPullvsE_simLC_fromIndirectHits_;
0232   MonitorElement* meCluTPullvsEta_simLC_fromIndirectHits_;
0233   MonitorElement* meCluRhoRes_simLC_fromIndirectHits_;
0234   MonitorElement* meCluPhiRes_simLC_fromIndirectHits_;
0235   MonitorElement* meCluLocalXRes_simLC_fromIndirectHits_;
0236 
0237   MonitorElement* meCluLocalYResZGlobPlus_simLC_fromIndirectHits_;
0238   MonitorElement* meCluLocalYResZGlobMinus_simLC_fromIndirectHits_;
0239 
0240   MonitorElement* meCluZRes_simLC_fromIndirectHits_;
0241   MonitorElement* meCluLocalXPull_simLC_fromIndirectHits_;
0242 
0243   MonitorElement* meCluLocalYPullZGlobPlus_simLC_fromIndirectHits_;
0244   MonitorElement* meCluLocalYPullZGlobMinus_simLC_fromIndirectHits_;
0245 
0246   MonitorElement* meCluSingCrystalLocalYRes_simLC_fromIndirectHits_;
0247   MonitorElement* meCluSingCrystalLocalYResZGlobPlus_simLC_fromIndirectHits_;
0248   MonitorElement* meCluSingCrystalLocalYResZGlobMinus_simLC_fromIndirectHits_;
0249 
0250   MonitorElement* meCluMultiCrystalLocalYRes_simLC_fromIndirectHits_;
0251   MonitorElement* meCluMultiCrystalLocalYResZGlobPlus_simLC_fromIndirectHits_;
0252   MonitorElement* meCluMultiCrystalLocalYResZGlobMinus_simLC_fromIndirectHits_;
0253 
0254   MonitorElement* meCluCentralLocalYRes_simLC_fromIndirectHits_;
0255   MonitorElement* meCluCentralLocalYResZGlobPlus_simLC_fromIndirectHits_;
0256   MonitorElement* meCluCentralLocalYResZGlobMinus_simLC_fromIndirectHits_;
0257 
0258   MonitorElement* meCluForwardLocalYRes_simLC_fromIndirectHits_;
0259   MonitorElement* meCluForwardPlusLocalYRes_simLC_fromIndirectHits_;
0260   MonitorElement* meCluForwardMinusLocalYRes_simLC_fromIndirectHits_;
0261 
0262   MonitorElement* meCluZPull_simLC_fromIndirectHits_;
0263   MonitorElement* meCluYXLocalSim_simLC_fromIndirectHits_;
0264 
0265   MonitorElement* meUnmatchedCluEnergy_;
0266 
0267   // --- UncalibratedRecHits histograms
0268 
0269   MonitorElement* meUncEneLVsX_;
0270   MonitorElement* meUncEneRVsX_;
0271   MonitorElement* meUncTimeLVsX_;
0272   MonitorElement* meUncTimeRVsX_;
0273 
0274   static constexpr int nBinsQ_ = 20;
0275   static constexpr float binWidthQ_ = 30.;
0276   static constexpr int nBinsQEta_ = 3;
0277   static constexpr float binsQEta_[nBinsQEta_ + 1] = {0., 0.65, 1.15, 1.55};
0278 
0279   MonitorElement* meTimeResQ_[nBinsQ_];
0280   MonitorElement* meTimeResQvsEta_[nBinsQ_][nBinsQEta_];
0281 
0282   static constexpr int nBinsEta_ = 31;
0283   static constexpr float binWidthEta_ = 0.05;
0284   static constexpr int nBinsEtaQ_ = 7;
0285   static constexpr float binsEtaQ_[nBinsEtaQ_ + 1] = {0., 30., 60., 90., 120., 150., 360., 600.};
0286 
0287   MonitorElement* meTimeResEta_[nBinsEta_];
0288   MonitorElement* meTimeResEtavsQ_[nBinsEta_][nBinsEtaQ_];
0289 };
0290 
0291 bool BtlLocalRecoValidation::isSameCluster(const FTLCluster& clu1, const FTLCluster& clu2) {
0292   return clu1.id() == clu2.id() && clu1.size() == clu2.size() && clu1.x() == clu2.x() && clu1.y() == clu2.y() &&
0293          clu1.time() == clu2.time();
0294 }
0295 
0296 // ------------ constructor and destructor --------------
0297 BtlLocalRecoValidation::BtlLocalRecoValidation(const edm::ParameterSet& iConfig)
0298     : folder_(iConfig.getParameter<std::string>("folder")),
0299       hitMinEnergy_(iConfig.getParameter<double>("HitMinimumEnergy")),
0300       optionalPlots_(iConfig.getParameter<bool>("optionalPlots")),
0301       uncalibRecHitsPlots_(iConfig.getParameter<bool>("UncalibRecHitsPlots")),
0302       hitMinAmplitude_(iConfig.getParameter<double>("HitMinimumAmplitude")),
0303       mtdgeoToken_(esConsumes<MTDGeometry, MTDDigiGeometryRecord>()),
0304       mtdtopoToken_(esConsumes<MTDTopology, MTDTopologyRcd>()),
0305       cpeToken_(esConsumes<MTDClusterParameterEstimator, MTDCPERecord>(edm::ESInputTag("", "MTDCPEBase"))) {
0306   btlRecHitsToken_ = consumes<FTLRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsTag"));
0307   btlUncalibRecHitsToken_ =
0308       consumes<FTLUncalibratedRecHitCollection>(iConfig.getParameter<edm::InputTag>("uncalibRecHitsTag"));
0309   btlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("simHitsTag"));
0310   btlRecCluToken_ = consumes<FTLClusterCollection>(iConfig.getParameter<edm::InputTag>("recCluTag"));
0311   mtdTrackingHitToken_ = consumes<MTDTrackingDetSetVector>(iConfig.getParameter<edm::InputTag>("trkHitTag"));
0312   r2sAssociationMapToken_ = consumes<MtdRecoClusterToSimLayerClusterAssociationMap>(
0313       iConfig.getParameter<edm::InputTag>("r2sAssociationMapTag"));
0314 }
0315 
0316 BtlLocalRecoValidation::~BtlLocalRecoValidation() {}
0317 
0318 // ------------ method called for each event  ------------
0319 void BtlLocalRecoValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0320   using namespace edm;
0321   using namespace std;
0322   using namespace geant_units::operators;
0323 
0324   auto geometryHandle = iSetup.getTransientHandle(mtdgeoToken_);
0325   const MTDGeometry* geom = geometryHandle.product();
0326 
0327   auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
0328   const MTDTopology* topology = topologyHandle.product();
0329 
0330   auto const& cpe = iSetup.getData(cpeToken_);
0331 
0332   auto btlRecHitsHandle = makeValid(iEvent.getHandle(btlRecHitsToken_));
0333   auto btlSimHitsHandle = makeValid(iEvent.getHandle(btlSimHitsToken_));
0334   auto btlRecCluHandle = makeValid(iEvent.getHandle(btlRecCluToken_));
0335   auto mtdTrkHitHandle = makeValid(iEvent.getHandle(mtdTrackingHitToken_));
0336   const auto& r2sAssociationMap = iEvent.get(r2sAssociationMapToken_);
0337   MixCollection<PSimHit> btlSimHits(btlSimHitsHandle.product());
0338 
0339 #ifdef EDM_ML_DEBUG
0340   for (const auto& hits : *mtdTrkHitHandle) {
0341     if (MTDDetId(hits.id()).mtdSubDetector() == MTDDetId::MTDType::BTL) {
0342       LogDebug("BtlLocalRecoValidation") << "MTD cluster DetId " << hits.id() << " # cluster " << hits.size();
0343       for (const auto& hit : hits) {
0344         LogDebug("BtlLocalRecoValidation")
0345             << "MTD_TRH: " << hit.localPosition().x() << "," << hit.localPosition().y() << " : "
0346             << hit.localPositionError().xx() << "," << hit.localPositionError().yy() << " : " << hit.time() << " : "
0347             << hit.timeError();
0348       }
0349     }
0350   }
0351 #endif
0352 
0353   // --- Loop over the BTL SIM hits
0354   std::unordered_map<uint32_t, MTDHit> m_btlSimHits;
0355   for (auto const& simHit : btlSimHits) {
0356     // --- Use only hits compatible with the in-time bunch-crossing
0357     if (simHit.tof() < 0 || simHit.tof() > 25.)
0358       continue;
0359 
0360     DetId id = simHit.detUnitId();
0361 
0362     auto simHitIt = m_btlSimHits.emplace(id.rawId(), MTDHit()).first;
0363 
0364     // --- Accumulate the energy (in MeV) of SIM hits in the same detector cell
0365     (simHitIt->second).energy += convertUnitsTo(0.001_MeV, simHit.energyLoss());
0366 
0367     // --- Get the time of the first SIM hit in the cell
0368     if ((simHitIt->second).time == 0 || simHit.tof() < (simHitIt->second).time) {
0369       (simHitIt->second).time = simHit.tof();
0370 
0371       auto hit_pos = simHit.localPosition();
0372       (simHitIt->second).x = hit_pos.x();
0373       (simHitIt->second).y = hit_pos.y();
0374       (simHitIt->second).z = hit_pos.z();
0375     }
0376 
0377   }  // simHit loop
0378 
0379   // --- Loop over the BTL RECO hits
0380   unsigned int n_reco_btl = 0;
0381   unsigned int n_reco_btl_nosimhit = 0;
0382   for (const auto& recHit : *btlRecHitsHandle) {
0383     LogTrace("BtlLocalRecoValidation") << "@RH detid " << recHit.id().rawId() << " r/c/X/dX " << recHit.row() << " "
0384                                        << recHit.column() << " " << recHit.position() << " " << recHit.positionError()
0385                                        << " E,T,dT " << recHit.energy() << " " << recHit.time() << " "
0386                                        << recHit.timeError();
0387 
0388     BTLDetId detId = recHit.id();
0389     DetId geoId = detId.geographicalId(MTDTopologyMode::crysLayoutFromTopoMode(topology->getMTDTopologyMode()));
0390     const MTDGeomDet* thedet = geom->idToDet(geoId);
0391     if (thedet == nullptr)
0392       throw cms::Exception("BtlLocalRecoValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
0393                                                      << detId.rawId() << ") is invalid!" << std::dec << std::endl;
0394     const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
0395     const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
0396 
0397     Local3DPoint local_point(0., 0., 0.);
0398     local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
0399     const auto& global_point = thedet->toGlobal(local_point);
0400 
0401     meHitEnergy_->Fill(recHit.energy());
0402     meHitEnergyRUSlice_[detId.runit() - 1]->Fill(recHit.energy());
0403     meHitLogEnergy_->Fill(log10(recHit.energy()));
0404     meHitTime_->Fill(recHit.time());
0405     meHitTimeError_->Fill(recHit.timeError());
0406     meHitLongPos_->Fill(recHit.position());
0407 
0408     meOccupancy_->Fill(global_point.z(), global_point.phi());
0409 
0410     if (optionalPlots_) {
0411       meLocalOccupancy_->Fill(local_point.x() + recHit.position(), local_point.y());
0412     }
0413     meHitXlocal_->Fill(local_point.x());
0414     meHitYlocal_->Fill(local_point.y());
0415     meHitZlocal_->Fill(local_point.z());
0416     meHitZ_->Fill(global_point.z());
0417     meHitPhi_->Fill(global_point.phi());
0418     meHitEta_->Fill(global_point.eta());
0419 
0420     meHitTvsE_->Fill(recHit.energy(), recHit.time());
0421     meHitEvsPhi_->Fill(global_point.phi(), recHit.energy());
0422     meHitEvsEta_->Fill(global_point.eta(), recHit.energy());
0423     meHitEvsZ_->Fill(global_point.z(), recHit.energy());
0424     meHitTvsPhi_->Fill(global_point.phi(), recHit.time());
0425     meHitTvsEta_->Fill(global_point.eta(), recHit.time());
0426     meHitTvsZ_->Fill(global_point.z(), recHit.time());
0427 
0428     // Resolution histograms
0429     LogDebug("BtlLocalRecoValidation") << "RecoHit DetId= " << detId.rawId()
0430                                        << " sim hits in id= " << m_btlSimHits.count(detId.rawId());
0431     if (m_btlSimHits.count(detId.rawId()) == 1 && m_btlSimHits[detId.rawId()].energy > hitMinEnergy_) {
0432       float longpos_res = recHit.position() - convertMmToCm(m_btlSimHits[detId.rawId()].x);
0433       float time_res = recHit.time() - m_btlSimHits[detId.rawId()].time;
0434       float energy_res = recHit.energy() - m_btlSimHits[detId.rawId()].energy;
0435 
0436       Local3DPoint local_point_sim(convertMmToCm(m_btlSimHits[detId.rawId()].x),
0437                                    convertMmToCm(m_btlSimHits[detId.rawId()].y),
0438                                    convertMmToCm(m_btlSimHits[detId.rawId()].z));
0439       local_point_sim =
0440           topo.pixelToModuleLocalPoint(local_point_sim, detId.row(topo.nrows()), detId.column(topo.nrows()));
0441       const auto& global_point_sim = thedet->toGlobal(local_point_sim);
0442 
0443       meTimeRes_->Fill(time_res);
0444       meTimeResVsE_->Fill(recHit.energy(), time_res);
0445       meEnergyRes_->Fill(energy_res);
0446       meEnergyRelResVsE_->Fill(recHit.energy(), energy_res / recHit.energy());
0447 
0448       meLongPosPull_->Fill(longpos_res / recHit.positionError());
0449       meLongPosPullvsEta_->Fill(std::abs(global_point_sim.eta()), longpos_res / recHit.positionError());
0450       meLongPosPullvsE_->Fill(m_btlSimHits[detId.rawId()].energy, longpos_res / recHit.positionError());
0451 
0452       meTPullvsEta_->Fill(std::abs(global_point_sim.eta()), time_res / recHit.timeError());
0453       meTPullvsE_->Fill(m_btlSimHits[detId.rawId()].energy, time_res / recHit.timeError());
0454     } else if (m_btlSimHits.count(detId.rawId()) == 0) {
0455       n_reco_btl_nosimhit++;
0456       LogDebug("BtlLocalRecoValidation") << "BTL rec hit with no corresponding sim hit in crystal, DetId= "
0457                                          << detId.rawId() << " geoId= " << geoId.rawId() << " ene= " << recHit.energy()
0458                                          << " time= " << recHit.time();
0459     }
0460 
0461     n_reco_btl++;
0462 
0463   }  // recHit loop
0464 
0465   if (n_reco_btl > 0) {
0466     meNhits_->Fill(std::log10(n_reco_btl));
0467   }
0468   if (n_reco_btl_nosimhit == 0) {
0469     meUnmatchedRecHit_->Fill(-1.5);
0470   } else {
0471     meUnmatchedRecHit_->Fill(std::log10(n_reco_btl_nosimhit));
0472   }
0473 
0474   // --- Loop over the BTL RECO clusters ---
0475   unsigned int n_clus_btl(0);
0476   for (const auto& DetSetClu : *btlRecCluHandle) {
0477     for (const auto& cluster : DetSetClu) {
0478       if (cluster.energy() < hitMinEnergy_)
0479         continue;
0480       BTLDetId cluId = cluster.id();
0481       DetId detIdObject(cluId);
0482       const auto& genericDet = geom->idToDetUnit(detIdObject);
0483       if (genericDet == nullptr) {
0484         throw cms::Exception("BtlLocalRecoValidation")
0485             << "GeographicalID: " << std::hex << cluId << " is invalid!" << std::dec << std::endl;
0486       }
0487       n_clus_btl++;
0488 
0489       const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(genericDet->topology());
0490       const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
0491 
0492       MTDClusterParameterEstimator::ReturnType tuple = cpe.getParameters(cluster, *genericDet);
0493 
0494       // --- Cluster position in the module reference frame
0495       LocalPoint local_point(std::get<0>(tuple));
0496       const auto& global_point = genericDet->toGlobal(local_point);
0497 
0498       meCluEnergy_->Fill(cluster.energy());
0499       meCluTime_->Fill(cluster.time());
0500       meCluTimeError_->Fill(cluster.timeError());
0501       meCluPhi_->Fill(global_point.phi());
0502       meCluEta_->Fill(global_point.eta());
0503       meCluZvsPhi_->Fill(global_point.z(), global_point.phi());
0504       meCluHits_->Fill(cluster.size());
0505 
0506       // --- Get the SIM hits associated to the cluster and calculate
0507       //     the cluster SIM energy, time and position
0508 
0509       double cluEneSIM = 0.;
0510       double cluTimeSIM = 0.;
0511       double cluLocXSIM = 0.;
0512       double cluLocYSIM = 0.;
0513       double cluLocZSIM = 0.;
0514 
0515       if (optionalPlots_) {
0516         for (int ihit = 0; ihit < cluster.size(); ++ihit) {
0517           int hit_row = cluster.minHitRow() + cluster.hitOffset()[ihit * 2];
0518           int hit_col = cluster.minHitCol() + cluster.hitOffset()[ihit * 2 + 1];
0519 
0520           // Match the RECO hit to the corresponding SIM hit
0521           for (const auto& recHit : *btlRecHitsHandle) {
0522             BTLDetId hitId(recHit.id().rawId());
0523 
0524             if (m_btlSimHits.count(hitId.rawId()) == 0)
0525               continue;
0526 
0527             // Check the hit position
0528             if (hitId.mtdSide() != cluId.mtdSide() || hitId.mtdRR() != cluId.mtdRR() || recHit.row() != hit_row ||
0529                 recHit.column() != hit_col)
0530               continue;
0531 
0532             // Check the hit energy and time
0533             if (recHit.energy() != cluster.hitENERGY()[ihit] || recHit.time() != cluster.hitTIME()[ihit])
0534               continue;
0535 
0536             // SIM hit's position in the module reference frame
0537             Local3DPoint local_point_sim(convertMmToCm(m_btlSimHits[recHit.id().rawId()].x),
0538                                          convertMmToCm(m_btlSimHits[recHit.id().rawId()].y),
0539                                          convertMmToCm(m_btlSimHits[recHit.id().rawId()].z));
0540             local_point_sim =
0541                 topo.pixelToModuleLocalPoint(local_point_sim, hitId.row(topo.nrows()), hitId.column(topo.nrows()));
0542 
0543             // Calculate the SIM cluster's position in the module reference frame
0544             cluLocXSIM += local_point_sim.x() * m_btlSimHits[recHit.id().rawId()].energy;
0545             cluLocYSIM += local_point_sim.y() * m_btlSimHits[recHit.id().rawId()].energy;
0546             cluLocZSIM += local_point_sim.z() * m_btlSimHits[recHit.id().rawId()].energy;
0547 
0548             // Calculate the SIM cluster energy and time
0549             cluEneSIM += m_btlSimHits[recHit.id().rawId()].energy;
0550             cluTimeSIM += m_btlSimHits[recHit.id().rawId()].time * m_btlSimHits[recHit.id().rawId()].energy;
0551 
0552             break;
0553 
0554           }  // recHit loop
0555 
0556         }  // ihit loop
0557       }
0558 
0559       // Find the MTDTrackingRecHit corresponding to the cluster
0560       const MTDTrackingRecHit* comp(nullptr);
0561       bool matchClu = false;
0562       const auto& trkHits = (*mtdTrkHitHandle)[detIdObject];
0563       for (const auto& trkHit : trkHits) {
0564         if (isSameCluster(trkHit.mtdCluster(), cluster)) {
0565           comp = trkHit.clone();
0566           matchClu = true;
0567           break;
0568         }
0569       }
0570       if (!matchClu) {
0571         edm::LogWarning("BtlLocalRecoValidation")
0572             << "No valid TrackingRecHit corresponding to cluster, detId = " << detIdObject.rawId();
0573       }
0574 
0575       // --- Fill the cluster resolution histograms
0576       if (optionalPlots_) {
0577         if (cluTimeSIM > 0. && cluEneSIM > 0.) {
0578           cluTimeSIM /= cluEneSIM;
0579 
0580           Local3DPoint cluLocalPosSIM(cluLocXSIM / cluEneSIM, cluLocYSIM / cluEneSIM, cluLocZSIM / cluEneSIM);
0581           const auto& cluGlobalPosSIM = genericDet->toGlobal(cluLocalPosSIM);
0582 
0583           float time_res = cluster.time() - cluTimeSIM;
0584           float energy_res = cluster.energy() - cluEneSIM;
0585           meCluTimeRes_->Fill(time_res);
0586           meCluEnergyRes_->Fill(energy_res);
0587 
0588           float rho_res = global_point.perp() - cluGlobalPosSIM.perp();
0589           float phi_res = global_point.phi() - cluGlobalPosSIM.phi();
0590 
0591           meCluRhoRes_->Fill(rho_res);
0592           meCluPhiRes_->Fill(phi_res);
0593 
0594           float xlocal_res = local_point.x() - cluLocalPosSIM.x();
0595           float ylocal_res = local_point.y() - cluLocalPosSIM.y();
0596 
0597           float z_res = global_point.z() - cluGlobalPosSIM.z();
0598 
0599           meCluZRes_->Fill(z_res);
0600 
0601           if (matchClu && comp != nullptr) {
0602             meCluLocalXRes_->Fill(xlocal_res);
0603 
0604             if (global_point.z() > 0) {
0605               meCluLocalYResZGlobPlus_->Fill(ylocal_res);
0606               meCluLocalYPullZGlobPlus_->Fill(ylocal_res / std::sqrt(comp->localPositionError().yy()));
0607             } else {
0608               meCluLocalYResZGlobMinus_->Fill(ylocal_res);
0609               meCluLocalYPullZGlobMinus_->Fill(ylocal_res / std::sqrt(comp->localPositionError().yy()));
0610             }
0611             if (cluster.size() == 1) {  // single-crystal clusters
0612               meCluSingCrystalLocalYRes_->Fill(ylocal_res);
0613               if (global_point.z() > 0) {
0614                 meCluSingCrystalLocalYResZGlobPlus_->Fill(ylocal_res);
0615               } else {
0616                 meCluSingCrystalLocalYResZGlobMinus_->Fill(ylocal_res);
0617               }
0618             }  // end of single-crystal clusters
0619             else {
0620               if (cluster.size() > 1) {  // multi-crystal clusters
0621                 meCluMultiCrystalLocalYRes_->Fill(ylocal_res);
0622                 if (global_point.z() > 0) {
0623                   meCluMultiCrystalLocalYResZGlobPlus_->Fill(ylocal_res);
0624                 } else {
0625                   meCluMultiCrystalLocalYResZGlobMinus_->Fill(ylocal_res);
0626                 }
0627               }
0628             }  // end of multi-crystal clusters
0629 
0630             if (abs(global_point.eta()) < 0.3) {
0631               meCluCentralLocalYRes_->Fill(ylocal_res);
0632               if (global_point.z() > 0) {
0633                 meCluCentralLocalYResZGlobPlus_->Fill(ylocal_res);
0634               } else {
0635                 meCluCentralLocalYResZGlobMinus_->Fill(ylocal_res);
0636               }
0637 
0638             } else {
0639               if (abs(global_point.eta()) > 1) {
0640                 meCluForwardLocalYRes_->Fill(ylocal_res);
0641                 if (global_point.z() > 0) {
0642                   meCluForwardPlusLocalYRes_->Fill(ylocal_res);
0643                 } else {
0644                   meCluForwardMinusLocalYRes_->Fill(ylocal_res);
0645                 }
0646               }
0647             }
0648 
0649             meCluYXLocal_->Fill(local_point.x(), local_point.y());
0650             meCluYXLocalSim_->Fill(cluLocalPosSIM.x(), cluLocalPosSIM.y());
0651 
0652             meCluLocalXPull_->Fill(xlocal_res / std::sqrt(comp->localPositionError().xx()));
0653             meCluZPull_->Fill(z_res / std::sqrt(comp->globalPositionError().czz()));
0654             meCluXLocalErr_->Fill(std::sqrt(comp->localPositionError().xx()));
0655             meCluYLocalErr_->Fill(std::sqrt(comp->localPositionError().yy()));
0656           }
0657 
0658           meCluEnergyvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), cluster.energy());
0659           meCluHitsvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), cluster.size());
0660 
0661           meCluTResvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), time_res);
0662           meCluTResvsE_->Fill(cluEneSIM, time_res);
0663 
0664           meCluTPullvsEta_->Fill(std::abs(cluGlobalPosSIM.eta()), time_res / cluster.timeError());
0665           meCluTPullvsE_->Fill(cluEneSIM, time_res / cluster.timeError());
0666 
0667         }  // if ( cluTimeSIM > 0. &&  cluEneSIM > 0. )
0668         else {
0669           meUnmatchedCluEnergy_->Fill(std::log10(cluster.energy()));
0670         }
0671       }
0672 
0673       // --- Fill the cluster resolution histograms using MtdSimLayerClusters as mtd truth
0674       edm::Ref<edmNew::DetSetVector<FTLCluster>, FTLCluster> clusterRef = edmNew::makeRefTo(btlRecCluHandle, &cluster);
0675       auto itp = r2sAssociationMap.equal_range(clusterRef);
0676       if (itp.first != itp.second) {
0677         std::vector<MtdSimLayerClusterRef> simClustersRefs =
0678             (*itp.first).second;  // the range of itp.first, itp.second should be always 1
0679         for (unsigned int i = 0; i < simClustersRefs.size(); i++) {
0680           auto simClusterRef = simClustersRefs[i];
0681 
0682           float simClusEnergy = convertUnitsTo(0.001_MeV, (*simClusterRef).simLCEnergy());  // GeV --> MeV
0683           float simClusTime = (*simClusterRef).simLCTime();
0684           LocalPoint simClusLocalPos = (*simClusterRef).simLCPos();
0685           const auto& simClusGlobalPos = genericDet->toGlobal(simClusLocalPos);
0686           unsigned int idOffset = (*simClusterRef).trackIdOffset();
0687 
0688           float time_res = cluster.time() - simClusTime;
0689           float energy_res = cluster.energy() - simClusEnergy;
0690           float rho_res = global_point.perp() - simClusGlobalPos.perp();
0691           float phi_res = global_point.phi() - simClusGlobalPos.phi();
0692           float z_res = global_point.z() - simClusGlobalPos.z();
0693           float xlocal_res = local_point.x() - simClusLocalPos.x();
0694           float ylocal_res = local_point.y() - simClusLocalPos.y();
0695 
0696           meCluTrackIdOffset_->Fill(float(idOffset));
0697 
0698           // -- Fill for direct hits
0699           if (idOffset == 0) {
0700             meCluTimeRes_simLC_->Fill(time_res);
0701             meCluEnergyRes_simLC_->Fill(energy_res);
0702             meCluRhoRes_simLC_->Fill(rho_res);
0703             meCluPhiRes_simLC_->Fill(phi_res);
0704             meCluZRes_simLC_->Fill(z_res);
0705 
0706             if (matchClu && comp != nullptr) {
0707               meCluLocalXRes_simLC_->Fill(xlocal_res);
0708 
0709               if (global_point.z() > 0) {
0710                 meCluLocalYResZGlobPlus_simLC_->Fill(ylocal_res);
0711                 meCluLocalYPullZGlobPlus_simLC_->Fill(ylocal_res / std::sqrt(comp->localPositionError().yy()));
0712               } else {
0713                 meCluLocalYResZGlobMinus_simLC_->Fill(ylocal_res);
0714                 meCluLocalYPullZGlobMinus_simLC_->Fill(ylocal_res / std::sqrt(comp->localPositionError().yy()));
0715               }
0716               if (optionalPlots_) {
0717                 if (cluster.size() == 1) {  // single-crystal clusters
0718                   meCluSingCrystalLocalYRes_simLC_->Fill(ylocal_res);
0719                   if (global_point.z() > 0) {
0720                     meCluSingCrystalLocalYResZGlobPlus_simLC_->Fill(ylocal_res);
0721                   } else {
0722                     meCluSingCrystalLocalYResZGlobMinus_simLC_->Fill(ylocal_res);
0723                   }
0724                 }  // end of single-crystal clusters
0725                 else {
0726                   if (cluster.size() > 1) {  // multi-crystal clusters
0727                     meCluMultiCrystalLocalYRes_simLC_->Fill(ylocal_res);
0728                     if (global_point.z() > 0) {
0729                       meCluMultiCrystalLocalYResZGlobPlus_simLC_->Fill(ylocal_res);
0730                     } else {
0731                       meCluMultiCrystalLocalYResZGlobMinus_simLC_->Fill(ylocal_res);
0732                     }
0733                   }
0734                 }  // end of multi-crystal clusters
0735 
0736                 if (abs(global_point.eta()) < 0.3) {
0737                   meCluCentralLocalYRes_simLC_->Fill(ylocal_res);
0738                   if (global_point.z() > 0) {
0739                     meCluCentralLocalYResZGlobPlus_simLC_->Fill(ylocal_res);
0740                   } else {
0741                     meCluCentralLocalYResZGlobMinus_simLC_->Fill(ylocal_res);
0742                   }
0743                 } else {
0744                   if (abs(global_point.eta()) > 1) {
0745                     meCluForwardLocalYRes_simLC_->Fill(ylocal_res);
0746                     if (global_point.z() > 0) {
0747                       meCluForwardPlusLocalYRes_simLC_->Fill(ylocal_res);
0748                     } else {
0749                       meCluForwardMinusLocalYRes_simLC_->Fill(ylocal_res);
0750                     }
0751                   }
0752                 }
0753               }  //end of optional plots
0754 
0755               meCluLocalXPull_simLC_->Fill(xlocal_res / std::sqrt(comp->localPositionError().xx()));
0756               meCluZPull_simLC_->Fill(z_res / std::sqrt(comp->globalPositionError().czz()));
0757             }
0758 
0759             meCluTResvsEta_simLC_->Fill(std::abs(simClusGlobalPos.eta()), time_res);
0760             meCluTResvsE_simLC_->Fill(simClusEnergy, time_res);
0761 
0762             meCluTPullvsEta_simLC_->Fill(std::abs(simClusGlobalPos.eta()), time_res / cluster.timeError());
0763             meCluTPullvsE_simLC_->Fill(simClusEnergy, time_res / cluster.timeError());
0764 
0765           }  // if idOffset == 0
0766           else {
0767             meCluTimeRes_simLC_fromIndirectHits_->Fill(time_res);
0768             meCluEnergyRes_simLC_fromIndirectHits_->Fill(energy_res);
0769             meCluRhoRes_simLC_fromIndirectHits_->Fill(rho_res);
0770             meCluPhiRes_simLC_fromIndirectHits_->Fill(phi_res);
0771             meCluZRes_simLC_fromIndirectHits_->Fill(z_res);
0772 
0773             if (matchClu && comp != nullptr) {
0774               meCluLocalXRes_simLC_fromIndirectHits_->Fill(xlocal_res);
0775 
0776               if (global_point.z() > 0) {
0777                 meCluLocalYResZGlobPlus_simLC_fromIndirectHits_->Fill(ylocal_res);
0778                 meCluLocalYPullZGlobPlus_simLC_fromIndirectHits_->Fill(ylocal_res /
0779                                                                        std::sqrt(comp->localPositionError().yy()));
0780               } else {
0781                 meCluLocalYResZGlobMinus_simLC_fromIndirectHits_->Fill(ylocal_res);
0782                 meCluLocalYPullZGlobMinus_simLC_fromIndirectHits_->Fill(ylocal_res /
0783                                                                         std::sqrt(comp->localPositionError().yy()));
0784               }
0785               if (optionalPlots_) {
0786                 if (cluster.size() == 1) {  // single-crystal clusters
0787                   meCluSingCrystalLocalYRes_simLC_fromIndirectHits_->Fill(ylocal_res);
0788                   if (global_point.z() > 0) {
0789                     meCluSingCrystalLocalYResZGlobPlus_simLC_fromIndirectHits_->Fill(ylocal_res);
0790                   } else {
0791                     meCluSingCrystalLocalYResZGlobMinus_simLC_fromIndirectHits_->Fill(ylocal_res);
0792                   }
0793                 }  // end of single-crystal clusters
0794                 else {
0795                   if (cluster.size() > 1) {  // multi-crystal clusters
0796                     meCluMultiCrystalLocalYRes_simLC_fromIndirectHits_->Fill(ylocal_res);
0797                     if (global_point.z() > 0) {
0798                       meCluMultiCrystalLocalYResZGlobPlus_simLC_fromIndirectHits_->Fill(ylocal_res);
0799                     } else {
0800                       meCluMultiCrystalLocalYResZGlobMinus_simLC_fromIndirectHits_->Fill(ylocal_res);
0801                     }
0802                   }
0803                 }  // end of multi-crystal clusters
0804 
0805                 if (abs(global_point.eta()) < 0.3) {
0806                   meCluCentralLocalYRes_simLC_fromIndirectHits_->Fill(ylocal_res);
0807                   if (global_point.z() > 0) {
0808                     meCluCentralLocalYResZGlobPlus_simLC_fromIndirectHits_->Fill(ylocal_res);
0809                   } else {
0810                     meCluCentralLocalYResZGlobMinus_simLC_fromIndirectHits_->Fill(ylocal_res);
0811                   }
0812                 } else {
0813                   if (abs(global_point.eta()) > 1) {
0814                     meCluForwardLocalYRes_simLC_fromIndirectHits_->Fill(ylocal_res);
0815                     if (global_point.z() > 0) {
0816                       meCluForwardPlusLocalYRes_simLC_fromIndirectHits_->Fill(ylocal_res);
0817                     } else {
0818                       meCluForwardMinusLocalYRes_simLC_fromIndirectHits_->Fill(ylocal_res);
0819                     }
0820                   }
0821                 }
0822               }  //end of optional plots
0823 
0824               meCluLocalXPull_simLC_fromIndirectHits_->Fill(xlocal_res / std::sqrt(comp->localPositionError().xx()));
0825               meCluZPull_simLC_fromIndirectHits_->Fill(z_res / std::sqrt(comp->globalPositionError().czz()));
0826             }
0827 
0828             meCluTResvsEta_simLC_fromIndirectHits_->Fill(std::abs(simClusGlobalPos.eta()), time_res);
0829             meCluTResvsE_simLC_fromIndirectHits_->Fill(simClusEnergy, time_res);
0830 
0831             meCluTPullvsEta_simLC_fromIndirectHits_->Fill(std::abs(simClusGlobalPos.eta()),
0832                                                           time_res / cluster.timeError());
0833             meCluTPullvsE_simLC_fromIndirectHits_->Fill(simClusEnergy, time_res / cluster.timeError());
0834           }
0835 
0836         }  // simLayerClusterRefs loop
0837       }
0838 
0839     }  // cluster loop
0840 
0841   }  // DetSetClu loop
0842 
0843   if (n_clus_btl > 0)
0844     meNclusters_->Fill(log10(n_clus_btl));
0845 
0846   // --- This is to count the number of processed events, needed in the harvesting step
0847   meNevents_->Fill(0.5);
0848 
0849   // --- Loop over the BTL Uncalibrated RECO hits
0850   if (optionalPlots_) {
0851     auto btlUncalibRecHitsHandle = makeValid(iEvent.getHandle(btlUncalibRecHitsToken_));
0852     for (const auto& uRecHit : *btlUncalibRecHitsHandle) {
0853       BTLDetId detId = uRecHit.id();
0854 
0855       LogTrace("BtlLocalRecoValidation") << "@URH detid " << detId.rawId() << " A " << uRecHit.amplitude().first << " "
0856                                          << uRecHit.amplitude().second << " T " << uRecHit.time().first << " "
0857                                          << uRecHit.time().second;
0858 
0859       // --- Skip UncalibratedRecHits not matched to SimHits
0860       if (m_btlSimHits.count(detId.rawId()) != 1)
0861         continue;
0862 
0863       // --- Combine the information from the left and right BTL cell sides
0864 
0865       float nHits = 0.;
0866       float hit_amplitude = 0.;
0867       float hit_time = 0.;
0868 
0869       // left side:
0870       if (uRecHit.amplitude().first > 0.) {
0871         hit_amplitude += uRecHit.amplitude().first;
0872         hit_time += uRecHit.time().first;
0873         nHits += 1.;
0874       }
0875       // right side:
0876       if (uRecHit.amplitude().second > 0.) {
0877         hit_amplitude += uRecHit.amplitude().second;
0878         hit_time += uRecHit.time().second;
0879         nHits += 1.;
0880       }
0881 
0882       hit_amplitude /= nHits;
0883       hit_time /= nHits;
0884 
0885       LogDebug("BtlLocalRecoValidation") << "#unc " << nHits << " A/T " << hit_amplitude << " " << hit_time;
0886       if (nHits == 0.) {
0887         edm::LogWarning("BtlLocalRecoValidation") << "Empty uncalibrated hit in DetId " << detId;
0888         continue;
0889       }
0890 
0891       if (hit_amplitude < hitMinAmplitude_)
0892         continue;
0893 
0894       // --- Fill the histograms
0895 
0896       meUncEneRVsX_->Fill(uRecHit.position(), uRecHit.amplitude().first - hit_amplitude);
0897       meUncEneLVsX_->Fill(uRecHit.position(), uRecHit.amplitude().second - hit_amplitude);
0898 
0899       meUncTimeRVsX_->Fill(uRecHit.position(), uRecHit.time().first - hit_time);
0900       meUncTimeLVsX_->Fill(uRecHit.position(), uRecHit.time().second - hit_time);
0901 
0902       if (uncalibRecHitsPlots_) {
0903         DetId geoId = detId.geographicalId(MTDTopologyMode::crysLayoutFromTopoMode(topology->getMTDTopologyMode()));
0904         const MTDGeomDet* thedet = geom->idToDet(geoId);
0905         if (thedet == nullptr)
0906           throw cms::Exception("BtlLocalRecoValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
0907                                                          << detId.rawId() << ") is invalid!" << std::dec << std::endl;
0908         const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
0909         const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
0910 
0911         Local3DPoint local_point(0., 0., 0.);
0912         local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
0913         const auto& global_point = thedet->toGlobal(local_point);
0914 
0915         float time_res = hit_time - m_btlSimHits[detId.rawId()].time;
0916 
0917         // amplitude histograms
0918 
0919         int qBin = (int)(hit_amplitude / binWidthQ_);
0920         if (qBin > nBinsQ_ - 1)
0921           qBin = nBinsQ_ - 1;
0922 
0923         meTimeResQ_[qBin]->Fill(time_res);
0924 
0925         int etaBin = 0;
0926         for (int ibin = 1; ibin < nBinsQEta_; ++ibin)
0927           if (fabs(global_point.eta()) >= binsQEta_[ibin] && fabs(global_point.eta()) < binsQEta_[ibin + 1])
0928             etaBin = ibin;
0929 
0930         meTimeResQvsEta_[qBin][etaBin]->Fill(time_res);
0931 
0932         // eta histograms
0933 
0934         etaBin = (int)(fabs(global_point.eta()) / binWidthEta_);
0935         if (etaBin > nBinsEta_ - 1)
0936           etaBin = nBinsEta_ - 1;
0937 
0938         meTimeResEta_[etaBin]->Fill(time_res);
0939 
0940         qBin = 0;
0941         for (int ibin = 1; ibin < nBinsEtaQ_; ++ibin)
0942           if (hit_amplitude >= binsEtaQ_[ibin] && hit_amplitude < binsEtaQ_[ibin + 1])
0943             qBin = ibin;
0944 
0945         meTimeResEtavsQ_[etaBin][qBin]->Fill(time_res);
0946       }
0947     }  // uRecHit loop}
0948   }
0949 }
0950 
0951 // ------------ method for histogram booking ------------
0952 void BtlLocalRecoValidation::bookHistograms(DQMStore::IBooker& ibook,
0953                                             edm::Run const& run,
0954                                             edm::EventSetup const& iSetup) {
0955   ibook.setCurrentFolder(folder_);
0956 
0957   // --- histograms booking
0958 
0959   meNevents_ = ibook.book1D("BtlNevents", "Number of events", 1, 0., 1.);
0960 
0961   meNhits_ = ibook.book1D("BtlNhits", "Number of BTL RECO hits;log_{10}(N_{RECO})", 100, 0., 5.25);
0962 
0963   meHitEnergy_ = ibook.book1D("BtlHitEnergy", "BTL RECO hits energy;E_{RECO} [MeV]", 100, 0., 20.);
0964   for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0965     std::string name_Energy = "BtlHitEnergyRUSlice" + std::to_string(ihistoRU);
0966     std::string title_Energy = "BTL RECO hits energy (RU " + std::to_string(ihistoRU) + ");E_{RECO} [MeV])";
0967     meHitEnergyRUSlice_[ihistoRU] = ibook.book1D(name_Energy, title_Energy, 100, 0., 20.);
0968   }
0969   meHitLogEnergy_ = ibook.book1D("BtlHitLogEnergy", "BTL RECO hits energy;log_{10}(E_{RECO} [MeV])", 16, -0.1, 1.5);
0970   meHitTime_ = ibook.book1D("BtlHitTime", "BTL RECO hits ToA;ToA_{RECO} [ns]", 100, 0., 25.);
0971   meHitTimeError_ = ibook.book1D("BtlHitTimeError", "BTL RECO hits ToA error;#sigma^{ToA}_{RECO} [ns]", 50, 0., 0.1);
0972   meOccupancy_ = ibook.book2D(
0973       "BtlOccupancy", "BTL RECO hits occupancy;Z_{RECO} [cm]; #phi_{RECO} [rad]", 65, -260., 260., 126, -3.2, 3.2);
0974   if (optionalPlots_) {
0975     meLocalOccupancy_ = ibook.book2D(
0976         "BtlLocalOccupancy", "BTL RECO hits local occupancy;X_{RECO} [cm]; Y_{RECO} [cm]", 100, 10., 10., 60, -3., 3.);
0977   }
0978   meHitXlocal_ = ibook.book1D("BtlHitXlocal", "BTL RECO local X;X_{RECO}^{LOC} [cm]", 100, -10., 10.);
0979   meHitYlocal_ = ibook.book1D("BtlHitYlocal", "BTL RECO local Y;Y_{RECO}^{LOC} [cm]", 60, -3, 3);
0980   meHitZlocal_ = ibook.book1D("BtlHitZlocal", "BTL RECO local z;z_{RECO}^{LOC} [cm]", 8, -0.4, 0.4);
0981   meHitZ_ = ibook.book1D("BtlHitZ", "BTL RECO hits Z;Z_{RECO} [cm]", 100, -260., 260.);
0982   meHitPhi_ = ibook.book1D("BtlHitPhi", "BTL RECO hits #phi;#phi_{RECO} [rad]", 126, -3.2, 3.2);
0983   meHitEta_ = ibook.book1D("BtlHitEta", "BTL RECO hits #eta;#eta_{RECO}", 100, -1.55, 1.55);
0984   meHitTvsE_ =
0985       ibook.bookProfile("BtlHitTvsE", "BTL RECO ToA vs energy;E_{RECO} [MeV];ToA_{RECO} [ns]", 50, 0., 20., 0., 100.);
0986   meHitEvsPhi_ = ibook.bookProfile(
0987       "BtlHitEvsPhi", "BTL RECO energy vs #phi;#phi_{RECO} [rad];E_{RECO} [MeV]", 50, -3.2, 3.2, 0., 100.);
0988   meHitEvsEta_ = ibook.bookProfile(
0989       "BtlHitEvsEta", "BTL RECO energy vs #eta;#eta_{RECO};E_{RECO} [MeV]", 50, -1.55, 1.55, 0., 100.);
0990   meHitEvsZ_ =
0991       ibook.bookProfile("BtlHitEvsZ", "BTL RECO energy vs Z;Z_{RECO} [cm];E_{RECO} [MeV]", 50, -260., 260., 0., 100.);
0992   meHitTvsPhi_ = ibook.bookProfile(
0993       "BtlHitTvsPhi", "BTL RECO ToA vs #phi;#phi_{RECO} [rad];ToA_{RECO} [ns]", 50, -3.2, 3.2, 0., 100.);
0994   meHitTvsEta_ =
0995       ibook.bookProfile("BtlHitTvsEta", "BTL RECO ToA vs #eta;#eta_{RECO};ToA_{RECO} [ns]", 50, -1.6, 1.6, 0., 100.);
0996   meHitTvsZ_ =
0997       ibook.bookProfile("BtlHitTvsZ", "BTL RECO ToA vs Z;Z_{RECO} [cm];ToA_{RECO} [ns]", 50, -260., 260., 0., 100.);
0998   meHitLongPos_ = ibook.book1D("BtlLongPos", "BTL RECO hits longitudinal position;long. pos._{RECO}", 50, -5, 5);
0999   meTimeRes_ = ibook.book1D("BtlTimeRes", "BTL time resolution;T_{RECO}-T_{SIM}", 100, -0.5, 0.5);
1000   meTimeResVsE_ = ibook.bookProfile(
1001       "BtlTimeResvsE", "BTL time resolution vs hit energy;E_{RECO} [MeV];T_{RECO}-T_{SIM}", 50, 0., 20., -0.5, 0.5, "S");
1002   meEnergyRes_ = ibook.book1D("BtlEnergyRes", "BTL energy resolution;E_{RECO}-E_{SIM}", 100, -0.5, 0.5);
1003   meEnergyRelResVsE_ = ibook.bookProfile("BtlEnergyRelResvsE",
1004                                          "BTL relative energy resolution vs hit energy;E_{RECO} [MeV];E_{RECO}-E_{SIM}",
1005                                          50,
1006                                          0.,
1007                                          20.,
1008                                          -0.15,
1009                                          0.15,
1010                                          "S");
1011   meLongPosPull_ = ibook.book1D("BtlLongPosPull",
1012                                 "BTL longitudinal position pull;X^{loc}_{RECO}-X^{loc}_{SIM}/#sigma_{xloc_{RECO}}",
1013                                 100,
1014                                 -5.,
1015                                 5.);
1016   meLongPosPullvsE_ = ibook.bookProfile(
1017       "BtlLongposPullvsE",
1018       "BTL longitudinal position pull vs E;E_{SIM} [MeV];X^{loc}_{RECO}-X^{loc}_{SIM}/#sigma_{xloc_{RECO}}",
1019       20,
1020       0.,
1021       20.,
1022       -5.,
1023       5.,
1024       "S");
1025   meLongPosPullvsEta_ = ibook.bookProfile(
1026       "BtlLongposPullvsEta",
1027       "BTL longitudinal position pull vs #eta;|#eta_{RECO}|;X^{loc}_{RECO}-X^{loc}_{SIM}/#sigma_{xloc_{RECO}}",
1028       32,
1029       0,
1030       1.55,
1031       -5.,
1032       5.,
1033       "S");
1034   meTPullvsE_ = ibook.bookProfile(
1035       "BtlTPullvsE", "BTL time pull vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}", 20, 0., 20., -5., 5., "S");
1036   meTPullvsEta_ = ibook.bookProfile("BtlTPullvsEta",
1037                                     "BTL time pull vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
1038                                     30,
1039                                     0,
1040                                     1.55,
1041                                     -5.,
1042                                     5.,
1043                                     "S");
1044 
1045   meUnmatchedRecHit_ = ibook.book1D(
1046       "UnmatchedRecHit", "log10(#BTL crystals with rechits but no simhit);log10(#BTL rechits)", 80, -2., 6.);
1047 
1048   meNclusters_ = ibook.book1D("BtlNclusters", "Number of BTL RECO clusters;log_{10}(N_{RECO})", 100, 0., 5.25);
1049   meCluTime_ = ibook.book1D("BtlCluTime", "BTL cluster time ToA;ToA [ns]", 250, 0, 25);
1050   meCluTimeError_ = ibook.book1D("BtlCluTimeError", "BTL cluster time error;#sigma_{t} [ns]", 100, 0, 0.1);
1051   meCluEnergy_ = ibook.book1D("BtlCluEnergy", "BTL cluster energy;E_{RECO} [MeV]", 100, 0, 20);
1052   meCluPhi_ = ibook.book1D("BtlCluPhi", "BTL cluster #phi;#phi_{RECO} [rad]", 144, -3.2, 3.2);
1053   meCluEta_ = ibook.book1D("BtlCluEta", "BTL cluster #eta;#eta_{RECO}", 100, -1.55, 1.55);
1054   meCluHits_ = ibook.book1D("BtlCluHitNumber", "BTL hits per cluster; Cluster size", 10, 0, 10);
1055   meCluZvsPhi_ = ibook.book2D(
1056       "BtlOccupancy", "BTL cluster Z vs #phi;Z_{RECO} [cm]; #phi_{RECO} [rad]", 144, -260., 260., 50, -3.2, 3.2);
1057 
1058   if (optionalPlots_) {
1059     meCluEnergyvsEta_ = ibook.bookProfile(
1060         "BtlCluEnergyVsEta", "BTL cluster energy vs #eta; |#eta_{RECO}|; E_{RECO} [cm]", 30, 0., 1.55, 0., 20., "S");
1061     meCluHitsvsEta_ = ibook.bookProfile(
1062         "BtlCluHitsVsEta", "BTL hits per cluster vs #eta; |#eta_{RECO}|;Cluster size", 30, 0., 1.55, 0., 10., "S");
1063 
1064     meCluTimeRes_ = ibook.book1D("BtlCluTimeRes", "BTL cluster time resolution;T_{RECO}-T_{SIM} [ns]", 100, -0.5, 0.5);
1065     meCluEnergyRes_ =
1066         ibook.book1D("BtlCluEnergyRes", "BTL cluster energy resolution;E_{RECO}-E_{SIM} [MeV]", 100, -0.5, 0.5);
1067     meCluTResvsE_ = ibook.bookProfile("BtlCluTResvsE",
1068                                       "BTL cluster time resolution vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM}) [ns]",
1069                                       20,
1070                                       0.,
1071                                       20.,
1072                                       -0.5,
1073                                       0.5,
1074                                       "S");
1075     meCluTResvsEta_ = ibook.bookProfile("BtlCluTResvsEta",
1076                                         "BTL cluster time resolution vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM}) [ns]",
1077                                         30,
1078                                         0,
1079                                         1.55,
1080                                         -0.5,
1081                                         0.5,
1082                                         "S");
1083     meCluTPullvsE_ = ibook.bookProfile("BtlCluTPullvsE",
1084                                        "BTL cluster time pull vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
1085                                        20,
1086                                        0.,
1087                                        20.,
1088                                        -5.,
1089                                        5.,
1090                                        "S");
1091     meCluTPullvsEta_ =
1092         ibook.bookProfile("BtlCluTPullvsEta",
1093                           "BTL cluster time pull vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
1094                           30,
1095                           0,
1096                           1.55,
1097                           -5.,
1098                           5.,
1099                           "S");
1100     meCluRhoRes_ =
1101         ibook.book1D("BtlCluRhoRes", "BTL cluster #rho resolution;#rho_{RECO}-#rho_{SIM} [cm]", 100, -0.5, 0.5);
1102     meCluPhiRes_ =
1103         ibook.book1D("BtlCluPhiRes", "BTL cluster #phi resolution;#phi_{RECO}-#phi_{SIM} [rad]", 100, -0.03, 0.03);
1104     meCluZRes_ = ibook.book1D("BtlCluZRes", "BTL cluster Z resolution;Z_{RECO}-Z_{SIM} [cm]", 100, -0.2, 0.2);
1105     meCluLocalXRes_ =
1106         ibook.book1D("BtlCluLocalXRes", "BTL cluster local X resolution;X_{RECO}-X_{SIM} [cm]", 100, -3.1, 3.1);
1107     meCluLocalYResZGlobPlus_ = ibook.book1D(
1108         "BtlCluLocalYResZGlobPlus", "BTL cluster local Y resolution (glob Z > 0);Y_{RECO}-Y_{SIM} [cm]", 100, -0.2, 0.2);
1109     meCluLocalYResZGlobMinus_ = ibook.book1D("BtlCluLocalYResZGlobMinus",
1110                                              "BTL cluster local Y resolution (glob Z < 0);Y_{RECO}-Y_{SIM} [cm]",
1111                                              100,
1112                                              -0.2,
1113                                              0.2);
1114     meCluSingCrystalLocalYRes_ =
1115         ibook.book1D("BtlCluSingCrystalLocalYRes",
1116                      "BTL cluster local Y resolution (single Crystal clusters);Y_{RECO}-Y_{SIM} [cm]",
1117                      100,
1118                      -0.2,
1119                      0.2);
1120     meCluSingCrystalLocalYResZGlobPlus_ =
1121         ibook.book1D("BtlCluSingCrystalLocalYResZGlobPlus",
1122                      "BTL cluster local Y resolution (single Crystal clusters, Z glob > 0);Y_{RECO}-Y_{SIM} [cm]",
1123                      100,
1124                      -0.2,
1125                      0.2);
1126     meCluSingCrystalLocalYResZGlobMinus_ =
1127         ibook.book1D("BtlCluSingCrystalLocalYResZGlobMinus",
1128                      "BTL cluster local Y resolution (single Crystal clusters, Z glob < 0);Y_{RECO}-Y_{SIM} [cm]",
1129                      100,
1130                      -0.2,
1131                      0.2);
1132     meCluMultiCrystalLocalYRes_ =
1133         ibook.book1D("BtlCluMultiCrystalLocalYRes",
1134                      "BTL cluster local Y resolution (Multi-Crystal clusters);Y_{RECO}-Y_{SIM} [cm]",
1135                      100,
1136                      -0.2,
1137                      0.2);
1138     meCluMultiCrystalLocalYResZGlobPlus_ =
1139         ibook.book1D("BtlCluMultiCrystalLocalYResZGlobPlus",
1140                      "BTL cluster local Y resolution (Multi-Crystal clusters, Z glob > 0);Y_{RECO}-Y_{SIM} [cm]",
1141                      100,
1142                      -0.2,
1143                      0.2);
1144     meCluMultiCrystalLocalYResZGlobMinus_ =
1145         ibook.book1D("BtlCluMultiCrystalLocalYResZGlobMinus",
1146                      "BTL cluster local Y resolution (Multi-Crystal clusters, Z glob < 0);Y_{RECO}-Y_{SIM} [cm]",
1147                      100,
1148                      -0.2,
1149                      0.2);
1150     meCluCentralLocalYRes_ = ibook.book1D("BtlCluCentralLocalYRes",
1151                                           "BTL cluster local Y resolution (central region);Y_{RECO}-Y_{SIM} [cm]",
1152                                           100,
1153                                           -0.2,
1154                                           0.2);
1155     meCluCentralLocalYResZGlobPlus_ =
1156         ibook.book1D("BtlCluCentralLocalYResZGlobPlus",
1157                      "BTL cluster local Y resolution (central region, Z glob > 0);Y_{RECO}-Y_{SIM} [cm]",
1158                      100,
1159                      -0.2,
1160                      0.2);
1161     meCluCentralLocalYResZGlobMinus_ =
1162         ibook.book1D("BtlCluCentralLocalYResZGlobMinus",
1163                      "BTL cluster local Y resolution (central region, Z glob < 0);Y_{RECO}-Y_{SIM} [cm]",
1164                      100,
1165                      -0.2,
1166                      0.2);
1167     meCluForwardLocalYRes_ = ibook.book1D("BtlCluForwardLocalYRes",
1168                                           "BTL cluster local Y resolution (forward region);Y_{RECO}-Y_{SIM} [cm]",
1169                                           100,
1170                                           -0.2,
1171                                           0.2);
1172     meCluForwardPlusLocalYRes_ =
1173         ibook.book1D("BtlCluForwardPlusLocalYRes",
1174                      "BTL cluster local Y resolution (forward region, Z glob > 0);Y_{RECO}-Y_{SIM} [cm]",
1175                      100,
1176                      -0.2,
1177                      0.2);
1178     meCluForwardMinusLocalYRes_ =
1179         ibook.book1D("BtlCluForwardMinusLocalYRes",
1180                      "BTL cluster local Y resolution (forward region, Z glob < 0);Y_{RECO}-Y_{SIM} [cm]",
1181                      100,
1182                      -0.2,
1183                      0.2);
1184 
1185     meCluLocalYPullZGlobPlus_ = ibook.book1D("BtlCluLocalYPullZGlobPlus",
1186                                              "BTL cluster local Y pull (glob Z > 0);Y_{RECO}-Y_{SIM}/sigmaY_[RECO]",
1187                                              100,
1188                                              -5.,
1189                                              5.);
1190     meCluLocalYPullZGlobMinus_ = ibook.book1D("BtlCluLocalYPullZGlobMinus",
1191                                               "BTL cluster local Y pull (glob Z < 0);Y_{RECO}-Y_{SIM}/sigmaY_[RECO]",
1192                                               100,
1193                                               -5.,
1194                                               5.);
1195 
1196     meCluLocalXPull_ =
1197         ibook.book1D("BtlCluLocalXPull", "BTL cluster local X pull;X_{RECO}-X_{SIM}/sigmaX_[RECO]", 100, -5., 5.);
1198 
1199     meCluZPull_ = ibook.book1D("BtlCluZPull", "BTL cluster Z pull;Z_{RECO}-Z_{SIM}/sigmaZ_[RECO]", 100, -5., 5.);
1200     meCluXLocalErr_ = ibook.book1D("BtlCluXLocalErr", "BTL cluster X local error;sigmaX_{RECO,loc} [cm]", 20, 0., 2.);
1201     meCluYLocalErr_ = ibook.book1D("BtlCluYLocalErr", "BTL cluster Y local error;sigmaY_{RECO,loc} [cm]", 20, 0., 0.4);
1202     meCluYXLocal_ = ibook.book2D("BtlCluYXLocal",
1203                                  "BTL cluster local Y vs X;X^{local}_{RECO} [cm];Y^{local}_{RECO} [cm]",
1204                                  200,
1205                                  -9.5,
1206                                  9.5,
1207                                  200,
1208                                  -2.8,
1209                                  2.8);
1210     meCluYXLocalSim_ = ibook.book2D("BtlCluYXLocalSim",
1211                                     "BTL cluster local Y vs X;X^{local}_{SIM} [cm];Y^{local}_{SIM} [cm]",
1212                                     200,
1213                                     -9.5,
1214                                     9.5,
1215                                     200,
1216                                     -2.8,
1217                                     2.8);
1218     meUnmatchedCluEnergy_ =
1219         ibook.book1D("BtlUnmatchedCluEnergy", "BTL unmatched cluster log10(energy);log10(E_{RECO} [MeV])", 5, -3, 2);
1220   }
1221 
1222   // with MtdSimLayerCluster as truth
1223 
1224   meCluTrackIdOffset_ =
1225       ibook.book1D("BtlCluTrackIdOffset", "BTL cluster category (trackId offset); trackId offset", 4, 0.0, 4.0);
1226   meCluTimeRes_simLC_ = ibook.book1D("BtlCluTimeRes_simLC",
1227                                      "BTL cluster time resolution (wrt MtdSimLayerClusters);T_{RECO}-T_{SIM} [ns]",
1228                                      100,
1229                                      -0.5,
1230                                      0.5);
1231   meCluEnergyRes_simLC_ = ibook.book1D("BtlCluEnergyRes_simLC",
1232                                        "BTL cluster energy resolution (wrt MtdSimLayerClusters);E_{RECO}-E_{SIM} [MeV]",
1233                                        100,
1234                                        -0.5,
1235                                        0.5);
1236   meCluTResvsE_simLC_ = ibook.bookProfile(
1237       "BtlCluTResvsE_simLC",
1238       "BTL cluster time resolution (wrt MtdSimLayerClusters) vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM}) [ns]",
1239       20,
1240       0.,
1241       20.,
1242       -0.5,
1243       0.5,
1244       "S");
1245   meCluTResvsEta_simLC_ = ibook.bookProfile(
1246       "BtlCluTResvsEta_simLC",
1247       "BTL cluster time resolution (wrt MtdSimLayerClusters) vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM}) [ns]",
1248       30,
1249       0,
1250       1.55,
1251       -0.5,
1252       0.5,
1253       "S");
1254   meCluTPullvsE_simLC_ = ibook.bookProfile(
1255       "BtlCluTPullvsE_simLC",
1256       "BTL cluster time pull (wrt MtdSimLayerClusters) vs E;E_{SIM} [MeV];(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
1257       20,
1258       0.,
1259       20.,
1260       -5.,
1261       5.,
1262       "S");
1263   meCluTPullvsEta_simLC_ = ibook.bookProfile(
1264       "BtlCluTPullvsEta_simLC",
1265       "BTL cluster time pull (wrt MtdSimLayerClusters) vs #eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
1266       30,
1267       0,
1268       1.55,
1269       -5.,
1270       5.,
1271       "S");
1272   meCluRhoRes_simLC_ = ibook.book1D("BtlCluRhoRes_simLC",
1273                                     "BTL cluster #rho resolution (wrt MtdSimLayerClusters);#rho_{RECO}-#rho_{SIM} [cm]",
1274                                     100,
1275                                     -0.5,
1276                                     0.5);
1277   meCluPhiRes_simLC_ =
1278       ibook.book1D("BtlCluPhiRes_simLC",
1279                    "BTL cluster #phi resolution (wrt MtdSimLayerClusters);#phi_{RECO}-#phi_{SIM} [rad]",
1280                    100,
1281                    -0.03,
1282                    0.03);
1283   meCluZRes_simLC_ = ibook.book1D(
1284       "BtlCluZRes_simLC", "BTL cluster Z resolution (wrt MtdSimLayerClusters);Z_{RECO}-Z_{SIM} [cm]", 100, -0.2, 0.2);
1285   meCluLocalXRes_simLC_ = ibook.book1D("BtlCluLocalXRes_simLC",
1286                                        "BTL cluster local X resolution (wrt MtdSimLayerClusters);X_{RECO}-X_{SIM} [cm]",
1287                                        100,
1288                                        -3.1,
1289                                        3.1);
1290   meCluLocalYResZGlobPlus_simLC_ =
1291       ibook.book1D("BtlCluLocalYResZGlobPlus_simLC",
1292                    "BTL cluster local Y resolution (wrt MtdSimLayerClusters, glob Z > 0);Y_{RECO}-Y_{SIM} [cm]",
1293                    100,
1294                    -0.2,
1295                    0.2);
1296   meCluLocalYResZGlobMinus_simLC_ =
1297       ibook.book1D("BtlCluLocalYResZGlobMinus_simLC",
1298                    "BTL cluster local Y resolution (wrt MtdSimLayerClusters, glob Z < 0);Y_{RECO}-Y_{SIM} [cm]",
1299                    100,
1300                    -0.2,
1301                    0.2);
1302   if (optionalPlots_) {
1303     meCluSingCrystalLocalYRes_simLC_ = ibook.book1D(
1304         "BtlCluSingCrystalLocalYRes_simLC",
1305         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, single Crystal clusters);Y_{RECO}-Y_{SIM} [cm]",
1306         100,
1307         -0.2,
1308         0.2);
1309     meCluSingCrystalLocalYResZGlobPlus_simLC_ =
1310         ibook.book1D("BtlCluSingCrystalLocalYResZGlobPlus_simLC",
1311                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, single Crystal clusters, Z glob > "
1312                      "0);Y_{RECO}-Y_{SIM} [cm]",
1313                      100,
1314                      -0.2,
1315                      0.2);
1316     meCluSingCrystalLocalYResZGlobMinus_simLC_ =
1317         ibook.book1D("BtlCluSingCrystalLocalYResZGlobMinus_simLC",
1318                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, single Crystal clusters, Z glob < "
1319                      "0);Y_{RECO}-Y_{SIM} [cm]",
1320                      100,
1321                      -0.2,
1322                      0.2);
1323     meCluMultiCrystalLocalYRes_simLC_ = ibook.book1D(
1324         "BtlCluMultiCrystalLocalYRes_simLC",
1325         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, Multi-Crystal clusters);Y_{RECO}-Y_{SIM} [cm]",
1326         100,
1327         -0.2,
1328         0.2);
1329     meCluMultiCrystalLocalYResZGlobPlus_simLC_ =
1330         ibook.book1D("BtlCluMultiCrystalLocalYResZGlobPlus_simLC",
1331                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, Multi-Crystal clusters, Z glob > "
1332                      "0);Y_{RECO}-Y_{SIM} [cm]",
1333                      100,
1334                      -0.2,
1335                      0.2);
1336     meCluMultiCrystalLocalYResZGlobMinus_simLC_ =
1337         ibook.book1D("BtlCluMultiCrystalLocalYResZGlobMinus_simLC",
1338                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, Multi-Crystal clusters, Z glob < "
1339                      "0);Y_{RECO}-Y_{SIM} [cm]",
1340                      100,
1341                      -0.2,
1342                      0.2);
1343     meCluCentralLocalYRes_simLC_ =
1344         ibook.book1D("BtlCluCentralLocalYRes_simLC",
1345                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, central region);Y_{RECO}-Y_{SIM} [cm]",
1346                      100,
1347                      -0.2,
1348                      0.2);
1349     meCluCentralLocalYResZGlobPlus_simLC_ = ibook.book1D(
1350         "BtlCluCentralLocalYResZGlobPlus_simLC",
1351         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, central region, Z glob > 0);Y_{RECO}-Y_{SIM} [cm]",
1352         100,
1353         -0.2,
1354         0.2);
1355     meCluCentralLocalYResZGlobMinus_simLC_ = ibook.book1D(
1356         "BtlCluCentralLocalYResZGlobMinus_simLC",
1357         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, central region, Z glob < 0);Y_{RECO}-Y_{SIM} [cm]",
1358         100,
1359         -0.2,
1360         0.2);
1361     meCluForwardLocalYRes_simLC_ =
1362         ibook.book1D("BtlCluForwardLocalYRes_simLC",
1363                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, forward region);Y_{RECO}-Y_{SIM} [cm]",
1364                      100,
1365                      -0.2,
1366                      0.2);
1367     meCluForwardPlusLocalYRes_simLC_ = ibook.book1D(
1368         "BtlCluForwardPlusLocalYRes_simLC",
1369         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, forward region, Z glob > 0);Y_{RECO}-Y_{SIM} [cm]",
1370         100,
1371         -0.2,
1372         0.2);
1373     meCluForwardMinusLocalYRes_simLC_ = ibook.book1D(
1374         "BtlCluForwardMinusLocalYRes_simLC",
1375         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, forward region, Z glob < 0);Y_{RECO}-Y_{SIM} [cm]",
1376         100,
1377         -0.2,
1378         0.2);
1379   }
1380   meCluLocalYPullZGlobPlus_simLC_ = ibook.book1D("BtlCluLocalYPullZGlobPlus_simLC",
1381                                                  "BTL cluster local Y pull (glob Z > 0);Y_{RECO}-Y_{SIM}/sigmaY_[RECO]",
1382                                                  100,
1383                                                  -5.,
1384                                                  5.);
1385   meCluLocalYPullZGlobMinus_simLC_ =
1386       ibook.book1D("BtlCluLocalYPullZGlobMinus_simLC",
1387                    "BTL cluster local Y pull (wrt MtdSimLayerClusters, glob Z < 0);Y_{RECO}-Y_{SIM}/sigmaY_[RECO]",
1388                    100,
1389                    -5.,
1390                    5.);
1391 
1392   meCluLocalXPull_simLC_ =
1393       ibook.book1D("BtlCluLocalXPull_simLC",
1394                    "BTL cluster local X pull (wrt MtdSimLayerClusters);X_{RECO}-X_{SIM}/sigmaX_[RECO]",
1395                    100,
1396                    -5.,
1397                    5.);
1398 
1399   meCluZPull_simLC_ = ibook.book1D(
1400       "BtlCluZPull_simLC", "BTL cluster Z pull (wrt MtdSimLayerClusters);Z_{RECO}-Z_{SIM}/sigmaZ_[RECO]", 100, -5., 5.);
1401   if (optionalPlots_) {
1402     meCluYXLocalSim_simLC_ =
1403         ibook.book2D("BtlCluYXLocalSim_simLC",
1404                      "BTL cluster local Y vs X (MtdSimLayerClusters);X^{local}_{SIM} [cm];Y^{local}_{SIM} [cm]",
1405                      200,
1406                      -9.5,
1407                      9.5,
1408                      200,
1409                      -2.8,
1410                      2.8);
1411   }
1412 
1413   ///
1414 
1415   meCluTimeRes_simLC_fromIndirectHits_ =
1416       ibook.book1D("BtlCluTimeRes_simLC_fromIndirectHits",
1417                    "BTL cluster time resolution (wrt MtdSimLayerClusters, non-direct hits);T_{RECO}-T_{SIM} [ns]",
1418                    100,
1419                    -0.5,
1420                    0.5);
1421   meCluEnergyRes_simLC_fromIndirectHits_ =
1422       ibook.book1D("BtlCluEnergyRes_simLC_fromIndirectHits",
1423                    "BTL cluster energy resolution (wrt MtdSimLayerClusters, non-direct hits);E_{RECO}-E_{SIM} [MeV]",
1424                    100,
1425                    -0.5,
1426                    0.5);
1427   meCluTResvsE_simLC_fromIndirectHits_ =
1428       ibook.bookProfile("BtlCluTResvsE_simLC_fromIndirectHits",
1429                         "BTL cluster time resolution (wrt MtdSimLayerClusters, non-direct hits) vs E;E_{SIM} "
1430                         "[MeV];(T_{RECO}-T_{SIM}) [ns]",
1431                         20,
1432                         0.,
1433                         20.,
1434                         -0.5,
1435                         0.5,
1436                         "S");
1437   meCluTResvsEta_simLC_fromIndirectHits_ =
1438       ibook.bookProfile("BtlCluTResvsEta_simLC_fromIndirectHits",
1439                         "BTL cluster time resolution (wrt MtdSimLayerClusters, non-direct hits) vs "
1440                         "#eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM}) [ns]",
1441                         30,
1442                         0,
1443                         1.55,
1444                         -0.5,
1445                         0.5,
1446                         "S");
1447   meCluTPullvsE_simLC_fromIndirectHits_ =
1448       ibook.bookProfile("BtlCluTPullvsE_simLC_fromIndirectHits",
1449                         "BTL cluster time pull (wrt MtdSimLayerClusters, non-direct hits) vs E;E_{SIM} "
1450                         "[MeV];(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
1451                         20,
1452                         0.,
1453                         20.,
1454                         -5.,
1455                         5.,
1456                         "S");
1457   meCluTPullvsEta_simLC_fromIndirectHits_ =
1458       ibook.bookProfile("BtlCluTPullvsEta_simLC_fromIndirectHits",
1459                         "BTL cluster time pull (wrt MtdSimLayerClusters, non-direct hits) vs "
1460                         "#eta;|#eta_{RECO}|;(T_{RECO}-T_{SIM})/#sigma_{T_{RECO}}",
1461                         30,
1462                         0,
1463                         1.55,
1464                         -5.,
1465                         5.,
1466                         "S");
1467   meCluRhoRes_simLC_fromIndirectHits_ =
1468       ibook.book1D("BtlCluRhoRes_simLC_fromIndirectHits",
1469                    "BTL cluster #rho resolution (wrt MtdSimLayerClusters, non-direct hits);#rho_{RECO}-#rho_{SIM} [cm]",
1470                    100,
1471                    -0.5,
1472                    0.5);
1473   meCluPhiRes_simLC_fromIndirectHits_ = ibook.book1D(
1474       "BtlCluPhiRes_simLC_fromIndirectHits",
1475       "BTL cluster #phi resolution (wrt MtdSimLayerClusters, non-direct hits);#phi_{RECO}-#phi_{SIM} [rad]",
1476       100,
1477       -0.03,
1478       0.03);
1479   meCluZRes_simLC_fromIndirectHits_ =
1480       ibook.book1D("BtlCluZRes_simLC_fromIndirectHits",
1481                    "BTL cluster Z resolution (wrt MtdSimLayerClusters, non-direct hits);Z_{RECO}-Z_{SIM} [cm]",
1482                    100,
1483                    -0.2,
1484                    0.2);
1485   meCluLocalXRes_simLC_fromIndirectHits_ =
1486       ibook.book1D("BtlCluLocalXRes_simLC_fromIndirectHits",
1487                    "BTL cluster local X resolution (wrt MtdSimLayerClusters, non-direct hits);X_{RECO}-X_{SIM} [cm]",
1488                    100,
1489                    -3.1,
1490                    3.1);
1491   meCluLocalYResZGlobPlus_simLC_fromIndirectHits_ = ibook.book1D(
1492       "BtlCluLocalYResZGlobPlus_simLC_fromIndirectHits",
1493       "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, glob Z > 0);Y_{RECO}-Y_{SIM} [cm]",
1494       100,
1495       -0.2,
1496       0.2);
1497   meCluLocalYResZGlobMinus_simLC_fromIndirectHits_ = ibook.book1D(
1498       "BtlCluLocalYResZGlobMinus_simLC_fromIndirectHits",
1499       "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, glob Z < 0);Y_{RECO}-Y_{SIM} [cm]",
1500       100,
1501       -0.2,
1502       0.2);
1503   if (optionalPlots_) {
1504     meCluSingCrystalLocalYRes_simLC_fromIndirectHits_ =
1505         ibook.book1D("BtlCluSingCrystalLocalYRes_simLC_fromIndirectHits",
1506                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, single Crystal "
1507                      "clusters);Y_{RECO}-Y_{SIM} [cm]",
1508                      100,
1509                      -0.2,
1510                      0.2);
1511     meCluSingCrystalLocalYResZGlobPlus_simLC_fromIndirectHits_ = ibook.book1D(
1512         "BtlCluSingCrystalLocalYResZGlobPlus_simLC_fromIndirectHits",
1513         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, single Crystal clusters, Z glob > "
1514         "0);Y_{RECO}-Y_{SIM} [cm]",
1515         100,
1516         -0.2,
1517         0.2);
1518     meCluSingCrystalLocalYResZGlobMinus_simLC_fromIndirectHits_ = ibook.book1D(
1519         "BtlCluSingCrystalLocalYResZGlobMinus_simLC_fromIndirectHits",
1520         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, single Crystal clusters, Z glob < "
1521         "0);Y_{RECO}-Y_{SIM} [cm]",
1522         100,
1523         -0.2,
1524         0.2);
1525     meCluMultiCrystalLocalYRes_simLC_fromIndirectHits_ =
1526         ibook.book1D("BtlCluMultiCrystalLocalYRes_simLC_fromIndirectHits",
1527                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, Multi-Crystal "
1528                      "clusters);Y_{RECO}-Y_{SIM} [cm]",
1529                      100,
1530                      -0.2,
1531                      0.2);
1532     meCluMultiCrystalLocalYResZGlobPlus_simLC_fromIndirectHits_ = ibook.book1D(
1533         "BtlCluMultiCrystalLocalYResZGlobPlus_simLC_fromIndirectHits",
1534         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, Multi-Crystal clusters, Z glob > "
1535         "0);Y_{RECO}-Y_{SIM} [cm]",
1536         100,
1537         -0.2,
1538         0.2);
1539     meCluMultiCrystalLocalYResZGlobMinus_simLC_fromIndirectHits_ = ibook.book1D(
1540         "BtlCluMultiCrystalLocalYResZGlobMinus_simLC_fromIndirectHits",
1541         "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, Multi-Crystal clusters, Z glob < "
1542         "0);Y_{RECO}-Y_{SIM} [cm]",
1543         100,
1544         -0.2,
1545         0.2);
1546     meCluCentralLocalYRes_simLC_fromIndirectHits_ =
1547         ibook.book1D("BtlCluCentralLocalYRes_simLC_fromIndirectHits",
1548                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, central "
1549                      "region);Y_{RECO}-Y_{SIM} [cm]",
1550                      100,
1551                      -0.2,
1552                      0.2);
1553     meCluCentralLocalYResZGlobPlus_simLC_fromIndirectHits_ =
1554         ibook.book1D("BtlCluCentralLocalYResZGlobPlus_simLC_fromIndirectHits",
1555                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, central region, Z glob "
1556                      "> 0);Y_{RECO}-Y_{SIM} [cm]",
1557                      100,
1558                      -0.2,
1559                      0.2);
1560     meCluCentralLocalYResZGlobMinus_simLC_fromIndirectHits_ =
1561         ibook.book1D("BtlCluCentralLocalYResZGlobMinus_simLC_fromIndirectHits",
1562                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, central region, Z glob "
1563                      "< 0);Y_{RECO}-Y_{SIM} [cm]",
1564                      100,
1565                      -0.2,
1566                      0.2);
1567     meCluForwardLocalYRes_simLC_fromIndirectHits_ =
1568         ibook.book1D("BtlCluForwardLocalYRes_simLC_fromIndirectHits",
1569                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, forward "
1570                      "region);Y_{RECO}-Y_{SIM} [cm]",
1571                      100,
1572                      -0.2,
1573                      0.2);
1574     meCluForwardPlusLocalYRes_simLC_fromIndirectHits_ =
1575         ibook.book1D("BtlCluForwardPlusLocalYRes_simLC_fromIndirectHits",
1576                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, forward region, Z glob "
1577                      "> 0);Y_{RECO}-Y_{SIM} [cm]",
1578                      100,
1579                      -0.2,
1580                      0.2);
1581     meCluForwardMinusLocalYRes_simLC_fromIndirectHits_ =
1582         ibook.book1D("BtlCluForwardMinusLocalYRes_simLC_fromIndirectHits",
1583                      "BTL cluster local Y resolution (wrt MtdSimLayerClusters, non-direct hits, forward region, Z glob "
1584                      "< 0);Y_{RECO}-Y_{SIM} [cm]",
1585                      100,
1586                      -0.2,
1587                      0.2);
1588   }
1589   meCluLocalYPullZGlobPlus_simLC_fromIndirectHits_ =
1590       ibook.book1D("BtlCluLocalYPullZGlobPlus_simLC_fromIndirectHits",
1591                    "BTL cluster local Y pull (glob Z > 0);Y_{RECO}-Y_{SIM}/sigmaY_[RECO]",
1592                    100,
1593                    -5.,
1594                    5.);
1595   meCluLocalYPullZGlobMinus_simLC_fromIndirectHits_ = ibook.book1D(
1596       "BtlCluLocalYPullZGlobMinus_simLC_fromIndirectHits",
1597       "BTL cluster local Y pull (wrt MtdSimLayerClusters, non-direct hits, glob Z < 0);Y_{RECO}-Y_{SIM}/sigmaY_[RECO]",
1598       100,
1599       -5.,
1600       5.);
1601 
1602   meCluLocalXPull_simLC_fromIndirectHits_ =
1603       ibook.book1D("BtlCluLocalXPull_simLC_fromIndirectHits",
1604                    "BTL cluster local X pull (wrt MtdSimLayerClusters, non-direct hits);X_{RECO}-X_{SIM}/sigmaX_[RECO]",
1605                    100,
1606                    -5.,
1607                    5.);
1608 
1609   meCluZPull_simLC_fromIndirectHits_ =
1610       ibook.book1D("BtlCluZPull_simLC_fromIndirectHits",
1611                    "BTL cluster Z pull (wrt MtdSimLayerClusters, non-direct hits);Z_{RECO}-Z_{SIM}/sigmaZ_[RECO]",
1612                    100,
1613                    -5.,
1614                    5.);
1615   if (optionalPlots_) {
1616     meCluYXLocalSim_simLC_fromIndirectHits_ =
1617         ibook.book2D("BtlCluYXLocalSim_simLC_fromIndirectHits",
1618                      "BTL cluster local Y vs X (MtdSimLayerClusters);X^{local}_{SIM} [cm];Y^{local}_{SIM} [cm]",
1619                      200,
1620                      -9.5,
1621                      9.5,
1622                      200,
1623                      -2.8,
1624                      2.8);
1625   }
1626 
1627   // --- UncalibratedRecHits histograms
1628 
1629   if (optionalPlots_) {
1630     meUncEneLVsX_ = ibook.bookProfile("BTLUncEneLVsX",
1631                                       "BTL uncalibrated hit amplitude left - average vs X;X [cm];Delta(Q_{left}) [pC]",
1632                                       20,
1633                                       -5.,
1634                                       5.,
1635                                       -640.,
1636                                       640.,
1637                                       "S");
1638     meUncEneRVsX_ =
1639         ibook.bookProfile("BTLUncEneRVsX",
1640                           "BTL uncalibrated hit amplitude right - average vs X;X [cm];Delta(Q_{right}) [pC]",
1641                           20,
1642                           -5.,
1643                           5.,
1644                           -640.,
1645                           640.,
1646                           "S");
1647 
1648     meUncTimeLVsX_ = ibook.bookProfile("BTLUncTimeLVsX",
1649                                        "BTL uncalibrated hit time left - average vs X;X [cm];Delta(T_{left}) [MeV]",
1650                                        20,
1651                                        -5.,
1652                                        5.,
1653                                        -25.,
1654                                        25.,
1655                                        "S");
1656     meUncTimeRVsX_ = ibook.bookProfile("BTLUncTimeRVsX",
1657                                        "BTL uncalibrated hit time right - average vs X;X [cm];Delta(T_{right}) [MeV]",
1658                                        20,
1659                                        -5.,
1660                                        5.,
1661                                        -25.,
1662                                        25.,
1663                                        "S");
1664     if (uncalibRecHitsPlots_) {
1665       for (unsigned int ihistoQ = 0; ihistoQ < nBinsQ_; ++ihistoQ) {
1666         std::string hname = Form("TimeResQ_%d", ihistoQ);
1667         std::string htitle = Form("BTL time resolution (Q bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoQ);
1668         meTimeResQ_[ihistoQ] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
1669 
1670         for (unsigned int ihistoEta = 0; ihistoEta < nBinsQEta_; ++ihistoEta) {
1671           hname = Form("TimeResQvsEta_%d_%d", ihistoQ, ihistoEta);
1672           htitle =
1673               Form("BTL time resolution (Q bin = %d, |#eta| bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoQ, ihistoEta);
1674           meTimeResQvsEta_[ihistoQ][ihistoEta] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
1675 
1676         }  // ihistoEta loop
1677 
1678       }  // ihistoQ loop
1679 
1680       for (unsigned int ihistoEta = 0; ihistoEta < nBinsEta_; ++ihistoEta) {
1681         std::string hname = Form("TimeResEta_%d", ihistoEta);
1682         std::string htitle = Form("BTL time resolution (|#eta| bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoEta);
1683         meTimeResEta_[ihistoEta] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
1684 
1685         for (unsigned int ihistoQ = 0; ihistoQ < nBinsEtaQ_; ++ihistoQ) {
1686           hname = Form("TimeResEtavsQ_%d_%d", ihistoEta, ihistoQ);
1687           htitle =
1688               Form("BTL time resolution (|#eta| bin = %d, Q bin = %d);T_{RECO} - T_{SIM} [ns]", ihistoEta, ihistoQ);
1689           meTimeResEtavsQ_[ihistoEta][ihistoQ] = ibook.book1D(hname, htitle, 200, -0.3, 0.7);
1690 
1691         }  // ihistoQ loop
1692 
1693       }  // ihistoEta loop
1694     }
1695   }
1696 }
1697 
1698 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
1699 void BtlLocalRecoValidation::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
1700   edm::ParameterSetDescription desc;
1701 
1702   desc.add<std::string>("folder", "MTD/BTL/LocalReco");
1703   desc.add<edm::InputTag>("recHitsTag", edm::InputTag("mtdRecHits", "FTLBarrel"));
1704   desc.add<edm::InputTag>("uncalibRecHitsTag", edm::InputTag("mtdUncalibratedRecHits", "FTLBarrel"));
1705   desc.add<edm::InputTag>("simHitsTag", edm::InputTag("mix", "g4SimHitsFastTimerHitsBarrel"));
1706   desc.add<edm::InputTag>("recCluTag", edm::InputTag("mtdClusters", "FTLBarrel"));
1707   desc.add<edm::InputTag>("trkHitTag", edm::InputTag("mtdTrackingRecHits"));
1708   desc.add<edm::InputTag>("r2sAssociationMapTag", edm::InputTag("mtdRecoClusterToSimLayerClusterAssociation"));
1709   desc.add<double>("HitMinimumEnergy", 1.);  // [MeV]
1710   desc.add<bool>("optionalPlots", false);
1711   desc.add<bool>("UncalibRecHitsPlots", false);
1712   desc.add<double>("HitMinimumAmplitude", 30.);  // [pC]
1713 
1714   descriptions.add("btlLocalRecoValid", desc);
1715 }
1716 
1717 DEFINE_FWK_MODULE(BtlLocalRecoValidation);