9 #include "Math/ProbFuncMathCore.h" 11 #include "fastjet/contrib/ConstituentSubtractor.hh" 22 useModulatedRho_(conf.getParameter<
bool>(
"useModulatedRho")),
23 minFlowChi2Prob_(conf.getParameter<double>(
"minFlowChi2Prob")),
24 maxFlowChi2Prob_(conf.getParameter<double>(
"maxFlowChi2Prob")) {
58 std::vector<fastjet::PseudoJet> tempJets = fastjet::sorted_by_pt(
fjClusterSeq_->inclusive_jets(
jetPtMin_));
72 unsigned int bkgVecSize =
etaRanges->size();
74 throw cms::Exception(
"WrongBkgInput") <<
"Producer needs vector with background estimates\n";
76 if (bkgVecSize != (rhoRanges->size() + 1) || bkgVecSize != (rhomRanges->size() + 1)) {
78 <<
"Size of etaRanges (" << bkgVecSize <<
") and rhoRanges (" << rhoRanges->size() <<
") and/or rhomRanges (" 79 << rhomRanges->size() <<
") vectors inconsistent\n";
93 for (fastjet::PseudoJet& ijet : tempJets) {
96 std::vector<fastjet::PseudoJet>
particles, ghosts;
97 fastjet::SelectorIsPureGhost().sift(ijet.constituents(), ghosts,
particles);
98 unsigned long nParticles =
particles.size();
103 for (fastjet::PseudoJet& ighost : ghosts) {
104 double rhoModulationFactor = 1.;
105 double ghostPhi = ighost.phi_std();
119 int32_t ghostPos = -1;
120 auto const ighostEta = ighost.eta();
121 if (ighostEta <= etaRanges->at(0) || bkgVecSize == 1) {
123 }
else if (ighostEta >=
etaRanges->at(bkgVecSize - 1)) {
124 ghostPos = rhoRanges->size() - 1;
126 for (
unsigned int ie = 0; ie < (bkgVecSize - 1); ie++) {
133 double pt = rhoRanges->at(ghostPos) * ighost.area() * rhoModulationFactor;
134 double mass_squared =
135 std::pow(rhoModulationFactor * rhomRanges->at(ghostPos) * ighost.area() +
pt, 2) -
std::pow(
pt, 2);
136 double mass = (mass_squared > 0) ?
sqrt(mass_squared) : 0;
137 ighost.reset_momentum_PtYPhiM(
pt, ighost.rap(), ighost.phi(),
mass);
143 fastjet::contrib::ConstituentSubtractor subtractor;
145 subtractor.set_max_distance(
147 subtractor.set_alpha(
149 subtractor.set_do_mass_subtraction();
150 subtractor.set_remove_all_zero_pt_particles(
true);
152 std::vector<fastjet::PseudoJet> subtracted_particles = subtractor.do_subtraction(
particles, ghosts);
155 fastjet::PseudoJet subtracted_jet =
join(subtracted_particles);
156 if (subtracted_jet.perp() > 0.)
157 fjJets_.push_back(subtracted_jet);
167 descCSJetProducer.
add<
string>(
"jetCollInstanceName",
"");
169 descCSJetProducer.
add<
bool>(
"sumRecHits",
false);
172 descriptions.
add(
"CSJetProducer", descCSJetProducer);
176 desc.add<
double>(
"csRParam", -1.);
177 desc.add<
double>(
"csAlpha", 2.);
178 desc.add<
bool>(
"useModulatedRho",
false);
179 desc.add<
double>(
"minFlowChi2Prob", 0.05);
180 desc.add<
double>(
"maxFlowChi2Prob", 0.95);
184 desc.add<
edm::InputTag>(
"rhoFlowFitParams", {
"hiFJRhoFlowModulationProducer",
"rhoFlowFitParams"});
188 const double phi,
const double eventPlane2,
const double eventPlane3,
const double par1,
const double par2) {
191 return 1. + 2. *
par1 *
cos(2. * (phi - eventPlane2)) +
par2 *
cos(3. * (phi - eventPlane3));
T getParameter(std::string const &) const
std::vector< fastjet::PseudoJet > fjJets_
double csAlpha_
for constituent subtraction : R parameter
bool useModulatedRho_
for HI constituent subtraction : alpha (power of pt in metric)
#define DEFINE_FWK_MODULE(type)
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< std::vector< double > > rhoToken_
edm::EDGetTokenT< std::vector< double > > rhoFlowFitParamsToken_
double maxFlowChi2Prob_
flowFit chi2/ndof minimum compatability requirement
std::vector< fastjet::PseudoJet > fjInputs_
ClusterSequencePtr fjClusterSeq_
Cos< T >::type cos(const T &t)
double getModulatedRhoFactor(double phi, double eventPlane2, double eventPlane3, double par1, double par2)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Namespace of DDCMS conversion namespace.
edm::EDGetTokenT< std::vector< double > > etaToken_
flowFit chi2/ndof minimum compatability requirement
static std::string join(char **cmd)
void runAlgorithm(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< std::vector< double > > rhomToken_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
AreaDefinitionPtr fjAreaDefinition_
static void fillDescriptionsFromVirtualJetProducer(edm::ParameterSetDescription &desc)
static void fillDescriptionsFromCSJetProducer(edm::ParameterSetDescription &desc)
Power< A, B >::type pow(const A &a, const B &b)
double minFlowChi2Prob_
flag to turn on/off flow-modulated rho and rhom
JetDefPtr fjJetDefinition_