Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "Validation/RecoVertex/interface/BeamSpotHistogramMaker.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "FWCore/ServiceRegistry/interface/Service.h"
0005 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0006 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0007 #include "TH2F.h"
0008 #include "TH1F.h"
0009 #include "TProfile.h"
0010 
0011 BeamSpotHistogramMaker::BeamSpotHistogramMaker(edm::ConsumesCollector&& iC)
0012     : _currdir(nullptr), _histoParameters(), _rhm(iC) {}
0013 
0014 BeamSpotHistogramMaker::BeamSpotHistogramMaker(const edm::ParameterSet& iConfig, edm::ConsumesCollector&& iC)
0015     : _currdir(nullptr),
0016       _histoParameters(iConfig.getUntrackedParameter<edm::ParameterSet>("histoParameters", edm::ParameterSet())),
0017       _rhm(iC) {}
0018 
0019 BeamSpotHistogramMaker::~BeamSpotHistogramMaker() { delete _currdir; }
0020 
0021 void BeamSpotHistogramMaker::book(const std::string dirname) {
0022   edm::Service<TFileService> tfserv;
0023   TFileDirectory* currdir;
0024 
0025   if (!dirname.empty()) {
0026     currdir = new TFileDirectory(tfserv->mkdir(dirname));
0027     _currdir = currdir;
0028   } else {
0029     _currdir = &(tfserv->tFileDirectory());
0030   }
0031 
0032   edm::LogInfo("HistogramBooking") << "BeamSpot histogram booking in directory " << dirname;
0033 
0034   _hbsxrun = _rhm.makeTH1F("bsxrun",
0035                            "BeamSpot X position",
0036                            _histoParameters.getUntrackedParameter<unsigned int>("nBinX", 200),
0037                            _histoParameters.getUntrackedParameter<double>("xMin", -1.),
0038                            _histoParameters.getUntrackedParameter<double>("xMax", 1.));
0039 
0040   _hbsyrun = _rhm.makeTH1F("bsyrun",
0041                            "BeamSpot Y position",
0042                            _histoParameters.getUntrackedParameter<unsigned int>("nBinY", 200),
0043                            _histoParameters.getUntrackedParameter<double>("yMin", -1.),
0044                            _histoParameters.getUntrackedParameter<double>("yMax", 1.));
0045 
0046   _hbszrun = _rhm.makeTH1F("bszrun",
0047                            "BeamSpot Z position",
0048                            _histoParameters.getUntrackedParameter<unsigned int>("nBinZ", 200),
0049                            _histoParameters.getUntrackedParameter<double>("zMin", -1.),
0050                            _histoParameters.getUntrackedParameter<double>("zMax", 1.));
0051 
0052   _hbssigmaxrun = _rhm.makeTH1F("bssigmaxrun",
0053                                 "BeamSpot sigmaX",
0054                                 _histoParameters.getUntrackedParameter<unsigned int>("nBinSigmaX", 200),
0055                                 _histoParameters.getUntrackedParameter<double>("sigmaXMin", 0.),
0056                                 _histoParameters.getUntrackedParameter<double>("sigmaXMax", 0.025));
0057 
0058   _hbssigmayrun = _rhm.makeTH1F("bssigmayrun",
0059                                 "BeamSpot sigmaY",
0060                                 _histoParameters.getUntrackedParameter<unsigned int>("nBinSigmaY", 200),
0061                                 _histoParameters.getUntrackedParameter<double>("sigmaYMin", 0.),
0062                                 _histoParameters.getUntrackedParameter<double>("sigmaYMax", 0.025));
0063 
0064   _hbssigmazrun = _rhm.makeTH1F("bssigmazrun",
0065                                 "BeamSpot sigmaZ",
0066                                 _histoParameters.getUntrackedParameter<unsigned int>("nBinSigmaZ", 200),
0067                                 _histoParameters.getUntrackedParameter<double>("sigmaZMin", 0.),
0068                                 _histoParameters.getUntrackedParameter<double>("sigmaZMax", 15.));
0069 
0070   _hbsxvsorbrun =
0071       _rhm.makeTProfile("bsxvsorbrun", "BeamSpot X position vs orbit number", 1600, 0.5, 1600. * 16384 + 0.5);
0072   _hbsyvsorbrun =
0073       _rhm.makeTProfile("bsyvsorbrun", "BeamSpot Y position vs orbit number", 1600, 0.5, 1600. * 16384 + 0.5);
0074   _hbszvsorbrun =
0075       _rhm.makeTProfile("bszvsorbrun", "BeamSpot Z position vs orbit number", 1600, 0.5, 1600. * 16384 + 0.5);
0076   _hbssigmaxvsorbrun =
0077       _rhm.makeTProfile("bssigmaxvsorbrun", "BeamSpot sigmaX vs orbit number", 1600, 0.5, 1600. * 16384 + 0.5);
0078   _hbssigmayvsorbrun =
0079       _rhm.makeTProfile("bssigmayvsorbrun", "BeamSpot sigmaY vs orbit number", 1600, 0.5, 1600. * 16384 + 0.5);
0080   _hbssigmazvsorbrun =
0081       _rhm.makeTProfile("bssigmazvsorbrun", "BeamSpot sigmaZ vs orbit number", 1600, 0.5, 1600. * 16384 + 0.5);
0082 }
0083 
0084 void BeamSpotHistogramMaker::beginRun(const unsigned int nrun) {
0085   char runname[100];
0086   sprintf(runname, "run_%d", nrun);
0087 
0088   TFileDirectory* currdir = _currdir;
0089   if (currdir == nullptr) {
0090     edm::Service<TFileService> tfserv;
0091     currdir = &(tfserv->tFileDirectory());
0092   }
0093 
0094   _rhm.beginRun(nrun, *currdir);
0095 
0096   (*_hbsxrun)->GetXaxis()->SetTitle("X [cm]");
0097   (*_hbsxrun)->GetYaxis()->SetTitle("Events");
0098   (*_hbsyrun)->GetXaxis()->SetTitle("Y [cm]");
0099   (*_hbsyrun)->GetYaxis()->SetTitle("Events");
0100   (*_hbszrun)->GetXaxis()->SetTitle("Z [cm]");
0101   (*_hbszrun)->GetYaxis()->SetTitle("Events");
0102   (*_hbssigmaxrun)->GetXaxis()->SetTitle("sigmaX [cm]");
0103   (*_hbssigmaxrun)->GetYaxis()->SetTitle("Events");
0104   (*_hbssigmayrun)->GetXaxis()->SetTitle("sigmaY [cm]");
0105   (*_hbssigmayrun)->GetYaxis()->SetTitle("Events");
0106   (*_hbssigmazrun)->GetXaxis()->SetTitle("sigmaZ [cm]");
0107   (*_hbssigmazrun)->GetYaxis()->SetTitle("Events");
0108 
0109   (*_hbsxvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
0110   (*_hbsxvsorbrun)->GetYaxis()->SetTitle("X [cm]");
0111   (*_hbsxvsorbrun)->SetCanExtend(TH1::kAllAxes);
0112   (*_hbsyvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
0113   (*_hbsyvsorbrun)->GetYaxis()->SetTitle("Y [cm]");
0114   (*_hbsyvsorbrun)->SetCanExtend(TH1::kAllAxes);
0115   (*_hbszvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
0116   (*_hbszvsorbrun)->GetYaxis()->SetTitle("Z [cm]");
0117   (*_hbszvsorbrun)->SetCanExtend(TH1::kAllAxes);
0118   (*_hbssigmaxvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
0119   (*_hbssigmaxvsorbrun)->GetYaxis()->SetTitle("sigmaX [cm]");
0120   (*_hbssigmaxvsorbrun)->SetCanExtend(TH1::kAllAxes);
0121   (*_hbssigmayvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
0122   (*_hbssigmayvsorbrun)->GetYaxis()->SetTitle("sigmaY [cm]");
0123   (*_hbssigmayvsorbrun)->SetCanExtend(TH1::kAllAxes);
0124   (*_hbssigmazvsorbrun)->GetXaxis()->SetTitle("time [orbit#]");
0125   (*_hbssigmazvsorbrun)->GetYaxis()->SetTitle("sigmaZ [cm]");
0126   (*_hbssigmazvsorbrun)->SetCanExtend(TH1::kAllAxes);
0127 }
0128 
0129 void BeamSpotHistogramMaker::fill(const unsigned int orbit, const reco::BeamSpot& bs) {
0130   if (_hbsxrun && *_hbsxrun)
0131     (*_hbsxrun)->Fill(bs.x0());
0132   if (_hbsxvsorbrun && *_hbsxvsorbrun)
0133     (*_hbsxvsorbrun)->Fill(orbit, bs.x0());
0134 
0135   if (_hbsyrun && *_hbsyrun)
0136     (*_hbsyrun)->Fill(bs.y0());
0137   if (_hbsyvsorbrun && *_hbsyvsorbrun)
0138     (*_hbsyvsorbrun)->Fill(orbit, bs.y0());
0139 
0140   if (_hbszrun && *_hbszrun)
0141     (*_hbszrun)->Fill(bs.z0());
0142   if (_hbszvsorbrun && *_hbszvsorbrun)
0143     (*_hbszvsorbrun)->Fill(orbit, bs.z0());
0144 
0145   if (_hbssigmaxrun && *_hbssigmaxrun)
0146     (*_hbssigmaxrun)->Fill(bs.BeamWidthX());
0147   if (_hbssigmayrun && *_hbssigmayrun)
0148     (*_hbssigmayrun)->Fill(bs.BeamWidthY());
0149   if (_hbssigmazrun && *_hbssigmazrun)
0150     (*_hbssigmazrun)->Fill(bs.sigmaZ());
0151   if (_hbssigmaxvsorbrun && *_hbssigmaxvsorbrun)
0152     (*_hbssigmaxvsorbrun)->Fill(orbit, bs.BeamWidthX());
0153   if (_hbssigmayvsorbrun && *_hbssigmayvsorbrun)
0154     (*_hbssigmayvsorbrun)->Fill(orbit, bs.BeamWidthY());
0155   if (_hbssigmazvsorbrun && *_hbssigmazvsorbrun)
0156     (*_hbssigmazvsorbrun)->Fill(orbit, bs.sigmaZ());
0157 }