CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/CalibTracker/SiStripLorentzAngle/plugins/SiStripCalibLorentzAngle.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include <string>
00003 #include <iostream>
00004 #include <fstream>
00005 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00006 #include "CalibTracker/SiStripLorentzAngle/interface/SiStripCalibLorentzAngle.h"
00007 #include "FWCore/Framework/interface/ESHandle.h"
00008 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00009 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00010 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00011 #include "FWCore/ServiceRegistry/interface/Service.h"
00012 #include "CLHEP/Random/RandGauss.h"
00013 #include "CondFormats/DataRecord/interface/SiStripLorentzAngleRcd.h"
00014 #include "DQM/SiStripCommon/interface/SiStripHistoId.h"
00015 #include "DQMServices/Core/interface/MonitorElement.h"
00016 #include "DQMServices/Core/interface/DQMStore.h"
00017 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00018 #include "DQM/SiStripCommon/interface/ExtractTObject.h"
00019 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00020 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00021 
00022 SiStripCalibLorentzAngle::SiStripCalibLorentzAngle(edm::ParameterSet const& conf) : ConditionDBWriter<SiStripLorentzAngle>(conf) , conf_(conf){}
00023 
00024 void SiStripCalibLorentzAngle::algoBeginJob(const edm::EventSetup& c){
00025 
00026   c.get<TrackerDigiGeometryRecord>().get(estracker);
00027   tracker=&(*estracker); 
00028   
00029   //get magnetic field and geometry from ES
00030   edm::ESHandle<MagneticField> magfield_;
00031   c.get<IdealMagneticFieldRecord>().get(magfield_);
00032 
00033   edm::ESHandle<SiStripLorentzAngle> SiStripLorentzAngle_;
00034   c.get<SiStripLorentzAngleRcd>().get(SiStripLorentzAngle_);
00035   detid_la= SiStripLorentzAngle_->getLorentzAngles();
00036   
00037   DQMStore* dbe_ = edm::Service<DQMStore>().operator->();
00038   
00039   std::string inputFile_ =conf_.getUntrackedParameter<std::string>("fileName", "LAProfiles.root");
00040   std::string LAreport_ =conf_.getUntrackedParameter<std::string>("LA_Report", "LA_Report.txt");
00041   std::string NoEntriesHisto_ =conf_.getUntrackedParameter<std::string>("NoEntriesHisto", "NoEntriesHisto.txt");
00042   std::string Dir_Name_ =conf_.getUntrackedParameter<std::string>("Dir_Name", "SiStrip");
00043   
00044   LayerDB = conf_.getUntrackedParameter<bool>("LayerDB", false);
00045   
00046   CalibByMC = conf_.getUntrackedParameter<bool>("CalibByMC", false);
00047   
00048   dbe_->open(inputFile_);
00049   
00050   // use SistripHistoId for producing histogram id (and title)
00051   SiStripHistoId hidmanager;
00052   
00053   edm::LogInfo("SiStripCalibLorentzAngle")<<"### DIR-NAME = "<<Dir_Name_;
00054   histolist= dbe_->getAllContents(Dir_Name_);
00055   std::vector<MonitorElement*>::iterator histo;
00056   
00057   hFile = new TFile (conf_.getUntrackedParameter<std::string>("out_fileName").c_str(), "RECREATE" );
00058     
00059   LorentzAngle_Plots = hFile->mkdir("LorentzAngle_Plots");
00060   Rootple = LorentzAngle_Plots->mkdir("Rootple");
00061   MuH = LorentzAngle_Plots->mkdir("MuH");
00062   TIB_MuH = MuH->mkdir("TIB_MuH");
00063   TOB_MuH = MuH->mkdir("TOB_MuH");
00064   MuH_vs_Phi = LorentzAngle_Plots->mkdir("MuH_vs_Phi");
00065   TIB_Phi = MuH_vs_Phi->mkdir("TIB_Phi");
00066   TOB_Phi = MuH_vs_Phi->mkdir("TOB_Phi");
00067   MuH_vs_Eta = LorentzAngle_Plots->mkdir("MuH_vs_Eta");
00068   TIB_Eta = MuH_vs_Eta->mkdir("TIB_Eta");
00069   TOB_Eta = MuH_vs_Eta->mkdir("TOB_Eta");
00070   FirstIT_GoodFit_Histos = LorentzAngle_Plots->mkdir("1IT_GoodFit_Histos");
00071   TIB_1IT_GoodFit = FirstIT_GoodFit_Histos->mkdir("TIB_1IT_GoodFit");
00072   TOB_1IT_GoodFit = FirstIT_GoodFit_Histos->mkdir("TOB_1IT_GoodFit");
00073   SecondIT_GoodFit_Histos = LorentzAngle_Plots->mkdir("2IT_GoodFit_Histos");
00074   TIB_2IT_GoodFit = SecondIT_GoodFit_Histos->mkdir("TIB_2IT_GoodFit");
00075   TOB_2IT_GoodFit = SecondIT_GoodFit_Histos->mkdir("TOB_2IT_GoodFit");
00076   SecondIT_BadFit_Histos = LorentzAngle_Plots->mkdir("2IT_BadFit_Histos");
00077   TIB_2IT_BadFit = SecondIT_BadFit_Histos->mkdir("TIB_2IT_BadFit");
00078   TOB_2IT_BadFit = SecondIT_BadFit_Histos->mkdir("TOB_2IT_BadFit");
00079     
00080   TH1Ds["LA_TIB"] = new TH1D("TanLAPerTesla TIB","TanLAPerTesla TIB",1000,-0.5,0.5);
00081   TH1Ds["LA_TIB"]->SetDirectory(MuH);
00082   TH1Ds["LA_TOB"] = new TH1D("TanLAPerTesla TOB","TanLAPerTesla TOB",1000,-0.5,0.5);
00083   TH1Ds["LA_TOB"]->SetDirectory(MuH);
00084   TH1Ds["LA_err_TIB"] = new TH1D("TanLAPerTesla Error TIB","TanLAPerTesla Error TIB",1000,0,1);
00085   TH1Ds["LA_err_TIB"]->SetDirectory(MuH);
00086   TH1Ds["LA_err_TOB"] = new TH1D("TanLAPerTesla Error TOB","TanLAPerTesla Error TOB",1000,0,1);
00087   TH1Ds["LA_err_TOB"]->SetDirectory(MuH);
00088   TH1Ds["LA_chi2norm_TIB"] = new TH1D("TanLAPerTesla Chi2norm TIB","TanLAPerTesla Chi2norm TIB",2000,0,10);
00089   TH1Ds["LA_chi2norm_TIB"]->SetDirectory(MuH);
00090   TH1Ds["LA_chi2norm_TOB"] = new TH1D("TanLAPerTesla Chi2norm TOB","TanLAPerTesla Chi2norm TOB",2000,0,10);
00091   TH1Ds["LA_chi2norm_TOB"]->SetDirectory(MuH);
00092   TH1Ds["MagneticField"] = new TH1D("MagneticField","MagneticField",500,0,5);
00093   TH1Ds["MagneticField"]->SetDirectory(MuH);
00094   
00095   TH2Ds["LA_TIB_graph"] = new TH2D("TanLAPerTesla TIB Layers","TanLAPerTesla TIB Layers",60,0,5,1000,-0.3,0.3);
00096   TH2Ds["LA_TIB_graph"]->SetDirectory(MuH);
00097   TH2Ds["LA_TIB_graph"]->SetNdivisions(6);
00098   TH2Ds["LA_TOB_graph"] = new TH2D("TanLAPerTesla TOB Layers","TanLAPerTesla TOB Layers",80,0,7,1000,-0.3,0.3);
00099   TH2Ds["LA_TOB_graph"]->SetDirectory(MuH);
00100   TH2Ds["LA_TOB_graph"]->SetNdivisions(8);
00101 
00102   TH1Ds["LA_TIB_1"] = new TH1D("TanLAPerTesla TIB1","TanLAPerTesla TIB1",2000,-0.5,0.5);
00103   TH1Ds["LA_TIB_1"]->SetDirectory(TIB_MuH);
00104   TH1Ds["LA_TIB_1_mono"] = new TH1D("TanLAPerTesla TIB1 MONO","TanLAPerTesla TIB1 MONO",2000,-0.5,0.5);
00105   TH1Ds["LA_TIB_1_mono"]->SetDirectory(TIB_MuH);
00106   TH1Ds["LA_TIB_1_stereo"] = new TH1D("TanLAPerTesla TIB1 STEREO","TanLAPerTesla TIB1 STEREO",2000,-0.5,0.5);
00107   TH1Ds["LA_TIB_1_stereo"]->SetDirectory(TIB_MuH);
00108   TH1Ds["LA_TIB_2"] = new TH1D("TanLAPerTesla TIB2","TanLAPerTesla TIB2",2000,-0.5,0.5);
00109   TH1Ds["LA_TIB_2"]->SetDirectory(TIB_MuH);
00110   TH1Ds["LA_TIB_2_mono"] = new TH1D("TanLAPerTesla TIB2 MONO","TanLAPerTesla TIB2 MONO",2000,-0.5,0.5);
00111   TH1Ds["LA_TIB_2_mono"]->SetDirectory(TIB_MuH);
00112   TH1Ds["LA_TIB_2_stereo"] = new TH1D("TanLAPerTesla TIB2 STEREO","TanLAPerTesla TIB2 STEREO",2000,-0.5,0.5);
00113   TH1Ds["LA_TIB_2_stereo"]->SetDirectory(TIB_MuH);
00114   TH1Ds["LA_TIB_3"] = new TH1D("TanLAPerTesla_TIB 3","TanLAPerTesla TIB3",2000,-0.5,0.5);
00115   TH1Ds["LA_TIB_3"]->SetDirectory(TIB_MuH);
00116   TH1Ds["LA_TIB_4"] = new TH1D("TanLAPerTesla_TIB 4","TanLAPerTesla TIB4",2000,-0.5,0.5);
00117   TH1Ds["LA_TIB_4"]->SetDirectory(TIB_MuH);
00118   
00119   TH1Ds["LA_TOB_1"] = new TH1D("TanLAPerTesla TOB1","TanLAPerTesla TOB1",2000,-0.5,0.5);
00120   TH1Ds["LA_TOB_1"]->SetDirectory(TOB_MuH);
00121   TH1Ds["LA_TOB_1_mono"] = new TH1D("TanLAPerTesla TOB1 MONO","TanLAPerTesla TOB1 MONO",2000,-0.5,0.5);
00122   TH1Ds["LA_TOB_1_mono"]->SetDirectory(TOB_MuH);
00123   TH1Ds["LA_TOB_1_stereo"] = new TH1D("TanLAPerTesla TOB1 STEREO","TanLAPerTesla TOB1 STEREO",2000,-0.5,0.5);
00124   TH1Ds["LA_TOB_1_stereo"]->SetDirectory(TOB_MuH);
00125   TH1Ds["LA_TOB_2"] = new TH1D("TanLAPerTesla TOB2","TanLAPerTesla TOB2",2000,-0.5,0.5);
00126   TH1Ds["LA_TOB_2"]->SetDirectory(TOB_MuH);
00127   TH1Ds["LA_TOB_2_mono"] = new TH1D("TanLAPerTesla TOB2 MONO","TanLAPerTesla TOB2 MONO",2000,-0.5,0.5);
00128   TH1Ds["LA_TOB_2_mono"]->SetDirectory(TOB_MuH);
00129   TH1Ds["LA_TOB_2_stereo"] = new TH1D("TanLAPerTesla TOB2 STEREO","TanLAPerTesla TOB2 STEREO",2000,-0.5,0.5);
00130   TH1Ds["LA_TOB_2_stereo"]->SetDirectory(TOB_MuH);
00131   TH1Ds["LA_TOB_3"] = new TH1D("TanLAPerTesla TOB3","TanLAPerTesla TOB3",2000,-0.5,0.5);
00132   TH1Ds["LA_TOB_3"]->SetDirectory(TOB_MuH);
00133   TH1Ds["LA_TOB_4"] = new TH1D("TanLAPerTesla TOB4","TanLAPerTesla TOB4",2000,-0.5,0.5);
00134   TH1Ds["LA_TOB_4"]->SetDirectory(TOB_MuH);
00135   TH1Ds["LA_TOB_5"] = new TH1D("TanLAPerTesla TOB5","TanLAPerTesla TOB5",2000,-0.5,0.5);
00136   TH1Ds["LA_TOB_5"]->SetDirectory(TOB_MuH);
00137   TH1Ds["LA_TOB_6"] = new TH1D("TanLAPerTesla TOB6","TanLAPerTesla TOB6",2000,-0.5,0.5); 
00138   TH1Ds["LA_TOB_6"]->SetDirectory(TOB_MuH);
00139     
00140   TH2Ds["LA_phi_TIB"] = new TH2D("TanLAPerTesla vs Phi TIB","TanLAPerTesla vs Phi TIB",800,-4,4,600,-0.3,0.3);
00141   TH2Ds["LA_phi_TIB"]->SetDirectory(MuH_vs_Phi);
00142   TH2Ds["LA_phi_TOB"] = new TH2D("TanLAPerTesla vs Phi TOB","TanLAPerTesla vs Phi TOB",800,-4,4,600,-0.3,0.3);
00143   TH2Ds["LA_phi_TOB"]->SetDirectory(MuH_vs_Phi);
00144   
00145   TH2Ds["LA_phi_TIB1"] = new TH2D("TanLAPerTesla vs Phi TIB1","TanLAPerTesla vs Phi TIB1",800,-4,4,600,-0.3,0.3);
00146   TH2Ds["LA_phi_TIB1"]->SetDirectory(TIB_Phi);
00147   TH2Ds["LA_phi_TIB1_mono"] = new TH2D("TanLAPerTesla vs Phi TIB1 MONO","TanLAPerTesla vs Phi TIB1 MONO",800,-4,4,600,-0.3,0.3);
00148   TH2Ds["LA_phi_TIB1_mono"]->SetDirectory(TIB_Phi);
00149   TH2Ds["LA_phi_TIB1_stereo"] = new TH2D("TanLAPerTesla vs Phi TIB1 STEREO","TanLAPerTesla vs Phi TIB1 STEREO",800,-4,4,600,-0.3,0.3);
00150   TH2Ds["LA_phi_TIB1_stereo"]->SetDirectory(TIB_Phi);
00151   TH2Ds["LA_phi_TIB2"] = new TH2D("TanLAPerTesla vs Phi TIB2","TanLAPerTesla vs Phi TIB2",800,-4,4,600,-0.3,0.3);
00152   TH2Ds["LA_phi_TIB2"]->SetDirectory(TIB_Phi);
00153   TH2Ds["LA_phi_TIB2_mono"] = new TH2D("TanLAPerTesla vs Phi TIB2 MONO","TanLAPerTesla vs Phi TIB2 MONO",800,-4,4,600,-0.3,0.3);
00154   TH2Ds["LA_phi_TIB2_mono"]->SetDirectory(TIB_Phi);
00155   TH2Ds["LA_phi_TIB2_stereo"] = new TH2D("TanLAPerTesla vs Phi TIB2 STEREO","TanLAPerTesla vs Phi TIB2 STEREO",800,-4,4,600,-0.3,0.3);
00156   TH2Ds["LA_phi_TIB2_stereo"]->SetDirectory(TIB_Phi);
00157   TH2Ds["LA_phi_TIB3"] = new TH2D("TanLAPerTesla vs Phi TIB3","TanLAPerTesla vs Phi TIB3",800,-4,4,600,-0.3,0.3);
00158   TH2Ds["LA_phi_TIB3"]->SetDirectory(TIB_Phi);
00159   TH2Ds["LA_phi_TIB4"] = new TH2D("TanLAPerTesla vs Phi TIB4","TanLAPerTesla vs Phi TIB4",800,-4,4,600,-0.3,0.3);
00160   TH2Ds["LA_phi_TIB4"]->SetDirectory(TIB_Phi);
00161   
00162   TH2Ds["LA_phi_TOB1"] = new TH2D("TanLAPerTesla vs Phi TOB1","TanLAPerTesla vs Phi TOB1",800,-4,4,600,-0.3,0.3);
00163   TH2Ds["LA_phi_TOB1"]->SetDirectory(TOB_Phi);
00164   TH2Ds["LA_phi_TOB1_mono"] = new TH2D("TanLAPerTesla vs Phi TOB1 MONO","TanLAPerTesla vs Phi TOB1 MONO",800,-4,4,600,-0.3,0.3);
00165   TH2Ds["LA_phi_TOB1_mono"]->SetDirectory(TOB_Phi);
00166   TH2Ds["LA_phi_TOB1_stereo"] = new TH2D("TanLAPerTesla vs Phi TOB1 STEREO","TanLAPerTesla vs Phi TOB1 STEREO",800,-4,4,600,-0.3,0.3);
00167   TH2Ds["LA_phi_TOB1_stereo"]->SetDirectory(TOB_Phi);
00168   TH2Ds["LA_phi_TOB2"] = new TH2D("TanLAPerTesla vs Phi TOB2","TanLAPerTesla vs Phi TOB2",800,-4,4,600,-0.3,0.3);
00169   TH2Ds["LA_phi_TOB2"]->SetDirectory(TOB_Phi);
00170   TH2Ds["LA_phi_TOB2_mono"] = new TH2D("TanLAPerTesla vs Phi TOB2 MONO","TanLAPerTesla vs Phi TOB2 MONO",800,-4,4,600,-0.3,0.3);
00171   TH2Ds["LA_phi_TOB2_mono"]->SetDirectory(TOB_Phi);
00172   TH2Ds["LA_phi_TOB2_stereo"] = new TH2D("TanLAPerTesla vs Phi TOB2 STEREO","TanLAPerTesla vs Phi TOB2 STEREO",800,-4,4,600,-0.3,0.3);
00173   TH2Ds["LA_phi_TOB2_stereo"]->SetDirectory(TOB_Phi);
00174   TH2Ds["LA_phi_TOB3"] = new TH2D("TanLAPerTesla vs Phi TOB3","TanLAPerTesla vs Phi TOB3",800,-4,4,600,-0.3,0.3);
00175   TH2Ds["LA_phi_TOB3"]->SetDirectory(TOB_Phi);
00176   TH2Ds["LA_phi_TOB4"] = new TH2D("TanLAPerTesla vs Phi TOB4","TanLAPerTesla vs Phi TOB4",800,-4,4,600,-0.3,0.3);
00177   TH2Ds["LA_phi_TOB4"]->SetDirectory(TOB_Phi);
00178   TH2Ds["LA_phi_TOB5"] = new TH2D("TanLAPerTesla vs Phi TOB5","TanLAPerTesla vs Phi TOB5",800,-4,4,600,-0.3,0.3);
00179   TH2Ds["LA_phi_TOB5"]->SetDirectory(TOB_Phi);
00180   TH2Ds["LA_phi_TOB6"] = new TH2D("TanLAPerTesla vs Phi TOB6","TanLAPerTesla vs Phi TOB6",800,-4,4,600,-0.3,0.3);
00181   TH2Ds["LA_phi_TOB6"]->SetDirectory(TOB_Phi);
00182   
00183   TH2Ds["LA_eta_TIB"] = new TH2D("TanLAPerTesla vs Eta TIB","TanLAPerTesla vs Eta TIB",800,-2.6,2.6,600,-0.3,0.3);
00184   TH2Ds["LA_eta_TIB"]->SetDirectory(MuH_vs_Eta);
00185   TH2Ds["LA_eta_TOB"] = new TH2D("TanLAPerTesla vs Eta TOB","TanLAPerTesla vs Eta TOB",800,-2.6,2.6,600,-0.3,0.3);
00186   TH2Ds["LA_eta_TOB"]->SetDirectory(MuH_vs_Eta);
00187   
00188   TH2Ds["LA_eta_TIB1"] = new TH2D("TanLAPerTesla vs Eta TIB1","TanLAPerTesla vs Eta TIB1",800,-2.6,2.6,600,-0.3,0.3);
00189   TH2Ds["LA_eta_TIB1"]->SetDirectory(TIB_Eta);
00190   TH2Ds["LA_eta_TIB1_mono"] = new TH2D("TanLAPerTesla vs Eta TIB1 MONO","TanLAPerTesla vs Eta TIB1 MONO",800,-2.6,2.6,600,-0.3,0.3);
00191   TH2Ds["LA_eta_TIB1_mono"]->SetDirectory(TIB_Eta);
00192   TH2Ds["LA_eta_TIB1_stereo"] = new TH2D("TanLAPerTesla vs Eta TIB1 STEREO","TanLAPerTesla vs Eta TIB1 STEREO",800,-2.6,2.6,600,-0.3,0.3);
00193   TH2Ds["LA_eta_TIB1_stereo"]->SetDirectory(TIB_Eta);
00194   TH2Ds["LA_eta_TIB2"] = new TH2D("TanLAPerTesla vs Eta TIB2","TanLAPerTesla vs Eta TIB2",800,-2.6,2.6,600,-0.3,0.3);
00195   TH2Ds["LA_eta_TIB2"]->SetDirectory(TIB_Eta);
00196   TH2Ds["LA_eta_TIB2_mono"] = new TH2D("TanLAPerTesla vs Eta TIB2 MONO","TanLAPerTesla vs Eta TIB2 MONO",800,-2.6,2.6,600,-0.3,0.3);
00197   TH2Ds["LA_eta_TIB2_mono"]->SetDirectory(TIB_Eta);
00198   TH2Ds["LA_eta_TIB2_stereo"] = new TH2D("TanLAPerTesla vs Eta TIB2 STEREO","TanLAPerTesla vs Eta TIB2 STEREO",800,-2.6,2.6,600,-0.3,0.3);
00199   TH2Ds["LA_eta_TIB2_stereo"]->SetDirectory(TIB_Eta);
00200   TH2Ds["LA_eta_TIB3"] = new TH2D("TanLAPerTesla vs Eta TIB3","TanLAPerTesla vs Eta TIB3",800,-2.6,2.6,600,-0.3,0.3);
00201   TH2Ds["LA_eta_TIB3"]->SetDirectory(TIB_Eta);
00202   TH2Ds["LA_eta_TIB4"] = new TH2D("TanLAPerTesla vs Eta TIB4","TanLAPerTesla vs Eta TIB4",800,-2.6,2.6,600,-0.3,0.3);
00203   TH2Ds["LA_eta_TIB4"]->SetDirectory(TIB_Eta);
00204   
00205   TH2Ds["LA_eta_TOB1"] = new TH2D("TanLAPerTesla vs Eta TOB1","TanLAPerTesla vs Eta TOB1",800,-2.6,2.6,600,-0.3,0.3);
00206   TH2Ds["LA_eta_TOB1"]->SetDirectory(TIB_Eta);
00207   TH2Ds["LA_eta_TOB1_mono"] = new TH2D("TanLAPerTesla vs Eta TOB1 MONO","TanLAPerTesla vs Eta TOB1 MONO",800,-2.6,2.6,600,-0.3,0.3);
00208   TH2Ds["LA_eta_TOB1_mono"]->SetDirectory(TIB_Eta);
00209   TH2Ds["LA_eta_TOB1_stereo"] = new TH2D("TanLAPerTesla vs Eta TOB1 STEREO","TanLAPerTesla vs Eta TOB1 STEREO",800,-2.6,2.6,600,-0.3,0.3);
00210   TH2Ds["LA_eta_TOB1_stereo"]->SetDirectory(TIB_Eta);
00211   TH2Ds["LA_eta_TOB2"] = new TH2D("TanLAPerTesla vs Eta TOB2","TanLAPerTesla vs Eta TOB2",800,-2.6,2.6,600,-0.3,0.3);
00212   TH2Ds["LA_eta_TOB2"]->SetDirectory(TIB_Eta);
00213   TH2Ds["LA_eta_TOB2_mono"] = new TH2D("TanLAPerTesla vs Eta TOB2 MONO","TanLAPerTesla vs Eta TOB2 MONO",800,-2.6,2.6,600,-0.3,0.3);
00214   TH2Ds["LA_eta_TOB2_mono"]->SetDirectory(TIB_Eta);
00215   TH2Ds["LA_eta_TOB2_stereo"] = new TH2D("TanLAPerTesla vs Eta TOB2 STEREO","TanLAPerTesla vs Eta TOB2 STEREO",800,-2.6,2.6,600,-0.3,0.3);
00216   TH2Ds["LA_eta_TOB2_stereo"]->SetDirectory(TIB_Eta);
00217   TH2Ds["LA_eta_TOB3"] = new TH2D("TanLAPerTesla vs Eta TOB3","TanLAPerTesla vs Eta TOB3",800,-2.6,2.6,600,-0.3,0.3);
00218   TH2Ds["LA_eta_TOB3"]->SetDirectory(TIB_Eta);
00219   TH2Ds["LA_eta_TOB4"] = new TH2D("TanLAPerTesla vs Eta TOB4","TanLAPerTesla vs Eta TOB4",800,-2.6,2.6,600,-0.3,0.3);
00220   TH2Ds["LA_eta_TOB4"]->SetDirectory(TIB_Eta);
00221   TH2Ds["LA_eta_TOB5"] = new TH2D("TanLAPerTesla vs Eta TOB5","TanLAPerTesla vs Eta TOB5",800,-2.6,2.6,600,-0.3,0.3);
00222   TH2Ds["LA_eta_TOB5"]->SetDirectory(TIB_Eta);
00223   TH2Ds["LA_eta_TOB6"] = new TH2D("TanLAPerTesla vs Eta TOB6","TanLAPerTesla vs Eta TOB6",800,-2.6,2.6,600,-0.3,0.3);
00224   TH2Ds["LA_eta_TOB6"]->SetDirectory(TIB_Eta);
00225   
00226   ModuleTree = new TTree("ModuleTree", "ModuleTree");
00227   ModuleTree->Branch("histoEntries", &histoEntries, "histoEntries/F");
00228   ModuleTree->Branch("globalX", &globalX, "globalX/F");
00229   ModuleTree->Branch("globalY", &globalY, "globalY/F");
00230   ModuleTree->Branch("globalZ", &globalZ, "globalZ/F");
00231   ModuleTree->Branch("gphi", &gphi, "gphi/F");
00232   ModuleTree->Branch("geta", &geta, "geta/F");
00233   ModuleTree->Branch("gR", &gR, "gR/F");
00234   ModuleTree->Branch("goodFit", &goodFit, "goodFit/I");
00235   ModuleTree->Branch("goodFit1IT", &goodFit1IT, "goodFit1IT/I");
00236   ModuleTree->Branch("badFit", &badFit, "badFit/I");
00237   ModuleTree->Branch("TIB", &TIB, "TIB/I");
00238   ModuleTree->Branch("TOB", &TOB, "TOB/I");
00239   ModuleTree->Branch("Layer", &Layer, "Layer/I");
00240   ModuleTree->Branch("MonoStereo", &MonoStereo, "MonoStereo/I");
00241   ModuleTree->Branch("theBfield", &theBfield, "theBfield/F");
00242   ModuleTree->Branch("muH", &muH, "muH/F");
00243   
00244   ModuleTree->SetDirectory(Rootple);
00245      
00246   int histocounter = 0;
00247   int NotEnoughEntries = 0;
00248   int ZeroEntries = 0;
00249   int GoodFit = 0;
00250   int FirstIT_goodfit = 0;
00251   int FirstIT_badfit = 0;
00252   int SecondIT_badfit = 0;
00253   int SecondIT_goodfit = 0;
00254   int no_mod_histo = 0;
00255   float chi2norm = 0;
00256   LocalPoint p =LocalPoint(0,0,0);
00257   
00258   double ModuleRangeMin=conf_.getParameter<double>("ModuleFitXMin");
00259   double ModuleRangeMax=conf_.getParameter<double>("ModuleFitXMax");
00260   double ModuleRangeMin2IT=conf_.getParameter<double>("ModuleFit2ITXMin");
00261   double ModuleRangeMax2IT=conf_.getParameter<double>("ModuleFit2ITXMax");
00262   double FitCuts_Entries=conf_.getParameter<double>("FitCuts_Entries");
00263   double FitCuts_p0=conf_.getParameter<double>("FitCuts_p0");
00264   double FitCuts_p1=conf_.getParameter<double>("FitCuts_p1");
00265   double FitCuts_p2=conf_.getParameter<double>("FitCuts_p2");
00266   double FitCuts_chi2=conf_.getParameter<double>("FitCuts_chi2");
00267   double FitCuts_ParErr_p0=conf_.getParameter<double>("FitCuts_ParErr_p0");
00268   double p0_guess=conf_.getParameter<double>("p0_guess");
00269   double p1_guess=conf_.getParameter<double>("p1_guess");
00270   double p2_guess=conf_.getParameter<double>("p2_guess");
00271   
00272   double TIB1calib = 1.;
00273   double TIB2calib = 1.;
00274   double TIB3calib = 1.;
00275   double TIB4calib = 1.;
00276   double TOB1calib = 1.;
00277   double TOB2calib = 1.;
00278   double TOB3calib = 1.;
00279   double TOB4calib = 1.;
00280   double TOB5calib = 1.;
00281   double TOB6calib = 1.;
00282   
00283   if(CalibByMC==true){
00284   //Calibration factors evaluated by using MC analysis
00285   TIB1calib=conf_.getParameter<double>("TIB1calib");
00286   TIB2calib=conf_.getParameter<double>("TIB2calib");
00287   TIB3calib=conf_.getParameter<double>("TIB3calib");
00288   TIB4calib=conf_.getParameter<double>("TIB4calib");
00289   TOB1calib=conf_.getParameter<double>("TOB1calib");
00290   TOB2calib=conf_.getParameter<double>("TOB2calib");
00291   TOB3calib=conf_.getParameter<double>("TOB3calib");
00292   TOB4calib=conf_.getParameter<double>("TOB4calib");
00293   TOB5calib=conf_.getParameter<double>("TOB5calib");
00294   TOB6calib=conf_.getParameter<double>("TOB6calib");
00295   }
00296   
00297   
00298   TF1 *fitfunc= new TF1("fitfunc","([4]/[3])*[1]*(TMath::Abs(x-[0]))+[2]",-1,1);
00299   TF1 *fitfunc2IT= new TF1("fitfunc2IT","([4]/[3])*[1]*(TMath::Abs(x-[0]))+[2]",-1,1);
00300  
00301   ofstream NoEntries;
00302   NoEntries.open(NoEntriesHisto_.c_str());
00303   ofstream Rep;
00304   Rep.open(LAreport_.c_str());
00305   
00306   gStyle->SetOptStat(1110);
00307   
00308   for(histo=histolist.begin();histo!=histolist.end();++histo){
00309   
00310   FitFunction = 0;
00311   FitFunction2IT = 0;
00312   bool Good2ITFit = false;
00313   bool ModuleHisto = true;
00314   
00315   histoEntries = -99;
00316   gphi=-99;
00317   geta=-99;
00318   gz = -99;
00319   gR=-1;
00320   globalX = -99;
00321   globalY = -99;
00322   globalZ = -99;
00323   goodFit = 0;
00324   goodFit1IT = 0;
00325   badFit = 0;
00326   muH = -1;
00327   TIB = 0;
00328   TOB = 0;
00329   MonoStereo = -1;
00330   
00331     uint32_t id=hidmanager.getComponentId((*histo)->getName());
00332     DetId detid(id);
00333     StripSubdetector subid(id);
00334     const GeomDetUnit * stripdet;
00335     MonoStereo = subid.stereo();
00336     
00337     if(!(stripdet=tracker->idToDetUnit(subid))){
00338     no_mod_histo++;
00339     ModuleHisto=false;
00340     edm::LogInfo("SiStripCalibLorentzAngle")<<"### NO MODULE HISTOGRAM";}
00341     
00342     if(stripdet!=0 && ModuleHisto==true){
00343     
00344       if(subid.subdetId() == int (StripSubdetector::TIB)){
00345       TIBDetId TIBid=TIBDetId(subid);
00346       Layer = TIBid.layer();
00347       TIB = 1;}     
00348       if(subid.subdetId() == int (StripSubdetector::TOB)){
00349       TOBDetId TOBid=TOBDetId(subid);
00350       Layer = TOBid.layer();
00351       TOB = 1;}
00352       
00353     //get module coordinates
00354     const GlobalPoint gposition = (stripdet->surface()).toGlobal(p);
00355     histoEntries = (*histo)->getEntries();
00356     globalX = gposition.x();
00357     globalY = gposition.y();
00358     globalZ = gposition.z();
00359     gphi = gposition.phi();
00360     geta = gposition.eta();
00361     gR = sqrt(pow(gposition.x(),2)+pow(gposition.y(),2));
00362     gz = gposition.z();
00363     
00364     //get magnetic field
00365     const StripGeomDetUnit* det = dynamic_cast<const StripGeomDetUnit*>(estracker->idToDetUnit(detid));
00366     if (det==0){
00367     edm::LogError("SiStripCalibLorentzAngle") << "[SiStripCalibLorentzAngle::getNewObject] the detID " << id << " doesn't seem to belong to Tracker" <<std::endl;       
00368     continue;
00369     }
00370     LocalVector lbfield=(det->surface()).toLocal(magfield_->inTesla(det->surface().position()));
00371     theBfield = lbfield.mag();
00372     theBfield = (theBfield > 0) ? theBfield : 0.00001; 
00373     TH1Ds["MagneticField"]->Fill(theBfield);    
00374     }
00375     if(stripdet==0)continue;
00376       
00377     if(((*histo)->getEntries()<=FitCuts_Entries)&&ModuleHisto==true){
00378     if(((*histo)->getEntries()==0)&&ModuleHisto==true){    
00379     NoEntries<<"NO ENTRIES MODULE, ID = "<<id<<std::endl;    
00380     edm::LogInfo("SiStripCalibLorentzAngle")<<"### HISTOGRAM WITH 0 ENTRIES => TYPE:"<<subid.subdetId();
00381     ZeroEntries++;
00382     }else{    
00383     edm::LogInfo("SiStripCalibLorentzAngle")<<"### HISTOGRAM WITH NR. ENTRIES <= ENTRIES_CUT => TYPE:"<<subid.subdetId();
00384     NotEnoughEntries++;}    
00385     }
00386     
00387       std::string name;
00388       if(TIB==1){
00389       name+="TIB";
00390       }else{
00391       name+="TOB";}
00392       std::stringstream LayerStream;
00393       LayerStream<<Layer;
00394       name+=LayerStream.str();    
00395       std::stringstream idnum;
00396       idnum<<id;
00397       name+="_Id_";
00398       name+=idnum.str();
00399       
00400     gStyle->SetOptFit(111);  
00401     
00402     //Extract TProfile from Monitor Element to ProfileMap
00403     Profiles[name.c_str()] = new TProfile;   
00404     TProfile* theProfile=ExtractTObject<TProfile>().extract(*histo);   
00405     theProfile->Copy(*Profiles[name.c_str()]);    
00406     Profiles[name.c_str()]->SetName(name.c_str());
00407   
00408     if(((*histo)->getEntries()>FitCuts_Entries) && ModuleHisto==true){   
00409       histocounter++;            
00410       if(TIB==1){
00411       edm::LogInfo("SiStripCalibLorentzAngle")<<"TIB layer = "<<Layer;}    
00412       if(TOB==1){
00413       edm::LogInfo("SiStripCalibLorentzAngle")<<"TOB layer = "<<Layer;}     
00414       edm::LogInfo("SiStripCalibLorentzAngle")<<"id: "<<id;    
00415             
00416       float thickness=stripdet->specificSurface().bounds().thickness();
00417       const StripTopology& topol=(StripTopology&)stripdet->topology();
00418       float pitch = topol.localPitch(p);
00419                   
00420       fitfunc->SetParameter(0, p0_guess);
00421       fitfunc->SetParameter(1, p1_guess);
00422       fitfunc->SetParameter(2, p2_guess);
00423       fitfunc->FixParameter(3, pitch);
00424       fitfunc->FixParameter(4, thickness);
00425       
00426       Profiles[name.c_str()]->Fit("fitfunc","E","",ModuleRangeMin, ModuleRangeMax);
00427             
00428       FitFunction = Profiles[name.c_str()]->GetFunction("fitfunc");
00429       chi2norm = FitFunction->GetChisquare()/FitFunction->GetNDF();
00430       
00431       if(FitFunction->GetParameter(0)>FitCuts_p0 || FitFunction->GetParameter(1)<FitCuts_p1 || FitFunction->GetParameter(2)<FitCuts_p2 || chi2norm>FitCuts_chi2 || FitFunction->GetParError(0)<FitCuts_ParErr_p0){
00432       
00433       FirstIT_badfit++;   
00434             
00435       fitfunc2IT->SetParameter(0, p0_guess);
00436       fitfunc2IT->SetParameter(1, p1_guess);
00437       fitfunc2IT->SetParameter(2, p2_guess);
00438       fitfunc2IT->FixParameter(3, pitch);
00439       fitfunc2IT->FixParameter(4, thickness);
00440       
00441       //2nd Iteration
00442       Profiles[name.c_str()]->Fit("fitfunc2IT","E","",ModuleRangeMin2IT, ModuleRangeMax2IT);
00443       
00444       FitFunction = Profiles[name.c_str()]->GetFunction("fitfunc2IT");
00445       chi2norm = FitFunction->GetChisquare()/FitFunction->GetNDF();
00446                 
00447       //2nd Iteration failed
00448       
00449       if(FitFunction->GetParameter(0)>FitCuts_p0 || FitFunction->GetParameter(1)<FitCuts_p1 || FitFunction->GetParameter(2)<FitCuts_p2 || chi2norm>FitCuts_chi2 || FitFunction->GetParError(0)<FitCuts_ParErr_p0){
00450       
00451       if(subid.subdetId() == int (StripSubdetector::TIB)){
00452       Profiles[name.c_str()]->SetDirectory(TIB_2IT_BadFit);
00453       }else{
00454       Profiles[name.c_str()]->SetDirectory(TOB_2IT_BadFit);} 
00455             
00456       SecondIT_badfit++;
00457       badFit=1;     
00458           
00459       }
00460       
00461       //2nd Iteration ok
00462      
00463       if(FitFunction->GetParameter(0)<FitCuts_p0 && FitFunction->GetParameter(1)>FitCuts_p1 && FitFunction->GetParameter(2)>FitCuts_p2 && chi2norm<FitCuts_chi2 && FitFunction->GetParError(0)>FitCuts_ParErr_p0){
00464       
00465       if(subid.subdetId() == int (StripSubdetector::TIB)){
00466       Profiles[name.c_str()]->SetDirectory(TIB_2IT_GoodFit);
00467       }else{
00468       Profiles[name.c_str()]->SetDirectory(TOB_2IT_GoodFit);} 
00469       
00470       SecondIT_goodfit++;      
00471       Good2ITFit = true;
00472       
00473       }
00474             
00475       }  
00476                        
00477       if(FitFunction->GetParameter(0)<FitCuts_p0 && FitFunction->GetParameter(1)>FitCuts_p1 && FitFunction->GetParameter(2)>FitCuts_p2 && chi2norm<FitCuts_chi2 && FitFunction->GetParError(0)>FitCuts_ParErr_p0){
00478       
00479         if(Good2ITFit==false){
00480         
00481         FirstIT_goodfit++;
00482         goodFit1IT = 1;
00483         
00484         if(subid.subdetId() == int (StripSubdetector::TIB)){
00485         Profiles[name.c_str()]->SetDirectory(TIB_1IT_GoodFit);
00486         }else{
00487         Profiles[name.c_str()]->SetDirectory(TOB_1IT_GoodFit);} 
00488         
00489         }
00490         
00491         GoodFit++;
00492         goodFit=1;
00493         
00494         LorentzAngle_Plots->cd();
00495             
00496         edm::LogInfo("SiStripCalibLorentzAngle")<<FitFunction->GetParameter(0);
00497         
00498         muH = -(FitFunction->GetParameter(0))/theBfield;
00499         
00500         if(TIB==1){
00501         if(Layer==1) muH = muH/TIB1calib;
00502         if(Layer==2) muH = muH/TIB2calib;
00503         if(Layer==3) muH = muH/TIB3calib;
00504         if(Layer==4) muH = muH/TIB4calib;
00505         }
00506         if(TOB==1){
00507         if(Layer==1) muH = muH/TOB1calib;
00508         if(Layer==2) muH = muH/TOB2calib;
00509         if(Layer==3) muH = muH/TOB3calib;
00510         if(Layer==4) muH = muH/TOB4calib;
00511         if(Layer==5) muH = muH/TOB5calib;
00512         if(Layer==6) muH = muH/TOB6calib;
00513         }
00514         
00515         detid_la[id]= muH;
00516         
00517         if(TIB==1){
00518         
00519         TH1Ds["LA_TIB"]->Fill(muH);
00520         TH1Ds["LA_err_TIB"]->Fill(FitFunction->GetParError(0)/theBfield);
00521         TH1Ds["LA_chi2norm_TIB"]->Fill(chi2norm);
00522         TH2Ds["LA_phi_TIB"]->Fill(gphi,muH);
00523         TH2Ds["LA_eta_TIB"]->Fill(geta,muH);
00524         TH2Ds["LA_TIB_graph"]->Fill(Layer,muH);
00525         
00526         if(Layer==1){
00527         TH1Ds["LA_TIB_1"]->Fill(muH);
00528         TH2Ds["LA_phi_TIB1"]->Fill(gphi,muH);
00529         TH2Ds["LA_eta_TIB1"]->Fill(geta,muH);
00530         if(MonoStereo==0){
00531         TH1Ds["LA_TIB_1_mono"]->Fill(muH);
00532         TH2Ds["LA_phi_TIB1_mono"]->Fill(gphi,muH);
00533         TH2Ds["LA_eta_TIB1_mono"]->Fill(geta,muH);}
00534         if(MonoStereo==1){
00535         TH1Ds["LA_TIB_1_stereo"]->Fill(muH);
00536         TH2Ds["LA_phi_TIB1_stereo"]->Fill(gphi,muH);
00537         TH2Ds["LA_eta_TIB1_stereo"]->Fill(geta,muH);}
00538         }
00539         
00540         if(Layer==2){
00541         TH1Ds["LA_TIB_2"]->Fill(muH);
00542         TH2Ds["LA_phi_TIB2"]->Fill(gphi,muH);
00543         TH2Ds["LA_eta_TIB2"]->Fill(geta,muH);
00544         if(MonoStereo==0){
00545         TH1Ds["LA_TIB_2_mono"]->Fill(muH);
00546         TH2Ds["LA_phi_TIB2_mono"]->Fill(gphi,muH);
00547         TH2Ds["LA_eta_TIB2_mono"]->Fill(geta,muH);}
00548         if(MonoStereo==1){
00549         TH1Ds["LA_TIB_2_stereo"]->Fill(muH);
00550         TH2Ds["LA_phi_TIB2_stereo"]->Fill(gphi,muH);
00551         TH2Ds["LA_eta_TIB2_stereo"]->Fill(geta,muH);}
00552         }
00553         
00554         if(Layer==3){
00555         TH1Ds["LA_TIB_3"]->Fill(muH);
00556         TH2Ds["LA_phi_TIB3"]->Fill(gphi,muH);
00557         TH2Ds["LA_eta_TIB3"]->Fill(geta,muH);
00558         }
00559         
00560         if(Layer==4){
00561         TH1Ds["LA_TIB_4"]->Fill(muH);
00562         TH2Ds["LA_phi_TIB4"]->Fill(gphi,muH);
00563         TH2Ds["LA_eta_TIB4"]->Fill(geta,muH);
00564         }
00565         }       
00566         
00567         if(TOB==1){
00568         
00569         TH1Ds["LA_TOB"]->Fill(muH);
00570         TH1Ds["LA_err_TOB"]->Fill(FitFunction->GetParError(0)/theBfield);
00571         TH1Ds["LA_chi2norm_TOB"]->Fill(chi2norm);
00572         TH2Ds["LA_phi_TOB"]->Fill(gphi,muH);
00573         TH2Ds["LA_eta_TOB"]->Fill(geta,muH);
00574         TH2Ds["LA_TOB_graph"]->Fill(Layer,muH);
00575         
00576         if(Layer==1){
00577         TH1Ds["LA_TOB_1"]->Fill(muH);
00578         TH2Ds["LA_phi_TOB1"]->Fill(gphi,muH);
00579         TH2Ds["LA_eta_TOB1"]->Fill(geta,muH);
00580         if(MonoStereo==0){
00581         TH1Ds["LA_TOB_1_mono"]->Fill(muH);
00582         TH2Ds["LA_phi_TOB1_mono"]->Fill(gphi,muH);
00583         TH2Ds["LA_eta_TOB1_mono"]->Fill(geta,muH);}
00584         if(MonoStereo==1){
00585         TH1Ds["LA_TOB_1_stereo"]->Fill(muH);
00586         TH2Ds["LA_phi_TOB1_stereo"]->Fill(gphi,muH);
00587         TH2Ds["LA_eta_TOB1_stereo"]->Fill(geta,muH);}
00588         }
00589         
00590         if(Layer==2){
00591         TH1Ds["LA_TOB_2"]->Fill(muH);
00592         TH2Ds["LA_phi_TOB2"]->Fill(gphi,muH);
00593         TH2Ds["LA_eta_TOB2"]->Fill(geta,muH);
00594         if(MonoStereo==0){
00595         TH1Ds["LA_TOB_2_mono"]->Fill(muH);
00596         TH2Ds["LA_phi_TOB2_mono"]->Fill(gphi,muH);
00597         TH2Ds["LA_eta_TOB2_mono"]->Fill(geta,muH);}
00598         if(MonoStereo==1){
00599         TH1Ds["LA_TOB_2_stereo"]->Fill(muH);
00600         TH2Ds["LA_phi_TOB2_stereo"]->Fill(gphi,muH);
00601         TH2Ds["LA_eta_TOB2_stereo"]->Fill(geta,muH);}
00602         }
00603         
00604         if(Layer==3){
00605         TH1Ds["LA_TOB_3"]->Fill(muH);
00606         TH2Ds["LA_phi_TOB3"]->Fill(gphi,muH);
00607         TH2Ds["LA_eta_TOB3"]->Fill(geta,muH);
00608         }
00609         
00610         if(Layer==4){
00611         TH1Ds["LA_TOB_4"]->Fill(muH);
00612         TH2Ds["LA_phi_TOB4"]->Fill(gphi,muH);
00613         TH2Ds["LA_eta_TOB4"]->Fill(geta,muH);
00614         }
00615         
00616         if(Layer==5){
00617         TH1Ds["LA_TOB_5"]->Fill(muH);
00618         TH2Ds["LA_phi_TOB5"]->Fill(gphi,muH);
00619         TH2Ds["LA_eta_TOB5"]->Fill(geta,muH);
00620         }
00621         
00622         if(Layer==6){
00623         TH1Ds["LA_TOB_6"]->Fill(muH);
00624         TH2Ds["LA_phi_TOB6"]->Fill(gphi,muH);
00625         TH2Ds["LA_eta_TOB6"]->Fill(geta,muH);
00626         }
00627         }
00628         
00629        }
00630     }
00631     
00632     ModuleTree->Fill();
00633      
00634   }
00635   
00636   double GaussFitRange=conf_.getParameter<double>("GaussFitRange");
00637   
00638   TH1Ds["LA_TIB_1"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00639   mean_TIB1 = TH1Ds["LA_TIB_1"]->GetFunction("gaus")->GetParameter(1);
00640   float err_mean_TIB1 = TH1Ds["LA_TIB_1"]->GetFunction("gaus")->GetParError(1);
00641   float rms_TIB1 = TH1Ds["LA_TIB_1"]->GetFunction("gaus")->GetParameter(2);
00642   TH1Ds["LA_TIB_2"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00643   mean_TIB2 = TH1Ds["LA_TIB_2"]->GetFunction("gaus")->GetParameter(1);
00644   float err_mean_TIB2 = TH1Ds["LA_TIB_2"]->GetFunction("gaus")->GetParError(1);
00645   float rms_TIB2 = TH1Ds["LA_TIB_2"]->GetFunction("gaus")->GetParameter(2);
00646   TH1Ds["LA_TIB_3"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00647   mean_TIB3 = TH1Ds["LA_TIB_3"]->GetFunction("gaus")->GetParameter(1);
00648   float err_mean_TIB3 = TH1Ds["LA_TIB_3"]->GetFunction("gaus")->GetParError(1);
00649   float rms_TIB3 = TH1Ds["LA_TIB_3"]->GetFunction("gaus")->GetParameter(2);
00650   TH1Ds["LA_TIB_4"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00651   mean_TIB4 = TH1Ds["LA_TIB_4"]->GetFunction("gaus")->GetParameter(1);
00652   float err_mean_TIB4 = TH1Ds["LA_TIB_4"]->GetFunction("gaus")->GetParError(1);
00653   float rms_TIB4 = TH1Ds["LA_TIB_4"]->GetFunction("gaus")->GetParameter(2);
00654 
00655   TH1Ds["LA_TOB_1"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00656   mean_TOB1 = TH1Ds["LA_TOB_1"]->GetFunction("gaus")->GetParameter(1);
00657   float err_mean_TOB1 = TH1Ds["LA_TOB_1"]->GetFunction("gaus")->GetParError(1);
00658   float rms_TOB1 = TH1Ds["LA_TOB_1"]->GetFunction("gaus")->GetParameter(2);
00659   TH1Ds["LA_TOB_2"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00660   mean_TOB2 = TH1Ds["LA_TOB_2"]->GetFunction("gaus")->GetParameter(1);
00661   float err_mean_TOB2 = TH1Ds["LA_TOB_2"]->GetFunction("gaus")->GetParError(1);
00662   float rms_TOB2 = TH1Ds["LA_TOB_2"]->GetFunction("gaus")->GetParameter(2);
00663   TH1Ds["LA_TOB_3"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00664   mean_TOB3 = TH1Ds["LA_TOB_3"]->GetFunction("gaus")->GetParameter(1);
00665   float err_mean_TOB3 = TH1Ds["LA_TOB_3"]->GetFunction("gaus")->GetParError(1);
00666   float rms_TOB3 = TH1Ds["LA_TOB_3"]->GetFunction("gaus")->GetParameter(2);
00667   TH1Ds["LA_TOB_4"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00668   mean_TOB4 = TH1Ds["LA_TOB_4"]->GetFunction("gaus")->GetParameter(1);
00669   float err_mean_TOB4 = TH1Ds["LA_TOB_4"]->GetFunction("gaus")->GetParError(1);
00670   float rms_TOB4 = TH1Ds["LA_TOB_4"]->GetFunction("gaus")->GetParameter(2);
00671   TH1Ds["LA_TOB_5"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00672   mean_TOB5 = TH1Ds["LA_TOB_5"]->GetFunction("gaus")->GetParameter(1);
00673   float err_mean_TOB5 = TH1Ds["LA_TOB_5"]->GetFunction("gaus")->GetParError(1);
00674   float rms_TOB5 = TH1Ds["LA_TOB_5"]->GetFunction("gaus")->GetParameter(2);
00675   TH1Ds["LA_TOB_6"]->Fit("gaus","","",-GaussFitRange,GaussFitRange);  
00676   mean_TOB6 = TH1Ds["LA_TOB_6"]->GetFunction("gaus")->GetParameter(1);
00677   float err_mean_TOB6 = TH1Ds["LA_TOB_6"]->GetFunction("gaus")->GetParError(1);
00678   float rms_TOB6 = TH1Ds["LA_TOB_6"]->GetFunction("gaus")->GetParameter(2);
00679 
00680 int nlayersTIB = 4;
00681 float TIBx[4]={1,2,3,4};
00682 float TIBex[4]={0,0,0,0};
00683 float TIBy[4]={mean_TIB1, mean_TIB2, mean_TIB3, mean_TIB4};
00684 float TIBey[4]={err_mean_TIB1, err_mean_TIB2, err_mean_TIB3, err_mean_TIB4};
00685 
00686 int nlayersTOB = 6;
00687 float TOBx[6]={1,2,3,4,5,6};
00688 float TOBex[6]={0,0,0,0,0,0};
00689 float TOBy[6]={mean_TOB1, mean_TOB2, mean_TOB3, mean_TOB4, mean_TOB5, mean_TOB6};
00690 float TOBey[6]={err_mean_TOB1, err_mean_TOB2, err_mean_TOB3, err_mean_TOB4, err_mean_TOB5, err_mean_TOB6};
00691 
00692 TIB_graph = new TGraphErrors(nlayersTIB,TIBx,TIBy,TIBex,TIBey);
00693 TOB_graph = new TGraphErrors(nlayersTOB,TOBx,TOBy,TOBex,TOBey);
00694 
00695 //TF1 *fit_TIB= new TF1("fit_TIB","[0]",0,4);
00696 //TF1 *fit_TOB= new TF1("fit_TOB","[0]",0,6);
00697 
00698 gStyle->SetOptFit(111);
00699 gStyle->SetOptStat(111);
00700 
00701 TIB_graph->SetTitle("TIB Layers #mu_{H}");
00702 TIB_graph->GetXaxis()->SetTitle("Layers");
00703 TIB_graph->GetXaxis()->SetNdivisions(4);
00704 TIB_graph->GetYaxis()->SetTitle("#mu_{H}");
00705 TIB_graph->SetMarkerStyle(20);
00706 TIB_graph->GetYaxis()->SetTitleOffset(1.3);
00707 TIB_graph->Fit("fit_TIB","E","",1,4);
00708 meanMobility_TIB = TIB_graph->GetFunction("fit_TIB")->GetParameter(0);
00709 
00710 TOB_graph->SetTitle("TOB Layers #mu_{H}");
00711 TOB_graph->GetXaxis()->SetTitle("Layers");
00712 TOB_graph->GetXaxis()->SetNdivisions(6);
00713 TOB_graph->GetYaxis()->SetTitle("#mu_{H}");
00714 TOB_graph->SetMarkerStyle(20);
00715 TOB_graph->GetYaxis()->SetTitleOffset(1.3);
00716 TOB_graph->Fit("fit_TOB","E","",1,6);
00717 meanMobility_TOB = TOB_graph->GetFunction("fit_TOB")->GetParameter(0);
00718 
00719   TIB_graph->Write("TIB_graph");
00720   TOB_graph->Write("TOB_graph");
00721 
00722   Rep<<"- NR.OF TIB AND TOB MODULES = 7932"<<std::endl<<std::endl<<std::endl;
00723   Rep<<"- NO MODULE HISTOS FOUND = "<<no_mod_histo<<std::endl<<std::endl;
00724   Rep<<"- NR.OF HISTOS WITH ENTRIES > "<<FitCuts_Entries<<" = "<<histocounter<<std::endl<<std::endl;
00725   Rep<<"- NR.OF HISTOS WITH ENTRIES <= "<<FitCuts_Entries<<" (!=0) = "<<NotEnoughEntries<<std::endl<<std::endl;
00726   Rep<<"- NR.OF HISTOS WITH 0 ENTRIES = "<<ZeroEntries<<std::endl<<std::endl<<std::endl;
00727   Rep<<"- NR.OF GOOD FIT (FIRST IT + SECOND IT GOOD FIT)= "<<GoodFit<<std::endl<<std::endl;
00728   Rep<<"- NR.OF FIRST IT GOOD FIT = "<<FirstIT_goodfit<<std::endl<<std::endl;
00729   Rep<<"- NR.OF SECOND IT GOOD FIT = "<<SecondIT_goodfit<<std::endl<<std::endl;
00730   Rep<<"- NR.OF FIRST IT BAD FIT = "<<FirstIT_badfit<<std::endl<<std::endl;
00731   Rep<<"- NR.OF SECOND IT BAD FIT = "<<SecondIT_badfit<<std::endl<<std::endl<<std::endl;
00732   
00733   Rep<<"--------------- Mean MuH values per Layer -------------------"<<std::endl<<std::endl<<std::endl;
00734   Rep<<"TIB1 = "<<mean_TIB1<<" +- "<<err_mean_TIB1<<" RMS = "<<rms_TIB1<<std::endl;
00735   Rep<<"TIB2 = "<<mean_TIB2<<" +- "<<err_mean_TIB2<<" RMS = "<<rms_TIB2<<std::endl;
00736   Rep<<"TIB3 = "<<mean_TIB3<<" +- "<<err_mean_TIB3<<" RMS = "<<rms_TIB3<<std::endl;
00737   Rep<<"TIB4 = "<<mean_TIB4<<" +- "<<err_mean_TIB4<<" RMS = "<<rms_TIB4<<std::endl;
00738   Rep<<"TOB1 = "<<mean_TOB1<<" +- "<<err_mean_TOB1<<" RMS = "<<rms_TOB1<<std::endl;
00739   Rep<<"TOB2 = "<<mean_TOB2<<" +- "<<err_mean_TOB2<<" RMS = "<<rms_TOB2<<std::endl;
00740   Rep<<"TOB3 = "<<mean_TOB3<<" +- "<<err_mean_TOB3<<" RMS = "<<rms_TOB3<<std::endl;
00741   Rep<<"TOB4 = "<<mean_TOB4<<" +- "<<err_mean_TOB4<<" RMS = "<<rms_TOB4<<std::endl;
00742   Rep<<"TOB5 = "<<mean_TOB5<<" +- "<<err_mean_TOB5<<" RMS = "<<rms_TOB5<<std::endl;
00743   Rep<<"TOB6 = "<<mean_TOB6<<" +- "<<err_mean_TOB6<<" RMS = "<<rms_TOB6<<std::endl<<std::endl;
00744   Rep<<"Mean Hall Mobility TIB = "<<meanMobility_TIB<<" +- "<<TIB_graph->GetFunction("fit_TIB")->GetParError(0)<<std::endl;
00745   Rep<<"Mean Hall Mobility TOB = "<<meanMobility_TOB<<" +- "<<TOB_graph->GetFunction("fit_TOB")->GetParError(0)<<std::endl<<std::endl<<std::endl;
00746         
00747   Rep.close();
00748   NoEntries.close(); 
00749   
00750 hFile->Write();
00751 hFile->Close();
00752     
00753 }
00754 
00755 // Virtual destructor needed.
00756 
00757 SiStripCalibLorentzAngle::~SiStripCalibLorentzAngle(){
00758  delete hFile;
00759 }
00760   
00761 
00762 // Analyzer: Functions that gets called by framework every event
00763    
00764 SiStripLorentzAngle* SiStripCalibLorentzAngle::getNewObject(){
00765 
00766   SiStripLorentzAngle* LorentzAngle = new SiStripLorentzAngle();
00767   
00768   if(!LayerDB){
00769   for(std::map<uint32_t, float>::iterator it = detid_la.begin(); it != detid_la.end(); it++){
00770     
00771     float langle=it->second;
00772     if ( ! LorentzAngle->putLorentzAngle(it->first,langle) )
00773       edm::LogError("SiStripCalibLorentzAngle")<<"[SiStripCalibLorentzAngle::analyze] detid already exists"<<std::endl;
00774   }
00775   }
00776   
00777   else{
00778   
00779     const TrackerGeometry::DetIdContainer& Id = estracker->detIds();
00780     TrackerGeometry::DetIdContainer::const_iterator Iditer; 
00781        
00782     for(Iditer=Id.begin();Iditer!=Id.end();Iditer++){
00783   
00784     StripSubdetector subid(Iditer->rawId());
00785     
00786     hallMobility = 0.;
00787     
00788     if(subid.subdetId() == int (StripSubdetector::TIB)){
00789     TIBDetId TIBid=TIBDetId(subid);
00790     if(TIBid.layer()==1){
00791     hallMobility=mean_TIB1;}
00792     if(TIBid.layer()==2){
00793     hallMobility=mean_TIB2;}
00794     if(TIBid.layer()==3){
00795     hallMobility=mean_TIB3;}
00796     if(TIBid.layer()==4){
00797     hallMobility=mean_TIB4;}
00798     if (!LorentzAngle->putLorentzAngle(Iditer->rawId(),hallMobility)) edm::LogError("SiStripLorentzAngleGenerator")<<" detid already exists"<<std::endl;
00799     }
00800     
00801     if(subid.subdetId() == int (StripSubdetector::TOB)){
00802     TOBDetId TOBid=TOBDetId(subid);
00803     if(TOBid.layer()==1){
00804     hallMobility=mean_TOB1;}
00805     if(TOBid.layer()==2){
00806     hallMobility=mean_TOB2;}
00807     if(TOBid.layer()==3){
00808     hallMobility=mean_TOB3;}
00809     if(TOBid.layer()==4){
00810     hallMobility=mean_TOB4;}
00811     if(TOBid.layer()==5){
00812     hallMobility=mean_TOB5;}
00813     if(TOBid.layer()==6){
00814     hallMobility=mean_TOB6;}
00815     if (!LorentzAngle->putLorentzAngle(Iditer->rawId(),hallMobility)) edm::LogError("SiStripLorentzAngleGenerator")<<" detid already exists"<<std::endl;
00816     } 
00817      
00818     if( subid.subdetId() == int(StripSubdetector::TID) ) {
00819     hallMobility=meanMobility_TIB;
00820     if (!LorentzAngle->putLorentzAngle(Iditer->rawId(),hallMobility)) edm::LogError("SiStripLorentzAngleGenerator")<<" detid already exists"<<std::endl;
00821     } 
00822     
00823     if( subid.subdetId() == int(StripSubdetector::TEC) ) {
00824     TECDetId TECid = TECDetId(subid);
00825     if(TECid.ringNumber()<5 ) {
00826     hallMobility=meanMobility_TIB;
00827     }else{
00828     hallMobility=meanMobility_TOB;
00829     }
00830     if (!LorentzAngle->putLorentzAngle(Iditer->rawId(),hallMobility)) edm::LogError("SiStripLorentzAngleGenerator")<<" detid already exists"<<std::endl;
00831     }
00832                  
00833   }
00834   }
00835   
00836   return LorentzAngle;
00837   
00838 }
00839 
00840 
00841 
00842 
00843 
00844   
00845   
00846   
00847 
00848