ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/extras/JetAlgorithmAnalyzer.cc
Revision: 1.1
Committed: Wed Jun 30 11:50:52 2010 UTC (14 years, 10 months ago) by yilmaz
Content type: text/plain
Branch: MAIN
CVS Tags: HiForest_V02_85, HiForest_V02_84, HiForest_V02_83, HiForest_V02_82, HiForest_V02_81, HiForest_V02_80, HiForest_V02_79, HiForest_V02_78, HiForest_V02_77, HiForest_V02_76, HiForest_V02_75, HiForest_V02_74, HiForest_V02_73, HiForest_V02_72, HiForest_V02_71, HiForest_V02_70, HiForest_V02_69, HiForest_V02_68, HiForest_V02_67, HiForest_V02_66, HiForest_V02_65, HiForest_V02_64, HiForest_V02_63, HiForest_V02_62, HiForest_V02_61, HiForest_V02_60, HiForest_V02_59, HiForest_V02_58, HiForest_V02_57, HiForest_V02_56, HiForest_V02_55, HiForest_V02_54, HiForest_V02_53, HiForest_V02_52, HiForest_V02_51, HiForest_V02_50, HiForest_V02_49, HiForest_V02_48, HiForest_V02_47, HiForest_V02_46, HiForest_V02_45, HiForest_V02_44, HiForest_V02_43, HiForest_V02_42, HiForest_V02_41, HiForest_V02_40, HiForest_V02_39, HiForest_V02_38, HiForest_V02_37, HiForest_V02_36, HiForest_V02_35, HiForest_V02_34, HiForest_V02_33, HiForest_V02_32, HiForest_V02_31, HiForest_V02_30, HiForest_V02_27, HiForest_V02_26, QM_2012, HiForest_V02_25, HiForest_V02_24, HiForest_V02_23, HiForest_V02_22, HiForest_V02_21, HiForest_V02_20, HiForest_V02_19, HiForest_V02_18, HiForest_V02_17, HiForest_V02_16, HiForest_V02_15, HiForest_V02_14, HiForest_V02_13, HiForest_V02_12, HiForest_V02_11, HiForest_V02_10, HiForest_V02_09, HiForest_V02_08, HiForest_V02_07, HiForest_V02_06, HiForest_V02_05, HiForest_V02_04, HiForest_V02_03, HiForest_V02_02, HiForest_V02_01, HiForest_V02_00, hi44X_02, hi413_03, hi441_1, hi441_0, hi413_11, hi413_10, hi413_09, hi413_08, hi413_07, hi413_06, hi413_05, hi413_04, hi413_02, hi39X_01, tag_d20110915, cmssw39x_base, cmssw39X_base, HEAD
Branch point for: branch_44x, cmssw39x_branch
Log Message:
backup code

File Contents

