30 #include "fastjet/ClusterSequence.hh"
31 #include "fastjet/ClusterSequenceArea.hh"
32 #include "fastjet/Selector.hh"
33 #include "fastjet/SISConePlugin.hh"
36 using namespace fastjet;
38 FASTJET_BEGIN_NAMESPACE
44 auto_ptr<JetDefinition::Plugin>
plugin;
46 auto_ptr<ClusterSequence>
cs;
53 valarray<double> mom(4);
54 for (
int j=0;
j<=3;
j++) {
60 psjet.set_user_index(
i);
68 for (
int i=0;
i<njets;
i++) {
69 for (
int j=0;
j<=3;
j++) {
83 double * f77jets,
int & njets,
int * whichjet,
84 const double & ghost_maxrap = 0.0,
85 const int & nrepeat = 0,
const double & ghost_area = 0.0) {
91 if ( ghost_maxrap == 0.0 ) {
96 GhostedAreaSpec area_spec(ghost_maxrap,nrepeat,ghost_area);
97 AreaDefinition area_def(active_area, area_spec);
101 jets = sorted_by_pt(
cs->inclusive_jets(ptmin));
108 for(
int ii=0; ii<
npart; ++ii) whichjet[ii]=0;
111 for (
int kk=0;
kk<njets; ++
kk) {
112 vector<PseudoJet> constit =
cs->constituents(
jets[
kk]);
113 for(
unsigned int ll=0; ll<constit.size(); ++ll)
114 whichjet[constit[ll].user_index()]=kk+1;
120 FASTJET_END_NAMESPACE
122 using namespace fastjet::fwrapper;
158 const double &
R,
const double &
f,
159 double * f77jets,
int & njets,
int * whichjet) {
162 plugin.reset(
new SISConePlugin(R,f));
207 const double &
R,
const double &
f,
208 const double & ghost_rapmax,
const int & nrepeat,
const double & ghost_area,
209 double * f77jets,
int & njets,
int * whichjet) {
212 plugin.reset(
new SISConePlugin(R,f));
259 const double &
R,
const double & ptjetmin,
261 double * f77jets,
int & njets,
int * whichjet) {
265 jet_def = JetDefinition(kt_algorithm, R);
266 }
else if (palg == 0.0) {
267 jet_def = JetDefinition(cambridge_algorithm, R);
268 }
else if (palg == -1.0) {
269 jet_def = JetDefinition(antikt_algorithm, R);
271 jet_def = JetDefinition(genkt_algorithm, R, palg);
322 const double &
R,
const double & palg,
323 const double & ghost_rapmax,
const int & nrepeat,
const double & ghost_area,
324 double * f77jets,
int & njets,
int * whichjet) {
328 jet_def = JetDefinition(kt_algorithm, R);
329 }
else if (palg == 0.0) {
330 jet_def = JetDefinition(cambridge_algorithm, R);
331 }
else if (palg == -1.0) {
332 jet_def = JetDefinition(antikt_algorithm, R);
334 jet_def = JetDefinition(genkt_algorithm, R, palg);
360 int * constituent_indices,
int & nconstituents) {
361 assert(
cs.get() != 0);
362 assert(ijet > 0 && ijet <=
int(
jets.size()));
364 vector<PseudoJet> constituents =
cs->constituents(
jets[ijet-1]);
366 nconstituents = constituents.size();
367 for (
int i = 0;
i < nconstituents;
i++) {
368 constituent_indices[
i] = constituents[
i].cluster_hist_index()+1;
388 assert(ijet > 0 && ijet <=
int(
jets.size()));
389 const ClusterSequenceAreaBase * csab =
390 dynamic_cast<const ClusterSequenceAreaBase *
>(
cs.get());
393 return csab->area(
jets[ijet-1]);
411 assert(
cs.get() != 0);
412 return cs->exclusive_dmerge(n);
428 assert(
cs.get() != 0);
429 return cs->exclusive_dmerge_max(n);
446 double &
rho,
double & sigma,
double & meanarea) {
447 const ClusterSequenceAreaBase * csab =
448 dynamic_cast<const ClusterSequenceAreaBase *
>(
cs.get());
451 Selector range = SelectorRapRange(rapmin,rapmax) && SelectorPhiRange(phimin,phimax);
452 bool use_area_4vector =
false;
453 csab->get_median_rho_and_sigma(range,use_area_4vector,rho,sigma,meanarea);
455 Error(
"Clustering with area is necessary in order to be able to evaluate rho.");
void fastjetsisconewitharea_(const double *p, const int &npart, const double &R, const double &f, const double &ghost_rapmax, const int &nrepeat, const double &ghost_area, double *f77jets, int &njets, int *whichjet)
auto_ptr< ClusterSequence > cs
void fastjetppgenktwitharea_(const double *p, const int &npart, const double &R, const double &palg, const double &ghost_rapmax, const int &nrepeat, const double &ghost_area, double *f77jets, int &njets, int *whichjet)
double fastjetdmergemax_(const int &n)
void fastjetconstituents_(const int &ijet, int *constituent_indices, int &nconstituents)
void fastjetppgenkt_(const double *p, const int &npart, const double &R, const double &ptjetmin, const double &palg, double *f77jets, int &njets, int *whichjet)
void transfer_input_particles(const double *p, const int &npart)
helper routine to transfer fortran input particles into
Functor that operates on <T>
vector< PseudoJet > input_particles
void fastjetglobalrhoandsigma_(const double &rapmin, const double &rapmax, const double &phimin, const double &phimax, double &rho, double &sigma, double &meanarea)
void transfer_jets(double *f77jets, int &njets)
helper routine to help transfer jets -> f77jets[4*ijet+0..3]
double fastjetarea_(const int &ijet)
void transfer_cluster_transfer(const double *p, const int &npart, const JetDefinition &jet_def, const double &ptmin, double *f77jets, int &njets, int *whichjet, const double &ghost_maxrap=0.0, const int &nrepeat=0, const double &ghost_area=0.0)
void fastjetsiscone_(const double *p, const int &npart, const double &R, const double &f, double *f77jets, int &njets, int *whichjet)
double fastjetdmerge_(const int &n)