|
|
Go to the documentation of this file.
29 #include "fastjet/ClusterSequence.hh"
30 #include "fastjet/ClusterSequenceArea.hh"
31 #include "fastjet/Selector.hh"
32 #include "fastjet/SISConePlugin.hh"
37 FASTJET_BEGIN_NAMESPACE
43 unique_ptr<JetDefinition::Plugin>
plugin;
45 unique_ptr<ClusterSequence>
cs;
52 valarray<double> mom(4);
53 for (
int j=0;
j<=3;
j++) {
59 psjet.set_user_index(
i);
68 for (
int j=0;
j<=3;
j++) {
82 double * f77jets,
int &
njets,
int * whichjet,
83 const double & ghost_maxrap = 0.0,
84 const int & nrepeat = 0,
const double & ghost_area = 0.0) {
90 if ( ghost_maxrap == 0.0 ) {
95 GhostedAreaSpec area_spec(ghost_maxrap,nrepeat,ghost_area);
96 AreaDefinition area_def(active_area, area_spec);
111 vector<PseudoJet> constit =
cs->constituents(
jets[
kk]);
112 for(
unsigned int ll=0; ll<constit.size(); ++ll)
113 whichjet[constit[ll].user_index()]=
kk+1;
119 FASTJET_END_NAMESPACE
121 using namespace fastjet::fwrapper;
157 const double &
R,
const double &
f,
158 double * f77jets,
int &
njets,
int * whichjet) {
161 plugin.reset(
new SISConePlugin(
R,
f));
206 const double &
R,
const double &
f,
207 const double & ghost_rapmax,
const int & nrepeat,
const double & ghost_area,
208 double * f77jets,
int &
njets,
int * whichjet) {
211 plugin.reset(
new SISConePlugin(
R,
f));
258 const double &
R,
const double & ptjetmin,
260 double * f77jets,
int &
njets,
int * whichjet) {
264 jet_def = JetDefinition(kt_algorithm,
R);
265 }
else if (palg == 0.0) {
266 jet_def = JetDefinition(cambridge_algorithm,
R);
267 }
else if (palg == -1.0) {
268 jet_def = JetDefinition(antikt_algorithm,
R);
270 jet_def = JetDefinition(genkt_algorithm,
R, palg);
321 const double &
R,
const double & palg,
322 const double & ghost_rapmax,
const int & nrepeat,
const double & ghost_area,
323 double * f77jets,
int &
njets,
int * whichjet) {
327 jet_def = JetDefinition(kt_algorithm,
R);
328 }
else if (palg == 0.0) {
329 jet_def = JetDefinition(cambridge_algorithm,
R);
330 }
else if (palg == -1.0) {
331 jet_def = JetDefinition(antikt_algorithm,
R);
333 jet_def = JetDefinition(genkt_algorithm,
R, palg);
359 int * constituent_indices,
int & nconstituents) {
363 vector<PseudoJet> constituents =
cs->constituents(
jets[ijet-1]);
365 nconstituents = constituents.size();
366 for (
int i = 0;
i < nconstituents;
i++) {
367 constituent_indices[
i] = constituents[
i].cluster_hist_index()+1;
388 const ClusterSequenceAreaBase * csab =
389 dynamic_cast<const ClusterSequenceAreaBase *>(
cs.get());
390 if (csab !=
nullptr) {
392 return csab->area(
jets[ijet-1]);
411 return cs->exclusive_dmerge(
n);
428 return cs->exclusive_dmerge_max(
n);
445 double & rho,
double & sigma,
double & meanarea) {
446 const ClusterSequenceAreaBase * csab =
447 dynamic_cast<const ClusterSequenceAreaBase *>(
cs.get());
448 if (csab !=
nullptr) {
451 bool use_area_4vector =
false;
452 csab->get_median_rho_and_sigma(
range,use_area_4vector,rho,sigma,meanarea);
454 Error(
"Clustering with area is necessary in order to be able to evaluate rho.");
double fastjetdmergemax_(const int &n)
Functor that operates on <T>
unique_ptr< ClusterSequence > cs
double fastjetdmerge_(const int &n)
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)
void fastjetsiscone_(const double *p, const int &npart, const double &R, const double &f, double *f77jets, int &njets, int *whichjet)
void fastjetglobalrhoandsigma_(const double &rapmin, const double &rapmax, const double &phimin, const double &phimax, double &rho, double &sigma, double &meanarea)
vector< PseudoJet > input_particles
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_jets(double *f77jets, int &njets)
helper routine to help transfer jets -> f77jets[4*ijet+0..3]
edm::ErrorSummaryEntry Error
unique_ptr< JetDefinition::Plugin > plugin
void transfer_input_particles(const double *p, const int &npart)
helper routine to transfer fortran input particles into
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 fastjetconstituents_(const int &ijet, int *constituent_indices, int &nconstituents)
double fastjetarea_(const int &ijet)
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)