Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:38

0001 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0002 #include "DQMServices/Core/interface/DQMStore.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ServiceRegistry/interface/Service.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 #include "DQM/SiTrackerPhase2/interface/TrackerPhase2HarvestingUtil.h"
0011 class Phase2ITRecHitHarvester : public DQMEDHarvester {
0012 public:
0013   explicit Phase2ITRecHitHarvester(const edm::ParameterSet&);
0014   ~Phase2ITRecHitHarvester() override;
0015   void dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) override;
0016 
0017   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0018 
0019 private:
0020   void gausFitslices(MonitorElement* srcME, MonitorElement* meanME, MonitorElement* sigmaME);
0021   void dofitsForLayer(const std::string& iFolder, DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter);
0022 
0023   // ----------member data ---------------------------
0024   const edm::ParameterSet config_;
0025   const std::string topFolder_;
0026   const unsigned int nbarrelLayers_;
0027   const unsigned int ndisk1Rings_;  //FOR IT epix//FOR OT TEDD1 rings
0028   const unsigned int ndisk2Rings_;  //FOR IT epix//FOR OT TEDD2 rings
0029   const unsigned int fitThreshold_;
0030   const std::string ecapdisk1Name_;  //FOR IT Epix//FOR OT TEDD_1
0031   const std::string ecapdisk2Name_;  //FOR IT Fpix//FOR OT TEDD_2
0032   const std::string histoPhiname_;
0033   const std::string deltaXvsEtaname_;
0034   const std::string deltaXvsPhiname_;
0035   const std::string deltaYvsEtaname_;
0036   const std::string deltaYvsPhiname_;
0037 };
0038 
0039 Phase2ITRecHitHarvester::Phase2ITRecHitHarvester(const edm::ParameterSet& iConfig)
0040     : config_(iConfig),
0041       topFolder_(iConfig.getParameter<std::string>("TopFolder")),
0042       nbarrelLayers_(iConfig.getParameter<uint32_t>("NbarrelLayers")),
0043       ndisk1Rings_(iConfig.getParameter<uint32_t>("NDisk1Rings")),
0044       ndisk2Rings_(iConfig.getParameter<uint32_t>("NDisk2Rings")),
0045       fitThreshold_(iConfig.getParameter<uint32_t>("NFitThreshold")),
0046       ecapdisk1Name_(iConfig.getParameter<std::string>("EcapDisk1Name")),
0047       ecapdisk2Name_(iConfig.getParameter<std::string>("EcapDisk2Name")),
0048       deltaXvsEtaname_(iConfig.getParameter<std::string>("ResidualXvsEta")),
0049       deltaXvsPhiname_(iConfig.getParameter<std::string>("ResidualXvsPhi")),
0050       deltaYvsEtaname_(iConfig.getParameter<std::string>("ResidualYvsEta")),
0051       deltaYvsPhiname_(iConfig.getParameter<std::string>("ResidualYvsPhi")) {}
0052 
0053 Phase2ITRecHitHarvester::~Phase2ITRecHitHarvester() {}
0054 
0055 void Phase2ITRecHitHarvester::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0056   for (unsigned int i = 1; i <= nbarrelLayers_; i++) {
0057     std::string iFolder = topFolder_ + "/Barrel/Layer" + std::to_string(i);
0058     dofitsForLayer(iFolder, ibooker, igetter);
0059   }
0060   for (unsigned int iSide = 1; iSide <= 2; iSide++) {
0061     const std::string ecapbasedisk1 =
0062         topFolder_ + "/EndCap_Side" + std::to_string(iSide) + "/" + ecapdisk1Name_ + "/Ring";
0063     const std::string ecapbasedisk2 =
0064         topFolder_ + "/EndCap_Side" + std::to_string(iSide) + "/" + ecapdisk2Name_ + "/Ring";
0065     //EPix or TEDD_1
0066     for (unsigned int epixr = 1; epixr <= ndisk1Rings_; epixr++) {
0067       std::string iFolder = ecapbasedisk1 + std::to_string(epixr);
0068       dofitsForLayer(iFolder, ibooker, igetter);
0069     }
0070     //FPix or TEDD_2
0071     for (unsigned int fpixr = 1; fpixr <= ndisk2Rings_; fpixr++) {
0072       std::string iFolder = ecapbasedisk2 + std::to_string(fpixr);
0073       dofitsForLayer(iFolder, ibooker, igetter);
0074     }
0075   }
0076 }
0077 //Function for Layer/Ring
0078 void Phase2ITRecHitHarvester::dofitsForLayer(const std::string& iFolder,
0079                                              DQMStore::IBooker& ibooker,
0080                                              DQMStore::IGetter& igetter) {
0081   MonitorElement* deltaX_eta = igetter.get(iFolder + "/" + deltaXvsEtaname_);
0082   MonitorElement* deltaX_phi = igetter.get(iFolder + "/" + deltaXvsPhiname_);
0083   MonitorElement* deltaY_eta = igetter.get(iFolder + "/" + deltaYvsEtaname_);
0084   MonitorElement* deltaY_phi = igetter.get(iFolder + "/" + deltaYvsPhiname_);
0085 
0086   std::string resFolder = iFolder + "/ResolutionFromFit/";
0087 
0088   ibooker.cd();
0089   ibooker.setCurrentFolder(resFolder);
0090   MonitorElement* sigmaX_eta =
0091       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("resXvseta"), ibooker);
0092   MonitorElement* meanX_eta =
0093       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("meanXvseta"), ibooker);
0094 
0095   gausFitslices(deltaX_eta, meanX_eta, sigmaX_eta);
0096 
0097   MonitorElement* sigmaY_eta =
0098       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("resYvseta"), ibooker);
0099   MonitorElement* meanY_eta =
0100       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("meanYvseta"), ibooker);
0101   gausFitslices(deltaY_eta, meanY_eta, sigmaY_eta);
0102 
0103   MonitorElement* sigmaX_phi =
0104       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("resXvsphi"), ibooker);
0105   MonitorElement* meanX_phi =
0106       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("meanXvsphi"), ibooker);
0107   gausFitslices(deltaX_phi, meanX_phi, sigmaX_phi);
0108 
0109   MonitorElement* sigmaY_phi =
0110       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("resYvsphi"), ibooker);
0111   MonitorElement* meanY_phi =
0112       phase2tkharvestutil::book1DFromPSet(config_.getParameter<edm::ParameterSet>("meanYvsphi"), ibooker);
0113   gausFitslices(deltaY_phi, meanY_phi, sigmaY_phi);
0114 }
0115 void Phase2ITRecHitHarvester::gausFitslices(MonitorElement* srcME, MonitorElement* meanME, MonitorElement* sigmaME) {
0116   TH2F* histo = srcME->getTH2F();
0117 
0118   // Fit slices projected along Y from bins in X
0119   const unsigned int cont_min = fitThreshold_;  //Minimum number of entries
0120   for (int i = 1; i <= srcME->getNbinsX(); i++) {
0121     TString iString(i);
0122     TH1* histoY = histo->ProjectionY(" ", i, i);
0123     const unsigned int cont = histoY->GetEntries();
0124 
0125     if (cont >= cont_min) {
0126       float minfit = histoY->GetMean() - histoY->GetRMS();
0127       float maxfit = histoY->GetMean() + histoY->GetRMS();
0128 
0129       std::unique_ptr<TF1> fitFcn{new TF1(TString("g") + histo->GetName() + iString, "gaus", minfit, maxfit)};
0130       histoY->Fit(fitFcn.get(), "QR0 SERIAL", "", minfit, maxfit);
0131 
0132       double* par = fitFcn->GetParameters();
0133       const double* err = fitFcn->GetParErrors();
0134 
0135       meanME->setBinContent(i, par[1]);
0136       meanME->setBinError(i, err[1]);
0137 
0138       sigmaME->setBinContent(i, par[2]);
0139       sigmaME->setBinError(i, err[2]);
0140 
0141       if (histoY)
0142         delete histoY;
0143     } else {
0144       if (histoY)
0145         delete histoY;
0146       continue;
0147     }
0148   }
0149 }
0150 
0151 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0152 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0153 void Phase2ITRecHitHarvester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0154   // phase2ITrechitHarvester
0155   edm::ParameterSetDescription desc;
0156   desc.add<std::string>("TopFolder", "TrackerPhase2ITRecHitV");
0157   desc.add<unsigned int>("NbarrelLayers", 4);
0158   desc.add<unsigned int>("NDisk1Rings", 5);
0159   desc.add<unsigned int>("NDisk2Rings", 4);
0160   desc.add<std::string>("EcapDisk1Name", "EPix");
0161   desc.add<std::string>("EcapDisk2Name", "FPix");
0162   desc.add<std::string>("ResidualXvsEta", "Delta_X_vs_Eta");
0163   desc.add<std::string>("ResidualXvsPhi", "Delta_X_vs_Phi");
0164   desc.add<std::string>("ResidualYvsEta", "Delta_Y_vs_Eta");
0165   desc.add<std::string>("ResidualYvsPhi", "Delta_Y_vs_Phi");
0166 
0167   edm::ParameterSetDescription psd0;
0168   psd0.add<std::string>("name", "resolutionXFitvseta");
0169   psd0.add<std::string>("title", ";|#eta|; X-Resolution from fit [#mum]");
0170   psd0.add<int>("NxBins", 41);
0171   psd0.add<double>("xmax", 4.1);
0172   psd0.add<double>("xmin", 0.);
0173   psd0.add<bool>("switch", true);
0174   desc.add<edm::ParameterSetDescription>("resXvseta", psd0);
0175 
0176   edm::ParameterSetDescription psd1;
0177   psd1.add<std::string>("name", "resolutionYFitvseta");
0178   psd1.add<std::string>("title", ";|#eta|; Y-Resolution from fit [#mum]");
0179   psd1.add<int>("NxBins", 41);
0180   psd1.add<double>("xmax", 4.1);
0181   psd1.add<double>("xmin", 0.);
0182   psd1.add<bool>("switch", true);
0183   desc.add<edm::ParameterSetDescription>("resYvseta", psd1);
0184 
0185   edm::ParameterSetDescription psd2;
0186   psd2.add<std::string>("name", "resolutionXFitvsphi");
0187   psd2.add<std::string>("title", ";#phi; X-Resolution from fit [#mum]");
0188   psd2.add<int>("NxBins", 36);
0189   psd2.add<double>("xmax", M_PI);
0190   psd2.add<double>("xmin", -M_PI);
0191   psd2.add<bool>("switch", true);
0192   desc.add<edm::ParameterSetDescription>("resXvsphi", psd2);
0193 
0194   edm::ParameterSetDescription psd3;
0195   psd3.add<std::string>("name", "resolutionYFitvsphi");
0196   psd3.add<std::string>("title", ";#phi; Y-Resolution from fit [#mum]");
0197   psd3.add<int>("NxBins", 36);
0198   psd3.add<double>("xmax", M_PI);
0199   psd3.add<double>("xmin", -M_PI);
0200   psd3.add<bool>("switch", true);
0201   desc.add<edm::ParameterSetDescription>("resYvsphi", psd3);
0202 
0203   edm::ParameterSetDescription psd4;
0204   psd4.add<std::string>("name", "meanXFitvseta");
0205   psd4.add<std::string>("title", ";|#eta|; Mean residual X from fit [#mum]");
0206   psd4.add<int>("NxBins", 41);
0207   psd4.add<double>("xmax", 4.1);
0208   psd4.add<double>("xmin", 0.);
0209   psd4.add<bool>("switch", true);
0210   desc.add<edm::ParameterSetDescription>("meanXvseta", psd4);
0211 
0212   edm::ParameterSetDescription psd5;
0213   psd5.add<std::string>("name", "meanYFitvseta");
0214   psd5.add<std::string>("title", ";|#eta|; Mean residual Y from fit [#mum]");
0215   psd5.add<int>("NxBins", 41);
0216   psd5.add<double>("xmax", 4.1);
0217   psd5.add<double>("xmin", 0.);
0218   psd5.add<bool>("switch", true);
0219   desc.add<edm::ParameterSetDescription>("meanYvseta", psd5);
0220 
0221   edm::ParameterSetDescription psd6;
0222   psd6.add<std::string>("name", "meanXFitvsphi");
0223   psd6.add<std::string>("title", ";#phi; Mean residual X from fit [#mum]");
0224   psd6.add<int>("NxBins", 36);
0225   psd6.add<double>("xmax", M_PI);
0226   psd6.add<double>("xmin", -M_PI);
0227   psd6.add<bool>("switch", true);
0228   desc.add<edm::ParameterSetDescription>("meanXvsphi", psd6);
0229 
0230   edm::ParameterSetDescription psd7;
0231   psd7.add<std::string>("name", "meanYFitvsphi");
0232   psd7.add<std::string>("title", ";#phi; Mean residual Y from fit [#mum]");
0233   psd7.add<int>("NxBins", 36);
0234   psd7.add<double>("xmax", M_PI);
0235   psd7.add<double>("xmin", -M_PI);
0236   psd7.add<bool>("switch", true);
0237   desc.add<edm::ParameterSetDescription>("meanYvsphi", psd7);
0238 
0239   desc.add<unsigned int>("NFitThreshold", 100);
0240 
0241   descriptions.add("Phase2ITRechitHarvester", desc);
0242   // or use the following to generate the label from the module's C++ type
0243   //descriptions.addWithDefaultLabel(desc);
0244 }
0245 
0246 //define this as a plug-in
0247 DEFINE_FWK_MODULE(Phase2ITRecHitHarvester);