Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef Validation_DTRecHits_Histograms_h
0002 #define Validation_DTRecHits_Histograms_h
0003 
0004 /** \class Histograms
0005  *  Collection of histograms for DT RecHit and Segment test.
0006  *
0007  *  \author S. Bolognesi and G. Cerminara - INFN Torino
0008  */
0009 
0010 #include <cmath>
0011 #include <iostream>
0012 #include <string>
0013 
0014 #include <TH1F.h>
0015 
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 
0018 //---------------------------------------------------------------------------------------
0019 /// Function to fill an efficiency histograms with binomial errors
0020 inline void divide(dqm::legacy::MonitorElement* eff,
0021                    const dqm::legacy::MonitorElement* numerator,
0022                    const dqm::legacy::MonitorElement* denominator) {
0023   TH1* effH = eff->getTH1();
0024   TH1* numH = numerator->getTH1();
0025   TH1* denH = denominator->getTH1();
0026   effH->Divide(numH, denH);
0027 
0028   // Set the error accordingly to binomial statistics
0029   int bins = effH->GetNbinsX();
0030   for (int bin = 1; bin <= bins; ++bin) {
0031     float den = denH->GetBinContent(bin);
0032     float eff = effH->GetBinContent(bin);
0033     float err = 0;
0034     if (den != 0) {
0035       err = sqrt(eff * (1 - eff) / den);
0036     }
0037     effH->SetBinError(bin, err);
0038   }
0039   return;
0040 }
0041 
0042 //---------------------------------------------------------------------------------------
0043 /// A set of histograms of residuals and pulls for 1D RecHits
0044 class HRes1DHit {
0045 public:
0046   typedef dqm::legacy::DQMStore DQMStore;
0047   typedef dqm::legacy::MonitorElement MonitorElement;
0048 
0049   HRes1DHit(const std::string& name, DQMStore::IBooker& booker, bool doall = true, bool local = true) {
0050     std::string pre = "1D_";
0051     pre += name;
0052     doall_ = doall;
0053     booker.setCurrentFolder("DT/1DRecHits/Res/");
0054     if (doall) {
0055       hDist = booker.book1D(pre + "_hDist", "1D RHit distance from wire", 100, 0, 2.5);
0056       hResVsAngle =
0057           booker.book2D(pre + "_hResVsAngle", "1D RHit residual vs impact angle", 100, -1.2, 1.2, 100, -0.2, 0.2);
0058       hResVsDistFE =
0059           booker.book2D(pre + "_hResVsDistFE", "1D RHit residual vs FE distance", 100, 0., 400., 150, -0.5, 0.5);
0060       booker.setCurrentFolder("DT/1DRecHits/Pull/");
0061       hPullVsPos = booker.book2D(pre + "_hPullVsPos", "1D RHit pull vs position", 100, 0, 2.5, 100, -5, 5);
0062       hPullVsAngle = booker.book2D(pre + "_hPullVsAngle", "1D RHit pull vs impact angle", 100, -1.2, 1.2, 100, -5, 5);
0063       hPullVsDistFE = booker.book2D(pre + "_hPullVsDistFE", "1D RHit pull vs FE distance", 100, 0., 400., 100, -5, 5);
0064     }
0065     booker.setCurrentFolder("DT/1DRecHits/Res/");
0066     hRes = booker.book1D(pre + "_hRes", "1D RHit residual", 300, -0.5, 0.5);
0067     hResSt[0] = booker.book1D(pre + "_hResMB1", "1D RHit residual", 300, -0.5, 0.5);
0068     hResSt[1] = booker.book1D(pre + "_hResMB2", "1D RHit residual", 300, -0.5, 0.5);
0069     hResSt[2] = booker.book1D(pre + "_hResMB3", "1D RHit residual", 300, -0.5, 0.5);
0070     hResSt[3] = booker.book1D(pre + "_hResMB4", "1D RHit residual", 300, -0.5, 0.5);
0071     hResVsEta = booker.book2D(pre + "_hResVsEta", "1D RHit residual vs eta", 50, -1.25, 1.25, 150, -0.5, 0.5);
0072     hResVsPhi = booker.book2D(pre + "_hResVsPhi", "1D RHit residual vs phi", 100, -3.2, 3.2, 150, -0.5, 0.5);
0073     hResVsPos = booker.book2D(pre + "_hResVsPos", "1D RHit residual vs position", 100, 0, 2.5, 150, -0.5, 0.5);
0074 
0075     booker.setCurrentFolder("DT/1DRecHits/Pull/");
0076     hPull = booker.book1D(pre + "_hPull", "1D RHit pull", 100, -5, 5);
0077     hPullSt[0] = booker.book1D(pre + "_hPullMB1", "1D RHit residual", 100, -5, 5);
0078     hPullSt[1] = booker.book1D(pre + "_hPullMB2", "1D RHit residual", 100, -5, 5);
0079     hPullSt[2] = booker.book1D(pre + "_hPullMB3", "1D RHit residual", 100, -5, 5);
0080     hPullSt[3] = booker.book1D(pre + "_hPullMB4", "1D RHit residual", 100, -5, 5);
0081   }
0082 
0083   void fill(float distSimHit,
0084             float thetaSimHit,
0085             float distFESimHit,
0086             float distRecHit,
0087             float etaSimHit,
0088             float phiSimHit,
0089             float errRecHit,
0090             int station) {
0091     // Reso, pull
0092     float res = distRecHit - distSimHit;
0093     if (doall_) {
0094       hDist->Fill(distRecHit);
0095       hResVsAngle->Fill(thetaSimHit, res);
0096       hResVsDistFE->Fill(distFESimHit, res);
0097     }
0098     hRes->Fill(res);
0099     hResSt[station - 1]->Fill(res);
0100     hResVsEta->Fill(etaSimHit, res);
0101     hResVsPhi->Fill(phiSimHit, res);
0102     hResVsPos->Fill(distSimHit, res);
0103     if (errRecHit != 0) {
0104       float pull = res / errRecHit;
0105       hPull->Fill(pull);
0106       hPullSt[station - 1]->Fill(pull);
0107       if (doall_) {
0108         hPullVsPos->Fill(distSimHit, pull);
0109         hPullVsAngle->Fill(thetaSimHit, pull);
0110         hPullVsDistFE->Fill(distFESimHit, pull);
0111       }
0112     } else {
0113       std::cout << "Error: RecHit error = 0" << std::endl;
0114     }
0115   }
0116 
0117 private:
0118   dqm::reco::MonitorElement* hDist;
0119   dqm::reco::MonitorElement* hRes;
0120   dqm::reco::MonitorElement* hResSt[4];
0121   dqm::reco::MonitorElement* hResVsEta;
0122   dqm::reco::MonitorElement* hResVsPhi;
0123   dqm::reco::MonitorElement* hResVsPos;
0124   dqm::reco::MonitorElement* hResVsAngle;
0125   dqm::reco::MonitorElement* hResVsDistFE;
0126 
0127   dqm::reco::MonitorElement* hPull;
0128   dqm::reco::MonitorElement* hPullSt[4];
0129   dqm::reco::MonitorElement* hPullVsPos;
0130   dqm::reco::MonitorElement* hPullVsAngle;
0131   dqm::reco::MonitorElement* hPullVsDistFE;
0132   bool doall_;
0133   std::string name_;
0134 };
0135 
0136 //---------------------------------------------------------------------------------------
0137 /// A set of histograms fo efficiency computation for 1D RecHits (producer)
0138 class HEff1DHit {
0139 public:
0140   typedef dqm::legacy::DQMStore DQMStore;
0141   typedef dqm::legacy::MonitorElement MonitorElement;
0142 
0143   HEff1DHit(const std::string& name, DQMStore::IBooker& booker) {
0144     std::string pre = "1D_";
0145     pre += name;
0146     name_ = pre;
0147     booker.setCurrentFolder("DT/1DRecHits/");
0148     hEtaMuSimHit = booker.book1D(pre + "_hEtaMuSimHit", "SimHit Eta distribution", 100, -1.5, 1.5);
0149     hEtaRecHit = booker.book1D(pre + "_hEtaRecHit", "SimHit Eta distribution with 1D RecHit", 100, -1.5, 1.5);
0150     hPhiMuSimHit = booker.book1D(pre + "_hPhiMuSimHit", "SimHit Phi distribution", 100, -M_PI, M_PI);
0151     hPhiRecHit = booker.book1D(pre + "_hPhiRecHit", "SimHit Phi distribution with 1D RecHit", 100, -M_PI, M_PI);
0152     hDistMuSimHit = booker.book1D(pre + "_hDistMuSimHit", "SimHit Distance from wire distribution", 100, 0, 2.5);
0153     hDistRecHit =
0154         booker.book1D(pre + "_hDistRecHit", "SimHit Distance from wire distribution with 1D RecHit", 100, 0, 2.5);
0155   }
0156 
0157   void fill(float distSimHit, float etaSimHit, float phiSimHit, bool fillRecHit) {
0158     hEtaMuSimHit->Fill(etaSimHit);
0159     hPhiMuSimHit->Fill(phiSimHit);
0160     hDistMuSimHit->Fill(distSimHit);
0161     if (fillRecHit) {
0162       hEtaRecHit->Fill(etaSimHit);
0163       hPhiRecHit->Fill(phiSimHit);
0164       hDistRecHit->Fill(distSimHit);
0165     }
0166   }
0167 
0168 private:
0169   dqm::reco::MonitorElement* hEtaMuSimHit;
0170   dqm::reco::MonitorElement* hEtaRecHit;
0171 
0172   dqm::reco::MonitorElement* hPhiMuSimHit;
0173   dqm::reco::MonitorElement* hPhiRecHit;
0174 
0175   dqm::reco::MonitorElement* hDistMuSimHit;
0176   dqm::reco::MonitorElement* hDistRecHit;
0177 
0178   std::string name_;
0179 };
0180 
0181 //---------------------------------------------------------------------------------------
0182 /// A set of histograms fo efficiency computation for 1D RecHits (harvesting)
0183 class HEff1DHitHarvest {
0184 public:
0185   typedef dqm::legacy::DQMStore DQMStore;
0186   typedef dqm::legacy::MonitorElement MonitorElement;
0187 
0188   HEff1DHitHarvest(const std::string& name, DQMStore::IBooker& booker, DQMStore::IGetter& getter) {
0189     std::string pre = "1D_";
0190     pre += name;
0191     name_ = pre;
0192     booker.setCurrentFolder("DT/1DRecHits/");
0193     hEffVsEta = booker.book1D(pre + "_hEffVsEta", "1D RecHit Efficiency as a function of Eta", 100, -1.5, 1.5);
0194     hEffVsPhi = booker.book1D(pre + "_hEffVsPhi", "1D RecHit Efficiency as a function of Phi", 100, -M_PI, M_PI);
0195     hEffVsDist = booker.book1D(pre + "_hEffVsDist", "1D RecHit Efficiency as a function of Dist", 100, 0, 2.5);
0196 
0197     computeEfficiency(getter);
0198   }
0199 
0200   void computeEfficiency(DQMStore::IGetter& getter) {
0201     std::string pre = "DT/1DRecHits/" + name_;
0202     divide(hEffVsEta, getter.get(pre + "_hEtaMuRecHit"), getter.get(pre + "_hEtaMuSimHit"));
0203     divide(hEffVsPhi, getter.get(pre + "_hPhiMuRecHit"), getter.get(pre + "_hPhiMuSimHit"));
0204     divide(hEffVsDist, getter.get(pre + "_hDistMuRecHit"), getter.get(pre + "_hDistMuSimHit"));
0205   }
0206 
0207 private:
0208   MonitorElement* hEffVsEta;
0209   MonitorElement* hEffVsPhi;
0210   MonitorElement* hEffVsDist;
0211 
0212   std::string name_;
0213 };
0214 
0215 //---------------------------------------------------------------------------------------
0216 // Histos of residuals for 2D rechits
0217 class HRes2DHit {
0218 public:
0219   typedef dqm::legacy::DQMStore DQMStore;
0220   typedef dqm::legacy::MonitorElement MonitorElement;
0221 
0222   HRes2DHit(const std::string& name, DQMStore::IBooker& booker, bool doall = true, bool local = true) {
0223     doall_ = doall;
0224     std::string pre = "2D_";
0225     pre += name;
0226     booker.setCurrentFolder("DT/2DSegments/Res/");
0227     if (doall) {
0228       hRecAngle = booker.book1D(pre + "_hRecAngle", "Distribution of Rec segment angles;angle (rad)", 100, -1.5, 1.5);
0229       hSimAngle =
0230           booker.book1D(pre + "_hSimAngle", "Distribution of segment angles from SimHits;angle (rad)", 100, -1.5, 1.5);
0231       hRecVsSimAngle =
0232           booker.book2D(pre + "_hRecVsSimAngle", "Rec angle vs sim angle;angle (rad)", 100, -1.5, 1.5, 100, -1.5, 1.5);
0233       hResAngleVsEta = booker.book2D(pre + "_hResAngleVsEta",
0234                                      "Residual on 2D segment angle vs Eta; #eta; res (rad)",
0235                                      100,
0236                                      -2.5,
0237                                      2.5,
0238                                      200,
0239                                      -0.2,
0240                                      0.2);
0241       hResAngleVsPhi = booker.book2D(pre + "_hResAngleVsPhi",
0242                                      "Residual on 2D segment angle vs Phi; #phi (rad);res (rad)",
0243                                      100,
0244                                      -3.2,
0245                                      3.2,
0246                                      150,
0247                                      -0.2,
0248                                      0.2);
0249       hResPosVsEta = booker.book2D(
0250           pre + "_hResPosVsEta", "Residual on 2D segment position vs Eta;#eta;res (cm)", 100, -2.5, 2.5, 150, -0.2, 0.2);
0251       hResPosVsPhi = booker.book2D(pre + "_hResPosVsPhi",
0252                                    "Residual on 2D segment position vs Phi;#phi (rad);res (cm)",
0253                                    100,
0254                                    -3.2,
0255                                    3.2,
0256                                    150,
0257                                    -0.2,
0258                                    0.2);
0259       hResPosVsResAngle = booker.book2D(pre + "_hResPosVsResAngle",
0260                                         "Residual on 2D segment position vs Residual on 2D "
0261                                         "segment angle;angle (rad);res (cm)",
0262                                         100,
0263                                         -0.3,
0264                                         0.3,
0265                                         150,
0266                                         -0.2,
0267                                         0.2);
0268     }
0269     hResAngle = booker.book1D(
0270         pre + "_hResAngle", "Residual on 2D segment angle;angle_{rec}-angle_{sim} (rad)", 50, -0.01, 0.01);
0271     hResPos = booker.book1D(
0272         pre + "_hResPos", "Residual on 2D segment position (x at SL center);x_{rec}-x_{sim} (cm)", 150, -0.2, 0.2);
0273 
0274     booker.setCurrentFolder("DT/2DSegments/Pull/");
0275     hPullAngle = booker.book1D(
0276         pre + "_hPullAngle", "Pull on 2D segment angle;(angle_{rec}-angle_{sim})/#sigma (rad)", 150, -5, 5);
0277     hPullPos = booker.book1D(pre + "_hPullPos",
0278                              "Pull on 2D segment position (x at SL "
0279                              "center);(x_{rec}-x_{sim} (cm))/#sigma",
0280                              150,
0281                              -5,
0282                              5);
0283   }
0284 
0285   void fill(float angleSimSegment,
0286             float angleRecSegment,
0287             float posSimSegment,
0288             float posRecSegment,
0289             float etaSimSegment,
0290             float phiSimSegment,
0291             float sigmaPos,
0292             float sigmaAngle) {
0293     float resAngle = angleRecSegment - angleSimSegment;
0294     hResAngle->Fill(resAngle);
0295     float resPos = posRecSegment - posSimSegment;
0296     hResPos->Fill(resPos);
0297     hPullAngle->Fill(resAngle / sigmaAngle);
0298     hPullPos->Fill(resPos / sigmaPos);
0299     if (doall_) {
0300       hRecAngle->Fill(angleRecSegment);
0301       hSimAngle->Fill(angleSimSegment);
0302       hRecVsSimAngle->Fill(angleSimSegment, angleRecSegment);
0303       hResAngleVsEta->Fill(etaSimSegment, resAngle);
0304       hResAngleVsPhi->Fill(phiSimSegment, resAngle);
0305       hResPosVsEta->Fill(etaSimSegment, resPos);
0306       hResPosVsPhi->Fill(phiSimSegment, resPos);
0307       hResPosVsResAngle->Fill(resAngle, resPos);
0308     }
0309   }
0310 
0311 private:
0312   dqm::reco::MonitorElement* hRecAngle;
0313   dqm::reco::MonitorElement* hSimAngle;
0314   dqm::reco::MonitorElement* hRecVsSimAngle;
0315   dqm::reco::MonitorElement* hResAngle;
0316   dqm::reco::MonitorElement* hResAngleVsEta;
0317   dqm::reco::MonitorElement* hResAngleVsPhi;
0318   dqm::reco::MonitorElement* hResPos;
0319   dqm::reco::MonitorElement* hResPosVsEta;
0320   dqm::reco::MonitorElement* hResPosVsPhi;
0321   dqm::reco::MonitorElement* hResPosVsResAngle;
0322   dqm::reco::MonitorElement* hPullAngle;
0323   dqm::reco::MonitorElement* hPullPos;
0324 
0325   std::string name_;
0326   bool doall_;
0327 };
0328 
0329 //---------------------------------------------------------------------------------------
0330 // Histos for 2D RecHit efficiency (producer)
0331 class HEff2DHit {
0332 public:
0333   typedef dqm::legacy::DQMStore DQMStore;
0334   typedef dqm::legacy::MonitorElement MonitorElement;
0335 
0336   HEff2DHit(const std::string& name, DQMStore::IBooker& booker) {
0337     std::string pre = "2D_";
0338     pre += name;
0339     name_ = pre;
0340     booker.setCurrentFolder("DT/2DSegments/");
0341     hEtaSimSegm = booker.book1D(pre + "_hEtaSimSegm", "Eta of SimHit segment", 100, -1.5, 1.5);
0342     hEtaRecHit =
0343         booker.book1D(pre + "_hEtaRecHit", "Eta distribution of SimHit segment with 2D RecHit", 100, -1.5, 1.5);
0344     hPhiSimSegm = booker.book1D(pre + "_hPhiSimSegm", "Phi of SimHit segment", 100, -M_PI, M_PI);
0345     hPhiRecHit =
0346         booker.book1D(pre + "_hPhiRecHit", "Phi distribution of SimHit segment with 2D RecHit", 100, -M_PI, M_PI);
0347     hPosSimSegm = booker.book1D(pre + "_hPosSimSegm", "Position in SL of SimHit segment (cm)", 100, -250, 250);
0348     hPosRecHit =
0349         booker.book1D(pre + "_hPosRecHit", "Position in SL of SimHit segment with 2D RecHit (cm)", 100, -250, 250);
0350     hAngleSimSegm = booker.book1D(pre + "_hAngleSimSegm", "Angle of SimHit segment (rad)", 100, -2, 2);
0351     hAngleRecHit = booker.book1D(pre + "_hAngleRecHit", "Angle of SimHit segment with 2D RecHit (rad)", 100, -2, 2);
0352   }
0353 
0354   void fill(float etaSimSegm, float phiSimSegm, float posSimSegm, float angleSimSegm, bool fillRecHit) {
0355     hEtaSimSegm->Fill(etaSimSegm);
0356     hPhiSimSegm->Fill(phiSimSegm);
0357     hPosSimSegm->Fill(posSimSegm);
0358     hAngleSimSegm->Fill(angleSimSegm);
0359 
0360     if (fillRecHit) {
0361       hEtaRecHit->Fill(etaSimSegm);
0362       hPhiRecHit->Fill(phiSimSegm);
0363       hPosRecHit->Fill(posSimSegm);
0364       hAngleRecHit->Fill(angleSimSegm);
0365     }
0366   }
0367 
0368 private:
0369   dqm::reco::MonitorElement* hEtaSimSegm;
0370   dqm::reco::MonitorElement* hEtaRecHit;
0371   dqm::reco::MonitorElement* hPhiSimSegm;
0372   dqm::reco::MonitorElement* hPhiRecHit;
0373   dqm::reco::MonitorElement* hPosSimSegm;
0374   dqm::reco::MonitorElement* hPosRecHit;
0375   dqm::reco::MonitorElement* hAngleSimSegm;
0376   dqm::reco::MonitorElement* hAngleRecHit;
0377 
0378   std::string name_;
0379 };
0380 
0381 //---------------------------------------------------------------------------------------
0382 // Histos for 2D RecHit efficiency (harvesting)
0383 class HEff2DHitHarvest {
0384 public:
0385   typedef dqm::legacy::DQMStore DQMStore;
0386   typedef dqm::legacy::MonitorElement MonitorElement;
0387 
0388   HEff2DHitHarvest(const std::string& name, DQMStore::IBooker& booker, DQMStore::IGetter& getter) {
0389     std::string pre = "2D_";
0390     pre += name;
0391     name_ = pre;
0392     booker.setCurrentFolder("DT/2DSegments/");
0393     hEffVsEta = booker.book1D(pre + "_hEffVsEta", "2D RecHit Efficiency as a function of Eta", 100, -1.5, 1.5);
0394     hEffVsPhi = booker.book1D(pre + "_hEffVsPhi", "2D RecHit Efficiency as a function of Phi", 100, -M_PI, M_PI);
0395     hEffVsPos =
0396         booker.book1D(pre + "_hEffVsPos", "2D RecHit Efficiency as a function of position in SL", 100, -250, 250);
0397     hEffVsAngle = booker.book1D(pre + "_hEffVsAngle", "2D RecHit Efficiency as a function of angle", 100, -2, 2);
0398 
0399     computeEfficiency(getter);
0400   }
0401 
0402   void computeEfficiency(DQMStore::IGetter& getter) {
0403     std::string pre = "DT/2DSegments/" + name_;
0404     divide(hEffVsEta, getter.get(pre + "_hEtaRecHit"), getter.get(pre + "_hEtaSimSegm"));
0405     divide(hEffVsPhi, getter.get(pre + "_hPhiRecHit"), getter.get(pre + "_hPhiSimSegm"));
0406     divide(hEffVsPos, getter.get(pre + "_hPosRecHit"), getter.get(pre + "_hPosSimSegm"));
0407     divide(hEffVsAngle, getter.get(pre + "_hAngleRecHit"), getter.get(pre + "_hAngleSimSegm"));
0408   }
0409 
0410 private:
0411   MonitorElement* hEffVsEta;
0412   MonitorElement* hEffVsPhi;
0413   MonitorElement* hEffVsPos;
0414   MonitorElement* hEffVsAngle;
0415 
0416   std::string name_;
0417 };
0418 
0419 //---------------------------------------------------------------------------------------
0420 // Histos of residuals for 4D rechits
0421 class HRes4DHit {
0422 public:
0423   typedef dqm::legacy::DQMStore DQMStore;
0424   typedef dqm::legacy::MonitorElement MonitorElement;
0425 
0426   HRes4DHit(const std::string& name, DQMStore::IBooker& booker, bool doall = true, bool local = true) : local_(local) {
0427     std::string pre = "4D_";
0428     pre += name;
0429     doall_ = doall;
0430 
0431     booker.setCurrentFolder("DT/4DSegments/Res/");
0432     if (doall) {
0433       hRecAlpha =
0434           booker.book1D(pre + "_hRecAlpha", "4D RecHit alpha (RPhi) distribution;#alpha^{x} (rad)", 100, -1.5, 1.5);
0435       hRecBeta = booker.book1D(pre + "_hRecBeta", "4D RecHit beta distribution:#alpha^{y} (rad)", 100, -1.5, 1.5);
0436 
0437       hSimAlpha = booker.book1D(
0438           pre + "_hSimAlpha", "4D segment from SimHit alpha (RPhi) distribution;i#alpha^{x} (rad)", 100, -1.5, 1.5);
0439       hSimBeta =
0440           booker.book1D(pre + "_hSimBeta", "4D segment from SimHit beta distribution;#alpha^{y} (rad)", 100, -1.5, 1.5);
0441       hRecVsSimAlpha = booker.book2D(pre + "_hRecVsSimAlpha",
0442                                      "4D segment rec alpha {v}s sim alpha (RPhi);#alpha^{x} (rad)",
0443                                      100,
0444                                      -1.5,
0445                                      1.5,
0446                                      100,
0447                                      -1.5,
0448                                      1.5);
0449       hRecVsSimBeta = booker.book2D(pre + "_hRecVsSimBeta",
0450                                     "4D segment rec beta vs sim beta (RZ);#alpha^{y} (rad)",
0451                                     100,
0452                                     -1.5,
0453                                     1.5,
0454                                     100,
0455                                     -1.5,
0456                                     1.5);
0457 
0458       hResAlphaVsEta = booker.book2D(pre + "_hResAlphaVsEta",
0459                                      "4D RecHit residual on #alpha_x direction vs "
0460                                      "eta;#eta;#alpha^{x}_{rec}-#alpha^{x}_{sim} (rad)",
0461                                      100,
0462                                      -2.5,
0463                                      2.5,
0464                                      100,
0465                                      -0.025,
0466                                      0.025);
0467       hResAlphaVsPhi = booker.book2D(pre + "_hResAlphaVsPhi",
0468                                      "4D RecHit residual on #alpha_x direction vs phi (rad);#phi "
0469                                      "(rad);#alpha^{x}_{rec}-#alpha^{x}_{sim} (rad)",
0470                                      100,
0471                                      -3.2,
0472                                      3.2,
0473                                      100,
0474                                      -0.025,
0475                                      0.025);
0476       hResBetaVsEta = booker.book2D(pre + "_hResBetaVsEta",
0477                                     "4D RecHit residual on beta direction vs "
0478                                     "eta;#eta;#alpha^{y}_{rec}-#alpha^{y}_{sim} (rad)",
0479                                     100,
0480                                     -2.5,
0481                                     2.5,
0482                                     200,
0483                                     -0.2,
0484                                     0.2);
0485       hResBetaVsPhi = booker.book2D(pre + "_hResBetaVsPhi",
0486                                     "4D RecHit residual on beta direction vs phi;#phi "
0487                                     "(rad);#alpha^{y}_{rec}-#alpha^{y}_{sim} (rad)",
0488                                     100,
0489                                     -3.2,
0490                                     3.2,
0491                                     200,
0492                                     -0.2,
0493                                     0.2);
0494 
0495       hResXVsEta = booker.book2D(pre + "_hResXVsEta",
0496                                  "4D RecHit residual on position (x) in "
0497                                  "chamber vs eta;#eta;x_{rec}-x_{sim} (cm)",
0498                                  100,
0499                                  -2.5,
0500                                  2.5,
0501                                  150,
0502                                  -0.3,
0503                                  0.3);
0504       hResXVsPhi = booker.book2D(pre + "_hResXVsPhi",
0505                                  "4D RecHit residual on position (x) in chamber vs "
0506                                  "phi;#phi (rad);x_{rec}-x_{sim} (cm)",
0507                                  100,
0508                                  -3.2,
0509                                  3.2,
0510                                  150,
0511                                  -0.3,
0512                                  0.3);
0513 
0514       hResYVsEta = booker.book2D(pre + "_hResYVsEta",
0515                                  "4D RecHit residual on position (y) in "
0516                                  "chamber vs eta;#eta;y_{rec}-y_{sim} (cm)",
0517                                  100,
0518                                  -2.5,
0519                                  2.5,
0520                                  150,
0521                                  -0.6,
0522                                  0.6);
0523       hResYVsPhi = booker.book2D(pre + "_hResYVsPhi",
0524                                  "4D RecHit residual on position (y) in chamber vs "
0525                                  "phi;#phi (rad);y_{rec}-y_{sim} (cm)",
0526                                  100,
0527                                  -3.2,
0528                                  3.2,
0529                                  150,
0530                                  -0.6,
0531                                  0.6);
0532 
0533       hResAlphaVsResBeta = booker.book2D(pre + "_hResAlphaVsResBeta",
0534                                          "4D RecHit residual on alpha vs residual on beta",
0535                                          200,
0536                                          -0.3,
0537                                          0.3,
0538                                          500,
0539                                          -0.15,
0540                                          0.15);
0541       hResXVsResY = booker.book2D(
0542           pre + "_hResXVsResY", "4D RecHit residual on X vs residual on Y", 150, -0.6, 0.6, 50, -0.3, 0.3);
0543       hResAlphaVsResX = booker.book2D(
0544           pre + "_hResAlphaVsResX", "4D RecHit residual on alpha vs residual on x", 150, -0.3, 0.3, 500, -0.15, 0.15);
0545 
0546       hResAlphaVsResY = booker.book2D(
0547           pre + "_hResAlphaVsResY", "4D RecHit residual on alpha vs residual on y", 150, -0.6, 0.6, 500, -0.15, 0.15);
0548 
0549       hRecBetaRZ = booker.book1D(pre + "_hRecBetaRZ", "4D RecHit beta distribution:#alpha^{y} (rad)", 100, -1.5, 1.5);
0550 
0551       hSimBetaRZ = booker.book1D(
0552           pre + "_hSimBetaRZ", "4D segment from SimHit beta distribution in RZ SL;#alpha^{y} (rad)", 100, -1.5, 1.5);
0553       hRecVsSimBetaRZ = booker.book2D(pre + "_hRecVsSimBetaRZ",
0554                                       "4D segment rec beta vs sim beta (RZ) in RZ SL;#alpha^{y} (rad)",
0555                                       100,
0556                                       -1.5,
0557                                       1.5,
0558                                       100,
0559                                       -1.5,
0560                                       1.5);
0561 
0562       hResBetaVsEtaRZ = booker.book2D(pre + "_hResBetaVsEtaRZ",
0563                                       "4D RecHit residual on beta direction vs eta;#eta in "
0564                                       "RZ SL;#alpha^{y}_{rec}-#alpha^{y}_{sim} (rad)",
0565                                       100,
0566                                       -2.5,
0567                                       2.5,
0568                                       200,
0569                                       -0.2,
0570                                       0.2);
0571       hResBetaVsPhiRZ = booker.book2D(pre + "_hResBetaVsPhiRZ",
0572                                       "4D RecHit residual on beta direction vs phi in RZ "
0573                                       "SL;#phi (rad);#alpha^{y}_{rec}-#alpha^{y}_{sim} (rad)",
0574                                       100,
0575                                       -3.2,
0576                                       3.2,
0577                                       200,
0578                                       -0.2,
0579                                       0.2);
0580       hResYVsEtaRZ = booker.book2D(pre + "_hResYVsEtaRZ",
0581                                    "4D RecHit residual on position (y) in chamber vs eta "
0582                                    "in RZ SL;#eta;y_{rec}-y_{sim} (cm)",
0583                                    100,
0584                                    -2.5,
0585                                    2.5,
0586                                    150,
0587                                    -0.6,
0588                                    0.6);
0589       hResYVsPhiRZ = booker.book2D(pre + "_hResYVsPhiRZ",
0590                                    "4D RecHit residual on position (y) in chamber vs phi "
0591                                    "in RZ SL;#phi (rad);y_{rec}-y_{sim} (cm)",
0592                                    100,
0593                                    -3.2,
0594                                    3.2,
0595                                    150,
0596                                    -0.6,
0597                                    0.6);
0598 
0599       booker.setCurrentFolder("DT/4DSegments/Pull/");
0600       hPullAlphaVsEta = booker.book2D(pre + "_hPullAlphaVsEta",
0601                                       "4D RecHit pull on #alpha_x direction vs "
0602                                       "eta;#eta;(#alpha^{x}_{rec}-#alpha^{x}_{sim})/#sigma",
0603                                       100,
0604                                       -2.5,
0605                                       2.5,
0606                                       100,
0607                                       -5,
0608                                       5);
0609       hPullAlphaVsPhi = booker.book2D(pre + "_hPullAlphaVsPhi",
0610                                       "4D RecHit pull on #alpha_x direction vs phi (rad);#phi "
0611                                       "(rad);(#alpha^{x}_{rec}-#alpha^{x}_{sim})/#sigma",
0612                                       100,
0613                                       -3.2,
0614                                       3.2,
0615                                       100,
0616                                       -5,
0617                                       5);
0618       hPullBetaVsEta = booker.book2D(pre + "_hPullBetaVsEta",
0619                                      "4D RecHit pull on beta direction vs "
0620                                      "eta;#eta;(#alpha^{y}_{rec}-#alpha^{y}_{sim})/#sigma",
0621                                      100,
0622                                      -2.5,
0623                                      2.5,
0624                                      200,
0625                                      -5,
0626                                      5);
0627       hPullBetaVsPhi = booker.book2D(pre + "_hPullBetaVsPhi",
0628                                      "4D RecHit pull on beta direction vs phi;#phi "
0629                                      "(rad);(#alpha^{y}_{rec}-#alpha^{y}_{sim})/#sigma",
0630                                      100,
0631                                      -3.2,
0632                                      3.2,
0633                                      200,
0634                                      -5,
0635                                      5);
0636       hPullXVsEta = booker.book2D(pre + "_hPullXVsEta",
0637                                   "4D RecHit pull on position (x) in chamber "
0638                                   "vs eta;#eta;(x_{rec}-x_{sim})#sigma",
0639                                   100,
0640                                   -2.5,
0641                                   2.5,
0642                                   150,
0643                                   -5,
0644                                   5);
0645       hPullXVsPhi = booker.book2D(pre + "_hPullXVsPhi",
0646                                   "4D RecHit pull on position (x) in chamber "
0647                                   "vs phi;#phi (rad);(x_{rec}-x_{sim})/#sigma",
0648                                   100,
0649                                   -3.2,
0650                                   3.2,
0651                                   150,
0652                                   -5,
0653                                   5);
0654       hPullYVsEta = booker.book2D(pre + "_hPullYVsEta",
0655                                   "4D RecHit pull on position (y) in chamber "
0656                                   "vs eta;#eta;(y_{rec}-y_{sim})/#sigma",
0657                                   100,
0658                                   -2.5,
0659                                   2.5,
0660                                   150,
0661                                   -5,
0662                                   5);
0663       hPullYVsPhi = booker.book2D(pre + "_hPullYVsPhi",
0664                                   "4D RecHit pull on position (y) in chamber "
0665                                   "vs phi;#phi (rad);(y_{rec}-y_{sim})/#sigma",
0666                                   100,
0667                                   -3.2,
0668                                   3.2,
0669                                   150,
0670                                   -5,
0671                                   5);
0672       hPullBetaVsEtaRZ = booker.book2D(pre + "_hPullBetaVsEtaRZ",
0673                                        "4D RecHit pull on beta direction vs eta;#eta in RZ "
0674                                        "SL;(#alpha^{y}_{rec}-#alpha^{y}_{sim})/#sigma",
0675                                        100,
0676                                        -2.5,
0677                                        2.5,
0678                                        200,
0679                                        -5,
0680                                        5);
0681       hPullBetaVsPhiRZ = booker.book2D(pre + "_hPullBetaVsPhiRZ",
0682                                        "4D RecHit pull on beta direction vs phi in RZ SL;#phi "
0683                                        "(rad);(#alpha^{y}_{rec}-#alpha^{y}_{sim})/#sigma",
0684                                        100,
0685                                        -3.2,
0686                                        3.2,
0687                                        200,
0688                                        -5,
0689                                        5);
0690       hPullYVsEtaRZ = booker.book2D(pre + "_hPullYVsEtaRZ",
0691                                     "4D RecHit pull on position (y) in chamber vs eta in "
0692                                     "RZ SL;#eta;(y_{rec}-y_{sim})/#sigma",
0693                                     100,
0694                                     -2.5,
0695                                     2.5,
0696                                     150,
0697                                     -5,
0698                                     5);
0699       hPullYVsPhiRZ = booker.book2D(pre + "_hPullYVsPhiRZ",
0700                                     "4D RecHit pull on position (y) in chamber vs phi in "
0701                                     "RZ SL;#phi (rad);(y_{rec}-y_{sim})/#sigma",
0702                                     100,
0703                                     -3.2,
0704                                     3.2,
0705                                     150,
0706                                     -5,
0707                                     5);
0708     }
0709     booker.setCurrentFolder("DT/4DSegments/Res/");
0710     hResAlpha = booker.book1D(pre + "_hResAlpha",
0711                               "4D RecHit residual on #alpha_x "
0712                               "direction;#alpha^{x}_{rec}-#alpha^{x}_{sim} (rad)",
0713                               200,
0714                               -0.015,
0715                               0.015);
0716 
0717     hResBeta = booker.book1D(pre + "_hResBeta",
0718                              "4D RecHit residual on beta "
0719                              "direction;#alpha^{y}_{rec}-#alpha^{y}_{sim} (rad)",
0720                              200,
0721                              -0.1,
0722                              0.1);
0723     hResX = booker.book1D(
0724         pre + "_hResX", "4D RecHit residual on position (x) in chamber;x_{rec}-x_{sim} (cm)", 150, -0.15, 0.15);
0725     hResY = booker.book1D(
0726         pre + "_hResY", "4D RecHit residual on position (y) in chamber;y_{rec}-y_{sim} (cm)", 150, -0.6, 0.6);
0727 
0728     // histo in rz SL reference frame.
0729     hResBetaRZ = booker.book1D(pre + "_hResBetaRZ",
0730                                "4D RecHit residual on beta direction in RZ "
0731                                "SL;#alpha^{y}_{rec}-#alpha^{y}_{sim} (rad)",
0732                                200,
0733                                -0.1,
0734                                0.1);
0735 
0736     hResYRZ = booker.book1D(pre + "_hResYRZ",
0737                             "4D RecHit residual on position (y) in chamber in "
0738                             "RZ SL;y_{rec}-y_{sim} (cm)",
0739                             150,
0740                             -0.15,
0741                             0.15);
0742 
0743     // Pulls
0744     booker.setCurrentFolder("DT/4DSegments/Pull/");
0745 
0746     hPullAlpha = booker.book1D(pre + "_hPullAlpha",
0747                                "4D RecHit pull on #alpha_x "
0748                                "direction;(#alpha^{x}_{rec}-#alpha^{x}_{sim})/#sigma",
0749                                200,
0750                                -5,
0751                                5);
0752     hPullBeta = booker.book1D(pre + "_hPullBeta",
0753                               "4D RecHit pull on beta "
0754                               "direction;(#alpha^{y}_{rec}-#alpha^{y}_{sim})/#sigma",
0755                               200,
0756                               -5,
0757                               5);
0758 
0759     hPullX =
0760         booker.book1D(pre + "_hPullX", "4D RecHit pull on position (x) in chamber;(x_{rec}-x_{sim})#sigma", 150, -5, 5);
0761 
0762     hPullY = booker.book1D(
0763         pre + "_hPullY", "4D RecHit pull on position (y) in chamber;(y_{rec}-y_{sim})/#sigma", 150, -5, 5);
0764 
0765     hPullBetaRZ = booker.book1D(pre + "_hPullBetaRZ",
0766                                 "4D RecHit pull on beta direction in RZ "
0767                                 "SL;(#alpha^{y}_{rec}-#alpha^{y}_{sim})/#sigma",
0768                                 200,
0769                                 -5,
0770                                 5);
0771 
0772     hPullYRZ = booker.book1D(pre + "_hPullYRZ",
0773                              "4D RecHit pull on position (y) in chamber in RZ "
0774                              "SL;(y_{rec}-y_{sim})/#sigma",
0775                              150,
0776                              -5,
0777                              5);
0778 
0779     // NHits, t0
0780     if (local_) {
0781       booker.setCurrentFolder("DT/4DSegments/");
0782       hHitMult = booker.book2D(pre + "_hNHits", "NHits", 12, 0, 12, 6, 0, 6);
0783       ht0 = booker.book2D(pre + "_ht0", "t0", 200, -25, 25, 200, -25, 25);
0784     }
0785   }
0786 
0787   void fill(float simDirectionAlpha,
0788             float recDirectionAlpha,
0789             float simDirectionBeta,
0790             float recDirectionBeta,
0791             float simX,
0792             float recX,
0793             float simY,
0794             float recY,
0795             float simEta,
0796             float simPhi,
0797             float recYRZ,
0798             float simYRZ,
0799             float recBetaRZ,
0800             float simBetaRZ,
0801             float sigmaAlpha,
0802             float sigmaBeta,
0803             float sigmaX,
0804             float sigmaY,
0805             float sigmaBetaRZ,
0806             float sigmaYRZ,
0807             int nHitsPhi,
0808             int nHitsTheta,
0809             float t0Phi,
0810             float t0Theta) {
0811     float resAlpha = recDirectionAlpha - simDirectionAlpha;
0812     hResAlpha->Fill(resAlpha);
0813     hPullAlpha->Fill(resAlpha / sigmaAlpha);
0814     float resBeta = recDirectionBeta - simDirectionBeta;
0815     hResBeta->Fill(resBeta);
0816     hPullBeta->Fill(resBeta / sigmaBeta);
0817     float resX = recX - simX;
0818     hResX->Fill(resX);
0819     hPullX->Fill(resX / sigmaX);
0820     float resY = recY - simY;
0821     hResY->Fill(resY);
0822     hPullY->Fill(resY / sigmaY);
0823 
0824     float resBetaRZ = recBetaRZ - simBetaRZ;
0825     hResBetaRZ->Fill(resBetaRZ);
0826     hPullBetaRZ->Fill(resBetaRZ / sigmaBetaRZ);
0827     float resYRZ = recYRZ - simYRZ;
0828     hResYRZ->Fill(resYRZ);
0829     hPullYRZ->Fill(resYRZ / sigmaYRZ);
0830     if (doall_) {
0831       hRecAlpha->Fill(recDirectionAlpha);
0832       hRecBeta->Fill(recDirectionBeta);
0833       hSimAlpha->Fill(simDirectionAlpha);
0834       hSimBeta->Fill(simDirectionBeta);
0835       hRecVsSimAlpha->Fill(simDirectionAlpha, recDirectionAlpha);
0836       hRecVsSimBeta->Fill(simDirectionBeta, recDirectionBeta);
0837       hResAlphaVsEta->Fill(simEta, resAlpha);
0838       hResAlphaVsPhi->Fill(simPhi, resAlpha);
0839       hPullAlphaVsEta->Fill(simEta, resAlpha / sigmaAlpha);
0840       hPullAlphaVsPhi->Fill(simPhi, resAlpha / sigmaAlpha);
0841       hResBetaVsEta->Fill(simEta, resBeta);
0842       hResBetaVsPhi->Fill(simPhi, resBeta);
0843       hPullBetaVsEta->Fill(simEta, resBeta / sigmaBeta);
0844       hPullBetaVsPhi->Fill(simPhi, resBeta / sigmaBeta);
0845       hResXVsEta->Fill(simEta, resX);
0846       hResXVsPhi->Fill(simPhi, resX);
0847       hPullXVsEta->Fill(simEta, resX / sigmaX);
0848       hPullXVsPhi->Fill(simPhi, resX / sigmaX);
0849       hResYVsEta->Fill(simEta, resY);
0850       hResYVsPhi->Fill(simPhi, resY);
0851       hPullYVsEta->Fill(simEta, resY / sigmaY);
0852       hPullYVsPhi->Fill(simPhi, resY / sigmaY);
0853       hResAlphaVsResBeta->Fill(resBeta, resAlpha);
0854       hResXVsResY->Fill(resY, resX);
0855       hResAlphaVsResX->Fill(resX, resAlpha);
0856       hResAlphaVsResY->Fill(resY, resAlpha);
0857 
0858       // RZ SuperLayer
0859       hRecBetaRZ->Fill(recBetaRZ);
0860       hSimBetaRZ->Fill(simBetaRZ);
0861       hRecVsSimBetaRZ->Fill(simBetaRZ, recBetaRZ);
0862       hResBetaVsEtaRZ->Fill(simEta, resBetaRZ);
0863       hResBetaVsPhiRZ->Fill(simPhi, resBetaRZ);
0864       hPullBetaVsEtaRZ->Fill(simEta, resBetaRZ / sigmaBetaRZ);
0865       hPullBetaVsPhiRZ->Fill(simPhi, resBetaRZ / sigmaBetaRZ);
0866       hResYVsEtaRZ->Fill(simEta, resYRZ);
0867       hResYVsPhiRZ->Fill(simPhi, resYRZ);
0868       hPullYVsEtaRZ->Fill(simEta, resYRZ / sigmaYRZ);
0869       hPullYVsPhiRZ->Fill(simPhi, resYRZ / sigmaYRZ);
0870     }
0871     if (local_) {
0872       hHitMult->Fill(nHitsPhi, nHitsTheta);
0873       ht0->Fill(t0Phi, t0Theta);
0874     }
0875   }
0876 
0877 private:
0878   dqm::reco::MonitorElement* hRecAlpha;
0879   dqm::reco::MonitorElement* hRecBeta;
0880   dqm::reco::MonitorElement* hSimAlpha;
0881   dqm::reco::MonitorElement* hSimBeta;
0882   dqm::reco::MonitorElement* hRecVsSimAlpha;
0883   dqm::reco::MonitorElement* hRecVsSimBeta;
0884   dqm::reco::MonitorElement* hResAlpha;
0885   dqm::reco::MonitorElement* hResAlphaVsEta;
0886   dqm::reco::MonitorElement* hResAlphaVsPhi;
0887   dqm::reco::MonitorElement* hResBeta;
0888   dqm::reco::MonitorElement* hResBetaVsEta;
0889   dqm::reco::MonitorElement* hResBetaVsPhi;
0890   dqm::reco::MonitorElement* hResX;
0891   dqm::reco::MonitorElement* hResXVsEta;
0892   dqm::reco::MonitorElement* hResXVsPhi;
0893   dqm::reco::MonitorElement* hResY;
0894   dqm::reco::MonitorElement* hResYVsEta;
0895   dqm::reco::MonitorElement* hResYVsPhi;
0896   dqm::reco::MonitorElement* hResAlphaVsResBeta;
0897   dqm::reco::MonitorElement* hResXVsResY;
0898   dqm::reco::MonitorElement* hResAlphaVsResX;
0899   dqm::reco::MonitorElement* hResAlphaVsResY;
0900   dqm::reco::MonitorElement* hPullAlpha;
0901   dqm::reco::MonitorElement* hPullAlphaVsEta;
0902   dqm::reco::MonitorElement* hPullAlphaVsPhi;
0903   dqm::reco::MonitorElement* hPullBeta;
0904   dqm::reco::MonitorElement* hPullBetaVsEta;
0905   dqm::reco::MonitorElement* hPullBetaVsPhi;
0906   dqm::reco::MonitorElement* hPullX;
0907   dqm::reco::MonitorElement* hPullXVsEta;
0908   dqm::reco::MonitorElement* hPullXVsPhi;
0909   dqm::reco::MonitorElement* hPullY;
0910   dqm::reco::MonitorElement* hPullYVsEta;
0911   dqm::reco::MonitorElement* hPullYVsPhi;
0912 
0913   // RZ SL
0914   dqm::reco::MonitorElement* hRecBetaRZ;
0915   dqm::reco::MonitorElement* hSimBetaRZ;
0916   dqm::reco::MonitorElement* hRecVsSimBetaRZ;
0917   dqm::reco::MonitorElement* hResBetaRZ;
0918   dqm::reco::MonitorElement* hResBetaVsEtaRZ;
0919   dqm::reco::MonitorElement* hResBetaVsPhiRZ;
0920   dqm::reco::MonitorElement* hResYRZ;
0921   dqm::reco::MonitorElement* hResYVsEtaRZ;
0922   dqm::reco::MonitorElement* hResYVsPhiRZ;
0923   dqm::reco::MonitorElement* hPullBetaRZ;
0924   dqm::reco::MonitorElement* hPullBetaVsEtaRZ;
0925   dqm::reco::MonitorElement* hPullBetaVsPhiRZ;
0926   dqm::reco::MonitorElement* hPullYRZ;
0927   dqm::reco::MonitorElement* hPullYVsEtaRZ;
0928   dqm::reco::MonitorElement* hPullYVsPhiRZ;
0929 
0930   dqm::reco::MonitorElement* hHitMult;
0931   dqm::reco::MonitorElement* ht0;
0932 
0933   bool doall_;
0934   bool local_;
0935   std::string name_;
0936 };
0937 
0938 //---------------------------------------------------------------------------------------
0939 /// A set of histograms for efficiency 4D RecHits (producer)
0940 class HEff4DHit {
0941 public:
0942   typedef dqm::legacy::DQMStore DQMStore;
0943   typedef dqm::legacy::MonitorElement MonitorElement;
0944 
0945   HEff4DHit(const std::string& name, DQMStore::IBooker& booker) {
0946     std::string pre = "4D_";
0947     pre += name;
0948     name_ = pre;
0949     booker.setCurrentFolder("DT/4DSegments/");
0950     hEtaSimSegm = booker.book1D(pre + "_hEtaSimSegm", "Eta of SimHit segment", 100, -1.5, 1.5);
0951     hEtaRecHit =
0952         booker.book1D(pre + "_hEtaRecHit", "Eta distribution of SimHit segment with 4D RecHit", 100, -1.5, 1.5);
0953 
0954     hPhiSimSegm = booker.book1D(pre + "_hPhiSimSegm", "Phi of SimHit segment", 100, -M_PI, M_PI);
0955     hPhiRecHit =
0956         booker.book1D(pre + "_hPhiRecHit", "Phi distribution of SimHit segment with 4D RecHit", 100, -M_PI, M_PI);
0957 
0958     hXSimSegm = booker.book1D(pre + "_hXSimSegm", "X position in Chamber of SimHit segment (cm)", 100, -200, 200);
0959     hXRecHit =
0960         booker.book1D(pre + "_hXRecHit", "X position in Chamber of SimHit segment with 4D RecHit (cm)", 100, -200, 200);
0961 
0962     hYSimSegm = booker.book1D(pre + "_hYSimSegm", "Y position in Chamber of SimHit segment (cm)", 100, -200, 200);
0963     hYRecHit =
0964         booker.book1D(pre + "_hYRecHit", "Y position in Chamber of SimHit segment with 4D RecHit (cm)", 100, -200, 200);
0965 
0966     hAlphaSimSegm = booker.book1D(pre + "_hAlphaSimSegm", "Alpha of SimHit segment (rad)", 100, -1.5, 1.5);
0967     hAlphaRecHit = booker.book1D(pre + "_hAlphaRecHit", "Alpha of SimHit segment with 4D RecHit (rad)", 100, -1.5, 1.5);
0968 
0969     hBetaSimSegm = booker.book1D(pre + "_hBetaSimSegm", "Beta of SimHit segment (rad)", 100, -2, 2);
0970     hBetaRecHit = booker.book1D(pre + "_hBetaRecHit", "Beta of SimHit segment with 4D RecHit (rad)", 100, -2, 2);
0971 
0972     hNSeg = booker.book1D(pre + "_hNSeg", "Number of rec segment per sim seg", 20, 0, 20);
0973   }
0974 
0975   void fill(float etaSimSegm,
0976             float phiSimSegm,
0977             float xSimSegm,
0978             float ySimSegm,
0979             float alphaSimSegm,
0980             float betaSimSegm,
0981             bool fillRecHit,
0982             int nSeg) {
0983     hEtaSimSegm->Fill(etaSimSegm);
0984     hPhiSimSegm->Fill(phiSimSegm);
0985     hXSimSegm->Fill(xSimSegm);
0986     hYSimSegm->Fill(ySimSegm);
0987     hAlphaSimSegm->Fill(alphaSimSegm);
0988     hBetaSimSegm->Fill(betaSimSegm);
0989     hNSeg->Fill(nSeg);
0990 
0991     if (fillRecHit) {
0992       hEtaRecHit->Fill(etaSimSegm);
0993       hPhiRecHit->Fill(phiSimSegm);
0994       hXRecHit->Fill(xSimSegm);
0995       hYRecHit->Fill(ySimSegm);
0996       hAlphaRecHit->Fill(alphaSimSegm);
0997       hBetaRecHit->Fill(betaSimSegm);
0998     }
0999   }
1000 
1001 private:
1002   dqm::reco::MonitorElement* hEtaSimSegm;
1003   dqm::reco::MonitorElement* hEtaRecHit;
1004   dqm::reco::MonitorElement* hPhiSimSegm;
1005   dqm::reco::MonitorElement* hPhiRecHit;
1006   dqm::reco::MonitorElement* hXSimSegm;
1007   dqm::reco::MonitorElement* hXRecHit;
1008   dqm::reco::MonitorElement* hYSimSegm;
1009   dqm::reco::MonitorElement* hYRecHit;
1010   dqm::reco::MonitorElement* hAlphaSimSegm;
1011   dqm::reco::MonitorElement* hAlphaRecHit;
1012   dqm::reco::MonitorElement* hBetaSimSegm;
1013   dqm::reco::MonitorElement* hBetaRecHit;
1014 
1015   dqm::reco::MonitorElement* hNSeg;
1016 
1017   std::string name_;
1018 };
1019 
1020 //---------------------------------------------------------------------------------------
1021 /// A set of histograms for efficiency 4D RecHits (harvesting)
1022 class HEff4DHitHarvest {
1023 public:
1024   typedef dqm::legacy::DQMStore DQMStore;
1025   typedef dqm::legacy::MonitorElement MonitorElement;
1026 
1027   HEff4DHitHarvest(const std::string& name, DQMStore::IBooker& booker, DQMStore::IGetter& getter) {
1028     std::string pre = "4D_";
1029     pre += name;
1030     name_ = pre;
1031     booker.setCurrentFolder("DT/4DSegments/");
1032     hEffVsEta = booker.book1D(pre + "_hEffVsEta", "4D RecHit Efficiency as a function of Eta", 100, -1.5, 1.5);
1033     hEffVsPhi = booker.book1D(pre + "_hEffVsPhi", "4D RecHit Efficiency as a function of Phi", 100, -M_PI, M_PI);
1034     hEffVsX =
1035         booker.book1D(pre + "_hEffVsX", "4D RecHit Efficiency as a function of x position in Chamber", 100, -200, 200);
1036     hEffVsY =
1037         booker.book1D(pre + "_hEffVsY", "4D RecHit Efficiency as a function of y position in Chamber", 100, -200, 200);
1038     hEffVsAlpha = booker.book1D(pre + "_hEffVsAlpha", "4D RecHit Efficiency as a function of alpha", 100, -1.5, 1.5);
1039     hEffVsBeta = booker.book1D(pre + "_hEffVsBeta", "4D RecHit Efficiency as a function of beta", 100, -2, 2);
1040 
1041     computeEfficiency(getter);
1042   }
1043 
1044   void computeEfficiency(DQMStore::IGetter& getter) {
1045     std::string pre = "DT/4DSegments/" + name_;
1046     divide(hEffVsEta, getter.get(pre + "_hEtaRecHit"), getter.get(pre + "_hEtaSimSegm"));
1047     divide(hEffVsPhi, getter.get(pre + "_hPhiRecHit"), getter.get(pre + "_hPhiSimSegm"));
1048     divide(hEffVsX, getter.get(pre + "_hXRecHit"), getter.get(pre + "_hXSimSegm"));
1049     divide(hEffVsY, getter.get(pre + "_hYRecHit"), getter.get(pre + "_hYSimSegm"));
1050     divide(hEffVsAlpha, getter.get(pre + "_hAlphaRecHit"), getter.get(pre + "_hAlphaSimSegm"));
1051     divide(hEffVsBeta, getter.get(pre + "_hBetaRecHit"), getter.get(pre + "_hBetaSimSegm"));
1052   }
1053 
1054 private:
1055   MonitorElement* hEffVsEta;
1056   MonitorElement* hEffVsPhi;
1057 
1058   MonitorElement* hEffVsX;
1059   MonitorElement* hEffVsY;
1060 
1061   MonitorElement* hEffVsAlpha;
1062   MonitorElement* hEffVsBeta;
1063 
1064   std::string name_;
1065 };
1066 
1067 #endif  // Validation_DTRecHits_Histograms_h