AliPhysics  775474e (775474e)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliFJWrapper.h
Go to the documentation of this file.
1 #ifndef AliFJWrapper_H
2 #define AliFJWrapper_H
3 
4 #if !defined(__CINT__)
5 
6 #include <vector>
7 #include <TString.h>
8 #include "AliLog.h"
9 #include "FJ_includes.h"
10 #include "AliJetShape.h"
11 
12 
14 {
15  public:
16  AliFJWrapper(const char *name, const char *title);
17  virtual ~AliFJWrapper();
18 
19  virtual void AddInputVector (Double_t px, Double_t py, Double_t pz, Double_t E, Int_t index = -99999);
20  virtual void AddInputVector (const fastjet::PseudoJet& vec, Int_t index = -99999);
21  virtual void AddInputVectors(const std::vector<fastjet::PseudoJet>& vecs, Int_t offsetIndex = -99999);
22  virtual void AddInputGhost (Double_t px, Double_t py, Double_t pz, Double_t E, Int_t index = -99999);
23  virtual const char *ClassName() const { return "AliFJWrapper"; }
24  virtual void Clear(const Option_t* /*opt*/ = "");
25  virtual void ClearMemory();
26  virtual void CopySettingsFrom (const AliFJWrapper& wrapper);
27  virtual void GetMedianAndSigma(Double_t& median, Double_t& sigma, Int_t remove = 0) const;
28  fastjet::ClusterSequenceArea* GetClusterSequence() const { return fClustSeq; }
29  fastjet::ClusterSequence* GetClusterSequenceSA() const { return fClustSeqSA; }
30  fastjet::ClusterSequenceActiveAreaExplicitGhosts* GetClusterSequenceGhosts() const { return fClustSeqActGhosts; }
31  const std::vector<fastjet::PseudoJet>& GetInputVectors() const { return fInputVectors; }
32  const std::vector<fastjet::PseudoJet>& GetInputGhosts() const { return fInputGhosts; }
33  const std::vector<fastjet::PseudoJet>& GetInclusiveJets() const { return fInclusiveJets; }
34  const std::vector<fastjet::PseudoJet>& GetFilteredJets() const { return fFilteredJets; }
35  std::vector<fastjet::PseudoJet> GetJetConstituents(UInt_t idx) const;
36  std::vector<fastjet::PseudoJet> GetFilteredJetConstituents(UInt_t idx) const;
38  const char* GetName() const { return fName; }
39  const char* GetTitle() const { return fTitle; }
40  Double_t GetJetArea (UInt_t idx) const;
41  fastjet::PseudoJet GetJetAreaVector (UInt_t idx) const;
43  fastjet::PseudoJet GetFilteredJetAreaVector(UInt_t idx) const;
45  virtual std::vector<double> GetSubtractedJetsPts(Double_t median_pt = -1, Bool_t sorted = kFALSE);
48  Double_t NSubjettiness(Int_t N, Int_t Algorithm, Double_t Radius, Double_t Beta, Int_t Option=0, Int_t Measure=0, Double_t Beta_SD=0, Double_t ZCut=0.1);
49  Double32_t NSubjettinessDerivativeSub(Int_t N, Int_t Algorithm, Double_t Radius, Double_t Beta, Double_t JetR, fastjet::PseudoJet jet, Int_t Option=0, Int_t Measure=0);
50 #ifdef FASTJET_VERSION
51  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetMass() const {return fGenSubtractorInfoJetMass ; }
52  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetAngularity() const {return fGenSubtractorInfoJetAngularity ; }
53  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetpTD() const {return fGenSubtractorInfoJetpTD ; }
54  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetCircularity() const {return fGenSubtractorInfoJetCircularity ; }
55  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetSigma2() const {return fGenSubtractorInfoJetSigma2 ; }
56  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetConstituent() const {return fGenSubtractorInfoJetConstituent ; }
57  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetLeSub() const {return fGenSubtractorInfoJetLeSub ; }
58  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJet1subjettiness_kt() const {return fGenSubtractorInfoJet1subjettiness_kt ; }
59  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJet2subjettiness_kt() const {return fGenSubtractorInfoJet2subjettiness_kt ; }
60  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJet3subjettiness_kt() const {return fGenSubtractorInfoJet3subjettiness_kt ; }
61  const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetOpeningAngle_kt() const {return fGenSubtractorInfoJetOpeningAngle_kt ; }
62  const std::vector<fastjet::PseudoJet> GetConstituentSubtrJets() const {return fConstituentSubtrJets ; }
63  const std::vector<fastjet::PseudoJet> GetGroomedJets() const {return fGroomedJets ; }
64  Int_t CreateGenSub(); // fastjet::contrib::GenericSubtractor
65  Int_t CreateConstituentSub(); // fastjet::contrib::ConstituentSubtractor
66  Int_t CreateSoftDrop();
67 #endif
68  virtual std::vector<double> GetGRNumerator() const { return fGRNumerator ; }
69  virtual std::vector<double> GetGRDenominator() const { return fGRDenominator ; }
70  virtual std::vector<double> GetGRNumeratorSub() const { return fGRNumeratorSub ; }
71  virtual std::vector<double> GetGRDenominatorSub() const { return fGRDenominatorSub ; }
72 
73  virtual void RemoveLastInputVector();
74 
75  virtual Int_t Run();
76  virtual Int_t Filter();
78  virtual Int_t DoGenericSubtractionGR(Int_t ijet);
90  virtual Int_t DoSoftDrop();
91 
92  void SetName(const char* name) { fName = name; }
93  void SetTitle(const char* title) { fTitle = title; }
94  void SetStrategy(const fastjet::Strategy &strat) { fStrategy = strat; }
95  void SetAlgorithm(const fastjet::JetAlgorithm &algor) { fAlgor = algor; }
96  void SetRecombScheme(const fastjet::RecombinationScheme &scheme) { fScheme = scheme; }
97  void SetAreaType(const fastjet::AreaType &atype) { fAreaType = atype; }
98  void SetNRepeats(Int_t nrepeat) { fNGhostRepeats = nrepeat; }
99  void SetGhostArea(Double_t gharea) { fGhostArea = gharea; }
100  void SetMaxRap(Double_t maxrap) { fMaxRap = maxrap; }
101  void SetR(Double_t r) { fR = r; }
102  void SetGridScatter(Double_t gridSc) { fGridScatter = gridSc; }
103  void SetKtScatter(Double_t ktSc) { fKtScatter = ktSc; }
104  void SetMeanGhostKt(Double_t meankt) { fMeanGhostKt = meankt; }
105  void SetPluginAlgor(Int_t plugin) { fPluginAlgor = plugin; }
106  void SetUseArea4Vector(Bool_t useA4v) { fUseArea4Vector = useA4v; }
107  void SetupAlgorithmfromOpt(const char *option);
108  void SetupAreaTypefromOpt(const char *option);
109  void SetupSchemefromOpt(const char *option);
110  void SetupStrategyfromOpt(const char *option);
112  void SetLegacyFJ();
113  void SetUseExternalBkg(Bool_t b, Double_t rho, Double_t rhom) { fUseExternalBkg = b; fRho = rho; fRhom = rhom;}
114  void SetRMaxAndStep(Double_t rmax, Double_t dr) {fRMax = rmax; fDRStep = dr; }
115  void SetRhoRhom (Double_t rho, Double_t rhom) { fUseExternalBkg = kTRUE; fRho = rho; fRhom = rhom;} // if using rho,rhom then fUseExternalBkg is true
116  void SetMinJetPt(Double_t MinPt) {fMinJetPt=MinPt;}
117 
118  protected:
121  std::vector<fastjet::PseudoJet> fInputVectors;
122  std::vector<fastjet::PseudoJet> fInputGhosts;
123  std::vector<fastjet::PseudoJet> fInclusiveJets;
124  std::vector<fastjet::PseudoJet> fFilteredJets;
125  std::vector<double> fSubtractedJetsPt;
126  std::vector<fastjet::PseudoJet> fConstituentSubtrJets;
127  std::vector<fastjet::PseudoJet> fGroomedJets;
128  fastjet::AreaDefinition *fAreaDef;
129  fastjet::VoronoiAreaSpec *fVorAreaSpec;
130  fastjet::GhostedAreaSpec *fGhostedAreaSpec;
131  fastjet::JetDefinition *fJetDef;
132  fastjet::JetDefinition::Plugin *fPlugin;
133 #ifndef FASTJET_VERSION
134  fastjet::RangeDefinition *fRange;
135 #else
136  fastjet::Selector *fRange;
137 #endif
138  fastjet::ClusterSequenceArea *fClustSeq;
139  fastjet::ClusterSequence *fClustSeqSA;
140  fastjet::ClusterSequenceActiveAreaExplicitGhosts *fClustSeqActGhosts;
141  fastjet::Strategy fStrategy;
142  fastjet::JetAlgorithm fAlgor;
143  fastjet::RecombinationScheme fScheme;
144  fastjet::AreaType fAreaType;
150  // no setters for the moment - used default values in the constructor
155  // extra parameters
158  // condition to stop the grooming (rejection of soft splitting) z > fZcut theta^fBeta
159  Double_t fZcut; // fZcut = 0.1
160  Double_t fBeta; // fBeta = 0
161 #ifdef FASTJET_VERSION
162  fastjet::JetMedianBackgroundEstimator *fBkrdEstimator;
163  //from contrib package
164  fastjet::contrib::GenericSubtractor *fGenSubtractor;
165  fastjet::contrib::ConstituentSubtractor *fConstituentSubtractor;
166  fastjet::contrib::SoftDrop *fSoftDrop;
167  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetMass;
168  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoGRNum;
169  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoGRDen;
170  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetAngularity;
171  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetpTD;
172  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetCircularity;
173  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetSigma2;
174  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetConstituent;
175  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetLeSub;
176  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJet1subjettiness_kt;
177  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJet2subjettiness_kt;
178  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJet3subjettiness_kt;
179  std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetOpeningAngle_kt;
180 #endif
184  Double_t fRho; // pT background density
185  Double_t fRhom; // mT background density
188  std::vector<double> fGRNumerator;
189  std::vector<double> fGRDenominator;
190  std::vector<double> fGRNumeratorSub;
191  std::vector<double> fGRDenominatorSub;
192 
193  virtual void SubtractBackground(const Double_t median_pt = -1);
194 
195  private:
196  AliFJWrapper();
197  AliFJWrapper(const AliFJWrapper& wrapper);
198  AliFJWrapper& operator = (const AliFJWrapper& wrapper);
199 };
200 #endif
201 #endif
202 
203 #ifdef AliFJWrapper_CXX
204 #undef AliFJWrapper_CXX
205 
206 #if defined __GNUC__
207 #pragma GCC system_header
208 #endif
209 
210 namespace fj = fastjet;
211 
212 //_________________________________________________________________________________________________
213 AliFJWrapper::AliFJWrapper(const char *name, const char *title)
214  :
215  fName (name)
216  , fTitle (title)
217  , fInputVectors ( )
218  , fInputGhosts ( )
219  , fInclusiveJets ( )
220  , fFilteredJets ( )
221  , fSubtractedJetsPt ( )
222  , fConstituentSubtrJets ( )
223  , fSoftDrop ( )
224  , fAreaDef (0)
225  , fVorAreaSpec (0)
226  , fGhostedAreaSpec (0)
227  , fJetDef (0)
228  , fPlugin (0)
229  , fRange (0)
230  , fClustSeq (0)
231  , fClustSeqSA (0)
232  , fClustSeqActGhosts (0)
233  , fStrategy (fj::Best)
234  , fAlgor (fj::kt_algorithm)
235  , fScheme (fj::BIpt_scheme)
236  , fAreaType (fj::active_area)
237  , fNGhostRepeats (1)
238  , fGhostArea (0.005)
239  , fMaxRap (1.)
240  , fR (0.4)
241  , fGridScatter (1.0)
242  , fKtScatter (0.1)
243  , fMeanGhostKt (1e-100)
244  , fPluginAlgor (0)
245  , fMedUsedForBgSub (0)
246  , fUseArea4Vector (kFALSE)
247  , fZcut(0.1)
248  , fBeta(0)
249 #ifdef FASTJET_VERSION
250  , fBkrdEstimator (0)
251  , fGenSubtractor (0)
252  , fConstituentSubtractor (0)
253  , fGenSubtractorInfoJetMass ( )
254  , fGenSubtractorInfoGRNum ( )
255  , fGenSubtractorInfoGRDen ( )
256  , fGenSubtractorInfoJetAngularity ( )
257  , fGenSubtractorInfoJetpTD ( )
258  , fGenSubtractorInfoJetCircularity( )
259  , fGenSubtractorInfoJetSigma2()
260  , fGenSubtractorInfoJetConstituent ( )
261  , fGenSubtractorInfoJetLeSub ( )
262  , fGenSubtractorInfoJet1subjettiness_kt ( )
263  , fGenSubtractorInfoJet2subjettiness_kt ( )
264  , fGenSubtractorInfoJet3subjettiness_kt ( )
265  , fGenSubtractorInfoJetOpeningAngle_kt ( )
266 #endif
267  , fDoFilterArea (false)
268  , fLegacyMode (false)
269  , fUseExternalBkg (false)
270  , fRho (0)
271  , fRhom (0)
272  , fRMax(2.)
273  , fDRStep(0.04)
274  , fGRNumerator()
275  , fGRDenominator()
276  , fGRNumeratorSub()
277  , fGRDenominatorSub()
278 {
279  // Constructor.
280 }
281 
282 //_________________________________________________________________________________________________
284 {
285  // Destructor.
286  ClearMemory();
287 }
288 
289 //_________________________________________________________________________________________________
291 {
292  // Destructor.
293  if (fAreaDef) { delete fAreaDef; fAreaDef = NULL; }
294  if (fVorAreaSpec) { delete fVorAreaSpec; fVorAreaSpec = NULL; }
295  if (fGhostedAreaSpec) { delete fGhostedAreaSpec; fGhostedAreaSpec = NULL; }
296  if (fJetDef) { delete fJetDef; fJetDef = NULL; }
297  if (fPlugin) { delete fPlugin; fPlugin = NULL; }
298  if (fRange) { delete fRange; fRange = NULL; }
299  if (fClustSeq) { delete fClustSeq; fClustSeq = NULL; }
300  if (fClustSeqSA) { delete fClustSeqSA; fClustSeqSA = NULL; }
302  #ifdef FASTJET_VERSION
303  if (fBkrdEstimator) { delete fBkrdEstimator; fBkrdEstimator = NULL; }
304  if (fGenSubtractor) { delete fGenSubtractor; fGenSubtractor = NULL; }
305  if (fConstituentSubtractor) { delete fConstituentSubtractor; fConstituentSubtractor = NULL; }
306  if (fSoftDrop) { delete fSoftDrop; fSoftDrop = NULL;}
307  #endif
308 }
309 
310 //_________________________________________________________________________________________________
311 void AliFJWrapper::CopySettingsFrom(const AliFJWrapper& wrapper)
312 {
313  // Copy some settings.
314  // You very often want to keep most of the settings
315  // but change only the algorithm or R - do it after call to this function
316 
317  fStrategy = wrapper.fStrategy;
318  fAlgor = wrapper.fAlgor;
319  fScheme = wrapper.fScheme;
320  fAreaType = wrapper.fAreaType;
321  fNGhostRepeats = wrapper.fNGhostRepeats;
322  fGhostArea = wrapper.fGhostArea;
323  fMaxRap = wrapper.fMaxRap;
324  fR = wrapper.fR;
325  fGridScatter = wrapper.fGridScatter;
326  fKtScatter = wrapper.fKtScatter;
327  fMeanGhostKt = wrapper.fMeanGhostKt;
328  fPluginAlgor = wrapper.fPluginAlgor;
330  fZcut = wrapper.fZcut;
331  fBeta = wrapper.fBeta;
332  fLegacyMode = wrapper.fLegacyMode;
334  fRho = wrapper.fRho;
335  fRhom = wrapper.fRhom;
336 }
337 
338 //_________________________________________________________________________________________________
339 void AliFJWrapper::Clear(const Option_t */*opt*/)
340 {
341  // Simply clear the input vectors.
342  // Make sure done on every event if the instance is reused
343  // Reset the median to zero.
344 
345  fInputVectors.clear();
346  fInputGhosts.clear();
347  fMedUsedForBgSub = 0;
348 
349  // for the moment brute force delete everything
350  ClearMemory();
351 }
352 
353 //_________________________________________________________________________________________________
355 {
356  // Remove last input vector
357 
358  fInputVectors.pop_back();
359 }
360 
361 //_________________________________________________________________________________________________
363 {
364  // Make the input pseudojet.
365 
366  fastjet::PseudoJet inVec(px, py, pz, E);
367 
368  // Salvatore Aiola: not sure why this was done...
369  //if (index > -99999) {
370  inVec.set_user_index(index);
371  //} else {
372  //inVec.set_user_index(fInputVectors.size());
373  //}
374 
375  // add to the fj container of input vectors
376  fInputVectors.push_back(inVec);
377 }
378 
379 //_________________________________________________________________________________________________
380 void AliFJWrapper::AddInputVector(const fj::PseudoJet& vec, Int_t index)
381 {
382  // Add an input pseudojet.
383 
384  fj::PseudoJet inVec = vec;
385 
386  // Salvatore Aiola: not sure why this was done...
388  inVec.set_user_index(index);
389  //} else {
390  //inVec.set_user_index(fInputVectors.size());
391  //}
392 
393  // add to the fj container of input vectors
394  fInputVectors.push_back(inVec);
395 }
396 
397 //_________________________________________________________________________________________________
398 void AliFJWrapper::AddInputVectors(const std::vector<fj::PseudoJet>& vecs, Int_t offsetIndex)
399 {
400  // Add the input from vector of pseudojets.
401 
402  for (UInt_t i = 0; i < vecs.size(); ++i) {
403  fj::PseudoJet inVec = vecs[i];
404  if (offsetIndex > -99999)
405  inVec.set_user_index(fInputVectors.size() + offsetIndex);
406  // add to the fj container of input vectors
407  fInputVectors.push_back(inVec);
408  }
409 }
410 
411 //_________________________________________________________________________________________________
413 {
414  // Make the input pseudojet.
415 
416  fastjet::PseudoJet inVec(px, py, pz, E);
417 
418  if (index > -99999) {
419  inVec.set_user_index(index);
420  } else {
421  inVec.set_user_index(fInputGhosts.size());
422  }
423 
424  // add to the fj container of input vectors
425  fInputGhosts.push_back(inVec);
426  if (!fDoFilterArea) fDoFilterArea = kTRUE;
427 }
428 
429 //_________________________________________________________________________________________________
431 {
432  // Get the jet area.
433 
434  Double_t retval = -1; // really wrong area..
435  if ( idx < fInclusiveJets.size() ) {
436  retval = fClustSeq->area(fInclusiveJets[idx]);
437  } else {
438  AliError(Form("[e] ::GetJetArea wrong index: %d",idx));
439  }
440  return retval;
441 }
442 
443 //_________________________________________________________________________________________________
445 {
446  // Get the filtered jet area.
447 
448  Double_t retval = -1; // really wrong area..
449  if (fDoFilterArea && fClustSeqActGhosts && (idx<fFilteredJets.size())) {
450  retval = fClustSeqActGhosts->area(fFilteredJets[idx]);
451  } else {
452  AliError(Form("[e] ::GetFilteredJetArea wrong index: %d",idx));
453  }
454  return retval;
455 }
456 
457 //_________________________________________________________________________________________________
458 fastjet::PseudoJet AliFJWrapper::GetJetAreaVector(UInt_t idx) const
459 {
460  // Get the jet area as vector.
461  fastjet::PseudoJet retval;
462  if ( idx < fInclusiveJets.size() ) {
463  retval = fClustSeq->area_4vector(fInclusiveJets[idx]);
464  } else {
465  AliError(Form("[e] ::GetJetArea wrong index: %d",idx));
466  }
467  return retval;
468 }
469 
470 //_________________________________________________________________________________________________
471 fastjet::PseudoJet AliFJWrapper::GetFilteredJetAreaVector(UInt_t idx) const
472 {
473  // Get the jet area as vector.
474  fastjet::PseudoJet retval;
475  if (fDoFilterArea && fClustSeqActGhosts && (idx<fFilteredJets.size())) {
476  retval = fClustSeqActGhosts->area_4vector(fFilteredJets[idx]);
477  } else {
478  AliError(Form("[e] ::GetFilteredJetArea wrong index: %d",idx));
479  }
480  return retval;
481 }
482 
483 //_________________________________________________________________________________________________
484 std::vector<double> AliFJWrapper::GetSubtractedJetsPts(Double_t median_pt, Bool_t sorted)
485 {
486  // Get subtracted jets pTs, returns vector.
487 
488  SubtractBackground(median_pt);
489 
490  if (kTRUE == sorted) {
491  std::sort(fSubtractedJetsPt.begin(), fSubtractedJetsPt.begin());
492  }
493  return fSubtractedJetsPt;
494 }
495 
496 //_________________________________________________________________________________________________
498 {
499  // Get subtracted jets pTs, returns Double_t.
500 
501  Double_t retval = -99999.; // really wrong pt..
502  if ( idx < fSubtractedJetsPt.size() ) {
503  retval = fSubtractedJetsPt[idx];
504  }
505  return retval;
506 }
507 
508 //_________________________________________________________________________________________________
509 std::vector<fastjet::PseudoJet>
511 {
512  // Get jets constituents.
513 
514  std::vector<fastjet::PseudoJet> retval;
515 
516  if ( idx < fInclusiveJets.size() ) {
517  retval = fClustSeq->constituents(fInclusiveJets[idx]);
518  } else {
519  AliError(Form("[e] ::GetJetConstituents wrong index: %d",idx));
520  }
521 
522  return retval;
523 }
524 
525 //_________________________________________________________________________________________________
526 std::vector<fastjet::PseudoJet>
528 {
529  // Get jets constituents.
530 
531  std::vector<fastjet::PseudoJet> retval;
532 
533  if ( idx < fFilteredJets.size() ) {
534  if (fClustSeqSA) retval = fClustSeqSA->constituents(fFilteredJets[idx]);
535  if (fClustSeqActGhosts) retval = fClustSeqActGhosts->constituents(fFilteredJets[idx]);
536  } else {
537  AliError(Form("[e] ::GetFilteredJetConstituents wrong index: %d",idx));
538  }
539 
540  return retval;
541 }
542 
543 //_________________________________________________________________________________________________
544 void AliFJWrapper::GetMedianAndSigma(Double_t &median, Double_t &sigma, Int_t remove) const
545 {
546  // Get the median and sigma from fastjet.
547  // User can also do it on his own because the cluster sequence is exposed (via a getter)
548 
549  if (!fClustSeq) {
550  AliError("[e] Run the jfinder first.");
551  return;
552  }
553 
554  Double_t mean_area = 0;
555  try {
556  if(0 == remove) {
557  fClustSeq->get_median_rho_and_sigma(*fRange, fUseArea4Vector, median, sigma, mean_area);
558  } else {
559  std::vector<fastjet::PseudoJet> input_jets = sorted_by_pt(fClustSeq->inclusive_jets());
560  input_jets.erase(input_jets.begin(), input_jets.begin() + remove);
561  fClustSeq->get_median_rho_and_sigma(input_jets, *fRange, fUseArea4Vector, median, sigma, mean_area);
562  input_jets.clear();
563  }
564  } catch (fj::Error) {
565  AliError(" [w] FJ Exception caught.");
566  median = -1.;
567  sigma = -1;
568  }
569 }
570 
571 //_________________________________________________________________________________________________
573 {
574  // Run the actual jet finder.
575 
576  if (fAreaType == fj::voronoi_area) {
577  // Rfact - check dependence - default is 1.
578  // NOTE: hardcoded variable!
579  fVorAreaSpec = new fj::VoronoiAreaSpec(1.);
580  fAreaDef = new fj::AreaDefinition(*fVorAreaSpec);
581  } else {
582  fGhostedAreaSpec = new fj::GhostedAreaSpec(fMaxRap,
584  fGhostArea,
585  fGridScatter,
586  fKtScatter,
587  fMeanGhostKt);
588 
589  fAreaDef = new fj::AreaDefinition(*fGhostedAreaSpec, fAreaType);
590  }
591 
592  // this is acceptable by fastjet:
593 #ifndef FASTJET_VERSION
594  fRange = new fj::RangeDefinition(fMaxRap - 0.95 * fR);
595 #else
596  fRange = new fj::Selector(fj::SelectorAbsRapMax(fMaxRap - 0.95 * fR));
597 #endif
598 
599  if (fAlgor == fj::plugin_algorithm) {
600  if (fPluginAlgor == 0) {
601  // SIS CONE ALGOR
602  // NOTE: hardcoded split parameter
603  Double_t overlap_threshold = 0.75; // NOTE: this actually splits a lot: thr/min(pt1,pt2)
604  fPlugin = new fj::SISConePlugin(fR,
605  overlap_threshold,
606  0, //search of stable cones - zero = until no more
607  1.0); // this should be seed effectively for proto jets
608  fJetDef = new fastjet::JetDefinition(fPlugin);
609  } else if (fPluginAlgor == 1) {
610  // CDF cone
611  // NOTE: hardcoded split parameter
612  Double_t overlap_threshold = 0.75; // NOTE: this actually splits a lot: thr/min(pt1,pt2)
613  fPlugin = new fj::CDFMidPointPlugin(fR,
614  overlap_threshold,
615  1.0, //search of stable cones - zero = until no more
616  1.0); // this should be seed effectively for proto jets
617  fJetDef = new fastjet::JetDefinition(fPlugin);
618  } else {
619  AliError("[e] Unrecognized plugin number!");
620  }
621  } else {
622  fJetDef = new fj::JetDefinition(fAlgor, fR, fScheme, fStrategy);
623  }
624 
625  try {
626  fClustSeq = new fj::ClusterSequenceArea(fInputVectors, *fJetDef, *fAreaDef);
627  } catch (fj::Error) {
628  AliError(" [w] FJ Exception caught.");
629  return -1;
630  }
631 
632  // FJ3 :: Define an JetMedianBackgroundEstimator just in case it will be used
633 #ifdef FASTJET_VERSION
634  fBkrdEstimator = new fj::JetMedianBackgroundEstimator(fj::SelectorAbsRapMax(fMaxRap));
635 #endif
636 
637  if (fLegacyMode) { SetLegacyFJ(); } // for FJ 2.x even if fLegacyMode is set, SetLegacyFJ is dummy
638 
639  // inclusive jets:
640  fInclusiveJets.clear();
641  fInclusiveJets = fClustSeq->inclusive_jets(0.0);
642 
643  return 0;
644 }
645 
646 //_________________________________________________________________________________________________
648 {
649 //
650 // AliFJWrapper::Filter
651 //
652 
653  fJetDef = new fj::JetDefinition(fAlgor, fR, fScheme, fStrategy);
654 
655  if (fDoFilterArea) {
656  if (fInputGhosts.size()>0) {
657  try {
658  fClustSeqActGhosts = new fj::ClusterSequenceActiveAreaExplicitGhosts(fInputVectors,
659  *fJetDef,
660  fInputGhosts,
661  fGhostArea);
662  } catch (fj::Error) {
663  AliError(" [w] FJ Exception caught.");
664  return -1;
665  }
666 
667  fFilteredJets.clear();
668  fFilteredJets = fClustSeqActGhosts->inclusive_jets(0.0);
669  } else {
670  return -1;
671  }
672  } else {
673  try {
674  fClustSeqSA = new fastjet::ClusterSequence(fInputVectors, *fJetDef);
675  } catch (fj::Error) {
676  AliError(" [w] FJ Exception caught.");
677  return -1;
678  }
679 
680  fFilteredJets.clear();
681  fFilteredJets = fClustSeqSA->inclusive_jets(0.0);
682  }
683 
684  return 0;
685 }
686 
687 //_________________________________________________________________________________________________
689 {
690  // This methods enable legacy behaviour (FastJet 2.x) when AliROOT is compiled with FastJet 3.x
691 #ifdef FASTJET_VERSION
692  std::cout << "WARNING! Setting FastJet in legacy mode" << std::endl;
693  if (fGhostedAreaSpec) { fGhostedAreaSpec->set_fj2_placement(kTRUE); }
694  if (fBkrdEstimator) {
695  fBkrdEstimator->set_provide_fj2_sigma(kTRUE);
696  fBkrdEstimator->set_use_area_4vector(kFALSE);
697  }
698 #endif
699 }
700 
701 //_________________________________________________________________________________________________
703 {
704  // Subtract the background (specify the value - see below the meaning).
705  // Negative argument means the bg will be determined with the current algorithm
706  // this is the default behavior. Zero is allowed
707  // Note: user may set the switch for area4vector based subtraction.
708 
709  Double_t median = 0;
710  Double_t sigma = 0;
711  Double_t mean_area = 0;
712 
713  // clear the subtracted jet pt's vector<double>
714  fSubtractedJetsPt.clear();
715 
716  // check what was specified (default is -1)
717  if (median_pt < 0) {
718  try {
719  fClustSeq->get_median_rho_and_sigma(*fRange, fUseArea4Vector, median, sigma, mean_area);
720  }
721 
722  catch (fj::Error) {
723  AliError(" [w] FJ Exception caught.");
724  median = -9999.;
725  sigma = -1;
726  fMedUsedForBgSub = median;
727  return;
728  }
729  fMedUsedForBgSub = median;
730  } else {
731  // we do not know the sigma in this case
732  sigma = -1;
733  if (0.0 == median_pt) {
734  // AliWarning(Form(" [w] Median specified for bg subtraction is ZERO: %f", median_pt )) << endl;
735  fMedUsedForBgSub = 0.;
736  } else {
737  fMedUsedForBgSub = median_pt;
738  }
739  }
740 
741  // subtract:
742  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
743  if ( fUseArea4Vector ) {
744  // subtract the background using the area4vector
745  fj::PseudoJet area4v = fClustSeq->area_4vector(fInclusiveJets[i]);
746  fj::PseudoJet jet_sub = fInclusiveJets[i] - area4v * fMedUsedForBgSub;
747  fSubtractedJetsPt.push_back(jet_sub.perp()); // here we put only the pt of the jet - note: this can be negative
748  } else {
749  // subtract the background using scalars
750  // fj::PseudoJet jet_sub = fInclusiveJets[i] - area * fMedUsedForBgSub_;
751  Double_t area = fClustSeq->area(fInclusiveJets[i]);
752  // standard subtraction
753  Double_t pt_sub = fInclusiveJets[i].perp() - fMedUsedForBgSub * area;
754  fSubtractedJetsPt.push_back(pt_sub); // here we put only the pt of the jet - note: this can be negative
755  }
756  }
757 }
758 
759 //_________________________________________________________________________________________________
761  //Do generic subtraction for jet mass
762 #ifdef FASTJET_VERSION
763  CreateGenSub();
764 
765  // Define jet shape
766  AliJetShapeMass shapeMass;
767 
768  // clear the generic subtractor info vector
769  fGenSubtractorInfoJetMass.clear();
770  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
771  fj::contrib::GenericSubtractorInfo info;
772  if(fInclusiveJets[i].perp()>1.e-4)
773  double subtracted_shape = (*fGenSubtractor)(shapeMass, fInclusiveJets[i], info);
774  fGenSubtractorInfoJetMass.push_back(info);
775  }
776 #endif
777  return 0;
778 }
779 
780 //_________________________________________________________________________________________________
782  //Do generic subtraction for jet mass
783 #ifdef FASTJET_VERSION
784  CreateGenSub();
785 
786  if(ijet>fInclusiveJets.size()) return 0;
787 
788  fGRNumerator.clear();
789  fGRDenominator.clear();
790  fGRNumeratorSub.clear();
791  fGRDenominatorSub.clear();
792 
793  // Define jet shape
794  for(Double_t r = 0.; r<fRMax; r+=fDRStep) {
795  AliJetShapeGRNum shapeGRNum(r,fDRStep);
796  AliJetShapeGRDen shapeGRDen(r,fDRStep);
797 
798  // clear the generic subtractor info vector
799  fGenSubtractorInfoGRNum.clear();
800  fGenSubtractorInfoGRDen.clear();
801  fj::contrib::GenericSubtractorInfo infoNum;
802  fj::contrib::GenericSubtractorInfo infoDen;
803  if(fInclusiveJets[ijet].perp()>1.e-4) {
804  double sub_num = (*fGenSubtractor)(shapeGRNum, fInclusiveJets[ijet], infoNum);
805  double sub_den = (*fGenSubtractor)(shapeGRDen, fInclusiveJets[ijet], infoDen);
806  }
807  fGenSubtractorInfoGRNum.push_back(infoNum);
808  fGenSubtractorInfoGRDen.push_back(infoDen);
809  fGRNumerator.push_back(infoNum.unsubtracted());
810  fGRDenominator.push_back(infoDen.unsubtracted());
811  fGRNumeratorSub.push_back(infoNum.second_order_subtracted());
812  fGRDenominatorSub.push_back(infoDen.second_order_subtracted());
813  }
814 #endif
815  return 0;
816 }
817 //_________________________________________________________________________________________________
819  //Do generic subtraction for jet mass
820 #ifdef FASTJET_VERSION
821  CreateGenSub();
822 
823  // Define jet shape
824  AliJetShapeAngularity shapeAngularity;
825 
826  // clear the generic subtractor info vector
827  fGenSubtractorInfoJetAngularity.clear();
828  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
829  fj::contrib::GenericSubtractorInfo infoAng;
830  if(fInclusiveJets[i].perp()>1.e-4)
831  double subtracted_shape = (*fGenSubtractor)(shapeAngularity, fInclusiveJets[i], infoAng);
832  fGenSubtractorInfoJetAngularity.push_back(infoAng);
833  }
834 #endif
835  return 0;
836 }
837 //_________________________________________________________________________________________________
839  //Do generic subtraction for jet mass
840 #ifdef FASTJET_VERSION
841  CreateGenSub();
842 
843  // Define jet shape
844  AliJetShapepTD shapepTD;
845 
846  // clear the generic subtractor info vector
847  fGenSubtractorInfoJetpTD.clear();
848  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
849  fj::contrib::GenericSubtractorInfo infopTD;
850  if(fInclusiveJets[i].perp()>1.e-4)
851  double subtracted_shape = (*fGenSubtractor)(shapepTD, fInclusiveJets[i], infopTD);
852  fGenSubtractorInfoJetpTD.push_back(infopTD);
853  }
854 #endif
855  return 0;
856 }
857 //_________________________________________________________________________________________________
859  //Do generic subtraction for jet mass
860 #ifdef FASTJET_VERSION
861  CreateGenSub();
862 
863  // Define jet shape
864  AliJetShapeCircularity shapecircularity;
865 
866  // clear the generic subtractor info vector
867  fGenSubtractorInfoJetCircularity.clear();
868  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
869  fj::contrib::GenericSubtractorInfo infoCirc;
870  if(fInclusiveJets[i].perp()>1.e-4)
871  double subtracted_shape = (*fGenSubtractor)(shapecircularity, fInclusiveJets[i], infoCirc);
872  fGenSubtractorInfoJetCircularity.push_back(infoCirc);
873  }
874 #endif
875  return 0;
876 }
877 //_________________________________________________________________________________________________
879  //Do generic subtraction for jet mass
880 #ifdef FASTJET_VERSION
881  CreateGenSub();
882 
883  // Define jet shape
884  AliJetShapeSigma2 shapesigma2;
885 
886  // clear the generic subtractor info vector
887  fGenSubtractorInfoJetSigma2.clear();
888  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
889  fj::contrib::GenericSubtractorInfo infoSigma;
890  if(fInclusiveJets[i].perp()>1.e-4)
891  double subtracted_shape = (*fGenSubtractor)(shapesigma2, fInclusiveJets[i], infoSigma);
892  fGenSubtractorInfoJetSigma2.push_back(infoSigma);
893  }
894 #endif
895  return 0;
896 }
897 //_________________________________________________________________________________________________
899  //Do generic subtraction for jet mass
900 #ifdef FASTJET_VERSION
901  CreateGenSub();
902 
903  // Define jet shape
904  AliJetShapeConstituent shapeconst;
905 
906  // clear the generic subtractor info vector
907  fGenSubtractorInfoJetConstituent.clear();
908  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
909  fj::contrib::GenericSubtractorInfo infoConst;
910  if(fInclusiveJets[i].perp()>1.e-4)
911  double subtracted_shape = (*fGenSubtractor)(shapeconst, fInclusiveJets[i], infoConst);
912  fGenSubtractorInfoJetConstituent.push_back(infoConst);
913  }
914 #endif
915  return 0;
916 }
917 
918 //_________________________________________________________________________________________________
920  //Do generic subtraction for jet mass
921 #ifdef FASTJET_VERSION
922  CreateGenSub();
923 
924  // Define jet shape
925  AliJetShapeLeSub shapeLeSub;
926 
927  // clear the generic subtractor info vector
928  fGenSubtractorInfoJetLeSub.clear();
929  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
930  fj::contrib::GenericSubtractorInfo infoLeSub;
931  if(fInclusiveJets[i].perp()>1.e-4)
932  double subtracted_shape = (*fGenSubtractor)(shapeLeSub, fInclusiveJets[i], infoLeSub);
933  fGenSubtractorInfoJetLeSub.push_back(infoLeSub);
934  }
935 #endif
936  return 0;
937 }
938 
939 //_________________________________________________________________________________________________
941  //Do generic subtraction for 1subjettiness
942 #ifdef FASTJET_VERSION
943  CreateGenSub();
944 
945  // Define jet shape
946  AliJetShape1subjettiness_kt shape1subjettiness_kt;
947 
948  // clear the generic subtractor info vector
949  fGenSubtractorInfoJet1subjettiness_kt.clear();
950  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
951  fj::contrib::GenericSubtractorInfo info1subjettiness_kt;
952  if(fInclusiveJets[i].perp()>1.e-4)
953  double subtracted_shape = (*fGenSubtractor)(shape1subjettiness_kt, fInclusiveJets[i], info1subjettiness_kt);
954  fGenSubtractorInfoJet1subjettiness_kt.push_back(info1subjettiness_kt);
955  }
956 #endif
957  return 0;
958 }
959 
960 //_________________________________________________________________________________________________
962  //Do generic subtraction for 2subjettiness
963 #ifdef FASTJET_VERSION
964  CreateGenSub();
965 
966  // Define jet shape
967  AliJetShape2subjettiness_kt shape2subjettiness_kt;
968 
969  // clear the generic subtractor info vector
970  fGenSubtractorInfoJet2subjettiness_kt.clear();
971  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
972  fj::contrib::GenericSubtractorInfo info2subjettiness_kt;
973  if(fInclusiveJets[i].perp()>1.e-4)
974  double subtracted_shape = (*fGenSubtractor)(shape2subjettiness_kt, fInclusiveJets[i], info2subjettiness_kt);
975  fGenSubtractorInfoJet2subjettiness_kt.push_back(info2subjettiness_kt);
976  }
977 #endif
978  return 0;
979 }
980 
981 //_________________________________________________________________________________________________
983  //Do generic subtraction for 3subjettiness
984 #ifdef FASTJET_VERSION
985  CreateGenSub();
986 
987  // Define jet shape
988  AliJetShape3subjettiness_kt shape3subjettiness_kt;
989 
990  // clear the generic subtractor info vector
991  fGenSubtractorInfoJet3subjettiness_kt.clear();
992  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
993  fj::contrib::GenericSubtractorInfo info3subjettiness_kt;
994  if(fInclusiveJets[i].perp()>1.e-4)
995  double subtracted_shape = (*fGenSubtractor)(shape3subjettiness_kt, fInclusiveJets[i], info3subjettiness_kt);
996  fGenSubtractorInfoJet3subjettiness_kt.push_back(info3subjettiness_kt);
997  }
998 #endif
999  return 0;
1000 }
1001 
1002 //_________________________________________________________________________________________________
1004  //Do generic subtraction for 2subjettiness axes opening angle
1005 #ifdef FASTJET_VERSION
1006  CreateGenSub();
1007 
1008  // Define jet shape
1009  AliJetShapeOpeningAngle_kt shapeOpeningAngle_kt;
1010 
1011  // clear the generic subtractor info vector
1012  fGenSubtractorInfoJetOpeningAngle_kt.clear();
1013  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
1014  fj::contrib::GenericSubtractorInfo infoOpeningAngle_kt;
1015  if(fInclusiveJets[i].perp()>1.e-4)
1016  double subtracted_shape = (*fGenSubtractor)(shapeOpeningAngle_kt, fInclusiveJets[i], infoOpeningAngle_kt);
1017  fGenSubtractorInfoJetOpeningAngle_kt.push_back(infoOpeningAngle_kt);
1018  }
1019 #endif
1020  return 0;
1021 }
1022 
1023 //_________________________________________________________________________________________________
1025  //Do constituent subtraction
1026 #ifdef FASTJET_VERSION
1027  CreateConstituentSub();
1028  // fConstituentSubtractor->set_alpha(/* double alpha */);
1029  // fConstituentSubtractor->set_max_deltaR(/* double max_deltaR */);
1030 
1031  //clear constituent subtracted jets
1032  fConstituentSubtrJets.clear();
1033  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
1034  fj::PseudoJet subtracted_jet(0.,0.,0.,0.);
1035  if(fInclusiveJets[i].perp()>0.)
1036  subtracted_jet = (*fConstituentSubtractor)(fInclusiveJets[i]);
1037  fConstituentSubtrJets.push_back(subtracted_jet);
1038  }
1039  if(fConstituentSubtractor) { delete fConstituentSubtractor; fConstituentSubtractor = NULL; }
1040 
1041 #endif
1042  return 0;
1043 }
1044 
1045 //_________________________________________________________________________________________________
1047  //Do grooming
1048 #ifdef FASTJET_VERSION
1049  CreateSoftDrop();
1050 
1051  //clear groomed jets
1052  fGroomedJets.clear();
1053  //fastjet::Subtractor fjsub (fBkrdEstimator);
1054  //fSoftDrop->set_subtractor(&fjsub);
1055  //fSoftDrop->set_input_jet_is_subtracted(false); //??
1056 
1057  for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
1058  fj::PseudoJet groomed_jet(0.,0.,0.,0.);
1059  if(fInclusiveJets[i].perp()>0.){
1060  groomed_jet = (*fSoftDrop)(fInclusiveJets[i]);
1061  groomed_jet.set_user_index(i); //index of the corresponding inclusve jet
1062  if(groomed_jet!=0) fGroomedJets.push_back(groomed_jet);
1063  }
1064 
1065  }
1066  if(fSoftDrop) { delete fSoftDrop; fSoftDrop = NULL; }
1067 
1068 #endif
1069  return 0;
1070 }
1071 
1072 //_________________________________________________________________________________________________
1073 Int_t AliFJWrapper::CreateSoftDrop() {
1074  //Do grooming
1075  #ifdef FASTJET_VERSION
1076  if (fSoftDrop) { delete fSoftDrop; } // protect against memory leaks
1077 
1078  fSoftDrop = new fj::contrib::SoftDrop(fBeta,fZcut);
1079 
1080  #endif
1081  return 0;
1082 }
1083 
1084 
1085 //_________________________________________________________________________________________________
1086 Int_t AliFJWrapper::CreateGenSub() {
1087  //Do generic subtraction for jet mass
1088  #ifdef FASTJET_VERSION
1089  if (fGenSubtractor) { delete fGenSubtractor; } // protect against memory leaks
1090 
1091  if (fUseExternalBkg)
1092  { fGenSubtractor = new fj::contrib::GenericSubtractor(fRho,fRhom); }
1093  else
1094  {
1095  fGenSubtractor = new fj::contrib::GenericSubtractor(fBkrdEstimator);
1096  #if FASTJET_VERSION_NUMBER >= 30100
1097  fGenSubtractor->set_common_bge_for_rho_and_rhom(); // see contrib 1.020 GenericSubtractor.hh line 62
1098  #endif
1099  }
1100 
1101  #endif
1102  return 0;
1103 }
1104 
1105 //_________________________________________________________________________________________________
1106 Int_t AliFJWrapper::CreateConstituentSub() {
1107  //Do generic subtraction for jet mass
1108  #ifdef FASTJET_VERSION
1109  if (fConstituentSubtractor) { delete fConstituentSubtractor; } // protect against memory leaks
1110 
1111  // see ConstituentSubtractor.hh signatures
1112  // ConstituentSubtractor(double rho, double rhom=0, double alpha=0, double maxDeltaR=-1)
1113  if (fUseExternalBkg) { fConstituentSubtractor = new fj::contrib::ConstituentSubtractor(fRho,fRhom); }
1114  else { fConstituentSubtractor = new fj::contrib::ConstituentSubtractor(fBkrdEstimator); }
1115 
1116  #endif
1117  return 0;
1118 }
1119 
1120 //_________________________________________________________________________________________________
1121 void AliFJWrapper::SetupAlgorithmfromOpt(const char *option)
1122 {
1123  // Setup algorithm from char.
1124 
1125  std::string opt(option);
1126 
1127  if (!opt.compare("kt")) fAlgor = fj::kt_algorithm;
1128  if (!opt.compare("antikt")) fAlgor = fj::antikt_algorithm;
1129  if (!opt.compare("cambridge")) fAlgor = fj::cambridge_algorithm;
1130  if (!opt.compare("genkt")) fAlgor = fj::genkt_algorithm;
1131  if (!opt.compare("cambridge_passive")) fAlgor = fj::cambridge_for_passive_algorithm;
1132  if (!opt.compare("genkt_passive")) fAlgor = fj::genkt_for_passive_algorithm;
1133  if (!opt.compare("ee_kt")) fAlgor = fj::ee_kt_algorithm;
1134  if (!opt.compare("ee_genkt")) fAlgor = fj::ee_genkt_algorithm;
1135  if (!opt.compare("plugin")) fAlgor = fj::plugin_algorithm;
1136 }
1137 
1138 //_________________________________________________________________________________________________
1139 void AliFJWrapper::SetupAreaTypefromOpt(const char *option)
1140 {
1141  // Setup area type from char.
1142 
1143  std::string opt(option);
1144 
1145  if (!opt.compare("active")) fAreaType = fj::active_area;
1146  if (!opt.compare("invalid")) fAreaType = fj::invalid_area;
1147  if (!opt.compare("active_area_explicit_ghosts")) fAreaType = fj::active_area_explicit_ghosts;
1148  if (!opt.compare("one_ghost_passive")) fAreaType = fj::one_ghost_passive_area;
1149  if (!opt.compare("passive")) fAreaType = fj::passive_area;
1150  if (!opt.compare("voronoi")) fAreaType = fj::voronoi_area;
1151 }
1152 
1153 //_________________________________________________________________________________________________
1154 void AliFJWrapper::SetupSchemefromOpt(const char *option)
1155 {
1156  //
1157  // setup scheme from char
1158  //
1159 
1160  std::string opt(option);
1161 
1162  if (!opt.compare("BIpt")) fScheme = fj::BIpt_scheme;
1163  if (!opt.compare("BIpt2")) fScheme = fj::BIpt2_scheme;
1164  if (!opt.compare("E")) fScheme = fj::E_scheme;
1165  if (!opt.compare("pt")) fScheme = fj::pt_scheme;
1166  if (!opt.compare("pt2")) fScheme = fj::pt2_scheme;
1167  if (!opt.compare("Et")) fScheme = fj::Et_scheme;
1168  if (!opt.compare("Et2")) fScheme = fj::Et2_scheme;
1169 }
1170 
1171 //_________________________________________________________________________________________________
1172 void AliFJWrapper::SetupStrategyfromOpt(const char *option)
1173 {
1174  // Setup strategy from char.
1175 
1176  std::string opt(option);
1177 
1178  if (!opt.compare("Best")) fStrategy = fj::Best;
1179  if (!opt.compare("N2MinHeapTiled")) fStrategy = fj::N2MinHeapTiled;
1180  if (!opt.compare("N2Tiled")) fStrategy = fj::N2Tiled;
1181  if (!opt.compare("N2PoorTiled")) fStrategy = fj::N2PoorTiled;
1182  if (!opt.compare("N2Plain")) fStrategy = fj::N2Plain;
1183  if (!opt.compare("N3Dumb")) fStrategy = fj::N3Dumb;
1184  if (!opt.compare("NlnN")) fStrategy = fj::NlnN;
1185  if (!opt.compare("NlnN3pi")) fStrategy = fj::NlnN3pi;
1186  if (!opt.compare("NlnN4pi")) fStrategy = fj::NlnN4pi;
1187  if (!opt.compare("NlnNCam4pi")) fStrategy = fj::NlnNCam4pi;
1188  if (!opt.compare("NlnNCam2pi2R")) fStrategy = fj::NlnNCam2pi2R;
1189  if (!opt.compare("NlnNCam")) fStrategy = fj::NlnNCam;
1190  if (!opt.compare("plugin")) fStrategy = fj::plugin_strategy;
1191 }
1192 
1193 Double_t AliFJWrapper::NSubjettiness(Int_t N, Int_t Algorithm, Double_t Radius, Double_t Beta, Int_t Option, Int_t Measure, Double_t Beta_SD, Double_t ZCut){
1194 
1195 
1196  //Option 0=Nsubjettiness result, 1=opening angle between axes in Eta-Phi plane, 2=Distance between axes in Eta-Phi plane
1197 
1198  fJetDef = new fj::JetDefinition(fAlgor, fR*2, fScheme, fStrategy ); //the *2 is becasue of a handful of jets that end up missing a track for some reason.
1199 
1200  try {
1201  fClustSeqSA = new fastjet::ClusterSequence(fInputVectors, *fJetDef);
1202  // ClustSeqSA = new fastjet::ClusterSequenceArea(fInputVectors, *fJetDef, *fAreaDef);
1203  } catch (fj::Error) {
1204  AliError(" [w] FJ Exception caught.");
1205  return -1;
1206  }
1207  fFilteredJets.clear();
1208  fFilteredJets = fClustSeqSA->inclusive_jets(fMinJetPt-0.1); //becasue this is < not <=
1209  Double_t Result=-1;
1210  Double_t Result_SoftDrop=-1;
1211  std::vector<fastjet::PseudoJet> SubJet_Axes; //this is actually not subjet axis...but just axis
1212  fj::PseudoJet SubJet1_Axis;
1213  fj::PseudoJet SubJet2_Axis;
1214  std::vector<fastjet::PseudoJet> SubJets;
1215  fj::PseudoJet SubJet1;
1216  fj::PseudoJet SubJet2;
1217  if (Algorithm==0){
1218  if (Measure==0){
1219  fj::contrib::Nsubjettiness nSub(N, fj::contrib::KT_Axes(), fj::contrib::NormalizedMeasure(Beta,fR*2));
1220  Result= nSub.result(fFilteredJets[0]);
1221  SubJet_Axes=nSub.currentAxes();
1222  SubJets=nSub.currentSubjets();
1223  }
1224  else if (Measure==1){
1225  fj::contrib::Nsubjettiness nSub(N, fj::contrib::KT_Axes(), fj::contrib::UnnormalizedMeasure(Beta));
1226  Result= nSub.result(fFilteredJets[0]);
1227  SubJet_Axes=nSub.currentAxes();
1228  SubJets=nSub.currentSubjets();
1229  }
1230  }
1231  else if (Algorithm==1) {
1232  fj::contrib::Nsubjettiness nSub(N, fj::contrib::CA_Axes(), fj::contrib::NormalizedMeasure(Beta,fR));
1233  Result= nSub.result(fFilteredJets[0]);
1234  SubJet_Axes=nSub.currentAxes();
1235  SubJets=nSub.currentSubjets();
1236  }
1237  else if (Algorithm==2){
1238  fj::contrib::Nsubjettiness nSub(N, fj::contrib::AntiKT_Axes(Radius), fj::contrib::NormalizedMeasure(Beta,fR));
1239  Result= nSub.result(fFilteredJets[0]);
1240  SubJet_Axes=nSub.currentAxes();
1241  SubJets=nSub.currentSubjets();
1242  }
1243  else if (Algorithm==3) {
1244  fj::contrib::Nsubjettiness nSub(N, fj::contrib::WTA_KT_Axes(), fj::contrib::NormalizedMeasure(Beta,fR));
1245  Result= nSub.result(fFilteredJets[0]);
1246  SubJet_Axes=nSub.currentAxes();
1247  SubJets=nSub.currentSubjets();
1248  }
1249  else if (Algorithm==4) {
1250  fj::contrib::Nsubjettiness nSub(N, fj::contrib::WTA_CA_Axes(), fj::contrib::NormalizedMeasure(Beta,fR));
1251  Result= nSub.result(fFilteredJets[0]);
1252  SubJet_Axes=nSub.currentAxes();
1253  SubJets=nSub.currentSubjets();
1254  }
1255  else if (Algorithm==5){
1256  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_KT_Axes(), fj::contrib::NormalizedMeasure(Beta,fR));
1257  Result= nSub.result(fFilteredJets[0]);
1258  SubJet_Axes=nSub.currentAxes();
1259  SubJets=nSub.currentSubjets();
1260  }
1261  else if (Algorithm==6){
1262  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_CA_Axes(), fj::contrib::NormalizedMeasure(Beta,fR));
1263  Result= nSub.result(fFilteredJets[0]);
1264  SubJet_Axes=nSub.currentAxes();
1265  SubJets=nSub.currentSubjets();
1266  }
1267  else if (Algorithm==7){
1268  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_AntiKT_Axes(Radius), fj::contrib::NormalizedMeasure(Beta,fR));
1269  Result= nSub.result(fFilteredJets[0]);
1270  SubJet_Axes=nSub.currentAxes();
1271  SubJets=nSub.currentSubjets();
1272  }
1273  else if (Algorithm==8){
1274  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_WTA_KT_Axes(), fj::contrib::NormalizedMeasure(Beta,fR));
1275  Result= nSub.result(fFilteredJets[0]);
1276  SubJet_Axes=nSub.currentAxes();
1277  SubJets=nSub.currentSubjets();
1278  }
1279  else if (Algorithm==9){
1280  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_WTA_CA_Axes(), fj::contrib::NormalizedMeasure(Beta,fR));
1281  Result= nSub.result(fFilteredJets[0]);
1282  SubJet_Axes=nSub.currentAxes();
1283  SubJets=nSub.currentSubjets();
1284  }
1285  else if (Algorithm==10){
1286  fj::contrib::Nsubjettiness nSub(N, fj::contrib::MultiPass_Axes(100), fj::contrib::NormalizedMeasure(Beta,fR));
1287  Result= nSub.result(fFilteredJets[0]);
1288  SubJet_Axes=nSub.currentAxes();
1289  SubJets=nSub.currentSubjets();
1290  }
1291 
1292 
1293 
1294  SubJet1_Axis=SubJet_Axes[0];
1295  Double_t SubJet1_Eta=SubJet1_Axis.pseudorapidity();
1296  Double_t SubJet2_Eta;
1297  Double_t SubJet1_Phi=SubJet1_Axis.phi();
1298  if(SubJet1_Phi < -1*TMath::Pi()) SubJet1_Phi += (2*TMath::Pi());
1299  else if (SubJet1_Phi > TMath::Pi()) SubJet1_Phi -= (2*TMath::Pi());
1300  Double_t SubJet2_Phi;
1301  Double_t DeltaPhi=-5;
1302  if (SubJet_Axes.size()>1){
1303  SubJet2_Axis=SubJet_Axes[1];
1304  SubJet2_Eta=SubJet2_Axis.pseudorapidity();
1305  SubJet2_Phi=SubJet2_Axis.phi();
1306  if(SubJet2_Phi < -1*TMath::Pi()) SubJet2_Phi += (2*TMath::Pi());
1307  else if (SubJet2_Phi > TMath::Pi()) SubJet2_Phi -= (2*TMath::Pi());
1308  DeltaPhi=SubJet1_Phi-SubJet2_Phi;
1309  if(DeltaPhi < -1*TMath::Pi()) DeltaPhi += (2*TMath::Pi());
1310  else if (DeltaPhi > TMath::Pi()) DeltaPhi -= (2*TMath::Pi());
1311  }
1312 
1313 
1314  SubJet1=SubJets[0];
1315  Double_t SubJet1Eta=SubJet1.pseudorapidity();
1316  Double_t SubJet2Eta;
1317  Double_t SubJet1Phi=SubJet1.phi();
1318  if(SubJet1Phi < -1*TMath::Pi()) SubJet1Phi += (2*TMath::Pi());
1319  else if (SubJet1Phi > TMath::Pi()) SubJet1Phi -= (2*TMath::Pi());
1320  Double_t SubJet2Phi;
1321  Double_t DeltaPhiSubJets=-5;
1322  if (SubJet_Axes.size()>1){
1323  SubJet2=SubJets[1];
1324  SubJet2Eta=SubJet2.pseudorapidity();
1325  SubJet2Phi=SubJet2.phi();
1326  if(SubJet2Phi < -1*TMath::Pi()) SubJet2Phi += (2*TMath::Pi());
1327  else if (SubJet2Phi > TMath::Pi()) SubJet2Phi -= (2*TMath::Pi());
1328  DeltaPhiSubJets=SubJet1Phi-SubJet2Phi;
1329  if(DeltaPhiSubJets < -1*TMath::Pi()) DeltaPhiSubJets += (2*TMath::Pi());
1330  else if (DeltaPhiSubJets > TMath::Pi()) DeltaPhiSubJets -= (2*TMath::Pi());
1331  }
1332 
1333  //Added for quality control of the DeltaR-Nsubjettiness variable (comparing Nsubjettiness and soft drop results)
1334  fj::contrib::SoftDrop Soft_Drop(Beta_SD,ZCut);
1335  Soft_Drop.set_tagging_mode(); //if the first two subjets fail the soft drop criteria a jet = 0 is returned
1336  fj::PseudoJet Soft_Dropped_Jet=Soft_Drop(fFilteredJets[0]);
1337  if (Soft_Dropped_Jet==0) Result_SoftDrop=-1;
1338  else{
1339  if (Option==3) Result_SoftDrop=Soft_Dropped_Jet.structure_of<fj::contrib::SoftDrop>().delta_R();
1340  if (Option==4) Result_SoftDrop=Soft_Dropped_Jet.structure_of<fj::contrib::SoftDrop>().symmetry();
1341  }
1342 
1343  if (Option==0) return Result;
1344  else if (Option==1 && SubJet_Axes.size()>1 && N==2) return TMath::Sqrt(TMath::Power(SubJet1_Eta-SubJet2_Eta,2)+TMath::Power(DeltaPhi,2));
1345  else if (Option==2 && SubJet_Axes.size()>1 && N==2) return TMath::Sqrt(TMath::Power(SubJet1_Eta-SubJet2_Eta,2)+TMath::Power(DeltaPhi,2));
1346  else if ((Option==3 || Option==4) && N==2) return Result_SoftDrop;
1347  else if (Option==5 && SubJets.size()>1 && N==2) return SubJet1.perp();
1348  else if (Option==6 && SubJets.size()>1 && N==2) return SubJet2.perp();
1349  else if (Option==7 && SubJets.size()>1 && N==2) return TMath::Sqrt(TMath::Power(SubJet1Eta-SubJet2Eta,2)+TMath::Power(DeltaPhiSubJets,2));
1350  else return -2;
1351 }
1352 
1353 
1354 
1355 Double32_t AliFJWrapper::NSubjettinessDerivativeSub(Int_t N, Int_t Algorithm, Double_t Radius, Double_t Beta, Double_t JetR, fastjet::PseudoJet jet, Int_t Option, Int_t Measure){ //For derivative subtraction
1356 
1357  Double_t Result=-1;
1358  std::vector<fastjet::PseudoJet> SubJet_Axes;
1359  fj::PseudoJet SubJet1_Axis;
1360  fj::PseudoJet SubJet2_Axis;
1361  if (Algorithm==0){
1362  fj::contrib::Nsubjettiness nSub(N, fj::contrib::KT_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1363  Result= nSub.result(jet);
1364  SubJet_Axes=nSub.currentAxes();
1365  }
1366  else if (Algorithm==1) {
1367  fj::contrib::Nsubjettiness nSub(N, fj::contrib::CA_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1368  Result= nSub.result(jet);
1369  SubJet_Axes=nSub.currentAxes();
1370  }
1371  else if (Algorithm==2){
1372  fj::contrib::Nsubjettiness nSub(N, fj::contrib::AntiKT_Axes(Radius), fj::contrib::NormalizedMeasure(Beta,JetR));
1373  Result= nSub.result(jet);
1374  SubJet_Axes=nSub.currentAxes();
1375  }
1376  else if (Algorithm==3) {
1377  fj::contrib::Nsubjettiness nSub(N, fj::contrib::WTA_KT_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1378  Result= nSub.result(jet);
1379  SubJet_Axes=nSub.currentAxes();
1380  }
1381  else if (Algorithm==4) {
1382  fj::contrib::Nsubjettiness nSub(N, fj::contrib::WTA_CA_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1383  Result= nSub.result(jet);
1384  SubJet_Axes=nSub.currentAxes();
1385  }
1386  else if (Algorithm==5){
1387  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_KT_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1388  Result= nSub.result(jet);
1389  SubJet_Axes=nSub.currentAxes();
1390  }
1391  else if (Algorithm==6){
1392  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_CA_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1393  Result= nSub.result(jet);
1394  SubJet_Axes=nSub.currentAxes();
1395  }
1396  else if (Algorithm==7){
1397  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_AntiKT_Axes(Radius), fj::contrib::NormalizedMeasure(Beta,JetR));
1398  Result= nSub.result(jet);
1399  SubJet_Axes=nSub.currentAxes();
1400  }
1401  else if (Algorithm==8){
1402  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_WTA_KT_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1403  Result= nSub.result(jet);
1404  SubJet_Axes=nSub.currentAxes();
1405  }
1406  else if (Algorithm==9){
1407  fj::contrib::Nsubjettiness nSub(N, fj::contrib::OnePass_WTA_CA_Axes(), fj::contrib::NormalizedMeasure(Beta,JetR));
1408  Result= nSub.result(jet);
1409  SubJet_Axes=nSub.currentAxes();
1410  }
1411  else if (Algorithm==10){
1412  fj::contrib::Nsubjettiness nSub(N, fj::contrib::MultiPass_Axes(100), fj::contrib::NormalizedMeasure(Beta,JetR));
1413  Result= nSub.result(jet);
1414  SubJet_Axes=nSub.currentAxes();
1415  }
1416 
1417  SubJet1_Axis=SubJet_Axes[0];
1418  Double_t SubJet1_Eta=SubJet1_Axis.pseudorapidity();
1419  Double_t SubJet2_Eta;
1420  Double_t SubJet1_Phi=SubJet1_Axis.phi();
1421  if(SubJet1_Phi < -1*TMath::Pi()) SubJet1_Phi += (2*TMath::Pi());
1422  else if (SubJet1_Phi > TMath::Pi()) SubJet1_Phi -= (2*TMath::Pi());
1423  Double_t SubJet2_Phi;
1424  Double_t DeltaPhi=-5;
1425  if (SubJet_Axes.size()>1){
1426  SubJet2_Axis=SubJet_Axes[1];
1427  SubJet2_Eta=SubJet2_Axis.pseudorapidity();
1428  SubJet2_Phi=SubJet2_Axis.phi();
1429  if(SubJet2_Phi < -1*TMath::Pi()) SubJet2_Phi += (2*TMath::Pi());
1430  else if (SubJet2_Phi > TMath::Pi()) SubJet2_Phi -= (2*TMath::Pi());
1431  DeltaPhi=SubJet1_Phi-SubJet2_Phi;
1432  if(DeltaPhi < -1*TMath::Pi()) DeltaPhi += (2*TMath::Pi());
1433  else if (DeltaPhi > TMath::Pi()) DeltaPhi -= (2*TMath::Pi());
1434  }
1435 
1436 
1437  if (Option==0) return Result;
1438  else if (Option==1 && SubJet_Axes.size()>1 && N==2) return TMath::Sqrt(TMath::Power(SubJet1_Eta-SubJet2_Eta,2)+TMath::Power(DeltaPhi,2));
1439  else if (Option==2 && SubJet_Axes.size()>1 && N==2) return TMath::Sqrt(TMath::Power(SubJet1_Eta-SubJet2_Eta,2)+TMath::Power(DeltaPhi,2));
1440  else return -2;
1441 
1442 }
1443 
1444 
1445 
1446 
1447 
1448 
1449 
1450 
1451 
1452 #endif
Double32_t NSubjettinessDerivativeSub(Int_t N, Int_t Algorithm, Double_t Radius, Double_t Beta, Double_t JetR, fastjet::PseudoJet jet, Int_t Option=0, Int_t Measure=0)
fastjet::PseudoJet GetJetAreaVector(UInt_t idx) const
std::vector< double > fGRNumeratorSub
Definition: AliFJWrapper.h:190
void SetRMaxAndStep(Double_t rmax, Double_t dr)
Definition: AliFJWrapper.h:114
void SetUseExternalBkg(Bool_t b, Double_t rho, Double_t rhom)
Definition: AliFJWrapper.h:113
Double_t GetFilteredJetArea(UInt_t idx) const
std::vector< fastjet::PseudoJet > fInputVectors
Definition: AliFJWrapper.h:121
Double_t fR
Definition: AliFJWrapper.h:148
virtual Int_t Run()
double Double_t
Definition: External.C:58
virtual Int_t Filter()
const char * title
Definition: MakeQAPdf.C:26
virtual std::vector< double > GetSubtractedJetsPts(Double_t median_pt=-1, Bool_t sorted=kFALSE)
fastjet::ClusterSequence * fClustSeqSA
Definition: AliFJWrapper.h:139
virtual Int_t DoSoftDrop()
Int_t fPluginAlgor
Definition: AliFJWrapper.h:154
fastjet::ClusterSequenceActiveAreaExplicitGhosts * fClustSeqActGhosts
Definition: AliFJWrapper.h:140
virtual Int_t DoGenericSubtractionJetOpeningAngle_kt()
Double_t fKtScatter
Definition: AliFJWrapper.h:152
Int_t fNGhostRepeats
Definition: AliFJWrapper.h:145
virtual Int_t DoGenericSubtractionJetpTD()
Double_t GetMedianUsedForBgSubtraction() const
Definition: AliFJWrapper.h:37
virtual void AddInputVectors(const std::vector< fastjet::PseudoJet > &vecs, Int_t offsetIndex=-99999)
std::vector< fastjet::PseudoJet > fInclusiveJets
Definition: AliFJWrapper.h:123
fastjet::JetAlgorithm fAlgor
Definition: AliFJWrapper.h:142
const char * GetTitle() const
Definition: AliFJWrapper.h:39
fastjet::RecombinationScheme fScheme
Definition: AliFJWrapper.h:143
Double_t fRMax
Definition: AliFJWrapper.h:186
virtual void CopySettingsFrom(const AliFJWrapper &wrapper)
void SetupAlgorithmfromOpt(const char *option)
virtual Int_t DoGenericSubtractionJetAngularity()
virtual Int_t DoGenericSubtractionGR(Int_t ijet)
virtual Int_t DoGenericSubtractionJetMass()
void SetMinJetPt(Double_t MinPt)
Definition: AliFJWrapper.h:116
void SetMaxRap(Double_t maxrap)
Definition: AliFJWrapper.h:100
Double_t GetJetSubtractedPt(UInt_t idx) const
Double_t fMinJetPt
Definition: AliFJWrapper.h:149
fastjet::PseudoJet GetFilteredJetAreaVector(UInt_t idx) const
void SetRecombScheme(const fastjet::RecombinationScheme &scheme)
Definition: AliFJWrapper.h:96
fastjet::RangeDefinition * fRange
Definition: AliFJWrapper.h:134
Double_t fZcut
Definition: AliFJWrapper.h:159
fastjet::ClusterSequenceActiveAreaExplicitGhosts * GetClusterSequenceGhosts() const
Definition: AliFJWrapper.h:30
fastjet::ClusterSequence * GetClusterSequenceSA() const
Definition: AliFJWrapper.h:29
void SetStrategy(const fastjet::Strategy &strat)
Definition: AliFJWrapper.h:94
std::vector< fastjet::PseudoJet > fFilteredJets
Definition: AliFJWrapper.h:124
void SetupAreaTypefromOpt(const char *option)
virtual void ClearMemory()
fastjet::AreaType fAreaType
Definition: AliFJWrapper.h:144
void SetR(Double_t r)
Definition: AliFJWrapper.h:101
Double_t fMaxRap
Definition: AliFJWrapper.h:147
fastjet::ClusterSequenceArea * GetClusterSequence() const
Definition: AliFJWrapper.h:28
Double_t * sigma
const std::vector< fastjet::PseudoJet > & GetInclusiveJets() const
Definition: AliFJWrapper.h:33
int Int_t
Definition: External.C:63
AliFJWrapper & operator=(const AliFJWrapper &wrapper)
unsigned int UInt_t
Definition: External.C:33
Double_t GetJetArea(UInt_t idx) const
std::vector< double > fSubtractedJetsPt
Definition: AliFJWrapper.h:125
fastjet::JetDefinition::Plugin * fPlugin
Definition: AliFJWrapper.h:132
void SetLegacyMode(Bool_t mode)
Definition: AliFJWrapper.h:111
Bool_t fDoFilterArea
Definition: AliFJWrapper.h:181
void SetAlgorithm(const fastjet::JetAlgorithm &algor)
Definition: AliFJWrapper.h:95
virtual Int_t DoGenericSubtractionJetCircularity()
void SetupSchemefromOpt(const char *option)
const std::vector< fastjet::PseudoJet > & GetInputGhosts() const
Definition: AliFJWrapper.h:32
std::vector< fastjet::PseudoJet > fGroomedJets
Definition: AliFJWrapper.h:127
void SetUseArea4Vector(Bool_t useA4v)
Definition: AliFJWrapper.h:106
fastjet::ClusterSequenceArea * fClustSeq
Definition: AliFJWrapper.h:138
void SetName(const char *name)
Definition: AliFJWrapper.h:92
Double_t fBeta
Definition: AliFJWrapper.h:160
Definition: Option.C:68
virtual void SubtractBackground(const Double_t median_pt=-1)
std::vector< fastjet::PseudoJet > fInputGhosts
Definition: AliFJWrapper.h:122
Double_t fGridScatter
Definition: AliFJWrapper.h:151
virtual void Clear(const Option_t *="")
Bool_t fLegacyMode
Definition: AliFJWrapper.h:182
fastjet::VoronoiAreaSpec * fVorAreaSpec
Definition: AliFJWrapper.h:129
Int_t mode
Definition: anaM.C:40
Double_t fDRStep
Definition: AliFJWrapper.h:187
std::vector< double > fGRDenominatorSub
Definition: AliFJWrapper.h:191
Bool_t GetLegacyMode()
Definition: AliFJWrapper.h:46
fastjet::GhostedAreaSpec * fGhostedAreaSpec
Definition: AliFJWrapper.h:130
Bool_t fUseArea4Vector
Definition: AliFJWrapper.h:157
void SetNRepeats(Int_t nrepeat)
Definition: AliFJWrapper.h:98
virtual std::vector< double > GetGRNumeratorSub() const
Definition: AliFJWrapper.h:70
virtual ~AliFJWrapper()
std::vector< double > fGRDenominator
Definition: AliFJWrapper.h:189
virtual Int_t DoGenericSubtractionJet1subjettiness_kt()
Double_t fGhostArea
Definition: AliFJWrapper.h:146
virtual void GetMedianAndSigma(Double_t &median, Double_t &sigma, Int_t remove=0) const
TString fName
Definition: AliFJWrapper.h:119
void SetGhostArea(Double_t gharea)
Definition: AliFJWrapper.h:99
Double_t fMedUsedForBgSub
Definition: AliFJWrapper.h:156
std::vector< fastjet::PseudoJet > GetJetConstituents(UInt_t idx) const
fastjet::JetDefinition * fJetDef
Definition: AliFJWrapper.h:131
virtual Int_t DoGenericSubtractionJetSigma2()
void SetupStrategyfromOpt(const char *option)
const std::vector< fastjet::PseudoJet > & GetInputVectors() const
Definition: AliFJWrapper.h:31
fastjet::AreaDefinition * fAreaDef
Definition: AliFJWrapper.h:128
virtual std::vector< double > GetGRDenominatorSub() const
Definition: AliFJWrapper.h:71
Bool_t fUseExternalBkg
Definition: AliFJWrapper.h:183
void SetTitle(const char *title)
Definition: AliFJWrapper.h:93
Double_t NSubjettiness(Int_t N, Int_t Algorithm, Double_t Radius, Double_t Beta, Int_t Option=0, Int_t Measure=0, Double_t Beta_SD=0, Double_t ZCut=0.1)
void SetPluginAlgor(Int_t plugin)
Definition: AliFJWrapper.h:105
virtual Int_t DoConstituentSubtraction()
Double_t fRhom
Definition: AliFJWrapper.h:185
virtual std::vector< double > GetGRDenominator() const
Definition: AliFJWrapper.h:69
void SetRhoRhom(Double_t rho, Double_t rhom)
Definition: AliFJWrapper.h:115
void SetGridScatter(Double_t gridSc)
Definition: AliFJWrapper.h:102
std::vector< fastjet::PseudoJet > fConstituentSubtrJets
Definition: AliFJWrapper.h:126
void SetMeanGhostKt(Double_t meankt)
Definition: AliFJWrapper.h:104
const char Option_t
Definition: External.C:48
virtual void AddInputVector(Double_t px, Double_t py, Double_t pz, Double_t E, Int_t index=-99999)
virtual void AddInputGhost(Double_t px, Double_t py, Double_t pz, Double_t E, Int_t index=-99999)
virtual Int_t DoGenericSubtractionJet3subjettiness_kt()
bool Bool_t
Definition: External.C:53
void SetLegacyFJ()
std::vector< fastjet::PseudoJet > GetFilteredJetConstituents(UInt_t idx) const
void SetKtScatter(Double_t ktSc)
Definition: AliFJWrapper.h:103
virtual std::vector< double > GetGRNumerator() const
Definition: AliFJWrapper.h:68
virtual const char * ClassName() const
Definition: AliFJWrapper.h:23
virtual void RemoveLastInputVector()
void SetAreaType(const fastjet::AreaType &atype)
Definition: AliFJWrapper.h:97
Double_t fRho
Definition: AliFJWrapper.h:184
Bool_t GetDoFilterArea()
Definition: AliFJWrapper.h:47
virtual Int_t DoGenericSubtractionJetLeSub()
virtual Int_t DoGenericSubtractionJetConstituent()
const char * GetName() const
Definition: AliFJWrapper.h:38
const std::vector< fastjet::PseudoJet > & GetFilteredJets() const
Definition: AliFJWrapper.h:34
virtual Int_t DoGenericSubtractionJet2subjettiness_kt()
std::vector< double > fGRNumerator
Definition: AliFJWrapper.h:188
Double_t fMeanGhostKt
Definition: AliFJWrapper.h:153
fastjet::Strategy fStrategy
Definition: AliFJWrapper.h:141
TString fTitle
Definition: AliFJWrapper.h:120