#include <memory>
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/HcalCalibObjects/interface/HOCalibVariables.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "CommonTools/UtilAlgos/interface/TFileService.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "TMath.h"
#include "TFile.h"
#include "TH1F.h"
#include "TH2F.h"
#include "TTree.h"
#include "TProfile.h"
#include "TPostScript.h"
#include "TCanvas.h"
#include "TF1.h"
#include "TStyle.h"
#include "TMinuit.h"
#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
Go to the source code of this file.
Classes | |
class | HOCalibAnalyzer |
Functions | |
void | fcnbg (Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t flag) |
void | fcnsg (Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t flag) |
Double_t | gausX (Double_t *x, Double_t *par) |
Double_t | langaufun (Double_t *x, Double_t *par) |
static const edm::ParameterSetDescriptionFillerPluginFactory::PMaker< edm::ParameterSetDescriptionFiller< HOCalibAnalyzer > > | s_filler__LINE__ ("HOCalibAnalyzer") |
static const edm::MakerPluginFactory::PMaker< edm::WorkerMaker< HOCalibAnalyzer > > | s_maker__LINE__ ("HOCalibAnalyzer") |
void | set_mean (double &x, bool mdigi) |
void | set_sigma (double &x, bool mdigi) |
Double_t | totalfunc (Double_t *x, Double_t *par) |
Variables | |
std::vector< float > | cro_ssg [netamx][nphimx+1] |
static const int | etamap [4][21] |
int | ietafit |
int | iphifit |
static const int | mapx0m [9][2] = {{17, 19}, {14, 15}, {13, 16}, {9, 11}, {0, 0}, {8, 12}, {5, 6}, {4, 7}, {1, 3}} |
static const int | mapx0p [9][2] = {{3, 1}, {7, 4}, {6, 5}, {12, 8}, {0, 0}, {11, 9}, {16, 13}, {15, 14}, {19, 17}} |
static const int | mapx1 [6][3] = {{1, 4, 8}, {12, 7, 3}, {5, 9, 13}, {11, 6, 2}, {16, 15, 14}, {19, 18, 17}} |
static const int | mapx2 [6][3] = {{1, 4, 8}, {12, 7, 3}, {5, 9, 13}, {11, 6, 2}, {16, 15, 14}, {-1, -1, -1}} |
static const int | nbgpr = 3 |
static const int | netamx = 30 |
static const int | nphimx = 72 |
static const int | npixlebt [21] = {0, 0, 0, 0, 0, 1, 2, 3, 0, 4, 0, 6, 7, 8, 9, 0, 11, 13, 14, 15, 0} |
static const int | npixleft [21] = {0, 0, 1, 2, 0, 4, 5, 6, 0, 8, 0, 0, 11, 0, 13, 14, 15, 0, 17, 18, 0} |
static const int | npixleup [21] = {0, 4, 5, 6, 8, 9, 0, 11, 0, 13, 0, 15, 16, 0, 17, 18, 19, 0, 0, 0, 0} |
static const int | npixribt [21] = {0, 0, 0, 0, 1, 2, 3, 0, 4, 5, 0, 7, 0, 9, 0, 11, 12, 14, 15, 16, 0} |
static const int | npixrigh [21] = {0, 2, 3, 0, 5, 6, 7, 0, 9, 0, 0, 12, 0, 14, 15, 16, 0, 18, 19, 0, 0} |
static const int | npixriup [21] = {0, 5, 6, 7, 9, 0, 11, 12, 13, 14, 0, 16, 0, 17, 18, 19, 0, 0, 0, 0, 0} |
static const int | nsgpr = 7 |
static const int | phimap [4][21] |
std::vector< float > | sig_reg [netamx][nphimx+1] |
void fcnbg | ( | Int_t & | npar, |
Double_t * | gin, | ||
Double_t & | f, | ||
Double_t * | par, | ||
Int_t | flag | ||
) |
Definition at line 166 of file HOCalibAnalyzer.cc.
References cro_ssg, MillePedeFileConverter_cfg::e, ietafit, iphifit, cmsBatch::log, and SiStripPI::max.
Referenced by HOCalibAnalyzer::endJob().
void fcnsg | ( | Int_t & | npar, |
Double_t * | gin, | ||
Double_t & | f, | ||
Double_t * | par, | ||
Int_t | flag | ||
) |
Definition at line 176 of file HOCalibAnalyzer.cc.
References ietafit, iphifit, cmsBatch::log, sig_reg, and totalfunc().
Referenced by HOCalibAnalyzer::endJob().
Double_t gausX | ( | Double_t * | x, |
Double_t * | par | ||
) |
Definition at line 112 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::endJob(), and totalfunc().
Double_t langaufun | ( | Double_t * | x, |
Double_t * | par | ||
) |
Definition at line 114 of file HOCalibAnalyzer.cc.
References np, SimDataFormats::CaloAnalysis::sc, and geometryCSVtoXML::xx.
Referenced by HOCalibAnalyzer::endJob(), and totalfunc().
|
static |
|
static |
void set_mean | ( | double & | x, |
bool | mdigi | ||
) |
Definition at line 186 of file HOCalibAnalyzer.cc.
References SiStripPI::max, and min().
Referenced by HOCalibAnalyzer::endJob().
void set_sigma | ( | double & | x, |
bool | mdigi | ||
) |
Definition at line 196 of file HOCalibAnalyzer.cc.
References SiStripPI::max, and min().
Referenced by HOCalibAnalyzer::endJob().
Double_t totalfunc | ( | Double_t * | x, |
Double_t * | par | ||
) |
Definition at line 164 of file HOCalibAnalyzer.cc.
References gausX(), and langaufun().
Referenced by HOCalibAnalyzer::endJob(), and fcnsg().
Definition at line 105 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze(), and fcnbg().
|
static |
Definition at line 79 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
int ietafit |
Definition at line 102 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::endJob(), fcnbg(), and fcnsg().
int iphifit |
Definition at line 103 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::endJob(), fcnbg(), and fcnsg().
|
static |
Definition at line 77 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 76 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 72 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 74 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 99 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::endJob().
|
static |
Definition at line 97 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::endJob(), HOCalibAnalyzer::getHOieta(), HOCalibAnalyzer::HOCalibAnalyzer(), and HOCalibAnalyzer::invert_HOieta().
|
static |
Definition at line 98 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze(), HOCalibAnalyzer::endJob(), and HOCalibAnalyzer::HOCalibAnalyzer().
|
static |
Definition at line 92 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 90 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 94 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 93 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 91 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 95 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
|
static |
Definition at line 100 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::endJob().
|
static |
Definition at line 84 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze().
Definition at line 104 of file HOCalibAnalyzer.cc.
Referenced by HOCalibAnalyzer::analyze(), HOCalibAnalyzer::endJob(), and fcnsg().