Page MenuHomeHEPForge

No OneTemporary

Index: contrib/contribs/Centauro/trunk/Centauro.cc
===================================================================
--- contrib/contribs/Centauro/trunk/Centauro.cc (revision 1264)
+++ contrib/contribs/Centauro/trunk/Centauro.cc (revision 1265)
@@ -1,118 +1,134 @@
+// This code is part of Fastjet contrib
+
+// It is free software; you can redistribute it and/or modify it under
+// the terms of the GNU General Public License as published by the
+// Free Software Foundation; either version 2 of the License, or (at
+// your option) any later version.
+//
+// It is distributed in the hope that it will be useful, but WITHOUT
+// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
+// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
+// License for more details.
+//
+// You should have received a copy of the GNU General Public License
+// along with this code. If not, see <http://www.gnu.org/licenses/>.
+//----------------------------------------------------------------------x
+
#include "Centauro.hh"
#include "fastjet/NNH.hh"
// strings and streams
#include <sstream>
#include <limits>
FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
namespace contrib{
//----------------------------------------------------------------------
/// class that contains the algorithm parameters R, photon energy and photon longitudinal momentum.
class CentauroInfo {
public:
CentauroInfo(double Ri, double gammaEi, double gammaPzi)
{ R_ = Ri; gammaE_ = gammaEi; gammaPz_ = gammaPzi;}
double gammaPz() { return gammaPz_; }
double gammaE() { return gammaE_; }
double R() { return R_; }
private:
double R_, gammaE_, gammaPz_;
};
class CentauroBriefJet {
public: //For definitions see https://arxiv.org/abs/2006.10751
// n = (1,0,0,1)
// nbar = (1,0,0,-1)
// P = Q/2x(1,0,0,1) //proton
// q = Q(0,0,0,-1) //virtual photon
// etabar = -2Q/(nbar*q)*pT/(n*p) , so in the Breit frame: // etabar = +2*pT/(n*p) = 2*pT/(E-pz)
// The distance is: (for f=x)
// dij = [(etabar_i - etabar_j)^{2} + 2*etabar_i*etabar_j(1-cos(phi_i - phi_j))]/R^{2}
void init(const PseudoJet & jet, CentauroInfo * info) {
R = info->R();
gammaE = info->gammaE();
gammaPz = info->gammaPz();
// photon 4-momentum is q = (gammaE, 0 , 0 , gammaPz);
double norm = 1.0/sqrt(jet.modp2());
// pseudo-jet information needed to calculate distance
nx = jet.px() * norm;
ny = jet.py() * norm;
nz = jet.pz() * norm;
pT = jet.perp();
phi = jet.phi();
if(gammaE!=0 and gammaPz!=0){ //gammaE and gammaPz passed, so not running in Breit frame
Q = sqrt(-1.0*(gammaE*gammaE-gammaPz*gammaPz));
etabar = -2.0*(Q/(gammaE+gammaPz))*(pT/(jet.E()-jet.pz()));
}
else{ //gammaE and gammaPz not passed, so assume that it is running in the Breit frame
etabar = +2.0*pT/(jet.E()-jet.pz());
}
// beam distance
diB = 1.0;
}
double distance(const CentauroBriefJet * jet) const {
double dij = pow(etabar - jet->etabar, 2.0) + 2*etabar*jet->etabar*(1-cos(phi- jet->phi));
dij = dij/pow(R,2.0);
return dij;
}
double beam_distance() const {
return diB;
}
double pT, phi, nx, ny, nz;
double etabar;
double diB;
double R, gammaE, gammaPz, Q;
};
std::string CentauroPlugin::description () const {
std::ostringstream desc;
desc << "Centauro plugin with R = " << R();
if(gammaE()==0 and gammaPz()==0){
desc << " gamma E and gamma Pz parameters were not given --> assume you are giving particles momenta in Breit frame";
}
return desc.str();
}
void CentauroPlugin::run_clustering(fastjet::ClusterSequence & cs) const {
int njets = cs.jets().size();
CentauroInfo vinfo(R(), gammaE(), gammaPz());
NNH<CentauroBriefJet,CentauroInfo> nnh(cs.jets(),&vinfo);
while (njets > 0) {
int i, j, k;
double dij = nnh.dij_min(i, j);
if (j >= 0) {
cs.plugin_record_ij_recombination(i, j, dij, k);
nnh.merge_jets(i, j, cs.jets()[k], k);
} else {
cs.plugin_record_iB_recombination(i, dij);
nnh.remove_jet(i);
}
njets--;
}
}
} // namespace contrib
FASTJET_END_NAMESPACE
Index: contrib/contribs/Centauro/trunk/Centauro.hh
===================================================================
--- contrib/contribs/Centauro/trunk/Centauro.hh (revision 1264)
+++ contrib/contribs/Centauro/trunk/Centauro.hh (revision 1265)
@@ -1,55 +1,71 @@
+// This file is part of FastJet contrib.
+//
+// It is free software; you can redistribute it and/or modify it under
+// the terms of the GNU General Public License as published by the
+// Free Software Foundation; either version 2 of the License, or (at
+// your option) any later version.
+//
+// It is distributed in the hope that it will be useful, but WITHOUT
+// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
+// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
+// License for more details.
+//
+// You should have received a copy of the GNU General Public License
+// along with this code. If not, see <http://www.gnu.org/licenses/>.
+//----------------------------------------------------------------------
+
#ifndef __FASTJET_CONTRIB_CENTAUROJETALGORITHM_HH__
#define __FASTJET_CONTRIB_CENTAUROJETALGORITHM_HH__
#include <fastjet/internal/base.hh>
#include "fastjet/JetDefinition.hh"
#include "fastjet/ClusterSequence.hh"
FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
namespace contrib{
class CentauroPlugin : public JetDefinition::Plugin {
public:
/// Constructor for the Centauro Plugin class.
/// Three floating point arguments are specified to set the parameters
/// the radius parameter R, and gammaE and gammaPz are the energy and pz of virtual photon
CentauroPlugin (double R, double gammaE, double gammaPz) : _R(R), _gammaE(gammaE), _gammaPz(gammaPz){}
/// if only one argument is passed, assume that it runs in Breit frame so gammaE and gammaPz info not required so set
// so they are set to zero
CentauroPlugin (double R) : _R(R), _gammaE(0), _gammaPz(0){}
/// copy constructor
CentauroPlugin (const CentauroPlugin & plugin) {
*this = plugin;
}
// the things that are required by base class
virtual std::string description () const;
virtual void run_clustering(ClusterSequence &) const;
virtual double R() const {return _R;}
virtual double gammaE() const {return _gammaE;}
virtual double gammaPz() const {return _gammaPz;}
/// avoid the warning whenever the user requests "exclusive" jets
/// from the cluster sequence
virtual bool exclusive_sequence_meaningful() const {return true;}
private:
double _R;
double _gammaE;
double _gammaPz;
};
} // namespace contrib
FASTJET_END_NAMESPACE
#endif // __FASTJET_CONTRIB_CENTAUROJETALGORITHM_HH__
Index: contrib/contribs/Centauro/trunk/example.cc
===================================================================
--- contrib/contribs/Centauro/trunk/example.cc (revision 1264)
+++ contrib/contribs/Centauro/trunk/example.cc (revision 1265)
@@ -1,63 +1,83 @@
+// run it with
+// ./example < single-epDIS-event.dat
+//----------------------------------------------------------------------
+//
+//----------------------------------------------------------------------
+// This file is part of FastJet contrib.
+//
+// It is free software; you can redistribute it and/or modify it under
+// the terms of the GNU General Public License as published by the
+// Free Software Foundation; either version 2 of the License, or (at
+// your option) any later version.
+//
+// It is distributed in the hope that it will be useful, but WITHOUT
+// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
+// or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
+// License for more details.
+//
+// You should have received a copy of the GNU General Public License
+// along with this code. If not, see <http://www.gnu.org/licenses/>.
+
#include <iostream>
#include <iomanip>
#include "fastjet/JetDefinition.hh"
#include "fastjet/PseudoJet.hh"
#include "fastjet/Selector.hh"
#include "Centauro.hh"
#include "fastjet/EECambridgePlugin.hh"
using namespace std;
using namespace fastjet;
void read_event(vector<PseudoJet> &full_event);
//----------------------------------------------------------------------
int main(){
// example for the creation and execution of the fastjet plugin
// for the Centauro jet algorithm
//
// read in input particles
//----------------------------------------------------------
vector<PseudoJet> full_event;
read_event(full_event);
// create the jet definition using the plugin mechanism
//----------------------------------------------------------
//The Centauro jet algorithm is described in arXiv:2006.10751, "Asymmetric jet clustering in deep-inelastic scattering", Miguel Arratia, Yiannis Makris, Duff Neill, Felix Ringer, Nobuo Sato.
fastjet::contrib::CentauroPlugin * centauro_plugin = new fastjet::contrib::CentauroPlugin(1.0);
fastjet::JetDefinition jet_def(centauro_plugin);
ClusterSequence clust_seq(full_event, jet_def);
vector<PseudoJet> jets = clust_seq.inclusive_jets(0);
vector<PseudoJet> sortedJets = sorted_by_E(jets);
cout << "jets in inclusive clustering " << endl;
for (unsigned int i=0; i<sortedJets.size(); i++){
const PseudoJet &jet = sortedJets[i];
vector<fastjet::PseudoJet> constituents = jet.constituents();
cout << " rap = " << jet.rap() << " e " << jet.e() << " n constituents " << constituents.size() << endl;
}
cout << endl;
return 0;
}
//------------------------------------------------------------------------
// read the event
void read_event(vector<PseudoJet> &full_event){
string line;
std::cout << "Particles that are input" << std::endl;
while (getline(cin, line)) { istringstream linestream(line);
// take substrings to avoid problems when there are extra "pollution"
// characters (e.g. line-feed).
double px,py,pz,E;
linestream >> px >> py >> pz >> E; PseudoJet particle(px,py,pz,E);
std::cout << px << " " << py << " " << pz << " " << E <<std::endl;
// push event onto back of full_event vector
full_event.push_back(particle); }
}

File Metadata

Mime Type
text/x-diff
Expires
Tue, Nov 19, 7:08 PM (1 d, 12 h)
Storage Engine
blob
Storage Format
Raw Data
Storage Handle
3805767
Default Alt Text
(11 KB)

Event Timeline