![]() |
![]() |
Public Member Functions | |
LHEAnalyzer (const edm::ParameterSet ¶ms) | |
virtual | ~LHEAnalyzer () |
Protected Member Functions | |
virtual void | analyze (const edm::Event &event, const edm::EventSetup &es) |
Private Member Functions | |
void | fillDJRSched (unsigned int min, unsigned int max) |
Private Attributes | |
unsigned int | binsDelta |
unsigned int | binsPt |
double | defaultDeltaCut |
double | defaultPtCut |
boost::ptr_vector< JetClustering > | deltaClustering |
std::vector< unsigned int > | djrSched |
std::vector< TH1 * > | histoDelta |
std::vector< TH1 * > | histoPt |
JetInput | jetInput |
double | maxDelta |
unsigned int | maxDJR |
double | maxEta |
double | maxPt |
double | minDelta |
unsigned int | minDJR |
double | minPt |
std::auto_ptr< JetClustering > | ptClustering |
double | ptFraction |
edm::InputTag | sourceLabel |
bool | useEt |
Definition at line 38 of file LHEAnalyzer.cc.
LHEAnalyzer::LHEAnalyzer | ( | const edm::ParameterSet & | params | ) | [explicit] |
Definition at line 77 of file LHEAnalyzer.cc.
References edm::ParameterSet::addParameter(), binsDelta, binsPt, defaultDeltaCut, defaultPtCut, deltaClustering, fillDJRSched(), edm::ParameterSet::getParameter(), h, histoDelta, histoPt, i, maxDelta, maxDJR, maxPt, minDelta, minDJR, minPt, ptClustering, ptFraction, ss, tmp, and useEt.
00077 : 00078 sourceLabel(params.getParameter<edm::InputTag>("src")), 00079 jetInput(params.getParameter<edm::ParameterSet>("jetInput")), 00080 defaultDeltaCut(params.getParameter<double>("defaultDeltaCut")), 00081 defaultPtCut(params.getParameter<double>("defaultPtCut")), 00082 maxEta(params.getParameter<double>("maxEta")), 00083 useEt(params.getParameter<bool>("useEt")), 00084 ptFraction(params.getUntrackedParameter<double>("ptFraction", 0.75)), 00085 binsDelta(params.getParameter<unsigned int>("binsDelta")), 00086 minDelta(params.getParameter<double>("minDelta")), 00087 maxDelta(params.getParameter<double>("maxDelta")), 00088 binsPt(params.getParameter<unsigned int>("binsPt")), 00089 minPt(params.getParameter<double>("minPt")), 00090 maxPt(params.getParameter<double>("maxPt")), 00091 minDJR(params.getParameter<unsigned int>("minDJR")), 00092 maxDJR(params.getParameter<unsigned int>("maxDJR")) 00093 { 00094 edm::ParameterSet jetClusPSet = 00095 params.getParameter<edm::ParameterSet>("jetClustering"); 00096 00097 for(unsigned int i = 0; i <= binsDelta; i++) { 00098 double deltaCut = 00099 minDelta + (maxDelta - minDelta) * i / binsDelta; 00100 jetClusPSet.addParameter("coneRadius", deltaCut); 00101 edm::ParameterSet tmp; 00102 tmp.addParameter("algorithm", jetClusPSet); 00103 deltaClustering.push_back( 00104 new JetClustering(tmp, defaultPtCut * ptFraction)); 00105 } 00106 00107 jetClusPSet.addParameter("coneRadius", defaultDeltaCut); 00108 edm::ParameterSet tmp; 00109 tmp.addParameter("algorithm", jetClusPSet); 00110 ptClustering.reset(new JetClustering(tmp, minPt * ptFraction)); 00111 00112 fillDJRSched(minDJR <= 0 ? 1 : minDJR, maxDJR - 1); 00113 00114 edm::Service<TFileService> fs; 00115 for(unsigned int i = minDJR; i < maxDJR; i++) { 00116 std::ostringstream ss, ss2; 00117 ss << (i + 1) << "#rightarrow" << i << " jets"; 00118 ss2 << i; 00119 TH1 *h = fs->make<TH1D>(("delta" + ss2.str()).c_str(), 00120 ("DJR " + ss.str()).c_str(), 00121 binsDelta, minDelta, maxDelta); 00122 h->SetXTitle("p_{T} [GeV/c^2]"); 00123 h->SetYTitle("#delta#sigma [mb]"); 00124 00125 if (i == 0) { 00126 h->Delete(); 00127 h = 0; 00128 } 00129 histoDelta.push_back(h); 00130 00131 std::string what = useEt ? "E" : "p"; 00132 h = fs->make<TH1D>(("pt" + ss2.str()).c_str(), 00133 ("DJR " + ss.str()).c_str(), binsPt, 00134 std::log10(minPt), std::log10(maxPt)); 00135 h->SetXTitle(("log_{10}(" + what + 00136 "_{T} [GeV/c^{2}])").c_str()); 00137 h->SetYTitle("#delta#sigma [mb]"); 00138 00139 histoPt.push_back(h); 00140 } 00141 }
LHEAnalyzer::~LHEAnalyzer | ( | ) | [virtual] |
void LHEAnalyzer::analyze | ( | const edm::Event & | event, | |
const edm::EventSetup & | es | |||
) | [protected, virtual] |
Implements edm::EDAnalyzer.
Definition at line 159 of file LHEAnalyzer.cc.
References binsDelta, defaultPtCut, deltaClustering, djrSched, first, MCTruth2::genEvent, h, histoDelta, histoPt, i, int, iter, j, configurableAnalysis::Jet, jetInput, pfTauBenchmarkGeneric_cfi::jets, maxDJR, minDJR, HLT_VtxMuL3::result, edm::second(), python::multivaluedict::sort(), sourceLabel, and useEt.
00160 { 00161 using boost::bind; 00162 typedef JetClustering::Jet Jet; 00163 00164 edm::Handle<edm::HepMCProduct> hepmc; 00165 event.getByLabel(sourceLabel, hepmc); 00166 00167 std::auto_ptr<HepMC::GenEvent> clonedEvent; 00168 const HepMC::GenEvent *genEvent = hepmc->GetEvent(); 00169 if (!genEvent->signal_process_vertex()) { 00170 clonedEvent.reset(new HepMC::GenEvent(*genEvent)); 00171 const HepMC::GenVertex *signalVertex = 00172 LHEEvent::findSignalVertex(clonedEvent.get()); 00173 clonedEvent->set_signal_process_vertex( 00174 const_cast<HepMC::GenVertex*>(signalVertex)); 00175 genEvent = clonedEvent.get(); 00176 } 00177 00178 JetInput::ParticleVector particles = jetInput(genEvent); 00179 00180 std::vector<Jet> ptJets = (*ptClustering)(particles); 00181 std::sort(ptJets.begin(), ptJets.end(), 00182 bind(std::greater<double>(), 00183 bind(useEt ? &Jet::et : &Jet::pt, _1), 00184 bind(useEt ? &Jet::et : &Jet::pt, _2))); 00185 00186 typedef std::pair<int, int> Pair; 00187 std::vector<Pair> deltaJets(maxDJR - minDJR + 1, 00188 Pair(-1, binsDelta + 1)); 00189 00190 for(std::vector<unsigned int>::const_iterator djr = djrSched.begin(); 00191 djr != djrSched.end(); ++djr) { 00192 //std::cout << "DJR schedule " << (*djr + 1) << " -> " << *djr << std::endl; 00193 int result = -1; 00194 for(;;) { 00195 //for(int i = minDJR; i <= maxDJR; i++) 00196 //std::cout << "+++ " << i << ": (" << deltaJets[i - minDJR].first << ", " << deltaJets[i - minDJR].second << ")" << std::endl; 00197 int upper = binsDelta + 1; 00198 for(int i = *djr; i >= (int)minDJR; i--) { 00199 if (deltaJets[i - minDJR].second <= 00200 (int)binsDelta) { 00201 upper = deltaJets[i - minDJR].second; 00202 break; 00203 } 00204 } 00205 int lower = -1; 00206 for(int i = *djr + 1; i <= (int)maxDJR; i++) { 00207 if (deltaJets[i - minDJR].first >= 0) { 00208 lower = deltaJets[i - minDJR].first; 00209 break; 00210 } 00211 } 00212 //std::cout << "\t" << lower << " < " << upper << std::endl; 00213 00214 result = (lower + upper + 2) / 2 - 1; 00215 if (result == lower) 00216 break; 00217 else if (result < lower) { 00218 result = -1; 00219 break; 00220 } 00221 00222 std::vector<Jet> jets = 00223 deltaClustering[result](particles); 00224 unsigned int nJets = 0; 00225 for(std::vector<Jet>::const_iterator iter = 00226 jets.begin(); iter != jets.end(); ++iter) 00227 if ((useEt ? iter->et() : iter->pt()) 00228 > defaultPtCut) 00229 nJets++; 00230 00231 //std::cout << "\t---(" << *djr << ")--> bin " << result << ": " << nJets << " jets" << std::endl; 00232 00233 if (nJets < minDJR) 00234 nJets = minDJR; 00235 else if (nJets > maxDJR) 00236 nJets = maxDJR; 00237 00238 for(int j = nJets; j >= (int)minDJR; j--) { 00239 if (deltaJets[j - minDJR].first < 0 || 00240 result > deltaJets[j - minDJR].first) 00241 deltaJets[j - minDJR].first = result; 00242 } 00243 for(int j = nJets; j <= (int)maxDJR; j++) { 00244 if (deltaJets[j - minDJR].second < 00245 (int)binsDelta || 00246 result < deltaJets[j - minDJR].second) 00247 deltaJets[j - minDJR].second = result; 00248 } 00249 } 00250 00251 //std::cout << "final " << *djr << ": " << result << std::endl; 00252 TH1 *h = histoDelta[*djr - minDJR]; 00253 h->Fill(h->GetBinCenter(result + 1)); 00254 00255 h = histoPt[*djr - minDJR]; 00256 if (*djr >= ptJets.size()) 00257 h->Fill(-999.0); 00258 else if (useEt) 00259 h->Fill(std::log10(ptJets[*djr].et())); 00260 else 00261 h->Fill(std::log10(ptJets[*djr].pt())); 00262 } 00263 00264 if (minDJR <= 0) { 00265 TH1 *h = histoPt[0]; 00266 if (minDJR >= ptJets.size()) 00267 h->Fill(-999.0); 00268 else if (useEt) 00269 h->Fill(std::log10(ptJets[minDJR].et())); 00270 else 00271 h->Fill(std::log10(ptJets[minDJR].pt())); 00272 } 00273 }
Definition at line 147 of file LHEAnalyzer.cc.
References djrSched, and middle.
Referenced by LHEAnalyzer().
00148 { 00149 unsigned int middle = (min + max) / 2; 00150 00151 djrSched.push_back(middle); 00152 00153 if (min < middle) 00154 fillDJRSched(min, middle - 1); 00155 if (middle < max) 00156 fillDJRSched(middle + 1, max); 00157 }
unsigned int LHEAnalyzer::binsDelta [private] |
unsigned int LHEAnalyzer::binsPt [private] |
double LHEAnalyzer::defaultDeltaCut [private] |
double LHEAnalyzer::defaultPtCut [private] |
boost::ptr_vector<JetClustering> LHEAnalyzer::deltaClustering [private] |
std::vector<unsigned int> LHEAnalyzer::djrSched [private] |
std::vector<TH1*> LHEAnalyzer::histoDelta [private] |
std::vector<TH1*> LHEAnalyzer::histoPt [private] |
JetInput LHEAnalyzer::jetInput [private] |
double LHEAnalyzer::maxDelta [private] |
unsigned int LHEAnalyzer::maxDJR [private] |
double LHEAnalyzer::maxEta [private] |
Definition at line 56 of file LHEAnalyzer.cc.
double LHEAnalyzer::maxPt [private] |
double LHEAnalyzer::minDelta [private] |
unsigned int LHEAnalyzer::minDJR [private] |
double LHEAnalyzer::minPt [private] |
std::auto_ptr<JetClustering> LHEAnalyzer::ptClustering [private] |
double LHEAnalyzer::ptFraction [private] |
edm::InputTag LHEAnalyzer::sourceLabel [private] |
bool LHEAnalyzer::useEt [private] |