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