Fork me on GitHub

Changeset 8269ee2 in git for modules


Ignore:
Timestamp:
Nov 29, 2017, 8:04:09 PM (7 years ago)
Author:
GitHub <noreply@…>
Branches:
ImprovedOutputFile, Timing, dual_readout, llp, master
Children:
1d2d44b, dda357d
Parents:
5107603 (diff), 67b86c2 (diff)
Note: this is a merge changeset, the changes displayed below correspond to the merge itself.
Use the (diff) links above to see all the changes relative to each parent.
git-author:
Michele Selvaggi <michele.selvaggi@…> (11/29/17 20:04:09)
git-committer:
GitHub <noreply@…> (11/29/17 20:04:09)
Message:

Merge pull request #45 from uschnoor/VLC_jetclustering

Thanks!

Location:
modules
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • modules/FastJetFinder.cc

    r5107603 r8269ee2  
    6666#include "fastjet/contribs/Nsubjettiness/ExtraRecombiners.hh"
    6767
     68#include "fastjet/contribs/ValenciaPlugin/ValenciaPlugin.hh"
     69
    6870#include "fastjet/tools/Filter.hh"
    6971#include "fastjet/tools/Pruner.hh"
     
    7981FastJetFinder::FastJetFinder() :
    8082  fPlugin(0), fRecomb(0), fAxesDef(0), fMeasureDef(0), fNjettinessPlugin(0),
    81   fDefinition(0), fAreaDefinition(0), fItInputArray(0)
     83  fDefinition(0), fAreaDefinition(0), fItInputArray(0), fValenciaPlugin(0)
    8284{
    8385
     
    118120  fJetPTMin = GetDouble("JetPTMin", 10.0);
    119121
     122 
    120123  //-- N(sub)jettiness parameters --
    121124
     
    125128  fRcutOff = GetDouble("RcutOff", 0.8); // used only if Njettiness is used as jet clustering algo (case 8)
    126129  fN = GetInt("N", 2);                  // used only if Njettiness is used as jet clustering algo (case 8)
    127      
     130
     131  //-- Exclusive clustering for e+e- collisions --
     132 
     133  fNJets = GetInt("NJets",2);
     134  fExclusiveClustering = GetBool("ExclusiveClustering", false);
     135
     136  //-- Valencia Linear Collider algorithm
     137  fGamma = GetDouble("Gamma", 1.0);
     138  //fBeta parameter see above
     139 
    128140  fMeasureDef = new NormalizedMeasure(fBeta, fParameterR);
    129141   
     
    236248      fDefinition = new JetDefinition(fNjettinessPlugin);
    237249      break;
     250  case 9:
     251      fValenciaPlugin = new ValenciaPlugin(fParameterR, fBeta, fGamma);
     252      fDefinition = new JetDefinition(fValenciaPlugin);
     253      break;
     254
    238255  }
    239256
     
    294311  if(fAxesDef) delete fAxesDef;
    295312  if(fMeasureDef) delete fMeasureDef;
     313  if(fValenciaPlugin) delete static_cast<JetDefinition::Plugin*>(fValenciaPlugin);
     314
    296315}
    297316
     
    357376
    358377  outputList.clear();
    359   outputList = sorted_by_pt(sequence->inclusive_jets(fJetPTMin));
    360 
     378
     379  if(fExclusiveClustering)
     380    {
     381      outputList = sorted_by_pt(sequence->exclusive_jets( fNJets ));
     382    }
     383  else
     384    {
     385      outputList = sorted_by_pt(sequence->inclusive_jets(fJetPTMin));
     386    }
    361387
    362388  // loop over all jets and export them
  • modules/FastJetFinder.h

    r5107603 r8269ee2  
    4141  namespace contrib {
    4242    class NjettinessPlugin;
     43    class ValenciaPlugin;
    4344    class AxesDefinition;
    4445    class MeasureDefinition;   
     
    6162  void *fPlugin; //!
    6263  void *fRecomb; //!
     64
     65  fastjet::contrib::AxesDefinition *fAxesDef;
     66  fastjet::contrib::MeasureDefinition *fMeasureDef;   
     67
    6368  fastjet::contrib::NjettinessPlugin *fNjettinessPlugin; //!
    64 
     69  fastjet::contrib::ValenciaPlugin *fValenciaPlugin; //!
    6570  fastjet::JetDefinition *fDefinition; //!
    6671
    6772  Int_t fJetAlgorithm;
    6873  Double_t fParameterR;
     74
    6975  Double_t fJetPTMin;
    7076  Double_t fConeRadius;
     
    7783  Double_t fOverlapThreshold;
    7884
     85  //-- Exclusive clustering for e+e- collisions --
     86
     87  Int_t fNJets;
     88  Bool_t fExclusiveClustering;
     89
     90  //-- Valencia Linear Collider algorithm
     91  Double_t fGamma;
     92 
    7993  //-- N (sub)jettiness parameters --
    8094
    8195  Bool_t fComputeNsubjettiness;
    82   fastjet::contrib::AxesDefinition *fAxesDef;
    83   fastjet::contrib::MeasureDefinition *fMeasureDef;   
    8496  Double_t fBeta;
    8597  Int_t fAxisMode;
Note: See TracChangeset for help on using the changeset viewer.