]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGJE/EMCALJetTasks/AliFJWrapper.h
add GR for jet derivative method
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliFJWrapper.h
CommitLineData
48d874ff 1#ifndef AliFJWrapper_H
2#define AliFJWrapper_H
3
4// $Id$
5
a9c8e539 6#if !defined(__CINT__)
48d874ff 7
8#include <vector>
a9c8e539 9#include <TString.h>
48d874ff 10#include "AliLog.h"
11#include "FJ_includes.h"
8082a80b 12#include "AliJetShape.h"
48d874ff 13
a9c8e539 14class AliFJWrapper
48d874ff 15{
16 public:
17 AliFJWrapper(const char *name, const char *title);
18 virtual ~AliFJWrapper();
a9c8e539 19
48d874ff 20 virtual void AddInputVector (Double_t px, Double_t py, Double_t pz, Double_t E, Int_t index = -99999);
21 virtual void AddInputVector (const fastjet::PseudoJet& vec, Int_t index = -99999);
22 virtual void AddInputVectors(const std::vector<fastjet::PseudoJet>& vecs, Int_t offsetIndex = -99999);
a9c8e539 23 virtual const char *ClassName() const { return "AliFJWrapper"; }
24 virtual void Clear(const Option_t* /*opt*/ = "");
48d874ff 25 virtual void CopySettingsFrom (const AliFJWrapper& wrapper);
26 virtual void GetMedianAndSigma(Double_t& median, Double_t& sigma, Int_t remove = 0) const;
8082a80b 27 fastjet::ClusterSequenceArea* GetClusterSequence() const { return fClustSeq; }
28 const std::vector<fastjet::PseudoJet> GetInputVectors() const { return fInputVectors; }
29 const std::vector<fastjet::PseudoJet> GetInclusiveJets() const { return fInclusiveJets; }
30 std::vector<fastjet::PseudoJet> GetJetConstituents(UInt_t idx) const;
31 Double_t GetMedianUsedForBgSubtraction() const { return fMedUsedForBgSub; }
32 const char* GetName() const { return fName; }
33 const char* GetTitle() const { return fTitle; }
34 Double_t GetJetArea (UInt_t idx) const;
35 fastjet::PseudoJet GetJetAreaVector (UInt_t idx) const;
36 Double_t GetJetSubtractedPt (UInt_t idx) const;
37 virtual std::vector<double> GetSubtractedJetsPts(Double_t median_pt = -1, Bool_t sorted = kFALSE);
38 Bool_t GetLegacyMode() { return fLegacyMode; }
93daf3f7 39#ifdef FASTJET_VERSION
8082a80b 40 const std::vector<fastjet::contrib::GenericSubtractorInfo> GetGenSubtractorInfoJetMass() const {return fGenSubtractorInfoJetMass;}
41 const std::vector<fastjet::PseudoJet> GetConstituentSubtrJets() const { return fConstituentSubtrJets; }
42#endif
48d874ff 43
44 virtual Int_t Run();
8082a80b 45 virtual Int_t DoGenericSubtractionJetMass();
46 virtual Int_t DoConstituentSubtraction();
48d874ff 47
48 void SetStrategy(const fastjet::Strategy &strat) { fStrategy = strat; }
49 void SetAlgorithm(const fastjet::JetAlgorithm &algor) { fAlgor = algor; }
50 void SetRecombScheme(const fastjet::RecombinationScheme &scheme) { fScheme = scheme; }
51 void SetAreaType(const fastjet::AreaType &atype) { fAreaType = atype; }
52 void SetNRepeats(Int_t nrepeat) { fNGhostRepeats = nrepeat; }
53 void SetGhostArea(Double_t gharea) { fGhostArea = gharea; }
54 void SetMaxRap(Double_t maxrap) { fMaxRap = maxrap; }
55 void SetR(Double_t r) { fR = r; }
56 void SetGridScatter(Double_t gridSc) { fGridScatter = gridSc; }
57 void SetKtScatter(Double_t ktSc) { fKtScatter = ktSc; }
58 void SetMeanGhostKt(Double_t meankt) { fMeanGhostKt = meankt; }
59 void SetPluginAlgor(Int_t plugin) { fPluginAlgor = plugin; }
60 void SetUseArea4Vector(Bool_t useA4v) { fUseArea4Vector = useA4v; }
61 void SetupAlgorithmfromOpt(const char *option);
62 void SetupAreaTypefromOpt(const char *option);
63 void SetupSchemefromOpt(const char *option);
64 void SetupStrategyfromOpt(const char *option);
5d0a5007 65 void SetLegacyMode (Bool_t mode) { fLegacyMode ^= mode; }
66 void SetLegacyFJ();
8082a80b 67 void SetUseExternalBkg(Bool_t b, Double_t rho, Double_t rhom) { fUseExternalBkg = b; fRho = rho; fRhom = rhom;}
48d874ff 68
69 protected:
8082a80b 70 TString fName; //!
71 TString fTitle; //!
72 std::vector<fastjet::PseudoJet> fInputVectors; //!
73 std::vector<fastjet::PseudoJet> fInclusiveJets; //!
74 std::vector<double> fSubtractedJetsPt; //!
75 std::vector<fastjet::PseudoJet> fConstituentSubtrJets; //!
76 fastjet::AreaDefinition *fAreaDef; //!
77 fastjet::VoronoiAreaSpec *fVorAreaSpec; //!
78 fastjet::GhostedAreaSpec *fGhostedAreaSpec; //!
79 fastjet::JetDefinition *fJetDef; //!
80 fastjet::JetDefinition::Plugin *fPlugin; //!
81 fastjet::RangeDefinition *fRange; //!
82 fastjet::ClusterSequenceArea *fClustSeq; //!
83 fastjet::Strategy fStrategy; //!
84 fastjet::JetAlgorithm fAlgor; //!
85 fastjet::RecombinationScheme fScheme; //!
86 fastjet::AreaType fAreaType; //!
87 Int_t fNGhostRepeats; //!
88 Double_t fGhostArea; //!
89 Double_t fMaxRap; //!
90 Double_t fR; //!
48d874ff 91 // no setters for the moment - used default values in the constructor
8082a80b 92 Double_t fGridScatter; //!
93 Double_t fKtScatter; //!
94 Double_t fMeanGhostKt; //!
95 Int_t fPluginAlgor; //!
48d874ff 96 // extra parameters
8082a80b 97 Double_t fMedUsedForBgSub; //!
98 Bool_t fUseArea4Vector; //!
5d0a5007 99#ifdef FASTJET_VERSION
8082a80b 100 fastjet::JetMedianBackgroundEstimator *fBkrdEstimator; //!
93daf3f7 101 //from contrib package
8082a80b 102 fastjet::contrib::GenericSubtractor *fGenSubtractor; //!
103 std::vector<fastjet::contrib::GenericSubtractorInfo> fGenSubtractorInfoJetMass; //!
5d0a5007 104#endif
8082a80b 105 Bool_t fLegacyMode; //!
106 Bool_t fUseExternalBkg; //!
107 Double_t fRho; // pT background density
108 Double_t fRhom; // mT background density
48d874ff 109
110 virtual void SubtractBackground(const Double_t median_pt = -1);
111
112 private:
113 AliFJWrapper();
114 AliFJWrapper(const AliFJWrapper& wrapper);
115 AliFJWrapper& operator = (const AliFJWrapper& wrapper);
116};
117#endif
118#endif
119
120#ifdef AliFJWrapper_CXX
121#undef AliFJWrapper_CXX
122
123#if defined __GNUC__
124#pragma GCC system_header
125#endif
126
127namespace fj = fastjet;
128
129//_________________________________________________________________________________________________
130AliFJWrapper::AliFJWrapper(const char *name, const char *title)
a9c8e539 131 :
132 fName (name)
133 , fTitle (title)
48d874ff 134 , fInputVectors ( )
135 , fInclusiveJets ( )
136 , fSubtractedJetsPt ( )
8082a80b 137 , fConstituentSubtrJets ( )
48d874ff 138 , fAreaDef (0)
139 , fVorAreaSpec (0)
140 , fGhostedAreaSpec (0)
141 , fJetDef (0)
142 , fPlugin (0)
143 , fRange (0)
144 , fClustSeq (0)
145 , fStrategy (fj::Best)
146 , fAlgor (fj::kt_algorithm)
147 , fScheme (fj::BIpt_scheme)
148 , fAreaType (fj::active_area)
149 , fNGhostRepeats (1)
6d75c7ee 150 , fGhostArea (0.005)
48d874ff 151 , fMaxRap (1.)
152 , fR (0.4)
153 , fGridScatter (1.0)
7dc91d51 154 , fKtScatter (0.1)
48d874ff 155 , fMeanGhostKt (1e-100)
156 , fPluginAlgor (0)
157 , fMedUsedForBgSub (0)
158 , fUseArea4Vector (kFALSE)
5d0a5007 159#ifdef FASTJET_VERSION
160 , fBkrdEstimator (0)
8082a80b 161 , fGenSubtractor (0)
162 , fGenSubtractorInfoJetMass ( )
5d0a5007 163#endif
164 , fLegacyMode (false)
8082a80b 165 , fUseExternalBkg (false)
166 , fRho (0)
167 , fRhom (0)
48d874ff 168{
169 // Constructor.
170}
171
172//_________________________________________________________________________________________________
173AliFJWrapper::~AliFJWrapper()
174{
175 // Destructor.
176
177 delete fAreaDef;
178 delete fVorAreaSpec;
179 delete fGhostedAreaSpec;
180 delete fJetDef;
181 delete fPlugin;
182 delete fRange;
5d0a5007 183 delete fClustSeq;
184#ifdef FASTJET_VERSION
8082a80b 185 if (fBkrdEstimator) delete fBkrdEstimator;
186 if (fGenSubtractor) delete fGenSubtractor;
5d0a5007 187#endif
48d874ff 188}
189
190//_________________________________________________________________________________________________
191void AliFJWrapper::CopySettingsFrom(const AliFJWrapper& wrapper)
192{
193 // Copy some settings.
194 // You very often want to keep most of the settings
195 // but change only the algorithm or R - do it after call to this function
196
8082a80b 197 fStrategy = wrapper.fStrategy;
198 fAlgor = wrapper.fAlgor;
199 fScheme = wrapper.fScheme;
200 fAreaType = wrapper.fAreaType;
201 fNGhostRepeats = wrapper.fNGhostRepeats;
202 fGhostArea = wrapper.fGhostArea;
203 fMaxRap = wrapper.fMaxRap;
204 fR = wrapper.fR;
205 fGridScatter = wrapper.fGridScatter;
206 fKtScatter = wrapper.fKtScatter;
207 fMeanGhostKt = wrapper.fMeanGhostKt;
208 fPluginAlgor = wrapper.fPluginAlgor;
209 fUseArea4Vector = wrapper.fUseArea4Vector;
210 fLegacyMode = wrapper.fLegacyMode;
211 fUseExternalBkg = wrapper.fUseExternalBkg;
212 fRho = wrapper.fRho;
213 fRhom = wrapper.fRhom;
48d874ff 214}
215
216//_________________________________________________________________________________________________
a9c8e539 217void AliFJWrapper::Clear(const Option_t */*opt*/)
48d874ff 218{
219 // Simply clear the input vectors.
220 // Make sure done on every event if the instance is reused
221 // Reset the median to zero.
222
223 fInputVectors.clear();
224 fMedUsedForBgSub = 0;
225
226 // for the moment brute force delete everything
227 delete fAreaDef; fAreaDef = 0;
228 delete fVorAreaSpec; fVorAreaSpec = 0;
229 delete fGhostedAreaSpec; fGhostedAreaSpec = 0;
230 delete fJetDef; fJetDef = 0;
231 delete fPlugin; fPlugin = 0;
232 delete fRange; fRange = 0;
233 delete fClustSeq; fClustSeq = 0;
5d0a5007 234#ifdef FASTJET_VERSION
8082a80b 235 if (fBkrdEstimator) delete fBkrdEstimator ; fBkrdEstimator = 0;
236 if (fGenSubtractor) delete fGenSubtractor ; fGenSubtractor = 0;
5d0a5007 237#endif
48d874ff 238}
239
240//_________________________________________________________________________________________________
241void AliFJWrapper::AddInputVector(Double_t px, Double_t py, Double_t pz, Double_t E, Int_t index)
242{
243 // Make the input pseudojet.
244
245 fastjet::PseudoJet inVec(px, py, pz, E);
246
247 if (index > -99999) {
248 inVec.set_user_index(index);
249 } else {
250 inVec.set_user_index(fInputVectors.size());
251 }
252
253 // add to the fj container of input vectors
254 fInputVectors.push_back(inVec);
255}
256
257//_________________________________________________________________________________________________
258void AliFJWrapper::AddInputVector(const fj::PseudoJet& vec, Int_t index)
259{
260 // Add an input pseudojet.
261
262 fj::PseudoJet inVec = vec;
263
264 if (index > -99999) {
265 inVec.set_user_index(index);
266 } else {
267 inVec.set_user_index(fInputVectors.size());
268 }
269
270 // add to the fj container of input vectors
271 fInputVectors.push_back(inVec);
272}
273
274//_________________________________________________________________________________________________
275void AliFJWrapper::AddInputVectors(const std::vector<fj::PseudoJet>& vecs, Int_t offsetIndex)
276{
277 // Add the input from vector of pseudojets.
278
279 for (UInt_t i = 0; i < vecs.size(); ++i) {
280 fj::PseudoJet inVec = vecs[i];
48d874ff 281 if (offsetIndex > -99999)
282 inVec.set_user_index(fInputVectors.size() + offsetIndex);
283 // add to the fj container of input vectors
48d874ff 284 fInputVectors.push_back(inVec);
285 }
286}
287
288//_________________________________________________________________________________________________
289Double_t AliFJWrapper::GetJetArea(UInt_t idx) const
290{
291 // Get the jet area.
292
293 Double_t retval = -1; // really wrong area..
294 if ( idx < fInclusiveJets.size() ) {
7dc91d51 295 retval = fClustSeq->area(fInclusiveJets[idx]);
48d874ff 296 } else {
297 AliError(Form("[e] ::GetJetArea wrong index: %d",idx));
298 }
299 return retval;
300}
301
4e22c11c 302//_________________________________________________________________________________________________
303fastjet::PseudoJet AliFJWrapper::GetJetAreaVector(UInt_t idx) const
304{
305 // Get the jet area as vector.
306 fastjet::PseudoJet retval;
307 if ( idx < fInclusiveJets.size() ) {
308 retval = fClustSeq->area_4vector(fInclusiveJets[idx]);
309 } else {
310 AliError(Form("[e] ::GetJetArea wrong index: %d",idx));
311 }
312 return retval;
313}
314
48d874ff 315//_________________________________________________________________________________________________
316std::vector<double> AliFJWrapper::GetSubtractedJetsPts(Double_t median_pt, Bool_t sorted)
317{
318 // Get subtracted jets pTs, returns vector.
319
320 SubtractBackground(median_pt);
321
322 if (kTRUE == sorted) {
323 std::sort(fSubtractedJetsPt.begin(), fSubtractedJetsPt.begin());
324 }
325 return fSubtractedJetsPt;
326}
327
328//_________________________________________________________________________________________________
329Double_t AliFJWrapper::GetJetSubtractedPt(UInt_t idx) const
330{
331 // Get subtracted jets pTs, returns Double_t.
332
333 Double_t retval = -99999.; // really wrong pt..
334 if ( idx < fSubtractedJetsPt.size() ) {
335 retval = fSubtractedJetsPt[idx];
336 }
337 return retval;
338}
339
340//_________________________________________________________________________________________________
341std::vector<fastjet::PseudoJet>
342AliFJWrapper::GetJetConstituents(UInt_t idx) const
343{
344 // Get jets constituents.
345
346 std::vector<fastjet::PseudoJet> retval;
347
348 if ( idx < fInclusiveJets.size() ) {
349 retval = fClustSeq->constituents(fInclusiveJets[idx]);
350 } else {
351 AliError(Form("[e] ::GetJetConstituents wrong index: %d",idx));
352 }
353
354 return retval;
355}
356
357//_________________________________________________________________________________________________
358void AliFJWrapper::GetMedianAndSigma(Double_t &median, Double_t &sigma, Int_t remove) const
359{
360 // Get the median and sigma from fastjet.
361 // User can also do it on his own because the cluster sequence is exposed (via a getter)
362
363 if (!fClustSeq) {
364 AliError("[e] Run the jfinder first.");
b0a20945 365 return;
48d874ff 366 }
367
368 Double_t mean_area = 0;
369 try {
370 if(0 == remove) {
371 fClustSeq->get_median_rho_and_sigma(*fRange, fUseArea4Vector, median, sigma, mean_area);
372 } else {
373 std::vector<fastjet::PseudoJet> input_jets = sorted_by_pt(fClustSeq->inclusive_jets());
374 input_jets.erase(input_jets.begin(), input_jets.begin() + remove);
5d0a5007 375 fClustSeq->get_median_rho_and_sigma(input_jets, *fRange, fUseArea4Vector, median, sigma, mean_area);
48d874ff 376 input_jets.clear();
377 }
378 } catch (fj::Error) {
379 AliError(" [w] FJ Exception caught.");
380 median = -1.;
381 sigma = -1;
382 }
383}
384
385//_________________________________________________________________________________________________
386Int_t AliFJWrapper::Run()
387{
388 // Run the actual jet finder.
389
390 if (fAreaType == fj::voronoi_area) {
391 // Rfact - check dependence - default is 1.
392 // NOTE: hardcoded variable!
393 fVorAreaSpec = new fj::VoronoiAreaSpec(1.);
394 fAreaDef = new fj::AreaDefinition(*fVorAreaSpec);
395 } else {
396 fGhostedAreaSpec = new fj::GhostedAreaSpec(fMaxRap,
397 fNGhostRepeats,
398 fGhostArea,
399 fGridScatter,
400 fKtScatter,
401 fMeanGhostKt);
402
403 fAreaDef = new fj::AreaDefinition(*fGhostedAreaSpec, fAreaType);
404 }
405
406 // this is acceptable by fastjet:
407 fRange = new fj::RangeDefinition(fMaxRap - 0.95 * fR);
408
409 if (fAlgor == fj::plugin_algorithm) {
410 if (fPluginAlgor == 0) {
411 // SIS CONE ALGOR
412 // NOTE: hardcoded split parameter
413 Double_t overlap_threshold = 0.75; // NOTE: this actually splits a lot: thr/min(pt1,pt2)
414 fPlugin = new fj::SISConePlugin(fR,
415 overlap_threshold,
416 0, //search of stable cones - zero = until no more
417 1.0); // this should be seed effectively for proto jets
418 fJetDef = new fastjet::JetDefinition(fPlugin);
c515ed6d 419 } else if (fPluginAlgor == 1) {
420 // CDF cone
421 // NOTE: hardcoded split parameter
422 Double_t overlap_threshold = 0.75; // NOTE: this actually splits a lot: thr/min(pt1,pt2)
423 fPlugin = new fj::CDFMidPointPlugin(fR,
424 overlap_threshold,
425 1.0, //search of stable cones - zero = until no more
426 1.0); // this should be seed effectively for proto jets
427 fJetDef = new fastjet::JetDefinition(fPlugin);
48d874ff 428 } else {
429 AliError("[e] Unrecognized plugin number!");
430 }
431 } else {
432 fJetDef = new fj::JetDefinition(fAlgor, fR, fScheme, fStrategy);
433 }
434
435 try {
436 fClustSeq = new fj::ClusterSequenceArea(fInputVectors, *fJetDef, *fAreaDef);
437 } catch (fj::Error) {
438 AliError(" [w] FJ Exception caught.");
439 return -1;
440 }
441
5d0a5007 442 // FJ3 :: Define an JetMedianBackgroundEstimator just in case it will be used
443#ifdef FASTJET_VERSION
8082a80b 444 fBkrdEstimator = new fj::JetMedianBackgroundEstimator(fj::SelectorAbsRapMax(fMaxRap));
5d0a5007 445#endif
446
447 if (fLegacyMode) { SetLegacyFJ(); } // for FJ 2.x even if fLegacyMode is set, SetLegacyFJ is dummy
448
48d874ff 449 // inclusive jets:
450 fInclusiveJets.clear();
7dc91d51 451 fInclusiveJets = fClustSeq->inclusive_jets(0.0);
48d874ff 452
453 return 0;
454}
455
5d0a5007 456//_________________________________________________________________________________________________
457void AliFJWrapper::SetLegacyFJ()
458{
459 // This methods enable legacy behaviour (FastJet 2.x) when AliROOT is compiled with FastJet 3.x
460#ifdef FASTJET_VERSION
461 std::cout << "WARNING! Setting FastJet in legacy mode" << std::endl;
462 if (fGhostedAreaSpec) { fGhostedAreaSpec->set_fj2_placement(kTRUE); }
463 if (fBkrdEstimator) {
464 fBkrdEstimator->set_provide_fj2_sigma(kTRUE);
465 fBkrdEstimator->set_use_area_4vector(kFALSE);
466 }
467#endif
468}
469
48d874ff 470//_________________________________________________________________________________________________
471void AliFJWrapper::SubtractBackground(Double_t median_pt)
472{
473 // Subtract the background (specify the value - see below the meaning).
474 // Negative argument means the bg will be determined with the current algorithm
475 // this is the default behavior. Zero is allowed
476 // Note: user may set the switch for area4vector based subtraction.
477
478 Double_t median = 0;
479 Double_t sigma = 0;
480 Double_t mean_area = 0;
481
482 // clear the subtracted jet pt's vector<double>
483 fSubtractedJetsPt.clear();
5d0a5007 484
48d874ff 485 // check what was specified (default is -1)
486 if (median_pt < 0) {
487 try {
5d0a5007 488 fClustSeq->get_median_rho_and_sigma(*fRange, fUseArea4Vector, median, sigma, mean_area);
48d874ff 489 }
5d0a5007 490
48d874ff 491 catch (fj::Error) {
492 AliError(" [w] FJ Exception caught.");
493 median = -9999.;
494 sigma = -1;
495 fMedUsedForBgSub = median;
496 return;
497 }
498 fMedUsedForBgSub = median;
499 } else {
500 // we do not know the sigma in this case
501 sigma = -1;
502 if (0.0 == median_pt) {
503 // AliWarning(Form(" [w] Median specified for bg subtraction is ZERO: %f", median_pt )) << endl;
504 fMedUsedForBgSub = 0.;
505 } else {
506 fMedUsedForBgSub = median_pt;
507 }
508 }
509
510 // subtract:
511 for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
5d0a5007 512 if ( fUseArea4Vector ) {
48d874ff 513 // subtract the background using the area4vector
514 fj::PseudoJet area4v = fClustSeq->area_4vector(fInclusiveJets[i]);
515 fj::PseudoJet jet_sub = fInclusiveJets[i] - area4v * fMedUsedForBgSub;
516 fSubtractedJetsPt.push_back(jet_sub.perp()); // here we put only the pt of the jet - note: this can be negative
517 } else {
518 // subtract the background using scalars
519 // fj::PseudoJet jet_sub = fInclusiveJets[i] - area * fMedUsedForBgSub_;
520 Double_t area = fClustSeq->area(fInclusiveJets[i]);
521 // standard subtraction
522 Double_t pt_sub = fInclusiveJets[i].perp() - fMedUsedForBgSub * area;
523 fSubtractedJetsPt.push_back(pt_sub); // here we put only the pt of the jet - note: this can be negative
524 }
525 }
526}
527
8082a80b 528//_________________________________________________________________________________________________
529Int_t AliFJWrapper::DoGenericSubtractionJetMass() {
530 //Do generic subtraction for jet mass
93daf3f7 531#ifdef FASTJET_VERSION
8082a80b 532 if(fUseExternalBkg) fGenSubtractor = new fj::contrib::GenericSubtractor(fRho,fRhom);
533 else fGenSubtractor = new fj::contrib::GenericSubtractor(fBkrdEstimator);
534
535 // Define jet shape
536 AliJetShapeMass shapeMass;
537
538 // clear the generic subtractor info vector
539 fGenSubtractorInfoJetMass.clear();
540 for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
541 fj::contrib::GenericSubtractorInfo info;
542 if(fInclusiveJets[i].perp()>0.)
543 double subtracted_shape = (*fGenSubtractor)(shapeMass, fInclusiveJets[i], info);
544 fGenSubtractorInfoJetMass.push_back(info);
545 }
546
547#endif
548 return 0;
549}
550
551//_________________________________________________________________________________________________
552Int_t AliFJWrapper::DoConstituentSubtraction() {
553 //Do constituent subtraction
93daf3f7 554#ifdef FASTJET_VERSION
8082a80b 555 fj::contrib::ConstituentSubtractor *subtractor;
8ed8c5bb 556 if(fUseExternalBkg)
8082a80b 557 subtractor = new fj::contrib::ConstituentSubtractor(fRho,fRhom,kFALSE,kTRUE);
8082a80b 558 else subtractor = new fj::contrib::ConstituentSubtractor(fBkrdEstimator);
559
560 //clear constituent subtracted jets
561 fConstituentSubtrJets.clear();
562 for (unsigned i = 0; i < fInclusiveJets.size(); i++) {
563 fj::PseudoJet subtracted_jet(0.,0.,0.,0.);
564 if(fInclusiveJets[i].perp()>0.)
565 subtracted_jet = (*subtractor)(fInclusiveJets[i]);
566 fConstituentSubtrJets.push_back(subtracted_jet);
567 }
568 if(subtractor) delete subtractor;
569
570#endif
571 return 0;
572}
573
48d874ff 574//_________________________________________________________________________________________________
575void AliFJWrapper::SetupAlgorithmfromOpt(const char *option)
576{
577 // Setup algorithm from char.
578
579 std::string opt(option);
580
581 if (!opt.compare("kt")) fAlgor = fj::kt_algorithm;
582 if (!opt.compare("antikt")) fAlgor = fj::antikt_algorithm;
583 if (!opt.compare("cambridge")) fAlgor = fj::cambridge_algorithm;
584 if (!opt.compare("genkt")) fAlgor = fj::genkt_algorithm;
585 if (!opt.compare("cambridge_passive")) fAlgor = fj::cambridge_for_passive_algorithm;
586 if (!opt.compare("genkt_passive")) fAlgor = fj::genkt_for_passive_algorithm;
587 if (!opt.compare("ee_kt")) fAlgor = fj::ee_kt_algorithm;
588 if (!opt.compare("ee_genkt")) fAlgor = fj::ee_genkt_algorithm;
589 if (!opt.compare("plugin")) fAlgor = fj::plugin_algorithm;
590}
591
592//_________________________________________________________________________________________________
593void AliFJWrapper::SetupAreaTypefromOpt(const char *option)
594{
595 // Setup area type from char.
596
597 std::string opt(option);
598
599 if (!opt.compare("active")) fAreaType = fj::active_area;
600 if (!opt.compare("invalid")) fAreaType = fj::invalid_area;
601 if (!opt.compare("active_area_explicit_ghosts")) fAreaType = fj::active_area_explicit_ghosts;
602 if (!opt.compare("one_ghost_passive")) fAreaType = fj::one_ghost_passive_area;
603 if (!opt.compare("passive")) fAreaType = fj::passive_area;
604 if (!opt.compare("voronoi")) fAreaType = fj::voronoi_area;
605}
606
607//_________________________________________________________________________________________________
608void AliFJWrapper::SetupSchemefromOpt(const char *option)
609{
610 //
611 // setup scheme from char
612 //
613
614 std::string opt(option);
615
616 if (!opt.compare("BIpt")) fScheme = fj::BIpt_scheme;
617 if (!opt.compare("BIpt2")) fScheme = fj::BIpt2_scheme;
618 if (!opt.compare("E")) fScheme = fj::E_scheme;
619 if (!opt.compare("pt")) fScheme = fj::pt_scheme;
620 if (!opt.compare("pt2")) fScheme = fj::pt2_scheme;
621 if (!opt.compare("Et")) fScheme = fj::Et_scheme;
622 if (!opt.compare("Et2")) fScheme = fj::Et2_scheme;
623}
624
625//_________________________________________________________________________________________________
626void AliFJWrapper::SetupStrategyfromOpt(const char *option)
627{
628 // Setup strategy from char.
629
630 std::string opt(option);
631
632 if (!opt.compare("Best")) fStrategy = fj::Best;
633 if (!opt.compare("N2MinHeapTiled")) fStrategy = fj::N2MinHeapTiled;
634 if (!opt.compare("N2Tiled")) fStrategy = fj::N2Tiled;
635 if (!opt.compare("N2PoorTiled")) fStrategy = fj::N2PoorTiled;
636 if (!opt.compare("N2Plain")) fStrategy = fj::N2Plain;
637 if (!opt.compare("N3Dumb")) fStrategy = fj::N3Dumb;
638 if (!opt.compare("NlnN")) fStrategy = fj::NlnN;
639 if (!opt.compare("NlnN3pi")) fStrategy = fj::NlnN3pi;
640 if (!opt.compare("NlnN4pi")) fStrategy = fj::NlnN4pi;
641 if (!opt.compare("NlnNCam4pi")) fStrategy = fj::NlnNCam4pi;
642 if (!opt.compare("NlnNCam2pi2R")) fStrategy = fj::NlnNCam2pi2R;
643 if (!opt.compare("NlnNCam")) fStrategy = fj::NlnNCam;
644 if (!opt.compare("plugin")) fStrategy = fj::plugin_strategy;
645}
646#endif