Go to the documentation of this file.00001 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00002
00003 #include "DataFormats/Common/interface/Handle.h"
00004
00005 #include "DataFormats/JetReco/interface/CaloJet.h"
00006 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00007
00008 #include "DataFormats/METReco/interface/CaloMET.h"
00009 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00010
00011 #include "FWCore/Framework/interface/MakerMacros.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013
00014 #include "FWCore/Utilities/interface/InputTag.h"
00015
00016 #include "DQM/DataScouting/interface/AlphaTVarProducer.h"
00017
00018 #include "TVector3.h"
00019
00020 #include <memory>
00021 #include <vector>
00022
00023
00024
00025
00026 AlphaTVarProducer::AlphaTVarProducer(const edm::ParameterSet& iConfig) :
00027 inputJetTag_ (iConfig.getParameter<edm::InputTag>("inputJetTag")){
00028
00029 produces<std::vector<double> >();
00030
00031 LogDebug("") << "Inputs: "
00032 << inputJetTag_.encode() << " ";
00033 }
00034
00035 AlphaTVarProducer::~AlphaTVarProducer()
00036 {
00037 }
00038
00039
00040 void
00041 AlphaTVarProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00042 {
00043 using namespace std;
00044 using namespace edm;
00045 using namespace reco;
00046
00047
00048
00049 edm::Handle<reco::CaloJetCollection> calojet_handle;
00050 iEvent.getByLabel(inputJetTag_,calojet_handle);
00051
00052 std::auto_ptr<std::vector<double> > result(new std::vector<double>);
00053
00054 if (calojet_handle.isValid()){
00055 std::vector<TLorentzVector> myJets;
00056 reco::CaloJetCollection::const_iterator jetIt;
00057 for(jetIt = calojet_handle->begin(); jetIt != calojet_handle->end(); jetIt++){
00058 TLorentzVector j; j.SetPtEtaPhiE(jetIt->pt(),jetIt->eta(), jetIt->phi(), jetIt->energy());
00059 myJets.push_back(j);
00060 }
00061
00062 double alphaT = CalcAlphaT(myJets);
00063 double HT = CalcHT(myJets);
00064
00065 result->push_back(alphaT);
00066 result->push_back(HT);
00067 }
00068 iEvent.put(result);
00069 }
00070
00071 double
00072 AlphaTVarProducer::CalcAlphaT(std::vector<TLorentzVector> jets){
00073 std::vector<double> ETs;
00074 TVector3 MHT = CalcMHT(jets);
00075 float HT = CalcHT(jets);
00076
00077 for(unsigned int i = 0; i < jets.size(); i++){
00078 if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) ETs.push_back(jets[i].Et());
00079
00080 }
00081 if(ETs.size() < 2.) return 0.0;
00082 if(ETs.size() > 16.) return 0.0;
00083 float DHT = deltaHt(ETs);
00084
00085 float AlphaT = alphaT(HT,DHT,MHT.Mag());
00086
00087 return AlphaT;
00088
00089 }
00090
00091 double AlphaTVarProducer::deltaHt(const std::vector<double>& ETs) {
00092 if(ETs.size() > 16.) return 9999999;
00093 std::vector<double> diff( 1<<(ETs.size()-1) , 0. );
00094 for(unsigned i=0; i < diff.size(); i++)
00095 for(unsigned j=0; j < ETs.size(); j++)
00096 diff[i] += ETs[j] * ( 1 - 2 * (int(i>>j)&1) ) ;
00097 std::vector<double>::const_iterator it;
00098 double min=9999999;
00099 for(it = diff.begin(); it !=diff.end(); it++) if(*it<min) min = *it;
00100 return min;
00101 }
00102
00103 double AlphaTVarProducer::alphaT(const double HT, const double DHT, const double MHT) {
00104 return 0.5 * ( HT - DHT ) / sqrt( HT*HT - MHT*MHT );
00105 }
00106
00107 double
00108 AlphaTVarProducer::CalcHT(std::vector<TLorentzVector> jets){
00109
00110 double HT = 0;
00111 for(unsigned int i = 0; i < jets.size(); i++){
00112 if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) HT += jets[i].Et();
00113 }
00114
00115 return HT;
00116
00117 }
00118
00119 double AlphaTVarProducer::CalcMHT(std::vector<TLorentzVector> jets){
00120 TVector3 MHT;
00121 for(unsigned int i = 0; i < jets.size(); i++){
00122 if(jets[i].Et() > 50. && fabs(jets[i].Eta()) < 2.5) MHT -= jets[i].Vect();
00123 }
00124 MHT.SetZ(0.0);
00125 return MHT.Mag();
00126 }
00127
00128 DEFINE_FWK_MODULE(AlphaTVarProducer);