00001 //STARTHEADER 00002 // $Id: fastjetfortran_madfks.cc,v 1.1 2013/02/08 20:15:52 spadhi Exp $ 00003 // 00004 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin P. Salam and Gregory Soyez 00005 // 00006 //---------------------------------------------------------------------- 00007 // This file is part of FastJet. 00008 // 00009 // FastJet is free software; you can redistribute it and/or modify 00010 // it under the terms of the GNU General Public License as published by 00011 // the Free Software Foundation; either version 2 of the License, or 00012 // (at your option) any later version. 00013 // 00014 // The algorithms that underlie FastJet have required considerable 00015 // development and are described in hep-ph/0512210. If you use 00016 // FastJet as part of work towards a scientific publication, please 00017 // include a citation to the FastJet paper. 00018 // 00019 // FastJet is distributed in the hope that it will be useful, 00020 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00021 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00022 // GNU General Public License for more details. 00023 // 00024 // You should have received a copy of the GNU General Public License 00025 // along with FastJet. If not, see <http://www.gnu.org/licenses/>. 00026 //---------------------------------------------------------------------- 00027 //ENDHEADER 00028 00029 #include <iostream> 00030 #include "fastjet/ClusterSequence.hh" 00031 #include "fastjet/ClusterSequenceArea.hh" 00032 #include "fastjet/Selector.hh" 00033 #include "fastjet/SISConePlugin.hh" 00034 00035 using namespace std; 00036 using namespace fastjet; 00037 00038 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 00039 00042 namespace fwrapper { 00043 vector<PseudoJet> input_particles, jets; 00044 auto_ptr<JetDefinition::Plugin> plugin; 00045 JetDefinition jet_def; 00046 auto_ptr<ClusterSequence> cs; 00047 00049 void transfer_input_particles(const double * p, const int & npart) { 00050 input_particles.resize(0); 00051 input_particles.reserve(npart); 00052 for (int i=0; i<npart; i++) { 00053 valarray<double> mom(4); // mom[0..3] 00054 for (int j=0;j<=3; j++) { 00055 mom[j] = *(p++); 00056 // RF-MZ: reorder the arguments because in madfks energy goes first; in fastjet last. 00057 // mom[(j+3) % 4] = *(p++); 00058 } 00059 PseudoJet psjet(mom); 00060 psjet.set_user_index(i); 00061 input_particles.push_back(psjet); 00062 } 00063 } 00064 00066 void transfer_jets(double * f77jets, int & njets) { 00067 njets = jets.size(); 00068 for (int i=0; i<njets; i++) { 00069 for (int j=0;j<=3; j++) { 00070 *f77jets = jets[i][j]; 00071 // RF-MZ: reorder the arguments because in madfks energy goes first; in fastjet last. 00072 //*f77jets = jets[i][(j+3) % 4]; 00073 f77jets++; 00074 } 00075 } 00076 } 00077 00080 void transfer_cluster_transfer(const double * p, const int & npart, 00081 const JetDefinition & jet_def, 00082 const double & ptmin, 00083 double * f77jets, int & njets, int * whichjet, 00084 const double & ghost_maxrap = 0.0, 00085 const int & nrepeat = 0, const double & ghost_area = 0.0) { 00086 00087 // transfer p[4*ipart+0..3] -> input_particles[i] 00088 transfer_input_particles(p, npart); 00089 00090 // perform the clustering 00091 if ( ghost_maxrap == 0.0 ) { 00092 // cluster without areas 00093 cs.reset(new ClusterSequence(input_particles,jet_def)); 00094 } else { 00095 // cluster with areas 00096 GhostedAreaSpec area_spec(ghost_maxrap,nrepeat,ghost_area); 00097 AreaDefinition area_def(active_area, area_spec); 00098 cs.reset(new ClusterSequenceArea(input_particles,jet_def,area_def)); 00099 } 00100 // extract jets (pt-ordered) 00101 jets = sorted_by_pt(cs->inclusive_jets(ptmin)); 00102 00103 // transfer jets -> f77jets[4*ijet+0..3] 00104 transfer_jets(f77jets, njets); 00105 00106 // Determine which parton/particle ended-up in which jet 00107 // set all jet entrie to zero first 00108 for(int ii=0; ii<npart; ++ii) whichjet[ii]=0; 00109 00110 // Loop over jets and find constituents 00111 for (int kk=0; kk<njets; ++kk) { 00112 vector<PseudoJet> constit = cs->constituents(jets[kk]); 00113 for(unsigned int ll=0; ll<constit.size(); ++ll) 00114 whichjet[constit[ll].user_index()]=kk+1; 00115 } 00116 00117 } 00118 00119 } 00120 FASTJET_END_NAMESPACE 00121 00122 using namespace fastjet::fwrapper; 00123 00124 00125 extern "C" { 00126 00129 // 00130 // Corresponds to the following Fortran subroutine 00131 // interface structure: 00132 // 00133 // SUBROUTINE FASTJETSISCONE(P,NPART,R,F,F77JETS,NJETS,WHICHJET) 00134 // DOUBLE PRECISION P(4,*), R, F, F77JETS(4,*) 00135 // INTEGER NPART, NJETS, WHICHJET(*) 00136 // 00137 // where on input 00138 // 00139 // P the input particle 4-momenta 00140 // NPART the number of input momenta 00141 // R the radius parameter 00142 // F the overlap threshold 00143 // 00144 // and on output 00145 // 00146 // F77JETS the output jet momenta (whose second dim should be >= NPART) 00147 // sorted in order of decreasing p_t. 00148 // NJETS the number of output jets 00149 // WHICHJET(i) the jet of parton/particle 'i' 00150 // 00151 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its 00152 // momenta as a matrix of the form P(4000,5), whereas this fortran 00153 // interface to fastjet expects them as P(4,NPART), i.e. you must take 00154 // the transpose of the Pythia array and drop the fifth component 00155 // (particle mass). 00156 // 00157 void fastjetsiscone_(const double * p, const int & npart, 00158 const double & R, const double & f, 00159 double * f77jets, int & njets, int * whichjet) { 00160 00161 // prepare jet def 00162 plugin.reset(new SISConePlugin(R,f)); 00163 jet_def = plugin.get(); 00164 00165 // do everything 00166 // transfer_cluster_transfer(p,npart,jet_def,f77jets,njets,whichjet); 00167 } 00168 00169 00170 00175 // 00176 // Corresponds to the following Fortran subroutine 00177 // interface structure: 00178 // 00179 // SUBROUTINE FASTJETSISCONEWITHAREA(P,NPART,R,F,GHMAXRAP,NREP,GHAREA,F77JETS,NJETS,WHICHJET) 00180 // DOUBLE PRECISION P(4,*), R, F, F77JETS(4,*), GHMAXRAP, GHAREA 00181 // INTEGER NPART, NJETS, NREP, WHICHJET(*) 00182 // 00183 // where on input 00184 // 00185 // P the input particle 4-momenta 00186 // NPART the number of input momenta 00187 // R the radius parameter 00188 // F the overlap threshold 00189 // GHMAXRAP the maximum (abs) rapidity covered by ghosts (FastJet default 6.0) 00190 // NREP the number of repetitions used to evaluate the area (FastJet default 1) 00191 // GHAREA the area of a single ghost (FastJet default 0.01) 00192 // 00193 // and on output 00194 // 00195 // F77JETS the output jet momenta (whose second dim should be >= NPART) 00196 // sorted in order of decreasing p_t. 00197 // NJETS the number of output jets 00198 // WHICHJET(i) the jet of parton/particle 'i' 00199 // 00200 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its 00201 // momenta as a matrix of the form P(4000,5), whereas this fortran 00202 // interface to fastjet expects them as P(4,NPART), i.e. you must take 00203 // the transpose of the Pythia array and drop the fifth component 00204 // (particle mass). 00205 // 00206 void fastjetsisconewitharea_(const double * p, const int & npart, 00207 const double & R, const double & f, 00208 const double & ghost_rapmax, const int & nrepeat, const double & ghost_area, 00209 double * f77jets, int & njets, int * whichjet) { 00210 00211 // prepare jet def 00212 plugin.reset(new SISConePlugin(R,f)); 00213 jet_def = plugin.get(); 00214 00215 // do everything 00216 // transfer_cluster_transfer(p,npart,jet_def,f77jets,njets,whichjet,ghost_rapmax,nrepeat,ghost_area); 00217 } 00218 00219 00220 00224 // 00225 // Corresponds to the following Fortran subroutine 00226 // interface structure: 00227 // 00228 // SUBROUTINE FASTJETPPGENKT(P,NPART,R,PALG,F77JETS,NJETS,WHICHJET) 00229 // DOUBLE PRECISION P(4,*), R, PALG, F, F77JETS(4,*) 00230 // INTEGER NPART, NJETS, WHICHJET(*) 00231 // 00232 // where on input 00233 // 00234 // P the input particle 4-momenta 00235 // NPART the number of input momenta 00236 // R the radius parameter 00237 // PALG the power for the generalised kt alg 00238 // (1.0=kt, 0.0=C/A, -1.0 = anti-kt) 00239 // 00240 // and on output 00241 // 00242 // F77JETS the output jet momenta (whose second dim should be >= NPART) 00243 // sorted in order of decreasing p_t. 00244 // NJETS the number of output jets 00245 // WHICHJET(i) the jet of parton/particle 'i' 00246 // 00247 // For the values of PALG that correspond to "standard" cases (1.0=kt, 00248 // 0.0=C/A, -1.0 = anti-kt) this routine actually calls the direct 00249 // implementation of those algorithms, whereas for other values of 00250 // PALG it calls the generalised kt implementation. 00251 // 00252 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its 00253 // momenta as a matrix of the form P(4000,5), whereas this fortran 00254 // interface to fastjet expects them as P(4,NPART), i.e. you must take 00255 // the transpose of the Pythia array and drop the fifth component 00256 // (particle mass). 00257 // 00258 void fastjetppgenkt_(const double * p, const int & npart, 00259 const double & R, const double & ptjetmin, 00260 const double & palg, 00261 double * f77jets, int & njets, int * whichjet) { 00262 00263 // prepare jet def 00264 if (palg == 1.0) { 00265 jet_def = JetDefinition(kt_algorithm, R); 00266 } else if (palg == 0.0) { 00267 jet_def = JetDefinition(cambridge_algorithm, R); 00268 } else if (palg == -1.0) { 00269 jet_def = JetDefinition(antikt_algorithm, R); 00270 } else { 00271 jet_def = JetDefinition(genkt_algorithm, R, palg); 00272 } 00273 00274 // do everything 00275 transfer_cluster_transfer(p,npart,jet_def,ptjetmin,f77jets,njets,whichjet); 00276 } 00277 00278 00284 // 00285 // Corresponds to the following Fortran subroutine 00286 // interface structure: 00287 // 00288 // SUBROUTINE FASTJETPPGENKTWITHAREA(P,NPART,R,PALG,GHMAXRAP,NREP,GHAREA,F77JETS,NJETS,WHICHJET) 00289 // DOUBLE PRECISION P(4,*), R, PALG, GHMAXRAP, GHAREA, F77JETS(4,*) 00290 // INTEGER NPART, NREP, NJETS, WHICHJET(*) 00291 // 00292 // where on input 00293 // 00294 // P the input particle 4-momenta 00295 // NPART the number of input momenta 00296 // R the radius parameter 00297 // PALG the power for the generalised kt alg 00298 // (1.0=kt, 0.0=C/A, -1.0 = anti-kt) 00299 // GHMAXRAP the maximum (abs) rapidity covered by ghosts (FastJet default 6.0) 00300 // NREP the number of repetitions used to evaluate the area (FastJet default 1) 00301 // GHAREA the area of a single ghost (FastJet default 0.01) 00302 // 00303 // and on output 00304 // 00305 // F77JETS the output jet momenta (whose second dim should be >= NPART) 00306 // sorted in order of decreasing p_t. 00307 // NJETS the number of output jets 00308 // WHICHJET(i) the jet of parton/particle 'i' 00309 // 00310 // For the values of PALG that correspond to "standard" cases (1.0=kt, 00311 // 0.0=C/A, -1.0 = anti-kt) this routine actually calls the direct 00312 // implementation of those algorithms, whereas for other values of 00313 // PALG it calls the generalised kt implementation. 00314 // 00315 // NOTE: if you are interfacing fastjet to Pythia 6, Pythia stores its 00316 // momenta as a matrix of the form P(4000,5), whereas this fortran 00317 // interface to fastjet expects them as P(4,NPART), i.e. you must take 00318 // the transpose of the Pythia array and drop the fifth component 00319 // (particle mass). 00320 // 00321 void fastjetppgenktwitharea_(const double * p, const int & npart, 00322 const double & R, const double & palg, 00323 const double & ghost_rapmax, const int & nrepeat, const double & ghost_area, 00324 double * f77jets, int & njets, int * whichjet) { 00325 00326 // prepare jet def 00327 if (palg == 1.0) { 00328 jet_def = JetDefinition(kt_algorithm, R); 00329 } else if (palg == 0.0) { 00330 jet_def = JetDefinition(cambridge_algorithm, R); 00331 } else if (palg == -1.0) { 00332 jet_def = JetDefinition(antikt_algorithm, R); 00333 } else { 00334 jet_def = JetDefinition(genkt_algorithm, R, palg); 00335 } 00336 00337 // do everything 00338 // transfer_cluster_transfer(p,npart,jet_def,f77jets,njets,whichjet,ghost_rapmax,nrepeat,ghost_area); 00339 } 00340 00341 00350 // 00351 // Corresponds to the following Fortran subroutine 00352 // interface structure: 00353 // 00354 // SUBROUTINE FASTJETCONSTITUENTS(IJET,CONSTITUENT_INDICES,NCONSTITUENTS) 00355 // INTEGER IJET 00356 // INTEGER CONSTITUENT_INDICES(*) 00357 // INTEGER nconstituents 00358 // 00359 void fastjetconstituents_(const int & ijet, 00360 int * constituent_indices, int & nconstituents) { 00361 assert(cs.get() != 0); 00362 assert(ijet > 0 && ijet <= int(jets.size())); 00363 00364 vector<PseudoJet> constituents = cs->constituents(jets[ijet-1]); 00365 00366 nconstituents = constituents.size(); 00367 for (int i = 0; i < nconstituents; i++) { 00368 constituent_indices[i] = constituents[i].cluster_hist_index()+1; 00369 } 00370 } 00371 00372 00379 // 00380 // Corresponds to the following Fortran subroutine 00381 // interface structure: 00382 // 00383 // FUNCTION FASTJETAREA(IJET) 00384 // DOUBLE PRECISION FASTJETAREA 00385 // INTEGER IJET 00386 // 00387 double fastjetarea_(const int & ijet) { 00388 assert(ijet > 0 && ijet <= int(jets.size())); 00389 const ClusterSequenceAreaBase * csab = 00390 dynamic_cast<const ClusterSequenceAreaBase *>(cs.get()); 00391 if (csab != 0) { 00392 // we have areas and can use csab to access all the area-related info 00393 return csab->area(jets[ijet-1]); 00394 } else { 00395 return 0.; 00396 // Error("No area information associated to this jet."); 00397 } 00398 } 00399 00400 00403 // 00404 // Corresponds to the following Fortran interface 00405 // 00406 // FUNCTION FASTJETDMERGE(N) 00407 // DOUBLE PRECISION FASTJETDMERGE 00408 // INTEGER N 00409 // 00410 double fastjetdmerge_(const int & n) { 00411 assert(cs.get() != 0); 00412 return cs->exclusive_dmerge(n); 00413 } 00414 00415 00420 // 00421 // Corresponds to the following Fortran interface 00422 // 00423 // FUNCTION FASTJETDMERGEMAX(N) 00424 // DOUBLE PRECISION FASTJETDMERGEMAX 00425 // INTEGER N 00426 // 00427 double fastjetdmergemax_(const int & n) { 00428 assert(cs.get() != 0); 00429 return cs->exclusive_dmerge_max(n); 00430 } 00431 00432 00437 // 00438 // Corresponds to the following Fortran interface 00439 // 00440 // SUBROUTINE FASTJETGLOBALRHOANDSIGMA(RAPMIN,RAPMAX,PHIMIN,PHIMAX,RHO,SIGMA,MEANAREA) 00441 // DOUBLE PRECISION RAPMIN,RAPMAX,PHIMIN,PHIMAX 00442 // DOUBLE PRECISION RHO,SIGMA,MEANAREA 00443 // 00444 void fastjetglobalrhoandsigma_(const double & rapmin, const double & rapmax, 00445 const double & phimin, const double & phimax, 00446 double & rho, double & sigma, double & meanarea) { 00447 const ClusterSequenceAreaBase * csab = 00448 dynamic_cast<const ClusterSequenceAreaBase *>(cs.get()); 00449 if (csab != 0) { 00450 // we have areas and can use csab to access all the area-related info 00451 Selector range = SelectorRapRange(rapmin,rapmax) && SelectorPhiRange(phimin,phimax); 00452 bool use_area_4vector = false; 00453 csab->get_median_rho_and_sigma(range,use_area_4vector,rho,sigma,meanarea); 00454 } else { 00455 Error("Clustering with area is necessary in order to be able to evaluate rho."); 00456 } 00457 } 00458 00459 00460 00461 00462 }