Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:32

0001 ////////////////////////////////////////////////////////////////////////////////
0002 //
0003 // FastjetJetProducer
0004 // ------------------
0005 //
0006 //            04/21/2009 Philipp Schieferdecker <philipp.schieferdecker@cern.ch>
0007 ////////////////////////////////////////////////////////////////////////////////
0008 
0009 #include "RecoJets/JetProducers/plugins/FastjetJetProducer.h"
0010 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0012 
0013 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0014 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0015 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0016 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
0017 #include "DataFormats/Common/interface/Handle.h"
0018 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0019 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0020 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0021 
0022 #include "FWCore/Framework/interface/Event.h"
0023 #include "FWCore/Framework/interface/EventSetup.h"
0024 #include "FWCore/Framework/interface/ESHandle.h"
0025 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027 
0028 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
0029 #include "Geometry/Records/interface/CaloGeometryRecord.h"
0030 
0031 #include "fastjet/SISConePlugin.hh"
0032 #include "fastjet/CMSIterativeConePlugin.hh"
0033 #include "fastjet/ATLASConePlugin.hh"
0034 #include "fastjet/CDFMidPointPlugin.hh"
0035 #include "fastjet/tools/Filter.hh"
0036 #include "fastjet/tools/Pruner.hh"
0037 #include "fastjet/tools/MassDropTagger.hh"
0038 #include "fastjet/contrib/SoftDrop.hh"
0039 #include "fastjet/tools/JetMedianBackgroundEstimator.hh"
0040 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"
0041 #include "fastjet/tools/Subtractor.hh"
0042 #include "fastjet/contrib/ConstituentSubtractor.hh"
0043 #include "RecoJets/JetAlgorithms/interface/CMSBoostedTauSeedingAlgorithm.h"
0044 
0045 #include <iostream>
0046 #include <memory>
0047 #include <algorithm>
0048 #include <limits>
0049 #include <cmath>
0050 //#include <fstream>
0051 
0052 using namespace std;
0053 using namespace edm;
0054 
0055 ////////////////////////////////////////////////////////////////////////////////
0056 // construction / destruction
0057 ////////////////////////////////////////////////////////////////////////////////
0058 
0059 //______________________________________________________________________________
0060 FastjetJetProducer::FastjetJetProducer(const edm::ParameterSet& iConfig) : VirtualJetProducer(iConfig) {
0061   useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
0062   useOnlyOnePV_ = iConfig.getParameter<bool>("UseOnlyOnePV");
0063   dzTrVtxMax_ = iConfig.getParameter<double>("DzTrVtxMax");
0064   dxyTrVtxMax_ = iConfig.getParameter<double>("DxyTrVtxMax");
0065   minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
0066   maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
0067 
0068   useMassDropTagger_ = iConfig.getParameter<bool>("useMassDropTagger");
0069   muCut_ = iConfig.getParameter<double>("muCut");
0070   yCut_ = iConfig.getParameter<double>("yCut");
0071 
0072   useFiltering_ = iConfig.getParameter<bool>("useFiltering");
0073   rFilt_ = iConfig.getParameter<double>("rFilt");
0074   nFilt_ = iConfig.getParameter<int>("nFilt");
0075   useDynamicFiltering_ = iConfig.getParameter<bool>("useDynamicFiltering");
0076   if (useDynamicFiltering_)
0077     rFiltDynamic_ = std::make_shared<DynamicRfilt>(rFilt_, rFiltFactor_);
0078   rFiltFactor_ = iConfig.getParameter<double>("rFiltFactor");
0079 
0080   useTrimming_ = iConfig.getParameter<bool>("useTrimming");
0081   trimPtFracMin_ = iConfig.getParameter<double>("trimPtFracMin");
0082 
0083   usePruning_ = iConfig.getParameter<bool>("usePruning");
0084   zCut_ = iConfig.getParameter<double>("zcut");
0085   RcutFactor_ = iConfig.getParameter<double>("rcut_factor");
0086   useKtPruning_ = iConfig.getParameter<bool>("useKtPruning");
0087 
0088   useCMSBoostedTauSeedingAlgorithm_ = iConfig.getParameter<bool>("useCMSBoostedTauSeedingAlgorithm");
0089   subjetPtMin_ = iConfig.getParameter<double>("subjetPtMin");
0090   muMin_ = iConfig.getParameter<double>("muMin");
0091   muMax_ = iConfig.getParameter<double>("muMax");
0092   yMin_ = iConfig.getParameter<double>("yMin");
0093   yMax_ = iConfig.getParameter<double>("yMax");
0094   dRMin_ = iConfig.getParameter<double>("dRMin");
0095   dRMax_ = iConfig.getParameter<double>("dRMax");
0096   maxDepth_ = iConfig.getParameter<int>("maxDepth");
0097 
0098   useConstituentSubtraction_ = iConfig.getParameter<bool>("useConstituentSubtraction");
0099   csRho_EtaMax_ = iConfig.getParameter<double>("csRho_EtaMax");
0100   csRParam_ = iConfig.getParameter<double>("csRParam");
0101 
0102   useSoftDrop_ = iConfig.getParameter<bool>("useSoftDrop");
0103   zCut_ = iConfig.getParameter<double>("zcut");
0104   beta_ = iConfig.getParameter<double>("beta");
0105   R0_ = iConfig.getParameter<double>("R0");
0106 
0107   correctShape_ = iConfig.getParameter<bool>("correctShape");
0108   gridMaxRapidity_ = iConfig.getParameter<double>("gridMaxRapidity");
0109   gridSpacing_ = iConfig.getParameter<double>("gridSpacing");
0110 
0111   input_chrefcand_token_ =
0112       consumes<edm::View<reco::RecoChargedRefCandidate>>(iConfig.getParameter<edm::InputTag>("src"));
0113 
0114   if (useFiltering_ || useTrimming_ || usePruning_ || useMassDropTagger_ || useCMSBoostedTauSeedingAlgorithm_ ||
0115       useConstituentSubtraction_ || useSoftDrop_ || correctShape_)
0116     useExplicitGhosts_ = true;
0117 
0118   ////// adding exceptions
0119 
0120   if ((useMassDropTagger_) && ((muCut_ == -1) || (yCut_ == -1)))
0121     throw cms::Exception("useMassDropTagger")
0122         << "Parameters muCut and/or yCut for Mass Drop are not defined." << std::endl;
0123 
0124   if ((useFiltering_) && ((rFilt_ == -1) || (nFilt_ == -1))) {
0125     throw cms::Exception("useFiltering") << "Parameters rFilt and/or nFilt for Filtering are not defined." << std::endl;
0126     if ((useDynamicFiltering_) && (rFiltFactor_ == -1))
0127       throw cms::Exception("useDynamicFiltering")
0128           << "Parameters rFiltFactor for DynamicFiltering is not defined." << std::endl;
0129   }
0130 
0131   if ((useTrimming_) && ((rFilt_ == -1) || (trimPtFracMin_ == -1)))
0132     throw cms::Exception("useTrimming") << "Parameters rFilt and/or trimPtFracMin for Trimming are not defined."
0133                                         << std::endl;
0134 
0135   if ((usePruning_) && ((zCut_ == -1) || (RcutFactor_ == -1) || (nFilt_ == -1)))
0136     throw cms::Exception("usePruning") << "Parameters zCut and/or RcutFactor and/or nFilt for Pruning are not defined."
0137                                        << std::endl;
0138 
0139   if ((useCMSBoostedTauSeedingAlgorithm_) &&
0140       ((subjetPtMin_ == -1) || (maxDepth_ == -1) || (muMin_ == -1) || (muMax_ == -1) || (yMin_ == -1) ||
0141        (yMax_ == -1) || (dRMin_ == -1) || (dRMax_ == -1)))
0142     throw cms::Exception("useCMSBoostedTauSeedingAlgorithm")
0143         << "Parameters subjetPtMin, muMin, muMax, yMin, yMax, dRmin, dRmax, maxDepth for CMSBoostedTauSeedingAlgorithm "
0144            "are not defined."
0145         << std::endl;
0146 
0147   if (useConstituentSubtraction_ && (fjAreaDefinition_.get() == nullptr))
0148     throw cms::Exception("AreaMustBeSet")
0149         << "Logic error. The area definition must be set if you use constituent subtraction." << std::endl;
0150 
0151   if ((useConstituentSubtraction_) && ((csRho_EtaMax_ == -1) || (csRParam_ == -1)))
0152     throw cms::Exception("useConstituentSubtraction")
0153         << "Parameters csRho_EtaMax and/or csRParam for ConstituentSubtraction are not defined." << std::endl;
0154 
0155   if (useSoftDrop_ && usePruning_)
0156     throw cms::Exception("PruningAndSoftDrop")
0157         << "Logic error. Soft drop is a generalized pruning, do not run them together." << std::endl;
0158 
0159   if ((useSoftDrop_) && ((zCut_ == -1) || (beta_ == -1) || (R0_ == -1)))
0160     throw cms::Exception("useSoftDrop") << "Parameters zCut and/or beta and/or R0 for SoftDrop are not defined."
0161                                         << std::endl;
0162 
0163   if ((correctShape_) && ((gridMaxRapidity_ == -1) || (gridSpacing_ == -1)))
0164     throw cms::Exception("correctShape")
0165         << "Parameters gridMaxRapidity and/or gridSpacing for SoftDrop are not defined." << std::endl;
0166 }
0167 
0168 //______________________________________________________________________________
0169 FastjetJetProducer::~FastjetJetProducer() {}
0170 
0171 ////////////////////////////////////////////////////////////////////////////////
0172 // implementation of member functions
0173 ////////////////////////////////////////////////////////////////////////////////
0174 
0175 void FastjetJetProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0176   // for everything but track jets
0177   if (!makeTrackJet(jetTypeE)) {
0178     // use the default production from one collection
0179     VirtualJetProducer::produce(iEvent, iSetup);
0180 
0181   } else {  // produce trackjets from tracks grouped per primary vertex
0182 
0183     produceTrackJets(iEvent, iSetup);
0184   }
0185 
0186   // fjClusterSeq_ retains quite a lot of memory - about 1 to 7Mb at 200 pileup
0187   // depending on the exact configuration; and there are 24 FastjetJetProducers in the
0188   // sequence so this adds up to about 60 Mb. It's allocated every time runAlgorithm
0189   // is called, so safe to delete here.
0190   fjClusterSeq_.reset();
0191 }
0192 
0193 void FastjetJetProducer::produceTrackJets(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0194   // read in the track candidates
0195   edm::Handle<edm::View<reco::RecoChargedRefCandidate>> inputsHandle;
0196   iEvent.getByToken(input_chrefcand_token_, inputsHandle);
0197 
0198   // make collection with pointers so we can play around with it
0199   std::vector<edm::Ptr<reco::RecoChargedRefCandidate>> allInputs;
0200   std::vector<edm::Ptr<reco::Candidate>> origInputs;
0201   for (size_t i = 0; i < inputsHandle->size(); ++i) {
0202     allInputs.push_back(inputsHandle->ptrAt(i));
0203     origInputs.push_back(inputsHandle->ptrAt(i));
0204   }
0205 
0206   // read in the PV collection
0207   edm::Handle<reco::VertexCollection> pvCollection;
0208   iEvent.getByToken(input_vertex_token_, pvCollection);
0209   // define the overall output jet container
0210   auto jets = std::make_unique<std::vector<reco::TrackJet>>();
0211 
0212   // loop over the good vertices, clustering for each vertex separately
0213   for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
0214     if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_)
0215       continue;
0216 
0217     // clear the intermediate containers
0218     inputs_.clear();
0219     fjInputs_.clear();
0220     fjJets_.clear();
0221 
0222     // if only vertex-associated tracks should be used
0223     if (useOnlyVertexTracks_) {
0224       // loop over the tracks associated to the vertex
0225       for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
0226         // whether a match was found in the track candidate input
0227         bool found = false;
0228         // loop over input track candidates
0229         for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate>>::iterator itIn = allInputs.begin();
0230              itIn != allInputs.end();
0231              ++itIn) {
0232           // match the input track candidate to the track from the vertex
0233           reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
0234           // check if the tracks match
0235           if ((*itIn)->track() == trref) {
0236             found = true;
0237             // add this track candidate to the input for clustering
0238             inputs_.push_back(*itIn);
0239             // erase the track candidate from the total list of input, so we don't reuse it later
0240             allInputs.erase(itIn);
0241             // found the candidate track corresponding to the vertex track, so stop the loop
0242             break;
0243           }  // end if match found
0244         }    // end loop over input tracks
0245         // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
0246         if (!found)
0247           edm::LogInfo("FastjetTrackJetProducer")
0248               << "Ignoring a track at vertex which is not in input track collection!";
0249       }  // end loop over tracks associated to vertex
0250       // if all inpt track candidates should be used
0251     } else {
0252       // loop over input track candidates
0253       for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate>>::iterator itIn = allInputs.begin();
0254            itIn != allInputs.end();
0255            ++itIn) {
0256         // check if the track is close enough to the vertex
0257         float dz = (*itIn)->track()->dz(itVtx->position());
0258         float dxy = (*itIn)->track()->dxy(itVtx->position());
0259         if (fabs(dz) > dzTrVtxMax_)
0260           continue;
0261         if (fabs(dxy) > dxyTrVtxMax_)
0262           continue;
0263         bool closervtx = false;
0264         // now loop over the good vertices a second time
0265         for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end();
0266              ++itVtx2) {
0267           if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_)
0268             continue;
0269           // and check this track is closer to any other vertex (if more than 1 vertex considered)
0270           if (!useOnlyOnePV_ && itVtx != itVtx2 && fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
0271             closervtx = true;
0272             break;  // 1 closer vertex makes the track already not matched, so break
0273           }
0274         }
0275         // don't add this track if another vertex is found closer
0276         if (closervtx)
0277           continue;
0278         // add this track candidate to the input for clustering
0279         inputs_.push_back(*itIn);
0280         // erase the track candidate from the total list of input, so we don't reuse it later
0281         allInputs.erase(itIn);
0282         // take a step back in the loop since we just erased
0283         --itIn;
0284       }
0285     }
0286 
0287     // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
0288     fjInputs_.reserve(inputs_.size());
0289     inputTowers();
0290     LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";
0291 
0292     // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
0293     runAlgorithm(iEvent, iSetup);
0294     LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";
0295 
0296     // convert our jets and add to the overall jet vector
0297     for (unsigned int ijet = 0; ijet < fjJets_.size(); ++ijet) {
0298       // get the constituents from fastjet
0299       std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
0300       // convert them to CandidatePtr vector
0301       std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
0302       // fill the trackjet
0303       reco::TrackJet jet;
0304       // write the specifics to the jet (simultaneously sets 4-vector, vertex).
0305       writeSpecific(
0306           jet,
0307           reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
0308           vertex_,
0309           constituents);
0310       jet.setJetArea(0);
0311       jet.setPileup(0);
0312       jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int)(itVtx - pvCollection->begin())));
0313       jet.setVertex(itVtx->position());
0314       jets->push_back(jet);
0315     }
0316 
0317     if (useOnlyOnePV_)
0318       break;  // stop vertex loop if only one vertex asked for
0319   }           // end loop over vertices
0320 
0321   // put the jets in the collection
0322   LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
0323   iEvent.put(std::move(jets));
0324 
0325   // Clear the work vectors so that memory is free for other modules.
0326   // Use the trick of swapping with an empty vector so that the memory
0327   // is actually given back rather than silently kept.
0328   decltype(fjInputs_)().swap(fjInputs_);
0329   decltype(fjJets_)().swap(fjJets_);
0330   decltype(inputs_)().swap(inputs_);
0331 }
0332 
0333 //______________________________________________________________________________
0334 void FastjetJetProducer::runAlgorithm(edm::Event& iEvent, edm::EventSetup const& iSetup) {
0335   // run algorithm
0336   /*
0337   fjInputs_.clear();
0338   double px, py , pz, E;
0339   string line;
0340   std::ifstream fin("dump3.txt");
0341   while (getline(fin, line)){
0342     if (line == "#END") break;
0343     if (line.substr(0,1) == "#") {continue;}
0344     istringstream istr(line);
0345     istr >> px >> py >> pz >> E;
0346     // create a fastjet::PseudoJet with these components and put it onto
0347     // back of the input_particles vector
0348     fastjet::PseudoJet j(px,py,pz,E);
0349     //if ( fabs(j.rap()) < inputEtaMax )
0350     fjInputs_.push_back(fastjet::PseudoJet(px,py,pz,E)); 
0351   }
0352   fin.close();
0353   */
0354 
0355   if (!doAreaFastjet_ && !doRhoFastjet_) {
0356     fjClusterSeq_ = std::make_shared<fastjet::ClusterSequence>(fjInputs_, *fjJetDefinition_);
0357   } else if (voronoiRfact_ <= 0) {
0358     fjClusterSeq_ =
0359         ClusterSequencePtr(new fastjet::ClusterSequenceArea(fjInputs_, *fjJetDefinition_, *fjAreaDefinition_));
0360   } else {
0361     fjClusterSeq_ = ClusterSequencePtr(
0362         new fastjet::ClusterSequenceVoronoiArea(fjInputs_, *fjJetDefinition_, fastjet::VoronoiAreaSpec(voronoiRfact_)));
0363   }
0364 
0365   if (!(useMassDropTagger_ || useCMSBoostedTauSeedingAlgorithm_ || useTrimming_ || useFiltering_ || usePruning_ ||
0366         useSoftDrop_ || useConstituentSubtraction_)) {
0367     fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
0368   } else {
0369     fjJets_.clear();
0370 
0371     transformer_coll transformers;
0372 
0373     std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
0374 
0375     unique_ptr<fastjet::JetMedianBackgroundEstimator> bge_rho;
0376     if (useConstituentSubtraction_) {
0377       fastjet::Selector rho_range = fastjet::SelectorAbsRapMax(csRho_EtaMax_);
0378       bge_rho = std::make_unique<fastjet::JetMedianBackgroundEstimator>(
0379           rho_range, fastjet::JetDefinition(fastjet::kt_algorithm, csRParam_), *fjAreaDefinition_);
0380       bge_rho->set_particles(fjInputs_);
0381       fastjet::contrib::ConstituentSubtractor* constituentSubtractor =
0382           new fastjet::contrib::ConstituentSubtractor(bge_rho.get());
0383 
0384       transformers.push_back(transformer_ptr(constituentSubtractor));
0385     };
0386     if (useMassDropTagger_) {
0387       fastjet::MassDropTagger* md_tagger = new fastjet::MassDropTagger(muCut_, yCut_);
0388       transformers.push_back(transformer_ptr(md_tagger));
0389     }
0390     if (useCMSBoostedTauSeedingAlgorithm_) {
0391       fastjet::contrib::CMSBoostedTauSeedingAlgorithm* tau_tagger = new fastjet::contrib::CMSBoostedTauSeedingAlgorithm(
0392           subjetPtMin_, muMin_, muMax_, yMin_, yMax_, dRMin_, dRMax_, maxDepth_, verbosity_);
0393       transformers.push_back(transformer_ptr(tau_tagger));
0394     }
0395     if (useTrimming_) {
0396       fastjet::Filter* trimmer = new fastjet::Filter(fastjet::JetDefinition(fastjet::kt_algorithm, rFilt_),
0397                                                      fastjet::SelectorPtFractionMin(trimPtFracMin_));
0398       transformers.push_back(transformer_ptr(trimmer));
0399     }
0400     if ((useFiltering_) && (!useDynamicFiltering_)) {
0401       fastjet::Filter* filter = new fastjet::Filter(fastjet::JetDefinition(fastjet::cambridge_algorithm, rFilt_),
0402                                                     fastjet::SelectorNHardest(nFilt_));
0403       transformers.push_back(transformer_ptr(filter));
0404     }
0405 
0406     if ((usePruning_) && (!useKtPruning_)) {
0407       fastjet::Pruner* pruner = new fastjet::Pruner(fastjet::cambridge_algorithm, zCut_, RcutFactor_);
0408       transformers.push_back(transformer_ptr(pruner));
0409     }
0410 
0411     if (useDynamicFiltering_) {
0412       fastjet::Filter* filter =
0413           new fastjet::Filter(fastjet::Filter(&*rFiltDynamic_, fastjet::SelectorNHardest(nFilt_)));
0414       transformers.push_back(transformer_ptr(filter));
0415     }
0416 
0417     if (useKtPruning_) {
0418       fastjet::Pruner* pruner = new fastjet::Pruner(fastjet::kt_algorithm, zCut_, RcutFactor_);
0419       transformers.push_back(transformer_ptr(pruner));
0420     }
0421 
0422     if (useSoftDrop_) {
0423       fastjet::contrib::SoftDrop* sd = new fastjet::contrib::SoftDrop(beta_, zCut_, R0_);
0424       transformers.push_back(transformer_ptr(sd));
0425     }
0426 
0427     unique_ptr<fastjet::Subtractor> subtractor;
0428     unique_ptr<fastjet::GridMedianBackgroundEstimator> bge_rho_grid;
0429     if (correctShape_) {
0430       bge_rho_grid = std::make_unique<fastjet::GridMedianBackgroundEstimator>(gridMaxRapidity_, gridSpacing_);
0431       bge_rho_grid->set_particles(fjInputs_);
0432       subtractor = std::make_unique<fastjet::Subtractor>(bge_rho_grid.get());
0433       subtractor->set_use_rho_m();
0434       //subtractor->use_common_bge_for_rho_and_rhom(true);
0435     }
0436 
0437     for (std::vector<fastjet::PseudoJet>::const_iterator ijet = tempJets.begin(), ijetEnd = tempJets.end();
0438          ijet != ijetEnd;
0439          ++ijet) {
0440       fastjet::PseudoJet transformedJet = *ijet;
0441       bool passed = true;
0442       for (transformer_coll::const_iterator itransf = transformers.begin(), itransfEnd = transformers.end();
0443            itransf != itransfEnd;
0444            ++itransf) {
0445         if (transformedJet != 0) {
0446           transformedJet = (**itransf)(transformedJet);
0447         } else {
0448           passed = false;
0449         }
0450       }
0451 
0452       if (correctShape_) {
0453         transformedJet = (*subtractor)(transformedJet);
0454       }
0455 
0456       if (passed) {
0457         fjJets_.push_back(transformedJet);
0458       }
0459     }
0460   }
0461 }
0462 
0463 //______________________________________________________________________________
0464 void FastjetJetProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0465   edm::ParameterSetDescription descFastjetJetProducer;
0466   ////// From FastjetJetProducer
0467   fillDescriptionsFromFastJetProducer(descFastjetJetProducer);
0468   ///// From VirtualJetProducer
0469   VirtualJetProducer::fillDescriptionsFromVirtualJetProducer(descFastjetJetProducer);
0470   ////
0471   descFastjetJetProducer.add<string>("jetCollInstanceName", "");
0472   descFastjetJetProducer.add<bool>("sumRecHits", false);
0473 
0474   /////////////////////
0475   descriptions.add("FastjetJetProducer", descFastjetJetProducer);
0476 }
0477 
0478 void FastjetJetProducer::fillDescriptionsFromFastJetProducer(edm::ParameterSetDescription& desc) {
0479   desc.add<bool>("useMassDropTagger", false);
0480   desc.add<bool>("useFiltering", false);
0481   desc.add<bool>("useDynamicFiltering", false);
0482   desc.add<bool>("useTrimming", false);
0483   desc.add<bool>("usePruning", false);
0484   desc.add<bool>("useCMSBoostedTauSeedingAlgorithm", false);
0485   desc.add<bool>("useKtPruning", false);
0486   desc.add<bool>("useConstituentSubtraction", false);
0487   desc.add<bool>("useSoftDrop", false);
0488   desc.add<bool>("correctShape", false);
0489   desc.add<bool>("UseOnlyVertexTracks", false);
0490   desc.add<bool>("UseOnlyOnePV", false);
0491   desc.add<double>("muCut", -1.0);
0492   desc.add<double>("yCut", -1.0);
0493   desc.add<double>("rFilt", -1.0);
0494   desc.add<double>("rFiltFactor", -1.0);
0495   desc.add<double>("trimPtFracMin", -1.0);
0496   desc.add<double>("zcut", -1.0);
0497   desc.add<double>("rcut_factor", -1.0);
0498   desc.add<double>("csRho_EtaMax", -1.0);
0499   desc.add<double>("csRParam", -1.0);
0500   desc.add<double>("beta", -1.0);
0501   desc.add<double>("R0", -1.0);
0502   desc.add<double>("gridMaxRapidity", -1.0);  // For fixed-grid rho
0503   desc.add<double>("gridSpacing", -1.0);      // For fixed-grid rho
0504   desc.add<double>("DzTrVtxMax", 999999.);
0505   desc.add<double>("DxyTrVtxMax", 999999.);
0506   desc.add<double>("MaxVtxZ", 15.0);
0507   desc.add<double>("subjetPtMin", -1.0);
0508   desc.add<double>("muMin", -1.0);
0509   desc.add<double>("muMax", -1.0);
0510   desc.add<double>("yMin", -1.0);
0511   desc.add<double>("yMax", -1.0);
0512   desc.add<double>("dRMin", -1.0);
0513   desc.add<double>("dRMax", -1.0);
0514   desc.add<int>("maxDepth", -1);
0515   desc.add<int>("nFilt", -1);
0516   desc.add<int>("MinVtxNdof", 5);
0517 }
0518 
0519 ////////////////////////////////////////////////////////////////////////////////
0520 // define as cmssw plugin
0521 ////////////////////////////////////////////////////////////////////////////////
0522 
0523 DEFINE_FWK_MODULE(FastjetJetProducer);