CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
fastjetfortran_madfks.cc
Go to the documentation of this file.
1 //STARTHEADER
2 //
3 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
4 //
5 //----------------------------------------------------------------------
6 // This file is part of FastJet.
7 //
8 // FastJet is free software; you can redistribute it and/or modify
9 // it under the terms of the GNU General Public License as published by
10 // the Free Software Foundation; either version 2 of the License, or
11 // (at your option) any later version.
12 //
13 // The algorithms that underlie FastJet have required considerable
14 // development and are described in hep-ph/0512210. If you use
15 // FastJet as part of work towards a scientific publication, please
16 // include a citation to the FastJet paper.
17 //
18 // FastJet is distributed in the hope that it will be useful,
19 // but WITHOUT ANY WARRANTY; without even the implied warranty of
20 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 // GNU General Public License for more details.
22 //
23 // You should have received a copy of the GNU General Public License
24 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
25 //----------------------------------------------------------------------
26 //ENDHEADER
27 
28 #include <iostream>
29 #include "fastjet/ClusterSequence.hh"
30 #include "fastjet/ClusterSequenceArea.hh"
31 #include "fastjet/Selector.hh"
32 #include "fastjet/SISConePlugin.hh"
33 
34 using namespace std;
35 using namespace fastjet;
36 
37 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
38 
41 namespace fwrapper {
42  vector<PseudoJet> input_particles, jets;
43  auto_ptr<JetDefinition::Plugin> plugin;
44  JetDefinition jet_def;
45  auto_ptr<ClusterSequence> cs;
46 
48  void transfer_input_particles(const double * p, const int & npart) {
49  input_particles.resize(0);
50  input_particles.reserve(npart);
51  for (int i=0; i<npart; i++) {
52  valarray<double> mom(4); // mom[0..3]
53  for (int j=0;j<=3; j++) {
54  mom[j] = *(p++);
55  // RF-MZ: reorder the arguments because in madfks energy goes first; in fastjet last.
56  // mom[(j+3) % 4] = *(p++);
57  }
58  PseudoJet psjet(mom);
59  psjet.set_user_index(i);
60  input_particles.push_back(psjet);
61  }
62  }
63 
65  void transfer_jets(double * f77jets, int & njets) {
66  njets = jets.size();
67  for (int i=0; i<njets; i++) {
68  for (int j=0;j<=3; j++) {
69  *f77jets = jets[i][j];
70  // RF-MZ: reorder the arguments because in madfks energy goes first; in fastjet last.
71  //*f77jets = jets[i][(j+3) % 4];
72  f77jets++;
73  }
74  }
75  }
76 
79  void transfer_cluster_transfer(const double * p, const int & npart,
80  const JetDefinition & jet_def,
81  const double & ptmin,
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) {
85 
86  // transfer p[4*ipart+0..3] -> input_particles[i]
87  transfer_input_particles(p, npart);
88 
89  // perform the clustering
90  if ( ghost_maxrap == 0.0 ) {
91  // cluster without areas
92  cs.reset(new ClusterSequence(input_particles,jet_def));
93  } else {
94  // cluster with areas
95  GhostedAreaSpec area_spec(ghost_maxrap,nrepeat,ghost_area);
96  AreaDefinition area_def(active_area, area_spec);
97  cs.reset(new ClusterSequenceArea(input_particles,jet_def,area_def));
98  }
99  // extract jets (pt-ordered)
100  jets = sorted_by_pt(cs->inclusive_jets(ptmin));
101 
102  // transfer jets -> f77jets[4*ijet+0..3]
103  transfer_jets(f77jets, njets);
104 
105  // Determine which parton/particle ended-up in which jet
106  // set all jet entrie to zero first
107  for(int ii=0; ii<npart; ++ii) whichjet[ii]=0;
108 
109  // Loop over jets and find constituents
110  for (int kk=0; kk<njets; ++kk) {
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;
114  }
115 
116  }
117 
118 }
119 FASTJET_END_NAMESPACE
120 
121 using namespace fastjet::fwrapper;
122 
123 
124 extern "C" {
125 
128 //
129 // Corresponds to the following Fortran subroutine
130 // interface structure:
131 //
132 // SUBROUTINE FASTJETSISCONE(P,NPART,R,F,F77JETS,NJETS,WHICHJET)
133 // DOUBLE PRECISION P(4,*), R, F, F77JETS(4,*)
134 // INTEGER NPART, NJETS, WHICHJET(*)
135 //
136 // where on input
137 //
138 // P the input particle 4-momenta
139 // NPART the number of input momenta
140 // R the radius parameter
141 // F the overlap threshold
142 //
143 // and on output
144 //
145 // F77JETS the output jet momenta (whose second dim should be >= NPART)
146 // sorted in order of decreasing p_t.
147 // NJETS the number of output jets
148 // WHICHJET(i) the jet of parton/particle 'i'
149 //
150 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its
151 // momenta as a matrix of the form P(4000,5), whereas this fortran
152 // interface to fastjet expects them as P(4,NPART), i.e. you must take
153 // the transpose of the Pythia array and drop the fifth component
154 // (particle mass).
155 //
156 void fastjetsiscone_(const double * p, const int & npart,
157  const double & R, const double & f,
158  double * f77jets, int & njets, int * whichjet) {
159 
160  // prepare jet def
161  plugin.reset(new SISConePlugin(R,f));
162  jet_def = plugin.get();
163 
164  // do everything
165 // transfer_cluster_transfer(p,npart,jet_def,f77jets,njets,whichjet);
166 }
167 
168 
169 
174 //
175 // Corresponds to the following Fortran subroutine
176 // interface structure:
177 //
178 // SUBROUTINE FASTJETSISCONEWITHAREA(P,NPART,R,F,GHMAXRAP,NREP,GHAREA,F77JETS,NJETS,WHICHJET)
179 // DOUBLE PRECISION P(4,*), R, F, F77JETS(4,*), GHMAXRAP, GHAREA
180 // INTEGER NPART, NJETS, NREP, WHICHJET(*)
181 //
182 // where on input
183 //
184 // P the input particle 4-momenta
185 // NPART the number of input momenta
186 // R the radius parameter
187 // F the overlap threshold
188 // GHMAXRAP the maximum (abs) rapidity covered by ghosts (FastJet default 6.0)
189 // NREP the number of repetitions used to evaluate the area (FastJet default 1)
190 // GHAREA the area of a single ghost (FastJet default 0.01)
191 //
192 // and on output
193 //
194 // F77JETS the output jet momenta (whose second dim should be >= NPART)
195 // sorted in order of decreasing p_t.
196 // NJETS the number of output jets
197 // WHICHJET(i) the jet of parton/particle 'i'
198 //
199 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its
200 // momenta as a matrix of the form P(4000,5), whereas this fortran
201 // interface to fastjet expects them as P(4,NPART), i.e. you must take
202 // the transpose of the Pythia array and drop the fifth component
203 // (particle mass).
204 //
205 void fastjetsisconewitharea_(const double * p, const int & npart,
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) {
209 
210  // prepare jet def
211  plugin.reset(new SISConePlugin(R,f));
212  jet_def = plugin.get();
213 
214  // do everything
215 // transfer_cluster_transfer(p,npart,jet_def,f77jets,njets,whichjet,ghost_rapmax,nrepeat,ghost_area);
216 }
217 
218 
219 
223 //
224 // Corresponds to the following Fortran subroutine
225 // interface structure:
226 //
227 // SUBROUTINE FASTJETPPGENKT(P,NPART,R,PALG,F77JETS,NJETS,WHICHJET)
228 // DOUBLE PRECISION P(4,*), R, PALG, F, F77JETS(4,*)
229 // INTEGER NPART, NJETS, WHICHJET(*)
230 //
231 // where on input
232 //
233 // P the input particle 4-momenta
234 // NPART the number of input momenta
235 // R the radius parameter
236 // PALG the power for the generalised kt alg
237 // (1.0=kt, 0.0=C/A, -1.0 = anti-kt)
238 //
239 // and on output
240 //
241 // F77JETS the output jet momenta (whose second dim should be >= NPART)
242 // sorted in order of decreasing p_t.
243 // NJETS the number of output jets
244 // WHICHJET(i) the jet of parton/particle 'i'
245 //
246 // For the values of PALG that correspond to "standard" cases (1.0=kt,
247 // 0.0=C/A, -1.0 = anti-kt) this routine actually calls the direct
248 // implementation of those algorithms, whereas for other values of
249 // PALG it calls the generalised kt implementation.
250 //
251 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its
252 // momenta as a matrix of the form P(4000,5), whereas this fortran
253 // interface to fastjet expects them as P(4,NPART), i.e. you must take
254 // the transpose of the Pythia array and drop the fifth component
255 // (particle mass).
256 //
257 void fastjetppgenkt_(const double * p, const int & npart,
258  const double & R, const double & ptjetmin,
259  const double & palg,
260  double * f77jets, int & njets, int * whichjet) {
261 
262  // prepare jet def
263  if (palg == 1.0) {
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);
269  } else {
270  jet_def = JetDefinition(genkt_algorithm, R, palg);
271  }
272 
273  // do everything
274  transfer_cluster_transfer(p,npart,jet_def,ptjetmin,f77jets,njets,whichjet);
275 }
276 
277 
283 //
284 // Corresponds to the following Fortran subroutine
285 // interface structure:
286 //
287 // SUBROUTINE FASTJETPPGENKTWITHAREA(P,NPART,R,PALG,GHMAXRAP,NREP,GHAREA,F77JETS,NJETS,WHICHJET)
288 // DOUBLE PRECISION P(4,*), R, PALG, GHMAXRAP, GHAREA, F77JETS(4,*)
289 // INTEGER NPART, NREP, NJETS, WHICHJET(*)
290 //
291 // where on input
292 //
293 // P the input particle 4-momenta
294 // NPART the number of input momenta
295 // R the radius parameter
296 // PALG the power for the generalised kt alg
297 // (1.0=kt, 0.0=C/A, -1.0 = anti-kt)
298 // GHMAXRAP the maximum (abs) rapidity covered by ghosts (FastJet default 6.0)
299 // NREP the number of repetitions used to evaluate the area (FastJet default 1)
300 // GHAREA the area of a single ghost (FastJet default 0.01)
301 //
302 // and on output
303 //
304 // F77JETS the output jet momenta (whose second dim should be >= NPART)
305 // sorted in order of decreasing p_t.
306 // NJETS the number of output jets
307 // WHICHJET(i) the jet of parton/particle 'i'
308 //
309 // For the values of PALG that correspond to "standard" cases (1.0=kt,
310 // 0.0=C/A, -1.0 = anti-kt) this routine actually calls the direct
311 // implementation of those algorithms, whereas for other values of
312 // PALG it calls the generalised kt implementation.
313 //
314 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its
315 // momenta as a matrix of the form P(4000,5), whereas this fortran
316 // interface to fastjet expects them as P(4,NPART), i.e. you must take
317 // the transpose of the Pythia array and drop the fifth component
318 // (particle mass).
319 //
320 void fastjetppgenktwitharea_(const double * p, const int & npart,
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) {
324 
325  // prepare jet def
326  if (palg == 1.0) {
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);
332  } else {
333  jet_def = JetDefinition(genkt_algorithm, R, palg);
334  }
335 
336  // do everything
337 // transfer_cluster_transfer(p,npart,jet_def,f77jets,njets,whichjet,ghost_rapmax,nrepeat,ghost_area);
338 }
339 
340 
349 //
350 // Corresponds to the following Fortran subroutine
351 // interface structure:
352 //
353 // SUBROUTINE FASTJETCONSTITUENTS(IJET,CONSTITUENT_INDICES,NCONSTITUENTS)
354 // INTEGER IJET
355 // INTEGER CONSTITUENT_INDICES(*)
356 // INTEGER nconstituents
357 //
358 void fastjetconstituents_(const int & ijet,
359  int * constituent_indices, int & nconstituents) {
360  assert(cs.get() != 0);
361  assert(ijet > 0 && ijet <= int(jets.size()));
362 
363  vector<PseudoJet> constituents = cs->constituents(jets[ijet-1]);
364 
365  nconstituents = constituents.size();
366  for (int i = 0; i < nconstituents; i++) {
367  constituent_indices[i] = constituents[i].cluster_hist_index()+1;
368  }
369 }
370 
371 
378 //
379 // Corresponds to the following Fortran subroutine
380 // interface structure:
381 //
382 // FUNCTION FASTJETAREA(IJET)
383 // DOUBLE PRECISION FASTJETAREA
384 // INTEGER IJET
385 //
386 double fastjetarea_(const int & ijet) {
387  assert(ijet > 0 && ijet <= int(jets.size()));
388  const ClusterSequenceAreaBase * csab =
389  dynamic_cast<const ClusterSequenceAreaBase *>(cs.get());
390  if (csab != 0) {
391  // we have areas and can use csab to access all the area-related info
392  return csab->area(jets[ijet-1]);
393  } else {
394  return 0.;
395 // Error("No area information associated to this jet.");
396  }
397 }
398 
399 
402 //
403 // Corresponds to the following Fortran interface
404 //
405 // FUNCTION FASTJETDMERGE(N)
406 // DOUBLE PRECISION FASTJETDMERGE
407 // INTEGER N
408 //
409 double fastjetdmerge_(const int & n) {
410  assert(cs.get() != 0);
411  return cs->exclusive_dmerge(n);
412 }
413 
414 
419 //
420 // Corresponds to the following Fortran interface
421 //
422 // FUNCTION FASTJETDMERGEMAX(N)
423 // DOUBLE PRECISION FASTJETDMERGEMAX
424 // INTEGER N
425 //
426 double fastjetdmergemax_(const int & n) {
427  assert(cs.get() != 0);
428  return cs->exclusive_dmerge_max(n);
429 }
430 
431 
436 //
437 // Corresponds to the following Fortran interface
438 //
439 // SUBROUTINE FASTJETGLOBALRHOANDSIGMA(RAPMIN,RAPMAX,PHIMIN,PHIMAX,RHO,SIGMA,MEANAREA)
440 // DOUBLE PRECISION RAPMIN,RAPMAX,PHIMIN,PHIMAX
441 // DOUBLE PRECISION RHO,SIGMA,MEANAREA
442 //
443 void fastjetglobalrhoandsigma_(const double & rapmin, const double & rapmax,
444  const double & phimin, const double & phimax,
445  double & rho, double & sigma, double & meanarea) {
446  const ClusterSequenceAreaBase * csab =
447  dynamic_cast<const ClusterSequenceAreaBase *>(cs.get());
448  if (csab != 0) {
449  // we have areas and can use csab to access all the area-related info
450  Selector range = SelectorRapRange(rapmin,rapmax) && SelectorPhiRange(phimin,phimax);
451  bool use_area_4vector = false;
452  csab->get_median_rho_and_sigma(range,use_area_4vector,rho,sigma,meanarea);
453  } else {
454  Error("Clustering with area is necessary in order to be able to evaluate rho.");
455  }
456 }
457 
458 
459 
460 
461 }
int i
Definition: DBlmapReader.cc:9
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
auto_ptr< JetDefinition::Plugin > plugin
JetDefinition jet_def
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)
assert(m_qm.get())
double npart
Definition: HydjetWrapper.h:49
double fastjetdmergemax_(const int &n)
int ii
Definition: cuy.py:588
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
vector< PseudoJet > jets
int j
Definition: DBlmapReader.cc:9
double f[11][100]
Functor that operates on &lt;T&gt;
Definition: Selector.h:24
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 -&gt; f77jets[4*ijet+0..3]
double fastjetarea_(const int &ijet)
double ptmin
Definition: HydjetWrapper.h:90
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)