Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:58

0001 // system includes
0002 #include <memory>
0003 #include <string>
0004 #include <iostream>
0005 #include <fstream>
0006 
0007 // user includes
0008 #include "CLHEP/Random/RandGauss.h"
0009 #include "CalibTracker/SiStripLorentzAngle/interface/SiStripCalibLorentzAngle.h"
0010 #include "DQM/SiStripCommon/interface/ExtractTObject.h"
0011 #include "DQM/SiStripCommon/interface/SiStripHistoId.h"
0012 #include "DQMServices/Core/interface/DQMStore.h"
0013 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0014 #include "FWCore/ServiceRegistry/interface/Service.h"
0015 #include "Geometry/CommonTopologies/interface/StripTopology.h"
0016 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0017 
0018 SiStripCalibLorentzAngle::SiStripCalibLorentzAngle(edm::ParameterSet const& conf)
0019     : ConditionDBWriter<SiStripLorentzAngle>(conf),
0020       tTopoToken_(esConsumes<edm::Transition::BeginRun>()),
0021       tkGeomToken_(esConsumes<edm::Transition::BeginRun>()),
0022       magFieldToken_(esConsumes<edm::Transition::BeginRun>()),
0023       lorentzAngleToken_(esConsumes<edm::Transition::BeginRun>()),
0024       conf_(conf) {}
0025 
0026 void SiStripCalibLorentzAngle::algoBeginJob(const edm::EventSetup& c) {
0027   tTopo_ = &c.getData(tTopoToken_);
0028   tkGeom_ = &c.getData(tkGeomToken_);
0029   const auto& magField = c.getData(magFieldToken_);
0030   detid_la = c.getData(lorentzAngleToken_).getLorentzAngles();
0031 
0032   DQMStore* dbe_ = edm::Service<DQMStore>().operator->();
0033 
0034   std::string inputFile_ = conf_.getUntrackedParameter<std::string>("fileName", "LAProfiles.root");
0035   std::string LAreport_ = conf_.getUntrackedParameter<std::string>("LA_Report", "LA_Report.txt");
0036   std::string NoEntriesHisto_ = conf_.getUntrackedParameter<std::string>("NoEntriesHisto", "NoEntriesHisto.txt");
0037   std::string Dir_Name_ = conf_.getUntrackedParameter<std::string>("Dir_Name", "SiStrip");
0038 
0039   LayerDB = conf_.getUntrackedParameter<bool>("LayerDB", false);
0040 
0041   CalibByMC = conf_.getUntrackedParameter<bool>("CalibByMC", false);
0042 
0043   dbe_->open(inputFile_);
0044 
0045   // use SistripHistoId for producing histogram id (and title)
0046   SiStripHistoId hidmanager;
0047 
0048   edm::LogInfo("SiStripCalibLorentzAngle") << "### DIR-NAME = " << Dir_Name_;
0049   histolist = dbe_->getAllContents(Dir_Name_);
0050   std::vector<MonitorElement*>::iterator histo;
0051 
0052   hFile = new TFile(conf_.getUntrackedParameter<std::string>("out_fileName").c_str(), "RECREATE");
0053 
0054   LorentzAngle_Plots = hFile->mkdir("LorentzAngle_Plots");
0055   Rootple = LorentzAngle_Plots->mkdir("Rootple");
0056   MuH = LorentzAngle_Plots->mkdir("MuH");
0057   TIB_MuH = MuH->mkdir("TIB_MuH");
0058   TOB_MuH = MuH->mkdir("TOB_MuH");
0059   MuH_vs_Phi = LorentzAngle_Plots->mkdir("MuH_vs_Phi");
0060   TIB_Phi = MuH_vs_Phi->mkdir("TIB_Phi");
0061   TOB_Phi = MuH_vs_Phi->mkdir("TOB_Phi");
0062   MuH_vs_Eta = LorentzAngle_Plots->mkdir("MuH_vs_Eta");
0063   TIB_Eta = MuH_vs_Eta->mkdir("TIB_Eta");
0064   TOB_Eta = MuH_vs_Eta->mkdir("TOB_Eta");
0065   FirstIT_GoodFit_Histos = LorentzAngle_Plots->mkdir("1IT_GoodFit_Histos");
0066   TIB_1IT_GoodFit = FirstIT_GoodFit_Histos->mkdir("TIB_1IT_GoodFit");
0067   TOB_1IT_GoodFit = FirstIT_GoodFit_Histos->mkdir("TOB_1IT_GoodFit");
0068   SecondIT_GoodFit_Histos = LorentzAngle_Plots->mkdir("2IT_GoodFit_Histos");
0069   TIB_2IT_GoodFit = SecondIT_GoodFit_Histos->mkdir("TIB_2IT_GoodFit");
0070   TOB_2IT_GoodFit = SecondIT_GoodFit_Histos->mkdir("TOB_2IT_GoodFit");
0071   SecondIT_BadFit_Histos = LorentzAngle_Plots->mkdir("2IT_BadFit_Histos");
0072   TIB_2IT_BadFit = SecondIT_BadFit_Histos->mkdir("TIB_2IT_BadFit");
0073   TOB_2IT_BadFit = SecondIT_BadFit_Histos->mkdir("TOB_2IT_BadFit");
0074 
0075   TH1Ds["LA_TIB"] = new TH1D("TanLAPerTesla TIB", "TanLAPerTesla TIB", 1000, -0.5, 0.5);
0076   TH1Ds["LA_TIB"]->SetDirectory(MuH);
0077   TH1Ds["LA_TOB"] = new TH1D("TanLAPerTesla TOB", "TanLAPerTesla TOB", 1000, -0.5, 0.5);
0078   TH1Ds["LA_TOB"]->SetDirectory(MuH);
0079   TH1Ds["LA_err_TIB"] = new TH1D("TanLAPerTesla Error TIB", "TanLAPerTesla Error TIB", 1000, 0, 1);
0080   TH1Ds["LA_err_TIB"]->SetDirectory(MuH);
0081   TH1Ds["LA_err_TOB"] = new TH1D("TanLAPerTesla Error TOB", "TanLAPerTesla Error TOB", 1000, 0, 1);
0082   TH1Ds["LA_err_TOB"]->SetDirectory(MuH);
0083   TH1Ds["LA_chi2norm_TIB"] = new TH1D("TanLAPerTesla Chi2norm TIB", "TanLAPerTesla Chi2norm TIB", 2000, 0, 10);
0084   TH1Ds["LA_chi2norm_TIB"]->SetDirectory(MuH);
0085   TH1Ds["LA_chi2norm_TOB"] = new TH1D("TanLAPerTesla Chi2norm TOB", "TanLAPerTesla Chi2norm TOB", 2000, 0, 10);
0086   TH1Ds["LA_chi2norm_TOB"]->SetDirectory(MuH);
0087   TH1Ds["MagneticField"] = new TH1D("MagneticField", "MagneticField", 500, 0, 5);
0088   TH1Ds["MagneticField"]->SetDirectory(MuH);
0089 
0090   TH2Ds["LA_TIB_graph"] = new TH2D("TanLAPerTesla TIB Layers", "TanLAPerTesla TIB Layers", 60, 0, 5, 1000, -0.3, 0.3);
0091   TH2Ds["LA_TIB_graph"]->SetDirectory(MuH);
0092   TH2Ds["LA_TIB_graph"]->SetNdivisions(6);
0093   TH2Ds["LA_TOB_graph"] = new TH2D("TanLAPerTesla TOB Layers", "TanLAPerTesla TOB Layers", 80, 0, 7, 1000, -0.3, 0.3);
0094   TH2Ds["LA_TOB_graph"]->SetDirectory(MuH);
0095   TH2Ds["LA_TOB_graph"]->SetNdivisions(8);
0096 
0097   TH1Ds["LA_TIB_1"] = new TH1D("TanLAPerTesla TIB1", "TanLAPerTesla TIB1", 2000, -0.5, 0.5);
0098   TH1Ds["LA_TIB_1"]->SetDirectory(TIB_MuH);
0099   TH1Ds["LA_TIB_1_mono"] = new TH1D("TanLAPerTesla TIB1 MONO", "TanLAPerTesla TIB1 MONO", 2000, -0.5, 0.5);
0100   TH1Ds["LA_TIB_1_mono"]->SetDirectory(TIB_MuH);
0101   TH1Ds["LA_TIB_1_stereo"] = new TH1D("TanLAPerTesla TIB1 STEREO", "TanLAPerTesla TIB1 STEREO", 2000, -0.5, 0.5);
0102   TH1Ds["LA_TIB_1_stereo"]->SetDirectory(TIB_MuH);
0103   TH1Ds["LA_TIB_2"] = new TH1D("TanLAPerTesla TIB2", "TanLAPerTesla TIB2", 2000, -0.5, 0.5);
0104   TH1Ds["LA_TIB_2"]->SetDirectory(TIB_MuH);
0105   TH1Ds["LA_TIB_2_mono"] = new TH1D("TanLAPerTesla TIB2 MONO", "TanLAPerTesla TIB2 MONO", 2000, -0.5, 0.5);
0106   TH1Ds["LA_TIB_2_mono"]->SetDirectory(TIB_MuH);
0107   TH1Ds["LA_TIB_2_stereo"] = new TH1D("TanLAPerTesla TIB2 STEREO", "TanLAPerTesla TIB2 STEREO", 2000, -0.5, 0.5);
0108   TH1Ds["LA_TIB_2_stereo"]->SetDirectory(TIB_MuH);
0109   TH1Ds["LA_TIB_3"] = new TH1D("TanLAPerTesla_TIB 3", "TanLAPerTesla TIB3", 2000, -0.5, 0.5);
0110   TH1Ds["LA_TIB_3"]->SetDirectory(TIB_MuH);
0111   TH1Ds["LA_TIB_4"] = new TH1D("TanLAPerTesla_TIB 4", "TanLAPerTesla TIB4", 2000, -0.5, 0.5);
0112   TH1Ds["LA_TIB_4"]->SetDirectory(TIB_MuH);
0113 
0114   TH1Ds["LA_TOB_1"] = new TH1D("TanLAPerTesla TOB1", "TanLAPerTesla TOB1", 2000, -0.5, 0.5);
0115   TH1Ds["LA_TOB_1"]->SetDirectory(TOB_MuH);
0116   TH1Ds["LA_TOB_1_mono"] = new TH1D("TanLAPerTesla TOB1 MONO", "TanLAPerTesla TOB1 MONO", 2000, -0.5, 0.5);
0117   TH1Ds["LA_TOB_1_mono"]->SetDirectory(TOB_MuH);
0118   TH1Ds["LA_TOB_1_stereo"] = new TH1D("TanLAPerTesla TOB1 STEREO", "TanLAPerTesla TOB1 STEREO", 2000, -0.5, 0.5);
0119   TH1Ds["LA_TOB_1_stereo"]->SetDirectory(TOB_MuH);
0120   TH1Ds["LA_TOB_2"] = new TH1D("TanLAPerTesla TOB2", "TanLAPerTesla TOB2", 2000, -0.5, 0.5);
0121   TH1Ds["LA_TOB_2"]->SetDirectory(TOB_MuH);
0122   TH1Ds["LA_TOB_2_mono"] = new TH1D("TanLAPerTesla TOB2 MONO", "TanLAPerTesla TOB2 MONO", 2000, -0.5, 0.5);
0123   TH1Ds["LA_TOB_2_mono"]->SetDirectory(TOB_MuH);
0124   TH1Ds["LA_TOB_2_stereo"] = new TH1D("TanLAPerTesla TOB2 STEREO", "TanLAPerTesla TOB2 STEREO", 2000, -0.5, 0.5);
0125   TH1Ds["LA_TOB_2_stereo"]->SetDirectory(TOB_MuH);
0126   TH1Ds["LA_TOB_3"] = new TH1D("TanLAPerTesla TOB3", "TanLAPerTesla TOB3", 2000, -0.5, 0.5);
0127   TH1Ds["LA_TOB_3"]->SetDirectory(TOB_MuH);
0128   TH1Ds["LA_TOB_4"] = new TH1D("TanLAPerTesla TOB4", "TanLAPerTesla TOB4", 2000, -0.5, 0.5);
0129   TH1Ds["LA_TOB_4"]->SetDirectory(TOB_MuH);
0130   TH1Ds["LA_TOB_5"] = new TH1D("TanLAPerTesla TOB5", "TanLAPerTesla TOB5", 2000, -0.5, 0.5);
0131   TH1Ds["LA_TOB_5"]->SetDirectory(TOB_MuH);
0132   TH1Ds["LA_TOB_6"] = new TH1D("TanLAPerTesla TOB6", "TanLAPerTesla TOB6", 2000, -0.5, 0.5);
0133   TH1Ds["LA_TOB_6"]->SetDirectory(TOB_MuH);
0134 
0135   TH2Ds["LA_phi_TIB"] = new TH2D("TanLAPerTesla vs Phi TIB", "TanLAPerTesla vs Phi TIB", 800, -4, 4, 600, -0.3, 0.3);
0136   TH2Ds["LA_phi_TIB"]->SetDirectory(MuH_vs_Phi);
0137   TH2Ds["LA_phi_TOB"] = new TH2D("TanLAPerTesla vs Phi TOB", "TanLAPerTesla vs Phi TOB", 800, -4, 4, 600, -0.3, 0.3);
0138   TH2Ds["LA_phi_TOB"]->SetDirectory(MuH_vs_Phi);
0139 
0140   TH2Ds["LA_phi_TIB1"] = new TH2D("TanLAPerTesla vs Phi TIB1", "TanLAPerTesla vs Phi TIB1", 800, -4, 4, 600, -0.3, 0.3);
0141   TH2Ds["LA_phi_TIB1"]->SetDirectory(TIB_Phi);
0142   TH2Ds["LA_phi_TIB1_mono"] =
0143       new TH2D("TanLAPerTesla vs Phi TIB1 MONO", "TanLAPerTesla vs Phi TIB1 MONO", 800, -4, 4, 600, -0.3, 0.3);
0144   TH2Ds["LA_phi_TIB1_mono"]->SetDirectory(TIB_Phi);
0145   TH2Ds["LA_phi_TIB1_stereo"] =
0146       new TH2D("TanLAPerTesla vs Phi TIB1 STEREO", "TanLAPerTesla vs Phi TIB1 STEREO", 800, -4, 4, 600, -0.3, 0.3);
0147   TH2Ds["LA_phi_TIB1_stereo"]->SetDirectory(TIB_Phi);
0148   TH2Ds["LA_phi_TIB2"] = new TH2D("TanLAPerTesla vs Phi TIB2", "TanLAPerTesla vs Phi TIB2", 800, -4, 4, 600, -0.3, 0.3);
0149   TH2Ds["LA_phi_TIB2"]->SetDirectory(TIB_Phi);
0150   TH2Ds["LA_phi_TIB2_mono"] =
0151       new TH2D("TanLAPerTesla vs Phi TIB2 MONO", "TanLAPerTesla vs Phi TIB2 MONO", 800, -4, 4, 600, -0.3, 0.3);
0152   TH2Ds["LA_phi_TIB2_mono"]->SetDirectory(TIB_Phi);
0153   TH2Ds["LA_phi_TIB2_stereo"] =
0154       new TH2D("TanLAPerTesla vs Phi TIB2 STEREO", "TanLAPerTesla vs Phi TIB2 STEREO", 800, -4, 4, 600, -0.3, 0.3);
0155   TH2Ds["LA_phi_TIB2_stereo"]->SetDirectory(TIB_Phi);
0156   TH2Ds["LA_phi_TIB3"] = new TH2D("TanLAPerTesla vs Phi TIB3", "TanLAPerTesla vs Phi TIB3", 800, -4, 4, 600, -0.3, 0.3);
0157   TH2Ds["LA_phi_TIB3"]->SetDirectory(TIB_Phi);
0158   TH2Ds["LA_phi_TIB4"] = new TH2D("TanLAPerTesla vs Phi TIB4", "TanLAPerTesla vs Phi TIB4", 800, -4, 4, 600, -0.3, 0.3);
0159   TH2Ds["LA_phi_TIB4"]->SetDirectory(TIB_Phi);
0160 
0161   TH2Ds["LA_phi_TOB1"] = new TH2D("TanLAPerTesla vs Phi TOB1", "TanLAPerTesla vs Phi TOB1", 800, -4, 4, 600, -0.3, 0.3);
0162   TH2Ds["LA_phi_TOB1"]->SetDirectory(TOB_Phi);
0163   TH2Ds["LA_phi_TOB1_mono"] =
0164       new TH2D("TanLAPerTesla vs Phi TOB1 MONO", "TanLAPerTesla vs Phi TOB1 MONO", 800, -4, 4, 600, -0.3, 0.3);
0165   TH2Ds["LA_phi_TOB1_mono"]->SetDirectory(TOB_Phi);
0166   TH2Ds["LA_phi_TOB1_stereo"] =
0167       new TH2D("TanLAPerTesla vs Phi TOB1 STEREO", "TanLAPerTesla vs Phi TOB1 STEREO", 800, -4, 4, 600, -0.3, 0.3);
0168   TH2Ds["LA_phi_TOB1_stereo"]->SetDirectory(TOB_Phi);
0169   TH2Ds["LA_phi_TOB2"] = new TH2D("TanLAPerTesla vs Phi TOB2", "TanLAPerTesla vs Phi TOB2", 800, -4, 4, 600, -0.3, 0.3);
0170   TH2Ds["LA_phi_TOB2"]->SetDirectory(TOB_Phi);
0171   TH2Ds["LA_phi_TOB2_mono"] =
0172       new TH2D("TanLAPerTesla vs Phi TOB2 MONO", "TanLAPerTesla vs Phi TOB2 MONO", 800, -4, 4, 600, -0.3, 0.3);
0173   TH2Ds["LA_phi_TOB2_mono"]->SetDirectory(TOB_Phi);
0174   TH2Ds["LA_phi_TOB2_stereo"] =
0175       new TH2D("TanLAPerTesla vs Phi TOB2 STEREO", "TanLAPerTesla vs Phi TOB2 STEREO", 800, -4, 4, 600, -0.3, 0.3);
0176   TH2Ds["LA_phi_TOB2_stereo"]->SetDirectory(TOB_Phi);
0177   TH2Ds["LA_phi_TOB3"] = new TH2D("TanLAPerTesla vs Phi TOB3", "TanLAPerTesla vs Phi TOB3", 800, -4, 4, 600, -0.3, 0.3);
0178   TH2Ds["LA_phi_TOB3"]->SetDirectory(TOB_Phi);
0179   TH2Ds["LA_phi_TOB4"] = new TH2D("TanLAPerTesla vs Phi TOB4", "TanLAPerTesla vs Phi TOB4", 800, -4, 4, 600, -0.3, 0.3);
0180   TH2Ds["LA_phi_TOB4"]->SetDirectory(TOB_Phi);
0181   TH2Ds["LA_phi_TOB5"] = new TH2D("TanLAPerTesla vs Phi TOB5", "TanLAPerTesla vs Phi TOB5", 800, -4, 4, 600, -0.3, 0.3);
0182   TH2Ds["LA_phi_TOB5"]->SetDirectory(TOB_Phi);
0183   TH2Ds["LA_phi_TOB6"] = new TH2D("TanLAPerTesla vs Phi TOB6", "TanLAPerTesla vs Phi TOB6", 800, -4, 4, 600, -0.3, 0.3);
0184   TH2Ds["LA_phi_TOB6"]->SetDirectory(TOB_Phi);
0185 
0186   TH2Ds["LA_eta_TIB"] =
0187       new TH2D("TanLAPerTesla vs Eta TIB", "TanLAPerTesla vs Eta TIB", 800, -2.6, 2.6, 600, -0.3, 0.3);
0188   TH2Ds["LA_eta_TIB"]->SetDirectory(MuH_vs_Eta);
0189   TH2Ds["LA_eta_TOB"] =
0190       new TH2D("TanLAPerTesla vs Eta TOB", "TanLAPerTesla vs Eta TOB", 800, -2.6, 2.6, 600, -0.3, 0.3);
0191   TH2Ds["LA_eta_TOB"]->SetDirectory(MuH_vs_Eta);
0192 
0193   TH2Ds["LA_eta_TIB1"] =
0194       new TH2D("TanLAPerTesla vs Eta TIB1", "TanLAPerTesla vs Eta TIB1", 800, -2.6, 2.6, 600, -0.3, 0.3);
0195   TH2Ds["LA_eta_TIB1"]->SetDirectory(TIB_Eta);
0196   TH2Ds["LA_eta_TIB1_mono"] =
0197       new TH2D("TanLAPerTesla vs Eta TIB1 MONO", "TanLAPerTesla vs Eta TIB1 MONO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0198   TH2Ds["LA_eta_TIB1_mono"]->SetDirectory(TIB_Eta);
0199   TH2Ds["LA_eta_TIB1_stereo"] =
0200       new TH2D("TanLAPerTesla vs Eta TIB1 STEREO", "TanLAPerTesla vs Eta TIB1 STEREO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0201   TH2Ds["LA_eta_TIB1_stereo"]->SetDirectory(TIB_Eta);
0202   TH2Ds["LA_eta_TIB2"] =
0203       new TH2D("TanLAPerTesla vs Eta TIB2", "TanLAPerTesla vs Eta TIB2", 800, -2.6, 2.6, 600, -0.3, 0.3);
0204   TH2Ds["LA_eta_TIB2"]->SetDirectory(TIB_Eta);
0205   TH2Ds["LA_eta_TIB2_mono"] =
0206       new TH2D("TanLAPerTesla vs Eta TIB2 MONO", "TanLAPerTesla vs Eta TIB2 MONO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0207   TH2Ds["LA_eta_TIB2_mono"]->SetDirectory(TIB_Eta);
0208   TH2Ds["LA_eta_TIB2_stereo"] =
0209       new TH2D("TanLAPerTesla vs Eta TIB2 STEREO", "TanLAPerTesla vs Eta TIB2 STEREO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0210   TH2Ds["LA_eta_TIB2_stereo"]->SetDirectory(TIB_Eta);
0211   TH2Ds["LA_eta_TIB3"] =
0212       new TH2D("TanLAPerTesla vs Eta TIB3", "TanLAPerTesla vs Eta TIB3", 800, -2.6, 2.6, 600, -0.3, 0.3);
0213   TH2Ds["LA_eta_TIB3"]->SetDirectory(TIB_Eta);
0214   TH2Ds["LA_eta_TIB4"] =
0215       new TH2D("TanLAPerTesla vs Eta TIB4", "TanLAPerTesla vs Eta TIB4", 800, -2.6, 2.6, 600, -0.3, 0.3);
0216   TH2Ds["LA_eta_TIB4"]->SetDirectory(TIB_Eta);
0217 
0218   TH2Ds["LA_eta_TOB1"] =
0219       new TH2D("TanLAPerTesla vs Eta TOB1", "TanLAPerTesla vs Eta TOB1", 800, -2.6, 2.6, 600, -0.3, 0.3);
0220   TH2Ds["LA_eta_TOB1"]->SetDirectory(TIB_Eta);
0221   TH2Ds["LA_eta_TOB1_mono"] =
0222       new TH2D("TanLAPerTesla vs Eta TOB1 MONO", "TanLAPerTesla vs Eta TOB1 MONO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0223   TH2Ds["LA_eta_TOB1_mono"]->SetDirectory(TIB_Eta);
0224   TH2Ds["LA_eta_TOB1_stereo"] =
0225       new TH2D("TanLAPerTesla vs Eta TOB1 STEREO", "TanLAPerTesla vs Eta TOB1 STEREO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0226   TH2Ds["LA_eta_TOB1_stereo"]->SetDirectory(TIB_Eta);
0227   TH2Ds["LA_eta_TOB2"] =
0228       new TH2D("TanLAPerTesla vs Eta TOB2", "TanLAPerTesla vs Eta TOB2", 800, -2.6, 2.6, 600, -0.3, 0.3);
0229   TH2Ds["LA_eta_TOB2"]->SetDirectory(TIB_Eta);
0230   TH2Ds["LA_eta_TOB2_mono"] =
0231       new TH2D("TanLAPerTesla vs Eta TOB2 MONO", "TanLAPerTesla vs Eta TOB2 MONO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0232   TH2Ds["LA_eta_TOB2_mono"]->SetDirectory(TIB_Eta);
0233   TH2Ds["LA_eta_TOB2_stereo"] =
0234       new TH2D("TanLAPerTesla vs Eta TOB2 STEREO", "TanLAPerTesla vs Eta TOB2 STEREO", 800, -2.6, 2.6, 600, -0.3, 0.3);
0235   TH2Ds["LA_eta_TOB2_stereo"]->SetDirectory(TIB_Eta);
0236   TH2Ds["LA_eta_TOB3"] =
0237       new TH2D("TanLAPerTesla vs Eta TOB3", "TanLAPerTesla vs Eta TOB3", 800, -2.6, 2.6, 600, -0.3, 0.3);
0238   TH2Ds["LA_eta_TOB3"]->SetDirectory(TIB_Eta);
0239   TH2Ds["LA_eta_TOB4"] =
0240       new TH2D("TanLAPerTesla vs Eta TOB4", "TanLAPerTesla vs Eta TOB4", 800, -2.6, 2.6, 600, -0.3, 0.3);
0241   TH2Ds["LA_eta_TOB4"]->SetDirectory(TIB_Eta);
0242   TH2Ds["LA_eta_TOB5"] =
0243       new TH2D("TanLAPerTesla vs Eta TOB5", "TanLAPerTesla vs Eta TOB5", 800, -2.6, 2.6, 600, -0.3, 0.3);
0244   TH2Ds["LA_eta_TOB5"]->SetDirectory(TIB_Eta);
0245   TH2Ds["LA_eta_TOB6"] =
0246       new TH2D("TanLAPerTesla vs Eta TOB6", "TanLAPerTesla vs Eta TOB6", 800, -2.6, 2.6, 600, -0.3, 0.3);
0247   TH2Ds["LA_eta_TOB6"]->SetDirectory(TIB_Eta);
0248 
0249   ModuleTree = new TTree("ModuleTree", "ModuleTree");
0250   ModuleTree->Branch("histoEntries", &histoEntries, "histoEntries/F");
0251   ModuleTree->Branch("globalX", &globalX, "globalX/F");
0252   ModuleTree->Branch("globalY", &globalY, "globalY/F");
0253   ModuleTree->Branch("globalZ", &globalZ, "globalZ/F");
0254   ModuleTree->Branch("gphi", &gphi, "gphi/F");
0255   ModuleTree->Branch("geta", &geta, "geta/F");
0256   ModuleTree->Branch("gR", &gR, "gR/F");
0257   ModuleTree->Branch("goodFit", &goodFit, "goodFit/I");
0258   ModuleTree->Branch("goodFit1IT", &goodFit1IT, "goodFit1IT/I");
0259   ModuleTree->Branch("badFit", &badFit, "badFit/I");
0260   ModuleTree->Branch("TIB", &TIB, "TIB/I");
0261   ModuleTree->Branch("TOB", &TOB, "TOB/I");
0262   ModuleTree->Branch("Layer", &Layer, "Layer/I");
0263   ModuleTree->Branch("MonoStereo", &MonoStereo, "MonoStereo/I");
0264   ModuleTree->Branch("theBfield", &theBfield, "theBfield/F");
0265   ModuleTree->Branch("muH", &muH, "muH/F");
0266 
0267   ModuleTree->SetDirectory(Rootple);
0268 
0269   int histocounter = 0;
0270   int NotEnoughEntries = 0;
0271   int ZeroEntries = 0;
0272   int GoodFit = 0;
0273   int FirstIT_goodfit = 0;
0274   int FirstIT_badfit = 0;
0275   int SecondIT_badfit = 0;
0276   int SecondIT_goodfit = 0;
0277   int no_mod_histo = 0;
0278   float chi2norm = 0;
0279   LocalPoint p = LocalPoint(0, 0, 0);
0280 
0281   double ModuleRangeMin = conf_.getParameter<double>("ModuleFitXMin");
0282   double ModuleRangeMax = conf_.getParameter<double>("ModuleFitXMax");
0283   double ModuleRangeMin2IT = conf_.getParameter<double>("ModuleFit2ITXMin");
0284   double ModuleRangeMax2IT = conf_.getParameter<double>("ModuleFit2ITXMax");
0285   double FitCuts_Entries = conf_.getParameter<double>("FitCuts_Entries");
0286   double FitCuts_p0 = conf_.getParameter<double>("FitCuts_p0");
0287   double FitCuts_p1 = conf_.getParameter<double>("FitCuts_p1");
0288   double FitCuts_p2 = conf_.getParameter<double>("FitCuts_p2");
0289   double FitCuts_chi2 = conf_.getParameter<double>("FitCuts_chi2");
0290   double FitCuts_ParErr_p0 = conf_.getParameter<double>("FitCuts_ParErr_p0");
0291   double p0_guess = conf_.getParameter<double>("p0_guess");
0292   double p1_guess = conf_.getParameter<double>("p1_guess");
0293   double p2_guess = conf_.getParameter<double>("p2_guess");
0294 
0295   double TIB1calib = 1.;
0296   double TIB2calib = 1.;
0297   double TIB3calib = 1.;
0298   double TIB4calib = 1.;
0299   double TOB1calib = 1.;
0300   double TOB2calib = 1.;
0301   double TOB3calib = 1.;
0302   double TOB4calib = 1.;
0303   double TOB5calib = 1.;
0304   double TOB6calib = 1.;
0305 
0306   if (CalibByMC == true) {
0307     //Calibration factors evaluated by using MC analysis
0308     TIB1calib = conf_.getParameter<double>("TIB1calib");
0309     TIB2calib = conf_.getParameter<double>("TIB2calib");
0310     TIB3calib = conf_.getParameter<double>("TIB3calib");
0311     TIB4calib = conf_.getParameter<double>("TIB4calib");
0312     TOB1calib = conf_.getParameter<double>("TOB1calib");
0313     TOB2calib = conf_.getParameter<double>("TOB2calib");
0314     TOB3calib = conf_.getParameter<double>("TOB3calib");
0315     TOB4calib = conf_.getParameter<double>("TOB4calib");
0316     TOB5calib = conf_.getParameter<double>("TOB5calib");
0317     TOB6calib = conf_.getParameter<double>("TOB6calib");
0318   }
0319 
0320   auto fitfunc = std::make_unique<TF1>("fitfunc", "([4]/[3])*[1]*(TMath::Abs(x-[0]))+[2]", -1, 1);
0321   auto fitfunc2IT = std::make_unique<TF1>("fitfunc2IT", "([4]/[3])*[1]*(TMath::Abs(x-[0]))+[2]", -1, 1);
0322 
0323   std::ofstream NoEntries;
0324   NoEntries.open(NoEntriesHisto_.c_str());
0325   std::ofstream Rep;
0326   Rep.open(LAreport_.c_str());
0327 
0328   gStyle->SetOptStat(1110);
0329 
0330   for (histo = histolist.begin(); histo != histolist.end(); ++histo) {
0331     FitFunction = nullptr;
0332     FitFunction2IT = nullptr;
0333     bool Good2ITFit = false;
0334     bool ModuleHisto = true;
0335 
0336     histoEntries = -99;
0337     gphi = -99;
0338     geta = -99;
0339     gz = -99;
0340     gR = -1;
0341     globalX = -99;
0342     globalY = -99;
0343     globalZ = -99;
0344     goodFit = 0;
0345     goodFit1IT = 0;
0346     badFit = 0;
0347     muH = -1;
0348     TIB = 0;
0349     TOB = 0;
0350     MonoStereo = -1;
0351 
0352     uint32_t id = hidmanager.getComponentId((*histo)->getName());
0353     DetId detid(id);
0354     StripSubdetector subid(id);
0355     const GeomDetUnit* stripdet;
0356     MonoStereo = subid.stereo();
0357 
0358     if (!(stripdet = tkGeom_->idToDetUnit(subid))) {
0359       no_mod_histo++;
0360       ModuleHisto = false;
0361       edm::LogInfo("SiStripCalibLorentzAngle") << "### NO MODULE HISTOGRAM";
0362     }
0363 
0364     if (stripdet != nullptr && ModuleHisto == true) {
0365       if (subid.subdetId() == int(StripSubdetector::TIB)) {
0366         Layer = tTopo_->tibLayer(detid);
0367         TIB = 1;
0368       }
0369       if (subid.subdetId() == int(StripSubdetector::TOB)) {
0370         Layer = tTopo_->tobLayer(detid);
0371         TOB = 1;
0372       }
0373 
0374       //get module coordinates
0375       const GlobalPoint gposition = (stripdet->surface()).toGlobal(p);
0376       histoEntries = (*histo)->getEntries();
0377       globalX = gposition.x();
0378       globalY = gposition.y();
0379       globalZ = gposition.z();
0380       gphi = gposition.phi();
0381       geta = gposition.eta();
0382       gR = sqrt(pow(gposition.x(), 2) + pow(gposition.y(), 2));
0383       gz = gposition.z();
0384 
0385       //get magnetic field
0386       const StripGeomDetUnit* det = dynamic_cast<const StripGeomDetUnit*>(tkGeom_->idToDetUnit(detid));
0387       if (det == nullptr) {
0388         edm::LogError("SiStripCalibLorentzAngle") << "[SiStripCalibLorentzAngle::getNewObject] the detID " << id
0389                                                   << " doesn't seem to belong to Tracker" << std::endl;
0390         continue;
0391       }
0392       LocalVector lbfield = (det->surface()).toLocal(magField.inTesla(det->surface().position()));
0393       theBfield = lbfield.mag();
0394       theBfield = (theBfield > 0) ? theBfield : 0.00001;
0395       TH1Ds["MagneticField"]->Fill(theBfield);
0396     }
0397     if (stripdet == nullptr)
0398       continue;
0399 
0400     if (((*histo)->getEntries() <= FitCuts_Entries) && ModuleHisto == true) {
0401       if (((*histo)->getEntries() == 0) && ModuleHisto == true) {
0402         NoEntries << "NO ENTRIES MODULE, ID = " << id << std::endl;
0403         edm::LogInfo("SiStripCalibLorentzAngle") << "### HISTOGRAM WITH 0 ENTRIES => TYPE:" << subid.subdetId();
0404         ZeroEntries++;
0405       } else {
0406         edm::LogInfo("SiStripCalibLorentzAngle")
0407             << "### HISTOGRAM WITH NR. ENTRIES <= ENTRIES_CUT => TYPE:" << subid.subdetId();
0408         NotEnoughEntries++;
0409       }
0410     }
0411 
0412     std::string name;
0413     if (TIB == 1) {
0414       name += "TIB";
0415     } else {
0416       name += "TOB";
0417     }
0418     std::stringstream LayerStream;
0419     LayerStream << Layer;
0420     name += LayerStream.str();
0421     std::stringstream idnum;
0422     idnum << id;
0423     name += "_Id_";
0424     name += idnum.str();
0425 
0426     gStyle->SetOptFit(111);
0427 
0428     //Extract TProfile from Monitor Element to ProfileMap
0429     Profiles[name] = new TProfile;
0430     TProfile* theProfile = ExtractTObject<TProfile>().extract(*histo);
0431     theProfile->Copy(*Profiles[name]);
0432     Profiles[name]->SetName(name.c_str());
0433 
0434     if (((*histo)->getEntries() > FitCuts_Entries) && ModuleHisto == true) {
0435       histocounter++;
0436       if (TIB == 1) {
0437         edm::LogInfo("SiStripCalibLorentzAngle") << "TIB layer = " << Layer;
0438       }
0439       if (TOB == 1) {
0440         edm::LogInfo("SiStripCalibLorentzAngle") << "TOB layer = " << Layer;
0441       }
0442       edm::LogInfo("SiStripCalibLorentzAngle") << "id: " << id;
0443 
0444       float thickness = stripdet->specificSurface().bounds().thickness();
0445       const StripTopology& topol = (const StripTopology&)stripdet->topology();
0446       float pitch = topol.localPitch(p);
0447 
0448       fitfunc->SetParameter(0, p0_guess);
0449       fitfunc->SetParameter(1, p1_guess);
0450       fitfunc->SetParameter(2, p2_guess);
0451       fitfunc->FixParameter(3, pitch);
0452       fitfunc->FixParameter(4, thickness);
0453 
0454       Profiles[name]->Fit(fitfunc.get(), "E", "", ModuleRangeMin, ModuleRangeMax);
0455 
0456       FitFunction = fitfunc.get();
0457       chi2norm = FitFunction->GetChisquare() / FitFunction->GetNDF();
0458 
0459       if (FitFunction->GetParameter(0) > FitCuts_p0 || FitFunction->GetParameter(1) < FitCuts_p1 ||
0460           FitFunction->GetParameter(2) < FitCuts_p2 || chi2norm > FitCuts_chi2 ||
0461           FitFunction->GetParError(0) < FitCuts_ParErr_p0) {
0462         FirstIT_badfit++;
0463 
0464         fitfunc2IT->SetParameter(0, p0_guess);
0465         fitfunc2IT->SetParameter(1, p1_guess);
0466         fitfunc2IT->SetParameter(2, p2_guess);
0467         fitfunc2IT->FixParameter(3, pitch);
0468         fitfunc2IT->FixParameter(4, thickness);
0469 
0470         //2nd Iteration
0471         Profiles[name]->Fit(fitfunc2IT.get(), "E", "", ModuleRangeMin2IT, ModuleRangeMax2IT);
0472 
0473         FitFunction = fitfunc2IT.get();
0474         chi2norm = FitFunction->GetChisquare() / FitFunction->GetNDF();
0475 
0476         //2nd Iteration failed
0477 
0478         if (FitFunction->GetParameter(0) > FitCuts_p0 || FitFunction->GetParameter(1) < FitCuts_p1 ||
0479             FitFunction->GetParameter(2) < FitCuts_p2 || chi2norm > FitCuts_chi2 ||
0480             FitFunction->GetParError(0) < FitCuts_ParErr_p0) {
0481           if (subid.subdetId() == int(StripSubdetector::TIB)) {
0482             Profiles[name]->SetDirectory(TIB_2IT_BadFit);
0483           } else {
0484             Profiles[name]->SetDirectory(TOB_2IT_BadFit);
0485           }
0486 
0487           SecondIT_badfit++;
0488           badFit = 1;
0489         }
0490 
0491         //2nd Iteration ok
0492 
0493         if (FitFunction->GetParameter(0) < FitCuts_p0 && FitFunction->GetParameter(1) > FitCuts_p1 &&
0494             FitFunction->GetParameter(2) > FitCuts_p2 && chi2norm < FitCuts_chi2 &&
0495             FitFunction->GetParError(0) > FitCuts_ParErr_p0) {
0496           if (subid.subdetId() == int(StripSubdetector::TIB)) {
0497             Profiles[name]->SetDirectory(TIB_2IT_GoodFit);
0498           } else {
0499             Profiles[name]->SetDirectory(TOB_2IT_GoodFit);
0500           }
0501 
0502           SecondIT_goodfit++;
0503           Good2ITFit = true;
0504         }
0505       }
0506 
0507       if (FitFunction->GetParameter(0) < FitCuts_p0 && FitFunction->GetParameter(1) > FitCuts_p1 &&
0508           FitFunction->GetParameter(2) > FitCuts_p2 && chi2norm < FitCuts_chi2 &&
0509           FitFunction->GetParError(0) > FitCuts_ParErr_p0) {
0510         if (Good2ITFit == false) {
0511           FirstIT_goodfit++;
0512           goodFit1IT = 1;
0513 
0514           if (subid.subdetId() == int(StripSubdetector::TIB)) {
0515             Profiles[name]->SetDirectory(TIB_1IT_GoodFit);
0516           } else {
0517             Profiles[name]->SetDirectory(TOB_1IT_GoodFit);
0518           }
0519         }
0520 
0521         GoodFit++;
0522         goodFit = 1;
0523 
0524         LorentzAngle_Plots->cd();
0525 
0526         edm::LogInfo("SiStripCalibLorentzAngle") << FitFunction->GetParameter(0);
0527 
0528         muH = -(FitFunction->GetParameter(0)) / theBfield;
0529 
0530         if (TIB == 1) {
0531           if (Layer == 1)
0532             muH = muH / TIB1calib;
0533           if (Layer == 2)
0534             muH = muH / TIB2calib;
0535           if (Layer == 3)
0536             muH = muH / TIB3calib;
0537           if (Layer == 4)
0538             muH = muH / TIB4calib;
0539         }
0540         if (TOB == 1) {
0541           if (Layer == 1)
0542             muH = muH / TOB1calib;
0543           if (Layer == 2)
0544             muH = muH / TOB2calib;
0545           if (Layer == 3)
0546             muH = muH / TOB3calib;
0547           if (Layer == 4)
0548             muH = muH / TOB4calib;
0549           if (Layer == 5)
0550             muH = muH / TOB5calib;
0551           if (Layer == 6)
0552             muH = muH / TOB6calib;
0553         }
0554 
0555         detid_la[id] = muH;
0556 
0557         if (TIB == 1) {
0558           TH1Ds["LA_TIB"]->Fill(muH);
0559           TH1Ds["LA_err_TIB"]->Fill(FitFunction->GetParError(0) / theBfield);
0560           TH1Ds["LA_chi2norm_TIB"]->Fill(chi2norm);
0561           TH2Ds["LA_phi_TIB"]->Fill(gphi, muH);
0562           TH2Ds["LA_eta_TIB"]->Fill(geta, muH);
0563           TH2Ds["LA_TIB_graph"]->Fill(Layer, muH);
0564 
0565           if (Layer == 1) {
0566             TH1Ds["LA_TIB_1"]->Fill(muH);
0567             TH2Ds["LA_phi_TIB1"]->Fill(gphi, muH);
0568             TH2Ds["LA_eta_TIB1"]->Fill(geta, muH);
0569             if (MonoStereo == 0) {
0570               TH1Ds["LA_TIB_1_mono"]->Fill(muH);
0571               TH2Ds["LA_phi_TIB1_mono"]->Fill(gphi, muH);
0572               TH2Ds["LA_eta_TIB1_mono"]->Fill(geta, muH);
0573             }
0574             if (MonoStereo == 1) {
0575               TH1Ds["LA_TIB_1_stereo"]->Fill(muH);
0576               TH2Ds["LA_phi_TIB1_stereo"]->Fill(gphi, muH);
0577               TH2Ds["LA_eta_TIB1_stereo"]->Fill(geta, muH);
0578             }
0579           }
0580 
0581           if (Layer == 2) {
0582             TH1Ds["LA_TIB_2"]->Fill(muH);
0583             TH2Ds["LA_phi_TIB2"]->Fill(gphi, muH);
0584             TH2Ds["LA_eta_TIB2"]->Fill(geta, muH);
0585             if (MonoStereo == 0) {
0586               TH1Ds["LA_TIB_2_mono"]->Fill(muH);
0587               TH2Ds["LA_phi_TIB2_mono"]->Fill(gphi, muH);
0588               TH2Ds["LA_eta_TIB2_mono"]->Fill(geta, muH);
0589             }
0590             if (MonoStereo == 1) {
0591               TH1Ds["LA_TIB_2_stereo"]->Fill(muH);
0592               TH2Ds["LA_phi_TIB2_stereo"]->Fill(gphi, muH);
0593               TH2Ds["LA_eta_TIB2_stereo"]->Fill(geta, muH);
0594             }
0595           }
0596 
0597           if (Layer == 3) {
0598             TH1Ds["LA_TIB_3"]->Fill(muH);
0599             TH2Ds["LA_phi_TIB3"]->Fill(gphi, muH);
0600             TH2Ds["LA_eta_TIB3"]->Fill(geta, muH);
0601           }
0602 
0603           if (Layer == 4) {
0604             TH1Ds["LA_TIB_4"]->Fill(muH);
0605             TH2Ds["LA_phi_TIB4"]->Fill(gphi, muH);
0606             TH2Ds["LA_eta_TIB4"]->Fill(geta, muH);
0607           }
0608         }
0609 
0610         if (TOB == 1) {
0611           TH1Ds["LA_TOB"]->Fill(muH);
0612           TH1Ds["LA_err_TOB"]->Fill(FitFunction->GetParError(0) / theBfield);
0613           TH1Ds["LA_chi2norm_TOB"]->Fill(chi2norm);
0614           TH2Ds["LA_phi_TOB"]->Fill(gphi, muH);
0615           TH2Ds["LA_eta_TOB"]->Fill(geta, muH);
0616           TH2Ds["LA_TOB_graph"]->Fill(Layer, muH);
0617 
0618           if (Layer == 1) {
0619             TH1Ds["LA_TOB_1"]->Fill(muH);
0620             TH2Ds["LA_phi_TOB1"]->Fill(gphi, muH);
0621             TH2Ds["LA_eta_TOB1"]->Fill(geta, muH);
0622             if (MonoStereo == 0) {
0623               TH1Ds["LA_TOB_1_mono"]->Fill(muH);
0624               TH2Ds["LA_phi_TOB1_mono"]->Fill(gphi, muH);
0625               TH2Ds["LA_eta_TOB1_mono"]->Fill(geta, muH);
0626             }
0627             if (MonoStereo == 1) {
0628               TH1Ds["LA_TOB_1_stereo"]->Fill(muH);
0629               TH2Ds["LA_phi_TOB1_stereo"]->Fill(gphi, muH);
0630               TH2Ds["LA_eta_TOB1_stereo"]->Fill(geta, muH);
0631             }
0632           }
0633 
0634           if (Layer == 2) {
0635             TH1Ds["LA_TOB_2"]->Fill(muH);
0636             TH2Ds["LA_phi_TOB2"]->Fill(gphi, muH);
0637             TH2Ds["LA_eta_TOB2"]->Fill(geta, muH);
0638             if (MonoStereo == 0) {
0639               TH1Ds["LA_TOB_2_mono"]->Fill(muH);
0640               TH2Ds["LA_phi_TOB2_mono"]->Fill(gphi, muH);
0641               TH2Ds["LA_eta_TOB2_mono"]->Fill(geta, muH);
0642             }
0643             if (MonoStereo == 1) {
0644               TH1Ds["LA_TOB_2_stereo"]->Fill(muH);
0645               TH2Ds["LA_phi_TOB2_stereo"]->Fill(gphi, muH);
0646               TH2Ds["LA_eta_TOB2_stereo"]->Fill(geta, muH);
0647             }
0648           }
0649 
0650           if (Layer == 3) {
0651             TH1Ds["LA_TOB_3"]->Fill(muH);
0652             TH2Ds["LA_phi_TOB3"]->Fill(gphi, muH);
0653             TH2Ds["LA_eta_TOB3"]->Fill(geta, muH);
0654           }
0655 
0656           if (Layer == 4) {
0657             TH1Ds["LA_TOB_4"]->Fill(muH);
0658             TH2Ds["LA_phi_TOB4"]->Fill(gphi, muH);
0659             TH2Ds["LA_eta_TOB4"]->Fill(geta, muH);
0660           }
0661 
0662           if (Layer == 5) {
0663             TH1Ds["LA_TOB_5"]->Fill(muH);
0664             TH2Ds["LA_phi_TOB5"]->Fill(gphi, muH);
0665             TH2Ds["LA_eta_TOB5"]->Fill(geta, muH);
0666           }
0667 
0668           if (Layer == 6) {
0669             TH1Ds["LA_TOB_6"]->Fill(muH);
0670             TH2Ds["LA_phi_TOB6"]->Fill(gphi, muH);
0671             TH2Ds["LA_eta_TOB6"]->Fill(geta, muH);
0672           }
0673         }
0674       }
0675     }
0676 
0677     ModuleTree->Fill();
0678   }
0679 
0680   double GaussFitRange = conf_.getParameter<double>("GaussFitRange");
0681   auto gaus = std::make_unique<TF1>("gaus", "gaus");
0682 
0683   TH1Ds["LA_TIB_1"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0684   mean_TIB1 = gaus->GetParameter(1);
0685   float err_mean_TIB1 = gaus->GetParError(1);
0686   float rms_TIB1 = gaus->GetParameter(2);
0687   TH1Ds["LA_TIB_2"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0688   mean_TIB2 = gaus->GetParameter(1);
0689   float err_mean_TIB2 = gaus->GetParError(1);
0690   float rms_TIB2 = gaus->GetParameter(2);
0691   TH1Ds["LA_TIB_3"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0692   mean_TIB3 = gaus->GetParameter(1);
0693   float err_mean_TIB3 = gaus->GetParError(1);
0694   float rms_TIB3 = gaus->GetParameter(2);
0695   TH1Ds["LA_TIB_4"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0696   mean_TIB4 = gaus->GetParameter(1);
0697   float err_mean_TIB4 = gaus->GetParError(1);
0698   float rms_TIB4 = gaus->GetParameter(2);
0699 
0700   TH1Ds["LA_TOB_1"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0701   mean_TOB1 = gaus->GetParameter(1);
0702   float err_mean_TOB1 = gaus->GetParError(1);
0703   float rms_TOB1 = gaus->GetParameter(2);
0704   TH1Ds["LA_TOB_2"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0705   mean_TOB2 = gaus->GetParameter(1);
0706   float err_mean_TOB2 = gaus->GetParError(1);
0707   float rms_TOB2 = gaus->GetParameter(2);
0708   TH1Ds["LA_TOB_3"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0709   mean_TOB3 = gaus->GetParameter(1);
0710   float err_mean_TOB3 = gaus->GetParError(1);
0711   float rms_TOB3 = gaus->GetParameter(2);
0712   TH1Ds["LA_TOB_4"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0713   mean_TOB4 = gaus->GetParameter(1);
0714   float err_mean_TOB4 = gaus->GetParError(1);
0715   float rms_TOB4 = gaus->GetParameter(2);
0716   TH1Ds["LA_TOB_5"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0717   mean_TOB5 = gaus->GetParameter(1);
0718   float err_mean_TOB5 = gaus->GetParError(1);
0719   float rms_TOB5 = gaus->GetParameter(2);
0720   TH1Ds["LA_TOB_6"]->Fit(gaus.get(), "", "", -GaussFitRange, GaussFitRange);
0721   mean_TOB6 = gaus->GetParameter(1);
0722   float err_mean_TOB6 = gaus->GetParError(1);
0723   float rms_TOB6 = gaus->GetParameter(2);
0724 
0725   int nlayersTIB = 4;
0726   float TIBx[4] = {1, 2, 3, 4};
0727   float TIBex[4] = {0, 0, 0, 0};
0728   float TIBy[4] = {mean_TIB1, mean_TIB2, mean_TIB3, mean_TIB4};
0729   float TIBey[4] = {err_mean_TIB1, err_mean_TIB2, err_mean_TIB3, err_mean_TIB4};
0730 
0731   int nlayersTOB = 6;
0732   float TOBx[6] = {1, 2, 3, 4, 5, 6};
0733   float TOBex[6] = {0, 0, 0, 0, 0, 0};
0734   float TOBy[6] = {mean_TOB1, mean_TOB2, mean_TOB3, mean_TOB4, mean_TOB5, mean_TOB6};
0735   float TOBey[6] = {err_mean_TOB1, err_mean_TOB2, err_mean_TOB3, err_mean_TOB4, err_mean_TOB5, err_mean_TOB6};
0736 
0737   TIB_graph = new TGraphErrors(nlayersTIB, TIBx, TIBy, TIBex, TIBey);
0738   TOB_graph = new TGraphErrors(nlayersTOB, TOBx, TOBy, TOBex, TOBey);
0739 
0740   //TF1 *fit_TIB= new TF1("fit_TIB","[0]",0,4);
0741   //TF1 *fit_TOB= new TF1("fit_TOB","[0]",0,6);
0742 
0743   gStyle->SetOptFit(111);
0744   gStyle->SetOptStat(111);
0745 
0746   TIB_graph->SetTitle("TIB Layers #mu_{H}");
0747   TIB_graph->GetXaxis()->SetTitle("Layers");
0748   TIB_graph->GetXaxis()->SetNdivisions(4);
0749   TIB_graph->GetYaxis()->SetTitle("#mu_{H}");
0750   TIB_graph->SetMarkerStyle(20);
0751   TIB_graph->GetYaxis()->SetTitleOffset(1.3);
0752   TIB_graph->Fit("fit_TIB", "E", "", 1, 4);
0753   meanMobility_TIB = TIB_graph->GetFunction("fit_TIB")->GetParameter(0);
0754 
0755   TOB_graph->SetTitle("TOB Layers #mu_{H}");
0756   TOB_graph->GetXaxis()->SetTitle("Layers");
0757   TOB_graph->GetXaxis()->SetNdivisions(6);
0758   TOB_graph->GetYaxis()->SetTitle("#mu_{H}");
0759   TOB_graph->SetMarkerStyle(20);
0760   TOB_graph->GetYaxis()->SetTitleOffset(1.3);
0761   TOB_graph->Fit("fit_TOB", "E", "", 1, 6);
0762   meanMobility_TOB = TOB_graph->GetFunction("fit_TOB")->GetParameter(0);
0763 
0764   TIB_graph->Write("TIB_graph");
0765   TOB_graph->Write("TOB_graph");
0766 
0767   Rep << "- NR.OF TIB AND TOB MODULES = 7932" << std::endl << std::endl << std::endl;
0768   Rep << "- NO MODULE HISTOS FOUND = " << no_mod_histo << std::endl << std::endl;
0769   Rep << "- NR.OF HISTOS WITH ENTRIES > " << FitCuts_Entries << " = " << histocounter << std::endl << std::endl;
0770   Rep << "- NR.OF HISTOS WITH ENTRIES <= " << FitCuts_Entries << " (!=0) = " << NotEnoughEntries << std::endl
0771       << std::endl;
0772   Rep << "- NR.OF HISTOS WITH 0 ENTRIES = " << ZeroEntries << std::endl << std::endl << std::endl;
0773   Rep << "- NR.OF GOOD FIT (FIRST IT + SECOND IT GOOD FIT)= " << GoodFit << std::endl << std::endl;
0774   Rep << "- NR.OF FIRST IT GOOD FIT = " << FirstIT_goodfit << std::endl << std::endl;
0775   Rep << "- NR.OF SECOND IT GOOD FIT = " << SecondIT_goodfit << std::endl << std::endl;
0776   Rep << "- NR.OF FIRST IT BAD FIT = " << FirstIT_badfit << std::endl << std::endl;
0777   Rep << "- NR.OF SECOND IT BAD FIT = " << SecondIT_badfit << std::endl << std::endl << std::endl;
0778 
0779   Rep << "--------------- Mean MuH values per Layer -------------------" << std::endl << std::endl << std::endl;
0780   Rep << "TIB1 = " << mean_TIB1 << " +- " << err_mean_TIB1 << " RMS = " << rms_TIB1 << std::endl;
0781   Rep << "TIB2 = " << mean_TIB2 << " +- " << err_mean_TIB2 << " RMS = " << rms_TIB2 << std::endl;
0782   Rep << "TIB3 = " << mean_TIB3 << " +- " << err_mean_TIB3 << " RMS = " << rms_TIB3 << std::endl;
0783   Rep << "TIB4 = " << mean_TIB4 << " +- " << err_mean_TIB4 << " RMS = " << rms_TIB4 << std::endl;
0784   Rep << "TOB1 = " << mean_TOB1 << " +- " << err_mean_TOB1 << " RMS = " << rms_TOB1 << std::endl;
0785   Rep << "TOB2 = " << mean_TOB2 << " +- " << err_mean_TOB2 << " RMS = " << rms_TOB2 << std::endl;
0786   Rep << "TOB3 = " << mean_TOB3 << " +- " << err_mean_TOB3 << " RMS = " << rms_TOB3 << std::endl;
0787   Rep << "TOB4 = " << mean_TOB4 << " +- " << err_mean_TOB4 << " RMS = " << rms_TOB4 << std::endl;
0788   Rep << "TOB5 = " << mean_TOB5 << " +- " << err_mean_TOB5 << " RMS = " << rms_TOB5 << std::endl;
0789   Rep << "TOB6 = " << mean_TOB6 << " +- " << err_mean_TOB6 << " RMS = " << rms_TOB6 << std::endl << std::endl;
0790   Rep << "Mean Hall Mobility TIB = " << meanMobility_TIB << " +- " << TIB_graph->GetFunction("fit_TIB")->GetParError(0)
0791       << std::endl;
0792   Rep << "Mean Hall Mobility TOB = " << meanMobility_TOB << " +- " << TOB_graph->GetFunction("fit_TOB")->GetParError(0)
0793       << std::endl
0794       << std::endl
0795       << std::endl;
0796 
0797   Rep.close();
0798   NoEntries.close();
0799 
0800   hFile->Write();
0801   hFile->Close();
0802 }
0803 
0804 // Virtual destructor needed.
0805 
0806 SiStripCalibLorentzAngle::~SiStripCalibLorentzAngle() { delete hFile; }
0807 
0808 // Analyzer: Functions that gets called by framework every event
0809 
0810 std::unique_ptr<SiStripLorentzAngle> SiStripCalibLorentzAngle::getNewObject() {
0811   auto LorentzAngle = std::make_unique<SiStripLorentzAngle>();
0812 
0813   if (!LayerDB) {
0814     for (std::map<uint32_t, float>::iterator it = detid_la.begin(); it != detid_la.end(); it++) {
0815       float langle = it->second;
0816       if (!LorentzAngle->putLorentzAngle(it->first, langle))
0817         edm::LogError("SiStripCalibLorentzAngle")
0818             << "[SiStripCalibLorentzAngle::analyze] detid already exists" << std::endl;
0819     }
0820   }
0821 
0822   else {
0823     const TrackerGeometry::DetIdContainer& Id = tkGeom_->detIds();
0824     TrackerGeometry::DetIdContainer::const_iterator Iditer;
0825 
0826     for (Iditer = Id.begin(); Iditer != Id.end(); Iditer++) {
0827       StripSubdetector subid(Iditer->rawId());
0828 
0829       hallMobility = 0.;
0830 
0831       if (subid.subdetId() == int(StripSubdetector::TIB)) {
0832         uint32_t tibLayer = tTopo_->tibLayer(*Iditer);
0833         if (tibLayer == 1) {
0834           hallMobility = mean_TIB1;
0835         }
0836         if (tibLayer == 2) {
0837           hallMobility = mean_TIB2;
0838         }
0839         if (tibLayer == 3) {
0840           hallMobility = mean_TIB3;
0841         }
0842         if (tibLayer == 4) {
0843           hallMobility = mean_TIB4;
0844         }
0845         if (!LorentzAngle->putLorentzAngle(Iditer->rawId(), hallMobility))
0846           edm::LogError("SiStripLorentzAngleGenerator") << " detid already exists" << std::endl;
0847       }
0848 
0849       if (subid.subdetId() == int(StripSubdetector::TOB)) {
0850         uint32_t tobLayer = tTopo_->tobLayer(*Iditer);
0851         if (tobLayer == 1) {
0852           hallMobility = mean_TOB1;
0853         }
0854         if (tobLayer == 2) {
0855           hallMobility = mean_TOB2;
0856         }
0857         if (tobLayer == 3) {
0858           hallMobility = mean_TOB3;
0859         }
0860         if (tobLayer == 4) {
0861           hallMobility = mean_TOB4;
0862         }
0863         if (tobLayer == 5) {
0864           hallMobility = mean_TOB5;
0865         }
0866         if (tobLayer == 6) {
0867           hallMobility = mean_TOB6;
0868         }
0869         if (!LorentzAngle->putLorentzAngle(Iditer->rawId(), hallMobility))
0870           edm::LogError("SiStripLorentzAngleGenerator") << " detid already exists" << std::endl;
0871       }
0872 
0873       if (subid.subdetId() == int(StripSubdetector::TID)) {
0874         hallMobility = meanMobility_TIB;
0875         if (!LorentzAngle->putLorentzAngle(Iditer->rawId(), hallMobility))
0876           edm::LogError("SiStripLorentzAngleGenerator") << " detid already exists" << std::endl;
0877       }
0878 
0879       if (subid.subdetId() == int(StripSubdetector::TEC)) {
0880         if (tTopo_->tecRing(subid) < 5) {
0881           hallMobility = meanMobility_TIB;
0882         } else {
0883           hallMobility = meanMobility_TOB;
0884         }
0885         if (!LorentzAngle->putLorentzAngle(Iditer->rawId(), hallMobility))
0886           edm::LogError("SiStripLorentzAngleGenerator") << " detid already exists" << std::endl;
0887       }
0888     }
0889   }
0890 
0891   return LorentzAngle;
0892 }