# User Rev Content
1 yilmaz 1.1 #ifndef RecoJets_JetProducers_plugins_JetAlgorithmAnalyzer_h
2     #define RecoJets_JetProducers_plugins_JetAlgorithmAnalyzer_h
3    
4     #include "FWCore/Framework/interface/EDProducer.h"
5     #include "FWCore/ParameterSet/interface/ParameterSet.h"
6     #include "DataFormats/Candidate/interface/Candidate.h"
7     #include "DataFormats/Candidate/interface/CandidateFwd.h"
8     #include "DataFormats/JetReco/interface/Jet.h"
9     #include "DataFormats/JetReco/interface/CaloJet.h"
10     #include "DataFormats/JetReco/interface/PFJet.h"
11     #include "DataFormats/JetReco/interface/BasicJet.h"
12     #include "DataFormats/JetReco/interface/GenJet.h"
13     #include "FWCore/Framework/interface/Event.h"
14     #include "FWCore/Framework/interface/EventSetup.h"
15     #include "FWCore/Framework/interface/ESHandle.h"
16     #include "FWCore/Utilities/interface/CodedException.h"
17     #include "FWCore/MessageLogger/interface/MessageLogger.h"
18     #include "FWCore/Framework/interface/MakerMacros.h"
19    
20     #include "RecoJets/JetProducers/interface/PileUpSubtractor.h"
21    
22     #include "fastjet/JetDefinition.hh"
23     #include "fastjet/ClusterSequence.hh"
24     #include "fastjet/ClusterSequenceArea.hh"
25     #include "fastjet/PseudoJet.hh"
26     #include "fastjet/ActiveAreaSpec.hh"
27     #include "fastjet/SISConePlugin.hh"
28     #include "fastjet/CMSIterativeConePlugin.hh"
29     #include "fastjet/ATLASConePlugin.hh"
30     #include "fastjet/CDFMidPointPlugin.hh"
31    
32     #include <vector>
33     #include <iostream>
34     #include <memory>
35     #include <algorithm>
36     #include <limits>
37     #include <cmath>
38     #include <boost/shared_ptr.hpp>
39    
40    
41     #include "DataFormats/Common/interface/View.h"
42     #include "DataFormats/Common/interface/Handle.h"
43     #include "DataFormats/VertexReco/interface/Vertex.h"
44     #include "DataFormats/VertexReco/interface/VertexFwd.h"
45     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
46     #include "DataFormats/JetReco/interface/GenJetCollection.h"
47     #include "DataFormats/JetReco/interface/PFJetCollection.h"
48     #include "DataFormats/JetReco/interface/BasicJetCollection.h"
49     #include "DataFormats/JetReco/interface/TrackJetCollection.h"
50     #include "DataFormats/Candidate/interface/CandidateFwd.h"
51     #include "DataFormats/Candidate/interface/LeafCandidate.h"
52     #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
53     #include "Geometry/Records/interface/CaloGeometryRecord.h"
54    
55     #include "DataFormats/Math/interface/Vector3D.h"
56    
57     #include "FWCore/ServiceRegistry/interface/Service.h"
58     #include "CommonTools/UtilAlgos/interface/TFileService.h"
59     #include "DataFormats/HeavyIonEvent/interface/Centrality.h"
60     #include "TNtuple.h"
61     #include "TH2D.h"
62     #include <vector>
63    
64     static const int nSteps = 8;
65     static const double PI = 3.141592653589;
66    
67     class JetAlgorithmAnalyzer : public edm::EDProducer{
68    
69     public:
70     //
71     // construction/destruction
72     //
73     explicit JetAlgorithmAnalyzer(const edm::ParameterSet& iConfig);
74     virtual ~JetAlgorithmAnalyzer();
75    
76     protected:
77    
78     //
79     // member functions
80     //
81    
82     virtual void produce( edm::Event& iEvent, const edm::EventSetup& iSetup );
83     virtual void runAlgorithm( edm::Event& iEvent, const edm::EventSetup& iSetup );
84     virtual void output( edm::Event & iEvent, edm::EventSetup const& iSetup );
85     template< typename T >
86     void writeBkgJets( edm::Event & iEvent, edm::EventSetup const& iSetup );
87    
88     void fillNtuple(int output, const std::vector<fastjet::PseudoJet>& jets, int step);
89     void fillTowerNtuple(const std::vector<fastjet::PseudoJet>& jets, int step);
90     void fillJetNtuple(const std::vector<fastjet::PseudoJet>& jets, int step);
91     void fillBkgNtuple(const PileUpSubtractor* subtractor, int step);
92    
93     // VirtualJetProducer Stuff
94     //
95     // typedefs & structs
96     //
97     struct JetType {
98     enum Type {
99     BasicJet,
100     GenJet,
101     CaloJet,
102     PFJet,
103     TrackJet,
104     LastJetType // no real type, technical
105     };
106     static const char *names[];
107     static Type byName(const std::string &name);
108     };
109    
110     JetType::Type jetTypeE;
111    
112     inline bool makeCaloJet(const JetType::Type &fTag) {
113     return fTag == JetType::CaloJet;
114     }
115     inline bool makePFJet(const JetType::Type &fTag) {
116     return fTag == JetType::PFJet;
117     }
118     inline bool makeGenJet(const JetType::Type &fTag) {
119     return fTag == JetType::GenJet;
120     }
121     inline bool makeTrackJet(const JetType::Type &fTag) {
122     return fTag == JetType::TrackJet;
123     }
124     inline bool makeBasicJet(const JetType::Type &fTag) {
125     return fTag == JetType::BasicJet;
126     }
127    
128     public:
129     // typedefs
130     typedef boost::shared_ptr<fastjet::ClusterSequence> ClusterSequencePtr;
131     typedef boost::shared_ptr<fastjet::JetDefinition::Plugin> PluginPtr;
132     typedef boost::shared_ptr<fastjet::JetDefinition> JetDefPtr;
133     typedef boost::shared_ptr<fastjet::ActiveAreaSpec> ActiveAreaSpecPtr;
134     typedef boost::shared_ptr<fastjet::RangeDefinition> RangeDefPtr;
135     std::string jetType() const { return jetType_; }
136    
137     private:
138    
139     // trackjet clustering parameters
140     bool useOnlyVertexTracks_;
141     bool useOnlyOnePV_;
142     float dzTrVtxMax_;
143    
144     double phi0_;
145     int nFill_;
146     float etaMax_;
147     int iev_;
148     bool avoidNegative_;
149     bool doAnalysis_;
150     bool doMC_;
151    
152     int centBin_;
153     const CentralityBins* cbins_;
154     TNtuple* ntTowers;
155     TNtuple* ntJetTowers;
156     TNtuple* ntTowersFromEvtPlane;
157     TNtuple* ntJets;
158     TNtuple* ntPU;
159     TNtuple* ntRandom;
160     TNtuple* ntuple;
161     std::vector<TH2D*> hTowers;
162     std::vector<TH2D*> hJetTowers;
163     std::vector<TH2D*> hTowersFromEvtPlane;
164     std::vector<TH2D*> hJets;
165     std::vector<TH1D*> hPU;
166     std::vector<TH1D*> hMean;
167     std::vector<TH1D*> hSigma;
168     std::vector<TH2D*> hRandom;
169     const CaloGeometry *geo;
170     edm::Service<TFileService> f;
171    
172     // VirtualJetProducer Stuff
173    
174     protected:
175     virtual void inputTowers();
176     virtual bool isAnomalousTower(reco::CandidatePtr input);
177     virtual void copyConstituents(const std::vector<fastjet::PseudoJet>&fjConstituents,
178     reco::Jet* jet);
179    
180     // This will run the actual algorithm. This method is pure virtual and
181     // has no default.
182     template< typename T >
183     void writeJets( edm::Event & iEvent, edm::EventSetup const& iSetup );
184    
185     // This method copies the constituents from the fjConstituents method
186     // to an output of CandidatePtr's.
187     virtual std::vector<reco::CandidatePtr>
188     getConstituents(const std::vector<fastjet::PseudoJet>&fjConstituents);
189    
190    
191     protected:
192     std::string moduleLabel_; // label for this module
193     edm::InputTag src_; // input constituent source
194     edm::InputTag srcPVs_; // primary vertex source
195     std::string jetType_; // type of jet (Calo,PF,Basic,Gen)
196     std::string jetAlgorithm_; // the jet algorithm to use
197     double rParam_; // the R parameter to use
198     double inputEtMin_; // minimum et of input constituents
199     double inputEMin_; // minimum e of input constituents
200     double jetPtMin_; // minimum jet pt
201     bool doPVCorrection_; // correct to primary vertex?
202    
203     // for restricting inputs due to processing time
204     bool restrictInputs_; // restrict inputs to first "maxInputs" inputs.
205     unsigned int maxInputs_; // maximum number of inputs.
206    
207     // for fastjet jet area calculation
208     bool doAreaFastjet_; // calculate area w/ fastjet?
209     // for fastjet rho calculation
210     bool doRhoFastjet_; // calculate rho w/ fastjet?
211    
212     // for pileup offset correction
213     bool doPUOffsetCorr_; // add the pileup calculation from offset correction?
214     std::string puSubtractorName_;
215     // anomalous cell cuts
216     unsigned int maxBadEcalCells_; // maximum number of bad ECAL cells
217     unsigned int maxRecoveredEcalCells_; // maximum number of recovered ECAL cells
218     unsigned int maxProblematicEcalCells_; // maximum number of problematic ECAL cells
219     unsigned int maxBadHcalCells_; // maximum number of bad HCAL cells
220     unsigned int maxRecoveredHcalCells_; // maximum number of recovered HCAL cells
221     unsigned int maxProblematicHcalCells_; // maximum number of problematic HCAL cells
222    
223     std::vector<edm::Ptr<reco::Candidate> > inputs_; // input candidates [View, PtrVector and CandColle
224     reco::Particle::Point vertex_; // Primary vertex
225     ClusterSequencePtr fjClusterSeq_; // fastjet cluster sequence
226     JetDefPtr fjJetDefinition_; // fastjet jet definition
227     PluginPtr fjPlugin_; // fastjet plugin
228     ActiveAreaSpecPtr fjActiveArea_; // fastjet active area definition
229     RangeDefPtr fjRangeDef_; // range definition
230     std::vector<fastjet::PseudoJet> fjInputs_; // fastjet inputs
231     std::vector<fastjet::PseudoJet> fjJets_; // fastjet jets
232    
233     std::string jetCollInstanceName_; // instance name for output jet collection
234     boost::shared_ptr<PileUpSubtractor> subtractor_;
235    
236     };
237    
238    
239     #endif
240     ////////////////////////////////////////////////////////////////////////////////
241     //
242     // JetAlgorithmAnalyzer
243     // ------------------
244     //
245     // 04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
246     ////////////////////////////////////////////////////////////////////////////////
247    
248     //#include "RecoJets/JetProducers/plugins/JetAlgorithmAnalyzer.h"
249     //#include "JetAlgorithmAnalyzer.h"
250    
251     #include "RecoJets/JetProducers/interface/JetSpecific.h"
252    
253     #include "FWCore/Framework/interface/Event.h"
254     #include "FWCore/Framework/interface/EventSetup.h"
255     #include "FWCore/Framework/interface/ESHandle.h"
256     #include "FWCore/Utilities/interface/CodedException.h"
257     #include "FWCore/MessageLogger/interface/MessageLogger.h"
258     #include "FWCore/Framework/interface/MakerMacros.h"
259     #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
260     #include "FWCore/ServiceRegistry/interface/Service.h"
261    
262     #include "DataFormats/Common/interface/Handle.h"
263     #include "DataFormats/VertexReco/interface/Vertex.h"
264     #include "DataFormats/VertexReco/interface/VertexFwd.h"
265     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
266     #include "DataFormats/JetReco/interface/GenJetCollection.h"
267     #include "DataFormats/JetReco/interface/PFJetCollection.h"
268     #include "DataFormats/JetReco/interface/BasicJetCollection.h"
269     #include "DataFormats/Candidate/interface/CandidateFwd.h"
270     #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
271     #include "DataFormats/Candidate/interface/LeafCandidate.h"
272     #include "DataFormats/Math/interface/deltaR.h"
273    
274     #include "SimDataFormats/HiGenData/interface/GenHIEvent.h"
275    
276     #include "DataFormats/CaloTowers/interface/CaloTower.h"
277     #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
278     #include "Geometry/Records/interface/CaloGeometryRecord.h"
279    
280     #include "fastjet/SISConePlugin.hh"
281     #include "fastjet/CMSIterativeConePlugin.hh"
282     #include "fastjet/ATLASConePlugin.hh"
283     #include "fastjet/CDFMidPointPlugin.hh"
284    
285     #include "CLHEP/Random/RandomEngine.h"
286    
287     #include <iostream>
288     #include <memory>
289     #include <algorithm>
290     #include <limits>
291     #include <cmath>
292    
293     using namespace std;
294     using namespace edm;
295    
296     static const double pi = 3.14159265358979;
297    
298     CLHEP::HepRandomEngine* randomEngine;
299    
300    
301     ////////////////////////////////////////////////////////////////////////////////
302     // construction / destruction
303     ////////////////////////////////////////////////////////////////////////////////
304    
305     //______________________________________________________________________________
306     JetAlgorithmAnalyzer::JetAlgorithmAnalyzer(const edm::ParameterSet& iConfig)
307     : phi0_(0),
308     nFill_(5),
309     etaMax_(3),
310     iev_(0),
311     geo(0),
312     cbins_(0),
313     moduleLabel_ (iConfig.getParameter<string> ("@module_label"))
314     , src_ (iConfig.getParameter<edm::InputTag>("src"))
315     , srcPVs_ (iConfig.getParameter<edm::InputTag>("srcPVs"))
316     , jetType_ (iConfig.getParameter<string> ("jetType"))
317     , jetAlgorithm_ (iConfig.getParameter<string> ("jetAlgorithm"))
318     , rParam_ (iConfig.getParameter<double> ("rParam"))
319     , inputEtMin_ (iConfig.getParameter<double> ("inputEtMin"))
320     , inputEMin_ (iConfig.getParameter<double> ("inputEMin"))
321     , jetPtMin_ (iConfig.getParameter<double> ("jetPtMin"))
322     , doPVCorrection_(iConfig.getParameter<bool> ("doPVCorrection"))
323     , restrictInputs_(false)
324     , maxInputs_(99999999)
325     , doAreaFastjet_ (iConfig.getParameter<bool> ("doAreaFastjet"))
326     , doRhoFastjet_ (iConfig.getParameter<bool> ("doRhoFastjet"))
327     , doPUOffsetCorr_(iConfig.getParameter<bool> ("doPUOffsetCorr"))
328     , maxBadEcalCells_ (iConfig.getParameter<unsigned>("maxBadEcalCells"))
329     , maxRecoveredEcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredEcalCells"))
330     , maxProblematicEcalCells_(iConfig.getParameter<unsigned>("maxProblematicEcalCells"))
331     , maxBadHcalCells_ (iConfig.getParameter<unsigned>("maxBadHcalCells"))
332     , maxRecoveredHcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredHcalCells"))
333     , maxProblematicHcalCells_(iConfig.getParameter<unsigned>("maxProblematicHcalCells"))
334     , jetCollInstanceName_ ("")
335     {
336     edm::Service<RandomNumberGenerator> rng;
337     randomEngine = &(rng->getEngine());
338    
339     doAnalysis_ = iConfig.getUntrackedParameter<bool>("doAnalysis",true);
340     doMC_ = iConfig.getUntrackedParameter<bool>("doMC",true);
341    
342     if(doAnalysis_) centBin_ = iConfig.getUntrackedParameter<int>("centrality",0);
343    
344     avoidNegative_ = iConfig.getParameter<bool>("avoidNegative");
345    
346     if ( iConfig.exists("UseOnlyVertexTracks") )
347     useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
348     else
349     useOnlyVertexTracks_ = false;
350    
351     if ( iConfig.exists("UseOnlyOnePV") )
352     useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
353     else
354     useOnlyOnePV_ = false;
355    
356     if ( iConfig.exists("DrTrVtxMax") )
357     dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
358     else
359     dzTrVtxMax_ = false;
360    
361     if(doAnalysis_){
362     ntTowers = f->make<TNtuple>("ntTowers","Algorithm Analysis Towers","eta:phi:et:step:event");
363     ntTowersFromEvtPlane = f->make<TNtuple>("ntTowersFromEvtPlane","Algorithm Analysis Towers","eta:phi:et:step:event");
364     ntJetTowers = f->make<TNtuple>("ntTowers","Algorithm Analysis Towers","eta:phi:et:step:event");
365    
366     ntJets = f->make<TNtuple>("ntJets","Algorithm Analysis Jets","eta:phi:et:step:event");
367     ntPU = f->make<TNtuple>("ntPU","Algorithm Analysis Background","eta:mean:sigma:step:event");
368     ntRandom = f->make<TNtuple>("ntRandom","Algorithm Analysis Background","eta:phi:et:pu:event");
369    
370     ntuple = f->make<TNtuple>("nt","debug","ieta:eta:iphi:phi:pt:em:had");
371    
372     for(int i = 0; i < nSteps; ++i){
373     hTowers.push_back(f->make<TH2D>(Form("hTowers_step%d",i),"",200,-5.5,5.5,200,-0.02,6.3));
374     hJets.push_back(f->make<TH2D>(Form("hJets_step%d",i),"",200,-5.5,5.5,200,-0.02,6.3));
375     hTowersFromEvtPlane.push_back(f->make<TH2D>(Form("hTowersFromEvtPlane_step%d",i),"",200,-5.5,5.5,200,-0.02,6.3));
376    
377     hJetTowers.push_back(f->make<TH2D>(Form("hJetTowers_step%d",i),"",200,-5.5,5.5,200,-0.02,6.3));
378    
379     hPU.push_back(f->make<TH1D>(Form("hPU_step%d",i),"",200,-5.5,5.5));
380     hMean.push_back(f->make<TH1D>(Form("hMean_step%d",i),"",200,-5.5,5.5));
381     hSigma.push_back(f->make<TH1D>(Form("hSigma_step%d",i),"",200,-5.5,5.5));
382    
383     }
384    
385     }
386    
387    
388     if (jetAlgorithm_=="IterativeCone") {
389     fjPlugin_ = PluginPtr(new fastjet::CMSIterativeConePlugin(rParam_,1.0));
390     fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
391     }
392     else if (jetAlgorithm_=="Kt")
393     fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::kt_algorithm,rParam_));
394     else if (jetAlgorithm_=="AntiKt")
395     fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm,rParam_) );
396     else
397     throw cms::Exception("Invalid jetAlgorithm")
398     <<"Jet algorithm for MyVirtualJetProducer is invalid, Abort!\n";
399    
400     jetTypeE=JetType::byName(jetType_);
401    
402     if ( iConfig.exists("jetCollInstanceName") ) {
403     jetCollInstanceName_ = iConfig.getParameter<string>("jetCollInstanceName");
404     }
405    
406     if ( doPUOffsetCorr_ ) {
407     if ( jetTypeE != JetType::CaloJet ) {
408     throw cms::Exception("InvalidInput") << "Can only offset correct jets of type CaloJet";
409     }
410     // subtractor_ = boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig));
411     puSubtractorName_ = iConfig.getParameter<string> ("subtractorName");
412     if(puSubtractorName_.empty()){
413     edm::LogWarning("VirtualJetProducer") << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
414     subtractor_ = boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig));
415     }else{
416     subtractor_ = boost::shared_ptr<PileUpSubtractor>(PileUpSubtractorFactory::get()->create( puSubtractorName_, iConfig) );
417     }
418     }
419    
420     // do fasjet area / rho calcluation? => accept corresponding parameters
421     if ( doAreaFastjet_ || doRhoFastjet_ ) {
422     // default Ghost_EtaMax should be 5
423     double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
424     // default Active_Area_Repeats 1
425     int activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
426     // default GhostArea 0.01
427     double ghostArea = iConfig.getParameter<double> ("GhostArea");
428     fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
429     activeAreaRepeats,
430     ghostArea));
431     fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(ghostEtaMax) );
432     }
433    
434     // restrict inputs to first "maxInputs" towers?
435     if ( iConfig.exists("restrictInputs") ) {
436     restrictInputs_ = iConfig.getParameter<bool>("restrictInputs");
437     maxInputs_ = iConfig.getParameter<unsigned int>("maxInputs");
438     }
439    
440    
441     string alias=iConfig.getUntrackedParameter<string>("alias",moduleLabel_);
442    
443     // make the "produces" statements
444     // produces<reco::CaloJetCollection>();
445     produces<reco::CaloJetCollection>("randomCones");
446     produces<std::vector<bool> >("directions");
447    
448     if(doRhoFastjet_){
449     produces<double>("rho");
450     produces<double>("sigma");
451     }
452    
453     }
454    
455    
456     //______________________________________________________________________________
457     JetAlgorithmAnalyzer::~JetAlgorithmAnalyzer()
458     {
459     }
460    
461     void JetAlgorithmAnalyzer::fillNtuple(int output, const std::vector<fastjet::PseudoJet>& jets, int step){
462     if(!doAnalysis_) return;
463    
464     TNtuple* nt;
465     TH2D* h;
466    
467     if(output == 1){
468     nt = ntJets;
469     h = hJets[step];
470     }
471     if(output == 0){
472     nt = ntTowers;
473     h = hTowers[step];
474     }
475     if(output == 2){
476     nt = ntTowersFromEvtPlane;
477     h = hTowersFromEvtPlane[step];
478     }
479     if(output == 3){
480     nt = ntJetTowers;
481     h = hJetTowers[step];
482     }
483    
484     double totet = 0;
485     int ntow = 0;
486     int nj = jets.size();
487    
488     cout<<"step : "<<step<<endl;
489     cout<<"Size of input : "<<nj<<endl;
490    
491     for(unsigned int i = 0; i < jets.size(); ++i){
492     const fastjet::PseudoJet& jet = jets[i];
493    
494     double pt = jet.perp();
495    
496     if(output != 1){
497     reco::CandidatePtr const & itow = inputs_[ jet.user_index() ];
498     pt = itow->et();
499     }
500    
501     double phi = jet.phi();
502     if(output == 2){
503     phi = phi - phi0_;
504     if(phi < 0) phi += 2*PI;
505     if(phi > 2*PI) phi -= 2*PI;
506     }
507    
508     double eta = jet.eta();
509    
510     if(eta > 0 && eta < 0.08){
511     // if(fabs(eta) < 1.){
512     totet += pt;
513     ntow++;
514     }
515    
516     nt->Fill(jet.eta(),phi,pt,step,iev_);
517     h->Fill(jet.eta(),phi,pt);
518     }
519    
520     cout<<"-----------------------------"<<endl;
521     cout<<"STEP = "<<step<<endl;
522     cout<<"Total ET = "<<totet<<endl;
523     cout<<"Towers counted = "<<ntow<<endl;
524     cout<<"Average tower ET = "<<totet/ntow<<endl;
525     cout<<"-----------------------------"<<endl;
526    
527     }
528    
529    
530     void JetAlgorithmAnalyzer::fillTowerNtuple(const std::vector<fastjet::PseudoJet>& jets, int step){
531     fillNtuple(0,jets,step);
532     fillNtuple(2,jets,step);
533     }
534    
535     void JetAlgorithmAnalyzer::fillJetNtuple(const std::vector<fastjet::PseudoJet>& jets, int step){
536     fillNtuple(1,jets,step);
537    
538     for(unsigned int i = 0; i < jets.size(); ++i){
539     const fastjet::PseudoJet& jet = jets[i];
540     std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(jet));
541    
542    
543    
544    
545     fillNtuple(3,fjConstituents,step);
546     }
547     }
548    
549     void JetAlgorithmAnalyzer::fillBkgNtuple(const PileUpSubtractor* subtractor, int step){
550     if(!doAnalysis_) return;
551     CaloTowerCollection col;
552     for(int ieta = -29; ieta < 30; ++ieta){
553     if(ieta == 0) continue;
554     CaloTowerDetId id(ieta,1);
555     const GlobalPoint& hitpoint = geo->getPosition(id);
556     cout<<"iETA "<<ieta<<endl;
557     double eta = hitpoint.eta();
558     cout<<"eta "<<eta<<endl;
559     math::PtEtaPhiMLorentzVector p4(1,eta,1.,1.);
560     GlobalPoint pos(1,1,1);
561     CaloTower c(id,1.,1.,1.,1,1, p4, pos,pos);
562     col.push_back(c);
563     reco::CandidatePtr ptr(&col,col.size() - 1);
564     double mean = subtractor->getMeanAtTower(ptr);
565     double sigma = subtractor->getSigmaAtTower(ptr);
566     double pu = subtractor->getPileUpAtTower(ptr);
567     ntPU->Fill(eta,mean,sigma,step,iev_);
568     hPU[step]->Fill(eta,pu);
569     hMean[step]->Fill(eta,mean);
570     hSigma[step]->Fill(eta,sigma);
571     }
572     }
573    
574     void JetAlgorithmAnalyzer::produce(edm::Event& iEvent,const edm::EventSetup& iSetup)
575     {
576    
577     // if(doAnalysis_) subtractor_->setDebug(1);
578     // else subtractor_->setDebug(0);
579    
580     if(!geo){
581     edm::ESHandle<CaloGeometry> pGeo;
582     iSetup.get<CaloGeometryRecord>().get(pGeo);
583     geo = pGeo.product();
584     }
585    
586     if(doAnalysis_ && centBin_ >= 0){
587     edm::Handle<reco::Centrality> cent;
588     iEvent.getByLabel(edm::InputTag("hiCentrality"),cent);
589     if(!cbins_) cbins_ = getCentralityBinsFromDB(iSetup);
590     int bin = cbins_->getBin(cent->EtHFhitSum());
591     if(bin != centBin_) return;
592     }
593    
594     LogDebug("VirtualJetProducer") << "Entered produce\n";
595     //determine signal vertex
596     vertex_=reco::Jet::Point(0,0,0);
597     if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
598     LogDebug("VirtualJetProducer") << "Adding PV info\n";
599     edm::Handle<reco::VertexCollection> pvCollection;
600     iEvent.getByLabel(srcPVs_,pvCollection);
601     if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
602     }
603    
604     // For Pileup subtraction using offset correction:
605     // set up geometry map
606     if ( doPUOffsetCorr_ ) {
607     subtractor_->setupGeometryMap(iEvent, iSetup);
608     }
609    
610     // clear data
611     LogDebug("VirtualJetProducer") << "Clear data\n";
612     fjInputs_.clear();
613     fjJets_.clear();
614     inputs_.clear();
615    
616     if(doAnalysis_ && doMC_){
617     Handle<GenHIEvent> hi;
618     iEvent.getByLabel("heavyIon",hi);
619     phi0_ = hi->evtPlane();
620    
621     double b = hi->b();
622     cout<<"===================================="<<endl;
623     cout<<"IMPACT PARAMETER OF THE EVENT : "<<b<<endl;
624     cout<<"===================================="<<endl;
625    
626    
627     }
628    
629     // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
630     edm::Handle<reco::CandidateView> inputsHandle;
631     iEvent.getByLabel(src_,inputsHandle);
632     for (size_t i = 0; i < inputsHandle->size(); ++i) {
633     inputs_.push_back(inputsHandle->ptrAt(i));
634     }
635     LogDebug("VirtualJetProducer") << "Got inputs\n";
636    
637     // Convert candidates to fastjet::PseudoJets.
638     // Also correct to Primary Vertex. Will modify fjInputs_
639     // and use inputs_
640     fjInputs_.reserve(inputs_.size());
641     inputTowers();
642     LogDebug("VirtualJetProducer") << "Inputted towers\n";
643    
644     fillTowerNtuple(fjInputs_,0);
645     fillBkgNtuple(subtractor_.get(),0);
646    
647     // For Pileup subtraction using offset correction:
648     // Subtract pedestal.
649    
650     if ( doPUOffsetCorr_ ) {
651     subtractor_->reset(inputs_,fjInputs_,fjJets_);
652     subtractor_->calculatePedestal(fjInputs_);
653    
654     fillTowerNtuple(fjInputs_,1);
655     fillBkgNtuple(subtractor_.get(),1);
656     subtractor_->subtractPedestal(fjInputs_);
657    
658     fillTowerNtuple(fjInputs_,2);
659     fillBkgNtuple(subtractor_.get(),2);
660    
661     LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
662     }
663    
664     // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
665     // This will use fjInputs_
666     runAlgorithm( iEvent, iSetup );
667    
668     fillTowerNtuple(fjInputs_,3);
669     fillBkgNtuple(subtractor_.get(),3);
670     fillJetNtuple(fjJets_,3);
671    
672     if ( doPUOffsetCorr_ ) {
673     subtractor_->setAlgorithm(fjClusterSeq_);
674     }
675    
676     LogDebug("VirtualJetProducer") << "Ran algorithm\n";
677    
678     // For Pileup subtraction using offset correction:
679     // Now we find jets and need to recalculate their energy,
680     // mark towers participated in jet,
681     // remove occupied towers from the list and recalculate mean and sigma
682     // put the initial towers collection to the jet,
683     // and subtract from initial towers in jet recalculated mean and sigma of towers
684     if ( doPUOffsetCorr_ ) {
685     vector<fastjet::PseudoJet> orphanInput;
686     subtractor_->calculateOrphanInput(orphanInput);
687     fillTowerNtuple(orphanInput,4);
688     fillBkgNtuple(subtractor_.get(),4);
689     fillJetNtuple(fjJets_,4);
690    
691     //only the id's of the orphan input are used, not their energy
692     subtractor_->calculatePedestal(orphanInput);
693     fillTowerNtuple(orphanInput,5);
694     fillBkgNtuple(subtractor_.get(),5);
695     fillJetNtuple(fjJets_,5);
696    
697     subtractor_->offsetCorrectJets();
698     fillTowerNtuple(orphanInput,6);
699     fillBkgNtuple(subtractor_.get(),6);
700     fillJetNtuple(fjJets_,6);
701    
702     }
703    
704     // Write the output jets.
705     // This will (by default) call the member function template
706     // "writeJets", but can be overridden.
707     // this will use inputs_
708     output( iEvent, iSetup );
709     fillBkgNtuple(subtractor_.get(),7);
710     fillJetNtuple(fjJets_,7);
711    
712     LogDebug("VirtualJetProducer") << "Wrote jets\n";
713    
714     ++iev_;
715    
716     doAnalysis_ = false;
717     return;
718     }
719    
720     void JetAlgorithmAnalyzer::output(edm::Event & iEvent, edm::EventSetup const& iSetup)
721     {
722     // Write jets and constitutents. Will use fjJets_, inputs_
723     // and fjClusterSeq_
724     switch( jetTypeE ) {
725     case JetType::CaloJet :
726     writeJets<reco::CaloJet>( iEvent, iSetup);
727     writeBkgJets<reco::CaloJet>( iEvent, iSetup);
728     break;
729     case JetType::PFJet :
730     writeJets<reco::PFJet>( iEvent, iSetup);
731     writeBkgJets<reco::PFJet>( iEvent, iSetup);
732     break;
733     case JetType::GenJet :
734     writeJets<reco::GenJet>( iEvent, iSetup);
735     writeBkgJets<reco::GenJet>( iEvent, iSetup);
736     break;
737     case JetType::TrackJet :
738     writeJets<reco::TrackJet>( iEvent, iSetup);
739     writeBkgJets<reco::TrackJet>( iEvent, iSetup);
740     break;
741     case JetType::BasicJet :
742     writeJets<reco::BasicJet>( iEvent, iSetup);
743     writeBkgJets<reco::BasicJet>( iEvent, iSetup);
744     break;
745     default:
746     edm::LogError("InvalidInput") << " invalid jet type in VirtualJetProducer\n";
747     break;
748     };
749    
750     }
751    
752     template< typename T >
753     void JetAlgorithmAnalyzer::writeJets( edm::Event & iEvent, edm::EventSetup const& iSetup )
754     {
755     // produce output jet collection
756    
757     using namespace reco;
758    
759     std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
760     jets->reserve(fjJets_.size());
761    
762     for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
763     // allocate this jet
764     T jet;
765     // get the fastjet jet
766     const fastjet::PseudoJet& fjJet = fjJets_[ijet];
767     // get the constituents from fastjet
768     std::vector<fastjet::PseudoJet> fjConstituents =
769     sorted_by_pt(fjClusterSeq_->constituents(fjJet));
770     // convert them to CandidatePtr vector
771     std::vector<CandidatePtr> constituents =
772     getConstituents(fjConstituents);
773    
774     // calcuate the jet area
775     double jetArea=0.0;
776     if ( doAreaFastjet_ ) {
777     fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
778     dynamic_cast<fastjet::ClusterSequenceArea const *>(&*fjClusterSeq_);
779     jetArea = clusterSequenceWithArea->area(fjJet);
780     }
781    
782     // write the specifics to the jet (simultaneously sets 4-vector, vertex).
783     // These are overridden functions that will call the appropriate
784     // specific allocator.
785     writeSpecific(jet,
786     Particle::LorentzVector(fjJet.px(),
787     fjJet.py(),
788     fjJet.pz(),
789     fjJet.E()),
790     vertex_,
791     constituents, iSetup);
792    
793     jet.setJetArea (jetArea);
794    
795     if(doPUOffsetCorr_){
796     jet.setPileup(subtractor_->getPileUpEnergy(ijet));
797     }else{
798     jet.setPileup (0.0);
799     }
800    
801     // add to the list
802     jets->push_back(jet);
803     }
804    
805     // put the jets in the collection
806     // iEvent.put(jets);
807    
808     // calculate rho (median pT per unit area, for PU&UE subtraction down the line
809     std::auto_ptr<double> rho(new double(0.0));
810     std::auto_ptr<double> sigma(new double(0.0));
811    
812     if (doRhoFastjet_) {
813     fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
814     dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
815     double mean_area = 0;
816     clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
817     iEvent.put(rho,"rho");
818     iEvent.put(sigma,"sigma");
819     }
820     }
821    
822     template< typename T >
823     void JetAlgorithmAnalyzer::writeBkgJets( edm::Event & iEvent, edm::EventSetup const& iSetup )
824     {
825     // produce output jet collection
826    
827     using namespace reco;
828    
829     std::vector<fastjet::PseudoJet> fjFakeJets_;
830     std::vector<std::vector<reco::CandidatePtr> > constituents_;
831     vector<double> phiRandom;
832     vector<double> etaRandom;
833     vector<double> et;
834     vector<double> pileUp;
835     std::auto_ptr<std::vector<bool> > directions(new std::vector<bool>());
836     directions->reserve(nFill_);
837    
838     phiRandom.reserve(nFill_);
839     etaRandom.reserve(nFill_);
840     et.reserve(nFill_);
841     pileUp.reserve(nFill_);
842    
843     fjFakeJets_.reserve(nFill_);
844     constituents_.reserve(nFill_);
845    
846     constituents_.reserve(nFill_);
847     for(int ijet = 0; ijet < nFill_; ++ijet){
848     vector<reco::CandidatePtr> vec;
849     constituents_.push_back(vec);
850     directions->push_back(true);
851     }
852    
853     for(int ijet = 0; ijet < nFill_; ++ijet){
854     phiRandom[ijet] = 2*pi*randomEngine->flat() - pi;
855     etaRandom[ijet] = 2*etaMax_*randomEngine->flat() - etaMax_;
856     et[ijet] = 0;
857     pileUp[ijet] = 0;
858     }
859    
860     for (vector<fastjet::PseudoJet>::const_iterator input_object = fjInputs_.begin (),
861     fjInputsEnd = fjInputs_.end();
862     input_object != fjInputsEnd; ++input_object) {
863    
864     const reco::CandidatePtr & tower=inputs_[input_object->user_index()];
865     const CaloTower* ctc = dynamic_cast<const CaloTower*>(tower.get());
866     int ieta = ctc->id().ieta();
867     int iphi = ctc->id().iphi();
868     CaloTowerDetId id(ieta,iphi);
869     const GlobalPoint& hitpoint = geo->getPosition(id);
870     double towEta = hitpoint.eta();
871     double towPhi = hitpoint.phi();
872    
873     for(int ir = 0; ir < nFill_; ++ir){
874     if(reco::deltaR(towEta,towPhi,etaRandom[ir],phiRandom[ir]) > rParam_) continue;
875    
876     constituents_[ir].push_back(tower);
877    
878     double towet = tower->et();
879     double putow = subtractor_->getPileUpAtTower(tower);
880     double etadd = towet - putow;
881     if(avoidNegative_ && etadd < 0.) etadd = 0;
882     et[ir] += etadd;
883     pileUp[ir] += towet - etadd;
884     }
885     }
886     cout<<"Start filling jets"<<endl;
887    
888     for(int ir = 0; ir < nFill_; ++ir){
889    
890     if(doAnalysis_) ntRandom->Fill(etaRandom[ir],phiRandom[ir],et[ir],pileUp[ir],iev_);
891    
892     if(et[ir] < 0){
893     // cout<<"Flipping vector"<<endl;
894     (*directions)[ir] = false;
895     et[ir] = -et[ir];
896     }else{
897     // cout<<"Keep vector same"<<endl;
898     (*directions)[ir] = true;
899     }
900     cout<<"Lorentz"<<endl;
901    
902     math::PtEtaPhiMLorentzVector p(et[ir],etaRandom[ir],phiRandom[ir],0);
903     fastjet::PseudoJet jet(p.px(),p.py(),p.pz(),p.energy());
904     fjFakeJets_.push_back(jet);
905     }
906    
907     std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
908     jets->reserve(fjFakeJets_.size());
909    
910     for (unsigned int ijet=0;ijet<fjFakeJets_.size();++ijet) {
911     // allocate this jet
912     T jet;
913     // get the fastjet jet
914     const fastjet::PseudoJet& fjJet = fjFakeJets_[ijet];
915    
916     // convert them to CandidatePtr vector
917     std::vector<CandidatePtr> constituents =
918     constituents_[ijet];
919    
920     writeSpecific(jet,
921     Particle::LorentzVector(fjJet.px(),
922     fjJet.py(),
923     fjJet.pz(),
924     fjJet.E()),
925     vertex_,
926     constituents, iSetup);
927    
928     // calcuate the jet area
929     double jetArea=0.0;
930     jet.setJetArea (jetArea);
931     if(doPUOffsetCorr_){
932     jet.setPileup(pileUp[ijet]);
933     }else{
934     jet.setPileup (0.0);
935     }
936    
937     // add to the list
938     jets->push_back(jet);
939     }
940    
941     // put the jets in the collection
942     iEvent.put(jets,"randomCones");
943     iEvent.put(directions,"directions");
944     }
945    
946     //void JetAlgorithmAnalyzer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup, std::vector<fastjet::PseudoJet>& input )
947     void JetAlgorithmAnalyzer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup)
948     {
949     if ( !doAreaFastjet_ && !doRhoFastjet_) {
950     fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
951     } else {
952     fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjActiveArea_ ) );
953     }
954     fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
955    
956     }
957    
958     void JetAlgorithmAnalyzer::inputTowers( )
959     {
960     std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
961     inEnd = inputs_.end(), i = inBegin;
962     for (; i != inEnd; ++i ) {
963     reco::CandidatePtr input = *i;
964     if (isnan(input->pt())) continue;
965     if (input->et() <inputEtMin_) continue;
966     if (input->energy()<inputEMin_) continue;
967     if (isAnomalousTower(input)) continue;
968    
969     //Check consistency of kinematics
970    
971     const CaloTower* ctc = dynamic_cast<const CaloTower*>(input.get());
972     if(ctc){
973     int ieta = ctc->id().ieta();
974     int iphi = ctc->id().iphi();
975    
976     ntuple->Fill(ieta, input->eta(), iphi, input->phi(),input->et(),ctc->emEt(),ctc->hadEt());
977    
978    
979    
980     if(abs(ieta) < 5){
981    
982     if(0){
983     math::RhoEtaPhiVector v(1.4,input->eta(),input->phi());
984     GlobalPoint point(v.x(),v.y(),v.z());
985     // const DetId d = geo->getClosestCell(point);
986     // HcalDetId hd(d);
987     HcalDetId hd(0);
988     if(hd.ieta() != ieta || hd.iphi() != iphi){
989     cout<<"Inconsistent kinematics!!! ET = "<<input->pt()<<endl;
990     cout<<"ieta candidate : "<<ieta<<" ieta detid : "<<hd.ieta()<<endl;
991     cout<<"iphi candidate : "<<iphi<<" iphi detid : "<<hd.iphi()<<endl;
992     }
993    
994     }
995    
996     if(0){
997     HcalDetId det(HcalBarrel,ieta,iphi,1);
998    
999     if(geo->present(det)){
1000     double eta = geo->getPosition(det).eta();
1001     double phi = geo->getPosition(det).phi();
1002    
1003     if(input->eta() != eta || input->phi() != phi){
1004     cout<<"Inconsistent kinematics!!! ET = "<<input->pt()<<endl;
1005     cout<<"eta candidate : "<<input->eta()<<" eta detid : "<<eta<<endl;
1006     cout<<"phi candidate : "<<input->phi()<<" phi detid : "<<phi<<endl;
1007     }
1008     }else{
1009     cout<<"DetId not present in the Calo Geometry : ieta = "<<ieta<<" iphi = "<<iphi<<endl;
1010     }
1011     }
1012     }
1013    
1014     }
1015    
1016     if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
1017     const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
1018     math::PtEtaPhiMLorentzVector ct(tower->p4(vertex_));
1019     fjInputs_.push_back(fastjet::PseudoJet(ct.px(),ct.py(),ct.pz(),ct.energy()));
1020     }
1021     else {
1022     fjInputs_.push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
1023     input->energy()));
1024     }
1025     fjInputs_.back().set_user_index(i - inBegin);
1026     }
1027     }
1028    
1029     //______________________________________________________________________________
1030     bool JetAlgorithmAnalyzer::isAnomalousTower(reco::CandidatePtr input)
1031     {
1032     if (!makeCaloJet(jetTypeE)) return false;
1033     const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
1034     if (0==tower) return false;
1035     if (tower->numBadEcalCells() >maxBadEcalCells_ ||
1036     tower->numRecoveredEcalCells() >maxRecoveredEcalCells_ ||
1037     tower->numProblematicEcalCells()>maxProblematicEcalCells_||
1038     tower->numBadHcalCells() >maxBadHcalCells_ ||
1039     tower->numRecoveredHcalCells() >maxRecoveredHcalCells_ ||
1040     tower->numProblematicHcalCells()>maxProblematicHcalCells_) return true;
1041     return false;
1042     }
1043     //______________________________________________________________________________
1044     void JetAlgorithmAnalyzer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents,
1045     reco::Jet* jet)
1046     {
1047     for (unsigned int i=0;i<fjConstituents.size();++i)
1048     jet->addDaughter(inputs_[fjConstituents[i].user_index()]);
1049     }
1050    
1051    
1052     //______________________________________________________________________________
1053     vector<reco::CandidatePtr>
1054     JetAlgorithmAnalyzer::getConstituents(const vector<fastjet::PseudoJet>&fjConstituents)
1055     {
1056     vector<reco::CandidatePtr> result;
1057     for (unsigned int i=0;i<fjConstituents.size();i++) {
1058     int index = fjConstituents[i].user_index();
1059     reco::CandidatePtr candidate = inputs_[index];
1060     result.push_back(candidate);
1061     }
1062     return result;
1063     }
1064    
1065     JetAlgorithmAnalyzer::JetType::Type
1066     JetAlgorithmAnalyzer::JetType::byName(const string &name)
1067     {
1068     const char **pos = std::find(names, names + LastJetType, name);
1069     if (pos == names + LastJetType) {
1070     std::string errorMessage="Requested jetType not supported: "+name+"\n";
1071     throw cms::Exception("Configuration",errorMessage);
1072     }
1073     return (Type)(pos-names);
1074     }
1075    
1076     const char *JetAlgorithmAnalyzer::JetType::names[] = {
1077     "BasicJet","GenJet","CaloJet","PFJet","TrackJet"
1078     };
1079    
1080    
1081     ////////////////////////////////////////////////////////////////////////////////
1082     // define as cmssw plugin
1083     ////////////////////////////////////////////////////////////////////////////////
1084    
1085     DEFINE_FWK_MODULE(JetAlgorithmAnalyzer);
1086