Fork me on GitHub

Changeset e4c3fef in git for modules/FastJetFinder.cc


Ignore:
Timestamp:
Apr 17, 2014, 10:54:17 AM (11 years ago)
Author:
pavel <pavel@…>
Branches:
ImprovedOutputFile, Timing, dual_readout, llp, master
Children:
b96d99b
Parents:
2e6a81b
Message:

several minor corrections

File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/FastJetFinder.cc

    r2e6a81b re4c3fef  
    4747#include "fastjet/plugins/CDFCones/fastjet/CDFJetCluPlugin.hh"
    4848
    49 #include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh" 
     49#include "fastjet/contribs/Nsubjettiness/Nsubjettiness.hh"
    5050#include "fastjet/contribs/Nsubjettiness/Njettiness.hh"
    5151#include "fastjet/contribs/Nsubjettiness/NjettinessPlugin.hh"
     
    7676void FastJetFinder::Init()
    7777{
    78  
     78
    7979  JetDefinition::Plugin *plugin = NULL;
    8080  JetDefinition::Recombiner *recomb = NULL;
    8181  NjettinessPlugin *njet_plugin = NULL;
    82  
     82
    8383  // read eta ranges
    8484
     
    110110
    111111  //-- N(sub)jettiness parameters --
    112  
     112
    113113  fComputeNsubjettiness = GetBool("ComputeNsubjettiness", false);
    114114  fBeta = GetDouble("Beta", 1.0);
    115115  fAxisMode = GetInt("AxisMode", 1);
    116   fRcutOff = GetDouble("RcutOff", 0.8); //used only if Njettiness is used as jet clustering algo (case 8)
    117   fN = GetInt("N", 2);                  //used only if Njettiness is used as jet clustering algo (case 8)
    118  
     116  fRcutOff = GetDouble("RcutOff", 0.8); // used only if Njettiness is used as jet clustering algo (case 8)
     117  fN = GetInt("N", 2);                  // used only if Njettiness is used as jet clustering algo (case 8)
     118
    119119  // ---  Jet Area Parameters ---
    120120  fAreaAlgorithm = GetInt("AreaAlgorithm", 0);
    121121  fComputeRho = GetBool("ComputeRho", false);
    122  
     122
    123123  // - ghost based areas -
    124124  fGhostEtaMax = GetDouble("GhostEtaMax", 5.0);
     
    128128  fPtScatter = GetDouble("PtScatter", 0.1);
    129129  fMeanGhostPt = GetDouble("MeanGhostPt", 1.0E-100);
    130  
     130
    131131  // - voronoi based areas -
    132132  fEffectiveRfact = GetDouble("EffectiveRfact", 1.0);
     
    188188      break;
    189189  }
    190  
    191190
    192191  fPlugin = plugin;
    193192  fRecomb = recomb;
    194193  fNjettinessPlugin = njet_plugin;
    195  
     194
    196195  ClusterSequence::print_banner();
    197196
     
    299298    inputList.clear();
    300299    inputList = sequence->constituents(*itOutputList);
    301      
     300
    302301    for(itInputList = inputList.begin(); itInputList != inputList.end(); ++itInputList)
    303302    {
     
    308307      if(deta > detaMax) detaMax = deta;
    309308      if(dphi > dphiMax) dphiMax = dphi;
    310      
     309
    311310      time += TMath::Sqrt(constituent->Momentum.E())*(constituent->Position.T());
    312311      weightTime += TMath::Sqrt(constituent->Momentum.E());
    313    
     312
    314313      candidate->AddCandidate(constituent);
    315314    }
    316    
     315
    317316    avTime = time/weightTime;
    318317
     
    325324
    326325    // --- compute N-subjettiness with N = 1,2,3,4,5 ----
    327    
     326
    328327    if(fComputeNsubjettiness)
    329328    {
    330329      Njettiness::AxesMode axisMode;
    331      
    332       if (fAxisMode == 1) axisMode = Njettiness::wta_kt_axes;
    333       if (fAxisMode == 2) axisMode = Njettiness::onepass_wta_kt_axes;
    334       if (fAxisMode == 3) axisMode = Njettiness::kt_axes;
    335       if (fAxisMode == 4) axisMode = Njettiness::onepass_kt_axes;
    336        
     330
     331      switch(fAxisMode)
     332      {
     333        default:
     334        case 1:
     335          axisMode = Njettiness::wta_kt_axes;
     336          break;
     337        case 2:
     338          axisMode = Njettiness::onepass_wta_kt_axes;
     339          break;
     340        case 3:
     341          axisMode = Njettiness::kt_axes;
     342          break;
     343        case 4:
     344          axisMode = Njettiness::onepass_kt_axes;
     345          break;
     346      }
     347
    337348      Njettiness::MeasureMode measureMode = Njettiness::unnormalized_measure;
    338      
     349
    339350      Nsubjettiness nSub1(1, axisMode, measureMode, fBeta);
    340351      Nsubjettiness nSub2(2, axisMode, measureMode, fBeta);
    341       Nsubjettiness nSub3(3, axisMode, measureMode, fBeta); 
    342       Nsubjettiness nSub4(4, axisMode, measureMode, fBeta); 
    343       Nsubjettiness nSub5(5, axisMode, measureMode, fBeta); 
    344    
    345       candidate -> Tau1 = nSub1(*itOutputList);
    346       candidate -> Tau2 = nSub2(*itOutputList);
    347       candidate -> Tau3 = nSub3(*itOutputList);
    348       candidate -> Tau4 = nSub4(*itOutputList);
    349       candidate -> Tau5 = nSub5(*itOutputList);
     352      Nsubjettiness nSub3(3, axisMode, measureMode, fBeta);
     353      Nsubjettiness nSub4(4, axisMode, measureMode, fBeta);
     354      Nsubjettiness nSub5(5, axisMode, measureMode, fBeta);
     355
     356      candidate->Tau[0] = nSub1(*itOutputList);
     357      candidate->Tau[1] = nSub2(*itOutputList);
     358      candidate->Tau[2] = nSub3(*itOutputList);
     359      candidate->Tau[3] = nSub4(*itOutputList);
     360      candidate->Tau[4] = nSub5(*itOutputList);
    350361    }
    351362
Note: See TracChangeset for help on using the changeset viewer.