Fork me on GitHub

Changeset 1383 in svn


Ignore:
Timestamp:
Apr 17, 2014, 12:28:09 PM (11 years ago)
Author:
Pavel Demin
Message:

fix namespace usage

Location:
trunk/modules
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • trunk/modules/FastJetFinder.cc

    r1372 r1383  
    6060
    6161FastJetFinder::FastJetFinder() :
    62   fPlugin(0), fDefinition(0), fAreaDefinition(0), fItInputArray(0), fRecomb(0), fNjettinessPlugin(0)
     62  fPlugin(0), fRecomb(0), fNjettinessPlugin(0), fDefinition(0), fAreaDefinition(0), fItInputArray(0)
    6363{
    6464
     
    7676void FastJetFinder::Init()
    7777{
    78 
    79   JetDefinition::Plugin *plugin = NULL;
    80   JetDefinition::Recombiner *recomb = NULL;
    81   NjettinessPlugin *njet_plugin = NULL;
     78  JetDefinition::Plugin *plugin = 0;
     79  JetDefinition::Recombiner *recomb = 0;
     80  NjettinessPlugin *njetPlugin = 0;
    8281
    8382  // read eta ranges
     
    135134  {
    136135    case 1:
    137       fAreaDefinition = new fastjet::AreaDefinition(active_area_explicit_ghosts, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
     136      fAreaDefinition = new AreaDefinition(active_area_explicit_ghosts, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
    138137      break;
    139138    case 2:
    140       fAreaDefinition = new fastjet::AreaDefinition(one_ghost_passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
     139      fAreaDefinition = new AreaDefinition(one_ghost_passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
    141140      break;
    142141    case 3:
    143       fAreaDefinition = new fastjet::AreaDefinition(passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
     142      fAreaDefinition = new AreaDefinition(passive_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
    144143      break;
    145144    case 4:
    146       fAreaDefinition = new fastjet::AreaDefinition(VoronoiAreaSpec(fEffectiveRfact));
     145      fAreaDefinition = new AreaDefinition(VoronoiAreaSpec(fEffectiveRfact));
    147146      break;
    148147    case 5:
    149       fAreaDefinition = new fastjet::AreaDefinition(active_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
     148      fAreaDefinition = new AreaDefinition(active_area, GhostedAreaSpec(fGhostEtaMax, fRepeat, fGhostArea, fGridScatter, fPtScatter, fMeanGhostPt));
    150149      break;
    151150    default:
     
    158157  {
    159158    case 1:
    160       plugin = new fastjet::CDFJetCluPlugin(fSeedThreshold, fConeRadius, fAdjacencyCut, fMaxIterations, fIratch, fOverlapThreshold);
    161       fDefinition = new fastjet::JetDefinition(plugin);
     159      plugin = new CDFJetCluPlugin(fSeedThreshold, fConeRadius, fAdjacencyCut, fMaxIterations, fIratch, fOverlapThreshold);
     160      fDefinition = new JetDefinition(plugin);
    162161      break;
    163162    case 2:
    164       plugin = new fastjet::CDFMidPointPlugin(fSeedThreshold, fConeRadius, fConeAreaFraction, fMaxPairSize, fMaxIterations, fOverlapThreshold);
    165       fDefinition = new fastjet::JetDefinition(plugin);
     163      plugin = new CDFMidPointPlugin(fSeedThreshold, fConeRadius, fConeAreaFraction, fMaxPairSize, fMaxIterations, fOverlapThreshold);
     164      fDefinition = new JetDefinition(plugin);
    166165      break;
    167166    case 3:
    168       plugin = new fastjet::SISConePlugin(fConeRadius, fOverlapThreshold, fMaxIterations, fJetPTMin);
    169       fDefinition = new fastjet::JetDefinition(plugin);
     167      plugin = new SISConePlugin(fConeRadius, fOverlapThreshold, fMaxIterations, fJetPTMin);
     168      fDefinition = new JetDefinition(plugin);
    170169      break;
    171170    case 4:
    172       fDefinition = new fastjet::JetDefinition(fastjet::kt_algorithm, fParameterR);
     171      fDefinition = new JetDefinition(kt_algorithm, fParameterR);
    173172      break;
    174173    case 5:
    175       fDefinition = new fastjet::JetDefinition(fastjet::cambridge_algorithm, fParameterR);
     174      fDefinition = new JetDefinition(cambridge_algorithm, fParameterR);
    176175      break;
    177176    default:
    178177    case 6:
    179       fDefinition = new fastjet::JetDefinition(fastjet::antikt_algorithm, fParameterR);
     178      fDefinition = new JetDefinition(antikt_algorithm, fParameterR);
    180179      break;
    181180    case 7:
    182       recomb = new fastjet::contrib::WinnerTakeAllRecombiner();
    183       fDefinition = new fastjet::JetDefinition(fastjet::antikt_algorithm, fParameterR, recomb, Best);
     181      recomb = new WinnerTakeAllRecombiner();
     182      fDefinition = new JetDefinition(antikt_algorithm, fParameterR, recomb, Best);
    184183      break;
    185184    case 8:
    186       njet_plugin = new fastjet::contrib::NjettinessPlugin(fN, Njettiness::wta_kt_axes, Njettiness::unnormalized_cutoff_measure, fBeta, fRcutOff);
    187       fDefinition = new fastjet::JetDefinition(njet_plugin);
     185      njetPlugin = new NjettinessPlugin(fN, Njettiness::wta_kt_axes, Njettiness::unnormalized_cutoff_measure, fBeta, fRcutOff);
     186      fDefinition = new JetDefinition(njetPlugin);
    188187      break;
    189188  }
     
    191190  fPlugin = plugin;
    192191  fRecomb = recomb;
    193   fNjettinessPlugin = njet_plugin;
     192  fNjettinessPlugin = njetPlugin;
    194193
    195194  ClusterSequence::print_banner();
  • trunk/modules/FastJetFinder.h

    r1372 r1383  
    2525  class AreaDefinition;
    2626  class Selector;
     27  namespace contrib {
     28    class NjettinessPlugin;
     29  }
    2730}
    2831
     
    4245  void *fPlugin; //!
    4346  void *fRecomb; //!
    44   void *fNjettinessPlugin; //!
     47  fastjet::contrib::NjettinessPlugin *fNjettinessPlugin; //!
    4548
    4649  fastjet::JetDefinition *fDefinition; //!
Note: See TracChangeset for help on using the changeset viewer.