CMS 3D CMS Logo

HcalCalibrator.cc
Go to the documentation of this file.
1 // -*- C++ -*-F
2 //---------------------------------------------------------------------------------
3 // Package: HcalCalibrator
4 // Class: HcalCalibrator
5 //
19 //
20 // Original Author: "Anton Anastassov"
21 // Created: Tue Sept 24 09:13:48 CDT 2008
22 //
23 //
24 //_________________________________________________________________________________
25 
26 // system include files
27 #include <memory>
28 #include <fstream>
29 #include <iostream>
30 #include <string>
31 
32 // user include files
33 
36 
45 
50 
52 public:
53  explicit HcalCalibrator(const edm::ParameterSet&);
54  ~HcalCalibrator() override = default;
55 
56  // Added for running the CaloTower creation algorithm
57 
58 private:
59  void beginJob() override;
60  void analyze(const edm::Event&, const edm::EventSetup&) override;
61  void endJob() override;
62 
65 
68  const double mMinTargetE;
69  const double mMaxTargetE;
70  const double mMinCellE;
71  const double mMinEOverP;
72  const double mMaxEOverP;
73  const double mMaxTrkEmE;
74 
75  const double mMaxEtThirdJet;
76  const double mMinDPhiDiJets;
77  const bool mSumDepths;
78  const bool mSumSmallDepths;
79  const bool mCombinePhi;
80  const int mHbClusterSize;
81  const int mHeClusterSize;
82 
83  const bool mUseConeClustering;
84  const double mMaxConeDist;
85 
86  const int mCalibAbsIEtaMax;
87  const int mCalibAbsIEtaMin;
88 
89  const double mMaxProbeJetEmFrac;
90  const double mMaxTagJetEmFrac;
91  const double mMaxTagJetAbsEta;
92  const double mMinTagJetEt;
93  const double mMinProbeJetAbsEta;
94 
96  const bool mApplyPhiSymCorFlag;
97 
100 
103 
106 
108 };
109 
110 // constructor
111 
113  : mInputFileList(conf.getUntrackedParameter<std::string>("inputFileList")),
114  // mOutputFile(conf.getUntrackedParameter<std::string>("outputFile")),
115  mCalibType(conf.getUntrackedParameter<std::string>("calibType")),
116  mCalibMethod(conf.getUntrackedParameter<std::string>("calibMethod")),
117  mMinTargetE(conf.getUntrackedParameter<double>("minTargetE")),
118  mMaxTargetE(conf.getUntrackedParameter<double>("maxTargetE")),
119  mMinCellE(conf.getUntrackedParameter<double>("minCellE")),
120  mMinEOverP(conf.getUntrackedParameter<double>("minEOverP")),
121  mMaxEOverP(conf.getUntrackedParameter<double>("maxEOverP")),
122  mMaxTrkEmE(conf.getUntrackedParameter<double>("maxTrkEmE")),
123  mMaxEtThirdJet(conf.getUntrackedParameter<double>("maxEtThirdJet")),
124  mMinDPhiDiJets(conf.getUntrackedParameter<double>("minDPhiDiJets")),
125  mSumDepths(conf.getUntrackedParameter<bool>("sumDepths")),
126  mSumSmallDepths(conf.getUntrackedParameter<bool>("sumSmallDepths")),
127  mCombinePhi(conf.getUntrackedParameter<bool>("combinePhi")),
128  mHbClusterSize(conf.getUntrackedParameter<int>("hbClusterSize")),
129  mHeClusterSize(conf.getUntrackedParameter<int>("heClusterSize")),
130 
131  mUseConeClustering(conf.getUntrackedParameter<bool>("useConeClustering")),
132  mMaxConeDist(conf.getUntrackedParameter<double>("maxConeDist")),
133 
134  mCalibAbsIEtaMax(conf.getUntrackedParameter<int>("calibAbsIEtaMax")),
135  mCalibAbsIEtaMin(conf.getUntrackedParameter<int>("calibAbsIEtaMin")),
136  mMaxProbeJetEmFrac(conf.getUntrackedParameter<double>("maxProbeJetEmFrac")),
137  mMaxTagJetEmFrac(conf.getUntrackedParameter<double>("maxTagJetEmFrac")),
138  mMaxTagJetAbsEta(conf.getUntrackedParameter<double>("maxTagJetAbsEta")),
139  mMinTagJetEt(conf.getUntrackedParameter<double>("minTagJetEt")),
140  mMinProbeJetAbsEta(conf.getUntrackedParameter<double>("minProbeJetAbsEta")),
141  mPhiSymCorFileName(conf.getUntrackedParameter<std::string>("phiSymCorFileName")),
142  mApplyPhiSymCorFlag(conf.getUntrackedParameter<bool>("applyPhiSymCorFlag")),
143  mOutputCorCoefFileName(conf.getUntrackedParameter<std::string>("outputCorCoefFileName")),
144  mHistoFileName(conf.getUntrackedParameter<std::string>("histoFileName")),
147 
148 // ------------ method called to for each event ------------
149 
153 }
154 
155 // ------------ method called once each job just before starting event loop ------------
156 
158 
159 // ------------ method called once each job just after ending the event loop ------------
160 
162  if (mCalibType != "DI_JET" && mCalibType != "ISO_TRACK") {
163  edm::LogVerbatim("HcalCalibrator") << "\n\nUnknown calibration type " << mCalibType;
164  edm::LogVerbatim("HcalCalibrator") << "Please select ISO_TRACK or DI_JET in the python file.";
165  return;
166  }
167 
168  if (mCalibMethod != "L3" && mCalibMethod != "MATRIX_INV_OF_ETA_AVE" && mCalibMethod != "L3_AND_MTRX_INV") {
169  edm::LogVerbatim("HcalCalibrator") << "\n\nUnknown calibration method " << mCalibMethod;
170  edm::LogVerbatim("HcalCalibrator")
171  << "Supported methods for IsoTrack calibration are: L3, MATRIX_INV_OF_ETA_AVE, L3_AND_MTRX_INV";
172  edm::LogVerbatim("HcalCalibrator") << "For DiJets the supported method is L3";
173  return;
174  }
175 
176  if (mCalibType == "DI_JET" && mCalibMethod != "L3") {
177  edm::LogVerbatim("HcalCalibrator")
178  << "\n\nDiJet calibration can use only the L3 method. Please change the python file.";
179  return;
180  }
181 
182  if (mCalibAbsIEtaMin < 1 || mCalibAbsIEtaMax > 41 || mCalibAbsIEtaMin > mCalibAbsIEtaMax) {
183  edm::LogVerbatim("HcalCalibrator")
184  << "\n\nInvalid ABS(iEta) calibration range. Check calibAbsIEtaMin and calibAbsIEtaMax in the python file.";
185  return;
186  }
187 
188  hcalCalib* calibrator = new hcalCalib();
189 
190  // set the parameters controlling the calibratoration
191 
192  calibrator->SetCalibType(mCalibType);
193  calibrator->SetCalibMethod(mCalibMethod);
194  calibrator->SetMinTargetE(mMinTargetE);
195  calibrator->SetMaxTargetE(mMaxTargetE);
196  calibrator->SetMaxEtThirdJet(mMaxEtThirdJet);
197  calibrator->SetMinDPhiDiJets(mMinDPhiDiJets);
198  calibrator->SetSumDepthsFlag(mSumDepths);
200  calibrator->SetCombinePhiFlag(mCombinePhi);
201  calibrator->SetMinCellE(mMinCellE);
202  calibrator->SetMinEOverP(mMinEOverP);
203  calibrator->SetMaxEOverP(mMaxEOverP);
204  calibrator->SetMaxTrkEmE(mMaxTrkEmE);
205  calibrator->SetHbClusterSize(mHbClusterSize);
206  calibrator->SetHeClusterSize(mHeClusterSize);
207 
209  calibrator->SetConeMaxDist(mMaxConeDist);
210 
216 
217  calibrator->SetMinTagJetEt(mMinTagJetEt);
218 
223 
224  calibrator->SetHistoFileName(mHistoFileName);
225 
227 
228  std::ifstream inputFileList; // contains list of input root files
229 
230  TString files = mInputFileList;
231  inputFileList.open(files.Data());
232 
233  std::vector<TString> inputFiles;
234  while (!inputFileList.eof()) {
235  TString fileName;
237  if (!fileName.BeginsWith("#") && !fileName.Contains(" ") && fileName != "")
238  inputFiles.push_back(fileName);
239  }
240  inputFileList.close();
241 
242  edm::LogVerbatim("HcalCalibrator") << "\nInput files for processing:";
243  for (std::vector<TString>::iterator it = inputFiles.begin(); it != inputFiles.end(); ++it) {
244  edm::LogVerbatim("HcalCalibrator") << "file: " << it->Data();
245  }
246 
247  TChain* fChain = new TChain("hcalCalibTree");
248 
249  for (std::vector<TString>::iterator f_it = inputFiles.begin(); f_it != inputFiles.end(); ++f_it) {
250  fChain->Add(f_it->Data());
251  }
252 
253  fChain->Process(calibrator);
254 
255  if (fChain)
256  delete fChain;
257  delete calibrator;
258 
259  return;
260 }
261 
262 //define this as a plug-in
void SetMaxTagJetAbsEta(Float_t e)
Definition: hcalCalib.h:208
void SetPhiSymCorFileName(const TString &filename)
Definition: hcalCalib.h:205
Log< level::Info, true > LogVerbatim
void SetHbClusterSize(Int_t i)
Definition: hcalCalib.h:194
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
void SetOutputCorCoefFileName(const TString &filename)
Definition: hcalCalib.h:211
const int mCalibAbsIEtaMin
void SetMinTagJetEt(Float_t e)
Definition: hcalCalib.h:209
const std::string mInputFileList
void beginJob() override
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const double mMaxEtThirdJet
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tok_geom_
void SetMaxEOverP(Float_t e)
Definition: hcalCalib.h:190
const std::string mOutputCorCoefFileName
void SetMaxProbeJetEmFrac(Float_t f)
Definition: hcalCalib.h:206
void SetSumSmallDepthsFlag(Bool_t b)
Definition: hcalCalib.h:186
const std::string mOutputFile
void SetConeMaxDist(Float_t d)
Definition: hcalCalib.h:198
const std::string mCalibMethod
void SetApplyPhiSymCorFlag(Bool_t b)
Definition: hcalCalib.h:204
const CaloGeometry * mTheCaloGeometry
const std::string mPhiSymCorFileName
void SetCalibType(const TString &s)
Definition: hcalCalib.h:192
void SetMinDPhiDiJets(Float_t dphi)
Definition: hcalCalib.h:203
const bool mApplyPhiSymCorFlag
void SetSumDepthsFlag(Bool_t b)
Definition: hcalCalib.h:185
void analyze(const edm::Event &, const edm::EventSetup &) override
void SetUseConeClustering(Bool_t b)
Definition: hcalCalib.h:197
const double mMaxProbeJetEmFrac
const HcalTopology * mTheHcalTopology
const bool mSumSmallDepths
const double mMinProbeJetAbsEta
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > tok_htopo_
~HcalCalibrator() override=default
const std::string mHistoFileName
void SetMinTargetE(Float_t e)
Definition: hcalCalib.h:183
const bool mSumDepths
void SetCalibMethod(const TString &s)
Definition: hcalCalib.h:193
const double mMinDPhiDiJets
void SetCalibAbsIEtaMax(Int_t i)
Definition: hcalCalib.h:200
const double mMinTargetE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const double mMinTagJetEt
void SetMaxTargetE(Float_t e)
Definition: hcalCalib.h:184
void SetMaxTrkEmE(Float_t e)
Definition: hcalCalib.h:191
const double mMaxTargetE
HcalCalibrator(const edm::ParameterSet &)
const double mMaxTagJetEmFrac
const std::string mCalibType
void SetCaloGeometry(const CaloGeometry *g, const HcalTopology *topo)
Definition: hcalCalib.h:214
const double mMaxTagJetAbsEta
void SetMaxTagJetEmFrac(Float_t f)
Definition: hcalCalib.h:207
void SetCalibAbsIEtaMin(Int_t i)
Definition: hcalCalib.h:201
void SetCombinePhiFlag(Bool_t b)
Definition: hcalCalib.h:187
const double mMaxEOverP
void SetHistoFileName(const TString &filename)
Definition: hcalCalib.h:212
const double mMinEOverP
const double mMinCellE
const bool mCombinePhi
void SetHeClusterSize(Int_t i)
Definition: hcalCalib.h:195
const bool mUseConeClustering
void SetMinCellE(Float_t e)
Definition: hcalCalib.h:188
void SetMaxEtThirdJet(Float_t et)
Definition: hcalCalib.h:202
void endJob() override
void SetMinProbeJetAbsEta(Float_t e)
Definition: hcalCalib.h:210
const double mMaxTrkEmE
void SetMinEOverP(Float_t e)
Definition: hcalCalib.h:189
const double mMaxConeDist
const int mHeClusterSize
const int mCalibAbsIEtaMax
const int mHbClusterSize