ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/modifiedFiles/VirtualJetProducer.cc
Revision: 1.2
Committed: Thu Sep 15 22:04:26 2011 UTC (13 years, 7 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, HEAD
Branch point for: branch_44x
Changes since 1.1: +1 -0 lines
Log Message:
muons

File Contents

# User Rev Content
1 frankma 1.1 ////////////////////////////////////////////////////////////////////////////////
2     //
3     // VirtualJetProducer
4     // ------------------
5     //
6     // 04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
7     ////////////////////////////////////////////////////////////////////////////////
8    
9     #include "RecoJets/JetProducers/plugins/VirtualJetProducer.h"
10     #include "RecoJets/JetProducers/interface/JetSpecific.h"
11    
12     #include "FWCore/Framework/interface/Event.h"
13     #include "FWCore/Framework/interface/EventSetup.h"
14     #include "FWCore/Framework/interface/ESHandle.h"
15     #include "FWCore/Utilities/interface/CodedException.h"
16     #include "FWCore/MessageLogger/interface/MessageLogger.h"
17     #include "FWCore/Framework/interface/MakerMacros.h"
18    
19     #include "DataFormats/Common/interface/View.h"
20     #include "DataFormats/Common/interface/Handle.h"
21     #include "DataFormats/VertexReco/interface/Vertex.h"
22     #include "DataFormats/VertexReco/interface/VertexFwd.h"
23     #include "DataFormats/JetReco/interface/CaloJetCollection.h"
24     #include "DataFormats/JetReco/interface/GenJetCollection.h"
25     #include "DataFormats/JetReco/interface/PFJetCollection.h"
26     #include "DataFormats/JetReco/interface/BasicJetCollection.h"
27     #include "DataFormats/JetReco/interface/TrackJetCollection.h"
28     #include "DataFormats/Candidate/interface/CandidateFwd.h"
29     #include "DataFormats/Candidate/interface/LeafCandidate.h"
30    
31     #include "fastjet/SISConePlugin.hh"
32     #include "fastjet/CMSIterativeConePlugin.hh"
33     #include "fastjet/ATLASConePlugin.hh"
34     #include "fastjet/CDFMidPointPlugin.hh"
35    
36     #include <iostream>
37     #include <memory>
38     #include <algorithm>
39     #include <limits>
40     #include <cmath>
41    
42     using namespace std;
43    
44    
45     namespace reco {
46     namespace helper {
47     struct GreaterByPtPseudoJet {
48     bool operator()( const fastjet::PseudoJet & t1, const fastjet::PseudoJet & t2 ) const {
49     return t1.perp() > t2.perp();
50     }
51     };
52    
53     }
54     }
55    
56     //______________________________________________________________________________
57     const char *VirtualJetProducer::JetType::names[] = {
58     "BasicJet","GenJet","CaloJet","PFJet","TrackJet"
59     };
60    
61    
62     //______________________________________________________________________________
63     VirtualJetProducer::JetType::Type
64     VirtualJetProducer::JetType::byName(const string &name)
65     {
66     const char **pos = std::find(names, names + LastJetType, name);
67     if (pos == names + LastJetType) {
68     std::string errorMessage="Requested jetType not supported: "+name+"\n";
69     throw cms::Exception("Configuration",errorMessage);
70     }
71     return (Type)(pos-names);
72     }
73    
74    
75     void VirtualJetProducer::makeProduces( std::string alias, std::string tag )
76     {
77     if (makeCaloJet(jetTypeE)) {
78     produces<reco::CaloJetCollection>(tag).setBranchAlias(alias);
79     }
80     else if (makePFJet(jetTypeE)) {
81     produces<reco::PFJetCollection>(tag).setBranchAlias(alias);
82     }
83     else if (makeGenJet(jetTypeE)) {
84     produces<reco::GenJetCollection>(tag).setBranchAlias(alias);
85     }
86     else if (makeTrackJet(jetTypeE)) {
87     produces<reco::TrackJetCollection>(tag).setBranchAlias(alias);
88     }
89     else if (makeBasicJet(jetTypeE)) {
90     produces<reco::BasicJetCollection>(tag).setBranchAlias(alias);
91     }
92     }
93    
94     ////////////////////////////////////////////////////////////////////////////////
95     // construction / destruction
96     ////////////////////////////////////////////////////////////////////////////////
97    
98     //______________________________________________________________________________
99     VirtualJetProducer::VirtualJetProducer(const edm::ParameterSet& iConfig)
100     : moduleLabel_ (iConfig.getParameter<string> ("@module_label"))
101     , src_ (iConfig.getParameter<edm::InputTag>("src"))
102     , srcPVs_ (iConfig.getParameter<edm::InputTag>("srcPVs"))
103     , jetType_ (iConfig.getParameter<string> ("jetType"))
104     , jetAlgorithm_ (iConfig.getParameter<string> ("jetAlgorithm"))
105     , rParam_ (iConfig.getParameter<double> ("rParam"))
106     , inputEtMin_ (iConfig.getParameter<double> ("inputEtMin"))
107     , inputEMin_ (iConfig.getParameter<double> ("inputEMin"))
108     , jetPtMin_ (iConfig.getParameter<double> ("jetPtMin"))
109 yilmaz 1.2 //, puPtMin_ (iConfig.getParameter<double> ("puPtMin"))
110 frankma 1.1 , doPVCorrection_(iConfig.getParameter<bool> ("doPVCorrection"))
111     , restrictInputs_(false)
112     , maxInputs_(99999999)
113     , doAreaFastjet_ (iConfig.getParameter<bool> ("doAreaFastjet"))
114     , doRhoFastjet_ (iConfig.getParameter<bool> ("doRhoFastjet"))
115     , doPUOffsetCorr_(iConfig.getParameter<bool> ("doPUOffsetCorr"))
116     , maxBadEcalCells_ (iConfig.getParameter<unsigned>("maxBadEcalCells"))
117     , maxRecoveredEcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredEcalCells"))
118     , maxProblematicEcalCells_(iConfig.getParameter<unsigned>("maxProblematicEcalCells"))
119     , maxBadHcalCells_ (iConfig.getParameter<unsigned>("maxBadHcalCells"))
120     , maxRecoveredHcalCells_ (iConfig.getParameter<unsigned>("maxRecoveredHcalCells"))
121     , maxProblematicHcalCells_(iConfig.getParameter<unsigned>("maxProblematicHcalCells"))
122     , jetCollInstanceName_ ("")
123     {
124    
125     //
126     // additional parameters to think about:
127     // - overlap threshold (set to 0.75 for the time being)
128     // - p parameter for generalized kT (set to -2 for the time being)
129     // - fastjet PU subtraction parameters (not yet considered)
130     //
131     if (jetAlgorithm_=="SISCone") {
132     fjPlugin_ = PluginPtr( new fastjet::SISConePlugin(rParam_,0.75,0,0.0,false,
133     fastjet::SISConePlugin::SM_pttilde) );
134     fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(&*fjPlugin_) );
135     }
136     else if (jetAlgorithm_=="IterativeCone") {
137     fjPlugin_ = PluginPtr(new fastjet::CMSIterativeConePlugin(rParam_,1.0));
138     fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
139     }
140     else if (jetAlgorithm_=="CDFMidPoint") {
141     fjPlugin_ = PluginPtr(new fastjet::CDFMidPointPlugin(rParam_,0.75));
142     fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
143     }
144     else if (jetAlgorithm_=="ATLASCone") {
145     fjPlugin_ = PluginPtr(new fastjet::ATLASConePlugin(rParam_));
146     fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(&*fjPlugin_));
147     }
148     else if (jetAlgorithm_=="Kt")
149     fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::kt_algorithm,rParam_));
150     else if (jetAlgorithm_=="CambridgeAachen")
151     fjJetDefinition_= JetDefPtr(new fastjet::JetDefinition(fastjet::cambridge_algorithm,
152     rParam_) );
153     else if (jetAlgorithm_=="AntiKt")
154     fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::antikt_algorithm,rParam_) );
155     else if (jetAlgorithm_=="GeneralizedKt")
156     fjJetDefinition_= JetDefPtr( new fastjet::JetDefinition(fastjet::genkt_algorithm,
157     rParam_,-2) );
158     else
159     throw cms::Exception("Invalid jetAlgorithm")
160     <<"Jet algorithm for VirtualJetProducer is invalid, Abort!\n";
161    
162     jetTypeE=JetType::byName(jetType_);
163    
164     if ( iConfig.exists("jetCollInstanceName") ) {
165     jetCollInstanceName_ = iConfig.getParameter<string>("jetCollInstanceName");
166     }
167    
168     if ( doPUOffsetCorr_ ) {
169     if ( jetTypeE != JetType::CaloJet ) {
170     //throw cms::Exception("InvalidInput") << "Can only offset correct jets of type CaloJet";
171     std::cout<<" Implementing background subtraction using a pseudo calo grid"<<std::endl;
172     }
173    
174     if(iConfig.exists("subtractorName")) puSubtractorName_ = iConfig.getParameter<string> ("subtractorName");
175     else puSubtractorName_ = std::string();
176    
177     if(puSubtractorName_.empty()){
178     edm::LogWarning("VirtualJetProducer") << "Pile Up correction on; however, pile up type is not specified. Using default... \n";
179     subtractor_ = boost::shared_ptr<PileUpSubtractor>(new PileUpSubtractor(iConfig));
180     }else{
181     subtractor_ = boost::shared_ptr<PileUpSubtractor>(PileUpSubtractorFactory::get()->create( puSubtractorName_, iConfig));
182     }
183     }
184    
185     // do fasjet area / rho calcluation? => accept corresponding parameters
186     if ( doAreaFastjet_ || doRhoFastjet_ ) {
187     // Eta range of jets to be considered for Rho calculation
188     // Should be at most (jet acceptance - jet radius)
189     double rhoEtaMax=iConfig.getParameter<double>("Rho_EtaMax");
190     // default Ghost_EtaMax should be 5
191     double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
192     // default Active_Area_Repeats 1
193     int activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
194     // default GhostArea 0.01
195     double ghostArea = iConfig.getParameter<double> ("GhostArea");
196     fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
197     activeAreaRepeats,
198     ghostArea));
199     fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(rhoEtaMax) );
200     }
201    
202     // restrict inputs to first "maxInputs" towers?
203     if ( iConfig.exists("restrictInputs") ) {
204     restrictInputs_ = iConfig.getParameter<bool>("restrictInputs");
205     maxInputs_ = iConfig.getParameter<unsigned int>("maxInputs");
206     }
207    
208    
209     string alias=iConfig.getUntrackedParameter<string>("alias",moduleLabel_);
210    
211     // make the "produces" statements
212     makeProduces( alias, jetCollInstanceName_ );
213    
214     doFastJetNonUniform_ = false;
215     if(iConfig.exists("doFastJetNonUniform")) doFastJetNonUniform_ = iConfig.getParameter<bool> ("doFastJetNonUniform");
216     if(doFastJetNonUniform_){
217     puCenters_ = iConfig.getParameter<std::vector<double> >("puCenters");
218     puWidth_ = iConfig.getParameter<double>("puWidth");
219     }
220     produces<std::vector<double> >("rhos");
221     produces<std::vector<double> >("sigmas");
222     produces<double>("rho");
223     produces<double>("sigma");
224    
225     }
226    
227     //______________________________________________________________________________
228     VirtualJetProducer::~VirtualJetProducer()
229     {
230     }
231    
232    
233     ////////////////////////////////////////////////////////////////////////////////
234     // implementation of member functions
235     ////////////////////////////////////////////////////////////////////////////////
236    
237     //______________________________________________________________________________
238     void VirtualJetProducer::produce(edm::Event& iEvent,const edm::EventSetup& iSetup)
239     {
240     LogDebug("VirtualJetProducer") << "Entered produce\n";
241     //determine signal vertex
242     vertex_=reco::Jet::Point(0,0,0);
243     if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
244     LogDebug("VirtualJetProducer") << "Adding PV info\n";
245     edm::Handle<reco::VertexCollection> pvCollection;
246     iEvent.getByLabel(srcPVs_,pvCollection);
247     if (pvCollection->size()>0) vertex_=pvCollection->begin()->position();
248     }
249    
250     // For Pileup subtraction using offset correction:
251     // set up geometry map
252     if ( doPUOffsetCorr_ ) {
253     subtractor_->setupGeometryMap(iEvent, iSetup);
254     }
255    
256     // clear data
257     LogDebug("VirtualJetProducer") << "Clear data\n";
258     fjInputs_.clear();
259     fjJets_.clear();
260     inputs_.clear();
261    
262     // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
263     edm::Handle<reco::CandidateView> inputsHandle;
264     iEvent.getByLabel(src_,inputsHandle);
265     for (size_t i = 0; i < inputsHandle->size(); ++i) {
266     inputs_.push_back(inputsHandle->ptrAt(i));
267     }
268     LogDebug("VirtualJetProducer") << "Got inputs\n";
269    
270     // Convert candidates to fastjet::PseudoJets.
271     // Also correct to Primary Vertex. Will modify fjInputs_
272     // and use inputs_
273     fjInputs_.reserve(inputs_.size());
274     inputTowers();
275     LogDebug("VirtualJetProducer") << "Inputted towers\n";
276    
277     // For Pileup subtraction using offset correction:
278     // Subtract pedestal.
279     if ( doPUOffsetCorr_ ) {
280     subtractor_->reset(inputs_,fjInputs_,fjJets_);
281     subtractor_->calculatePedestal(fjInputs_);
282     subtractor_->subtractPedestal(fjInputs_);
283     LogDebug("VirtualJetProducer") << "Subtracted pedestal\n";
284     }
285    
286     // Run algorithm. Will modify fjJets_ and allocate fjClusterSeq_.
287     // This will use fjInputs_
288     runAlgorithm( iEvent, iSetup );
289     if ( doPUOffsetCorr_ ) {
290     subtractor_->setAlgorithm(fjClusterSeq_);
291     }
292    
293     LogDebug("VirtualJetProducer") << "Ran algorithm\n";
294    
295     // For Pileup subtraction using offset correction:
296     // Now we find jets and need to recalculate their energy,
297     // mark towers participated in jet,
298     // remove occupied towers from the list and recalculate mean and sigma
299     // put the initial towers collection to the jet,
300     // and subtract from initial towers in jet recalculated mean and sigma of towers
301     if ( doPUOffsetCorr_ ) {
302     vector<fastjet::PseudoJet> orphanInput;
303     subtractor_->calculateOrphanInput(orphanInput);
304     subtractor_->calculatePedestal(orphanInput);
305     subtractor_->offsetCorrectJets();
306     }
307    
308     // Write the output jets.
309     // This will (by default) call the member function template
310     // "writeJets", but can be overridden.
311     // this will use inputs_
312     output( iEvent, iSetup );
313     LogDebug("VirtualJetProducer") << "Wrote jets\n";
314    
315     return;
316     }
317    
318     //______________________________________________________________________________
319    
320     void VirtualJetProducer::inputTowers( )
321     {
322     std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
323     inEnd = inputs_.end(), i = inBegin;
324     for (; i != inEnd; ++i ) {
325     reco::CandidatePtr input = *i;
326     if (isnan(input->pt())) continue;
327     if (input->et() <inputEtMin_) continue;
328     if (input->energy()<inputEMin_) continue;
329     if (isAnomalousTower(input)) continue;
330     if (makeCaloJet(jetTypeE)&&doPVCorrection_) {
331     const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
332     math::PtEtaPhiMLorentzVector ct(tower->p4(vertex_));
333     fjInputs_.push_back(fastjet::PseudoJet(ct.px(),ct.py(),ct.pz(),ct.energy()));
334     }
335     else {
336     fjInputs_.push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
337     input->energy()));
338     }
339     fjInputs_.back().set_user_index(i - inBegin);
340     }
341    
342     if ( restrictInputs_ && fjInputs_.size() > maxInputs_ ) {
343     reco::helper::GreaterByPtPseudoJet pTComparator;
344     std::sort(fjInputs_.begin(), fjInputs_.end(), pTComparator);
345     fjInputs_.resize(maxInputs_);
346     edm::LogWarning("JetRecoTooManyEntries") << "Too many inputs in the event, limiting to first " << maxInputs_ << ". Output is suspect.";
347     }
348     }
349    
350     //______________________________________________________________________________
351     bool VirtualJetProducer::isAnomalousTower(reco::CandidatePtr input)
352     {
353     if (!makeCaloJet(jetTypeE)) return false;
354     const CaloTower* tower=dynamic_cast<const CaloTower*>(input.get());
355     if (0==tower) return false;
356     if (tower->numBadEcalCells() >maxBadEcalCells_ ||
357     tower->numRecoveredEcalCells() >maxRecoveredEcalCells_ ||
358     tower->numProblematicEcalCells()>maxProblematicEcalCells_||
359     tower->numBadHcalCells() >maxBadHcalCells_ ||
360     tower->numRecoveredHcalCells() >maxRecoveredHcalCells_ ||
361     tower->numProblematicHcalCells()>maxProblematicHcalCells_) return true;
362     return false;
363     }
364    
365    
366     //------------------------------------------------------------------------------
367     // This is pure virtual.
368     //______________________________________________________________________________
369     // void VirtualJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup,
370     // std::vector<edm::Ptr<reco::Candidate> > const & inputs_);
371    
372     //______________________________________________________________________________
373     void VirtualJetProducer::copyConstituents(const vector<fastjet::PseudoJet>& fjConstituents,
374     reco::Jet* jet)
375     {
376     for (unsigned int i=0;i<fjConstituents.size();++i)
377     jet->addDaughter(inputs_[fjConstituents[i].user_index()]);
378     }
379    
380    
381     //______________________________________________________________________________
382     vector<reco::CandidatePtr>
383     VirtualJetProducer::getConstituents(const vector<fastjet::PseudoJet>&fjConstituents)
384     {
385     vector<reco::CandidatePtr> result;
386     for (unsigned int i=0;i<fjConstituents.size();i++) {
387     int index = fjConstituents[i].user_index();
388     reco::CandidatePtr candidate = inputs_[index];
389     result.push_back(candidate);
390     }
391     return result;
392     }
393    
394     void VirtualJetProducer::output(edm::Event & iEvent, edm::EventSetup const& iSetup)
395     {
396     // Write jets and constitutents. Will use fjJets_, inputs_
397     // and fjClusterSeq_
398     switch( jetTypeE ) {
399     case JetType::CaloJet :
400     writeJets<reco::CaloJet>( iEvent, iSetup);
401     break;
402     case JetType::PFJet :
403     writeJets<reco::PFJet>( iEvent, iSetup);
404     break;
405     case JetType::GenJet :
406     writeJets<reco::GenJet>( iEvent, iSetup);
407     break;
408     case JetType::TrackJet :
409     writeJets<reco::TrackJet>( iEvent, iSetup);
410     break;
411     case JetType::BasicJet :
412     writeJets<reco::BasicJet>( iEvent, iSetup);
413     break;
414     default:
415     throw cms::Exception("InvalidInput") << "invalid jet type in VirtualJetProducer\n";
416     break;
417     };
418    
419     }
420    
421     template< typename T >
422     void VirtualJetProducer::writeJets( edm::Event & iEvent, edm::EventSetup const& iSetup )
423     {
424     // produce output jet collection
425    
426     using namespace reco;
427    
428     std::auto_ptr<std::vector<T> > jets(new std::vector<T>() );
429     jets->reserve(fjJets_.size());
430    
431     for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
432     // allocate this jet
433     T jet;
434     // get the fastjet jet
435     const fastjet::PseudoJet& fjJet = fjJets_[ijet];
436     // get the constituents from fastjet
437     std::vector<fastjet::PseudoJet> fjConstituents =
438     sorted_by_pt(fjClusterSeq_->constituents(fjJet));
439     // convert them to CandidatePtr vector
440     std::vector<CandidatePtr> constituents =
441     getConstituents(fjConstituents);
442    
443     // calcuate the jet area
444     double jetArea=0.0;
445     if ( doAreaFastjet_ ) {
446     fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
447     dynamic_cast<fastjet::ClusterSequenceArea const *>(&*fjClusterSeq_);
448     jetArea = clusterSequenceWithArea->area(fjJet);
449     }
450    
451     // write the specifics to the jet (simultaneously sets 4-vector, vertex).
452     // These are overridden functions that will call the appropriate
453     // specific allocator.
454     writeSpecific(jet,
455     Particle::LorentzVector(fjJet.px(),
456     fjJet.py(),
457     fjJet.pz(),
458     fjJet.E()),
459     vertex_,
460     constituents, iSetup);
461    
462     jet.setJetArea (jetArea);
463    
464     if(doPUOffsetCorr_){
465     jet.setPileup(subtractor_->getPileUpEnergy(ijet));
466     }else{
467     jet.setPileup (0.0);
468     }
469    
470     // add to the list
471     jets->push_back(jet);
472     }
473    
474     // put the jets in the collection
475     iEvent.put(jets);
476    
477     if (doRhoFastjet_) {
478     if(doFastJetNonUniform_){
479     std::auto_ptr<std::vector<double> > rhos(new std::vector<double>);
480     std::auto_ptr<std::vector<double> > sigmas(new std::vector<double>);
481     int nEta = puCenters_.size();
482     rhos->reserve(nEta);
483     sigmas->reserve(nEta);
484     fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
485     dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
486    
487     for(int ie = 0; ie < nEta; ++ie){
488     double eta = puCenters_[ie];
489     double rho = 0;
490     double sigma = 0;
491     double etamin=eta-puWidth_;
492     double etamax=eta+puWidth_;
493     fastjet::RangeDefinition range_rho(etamin,etamax);
494     clusterSequenceWithArea->get_median_rho_and_sigma(range_rho, true, rho, sigma);
495     rhos->push_back(rho);
496     sigmas->push_back(sigma);
497     }
498     iEvent.put(rhos,"rhos");
499     iEvent.put(sigmas,"sigmas");
500     }else{
501     std::auto_ptr<double> rho(new double(0.0));
502     std::auto_ptr<double> sigma(new double(0.0));
503     double mean_area = 0;
504    
505     fastjet::ClusterSequenceArea const * clusterSequenceWithArea =
506     dynamic_cast<fastjet::ClusterSequenceArea const *> ( &*fjClusterSeq_ );
507     clusterSequenceWithArea->get_median_rho_and_sigma(*fjRangeDef_,false,*rho,*sigma,mean_area);
508     iEvent.put(rho,"rho");
509     iEvent.put(sigma,"sigma");
510     }
511    
512     }
513     }
514    
515    
516