Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:43

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