7 #include "tbb/task_arena.h"
23 #include "fastjet/ClusterSequence.hh"
26 using namespace fastjet;
28 template <
typename TILES>
33 antikt_radius_(conf.getParameter<double>(
"antikt_radius")),
34 minNumLayerCluster_(conf.getParameter<int>(
"minNumLayerCluster")),
35 eidInputName_(conf.getParameter<std::
string>(
"eid_input_name")),
36 eidOutputNameEnergy_(conf.getParameter<std::
string>(
"eid_output_name_energy")),
37 eidOutputNameId_(conf.getParameter<std::
string>(
"eid_output_name_id")),
38 eidMinClusterEnergy_(conf.getParameter<double>(
"eid_min_cluster_energy")),
39 eidNLayers_(conf.getParameter<int>(
"eid_n_layers")),
40 eidNClusters_(conf.getParameter<int>(
"eid_n_clusters")){};
42 template <
typename TILES>
44 std::vector<ticl::Trackster> &
result) {
47 <<
"Creating FastJet with " << fjInputs.size() <<
" LayerClusters in input";
49 fastjet::ClusterSequence sequence(fjInputs, JetDefinition(antikt_algorithm, antikt_radius_));
50 auto jets = fastjet::sorted_by_pt(sequence.inclusive_jets(0));
52 edm::LogVerbatim(
"PatternRecogntionbyFastJet") <<
"FastJet produced " <<
jets.size() <<
" jets/trackster";
55 auto trackster_idx = result.size();
56 result.resize(trackster_idx +
jets.size());
57 for (
const auto &pj :
jets) {
58 if (pj.constituents().size() >
static_cast<unsigned int>(minNumLayerCluster_)) {
59 for (
const auto &component : pj.constituents()) {
60 result[trackster_idx].vertices().push_back(component.user_index());
61 result[trackster_idx].vertex_multiplicity().push_back(1);
64 <<
"Jet has " << pj.constituents().size() <<
" components that are stored in trackster " << trackster_idx;
71 <<
"Jet with " << pj.constituents().size() <<
" constituents discarded since too small wrt "
72 << minNumLayerCluster_;
79 template <
typename TILES>
82 std::vector<Trackster> &
result,
83 std::unordered_map<
int, std::vector<int>> &seedToTracksterAssociation) {
90 rhtools_.setGeometry(geom);
97 auto lastLayerPerSide =
static_cast<unsigned int>(rhtools_.lastLayer(isHFnose)) - 1;
98 unsigned int maxLayer = 2 * lastLayerPerSide - 1;
99 std::vector<fastjet::PseudoJet> fjInputs;
101 for (
unsigned int currentLayer = 0; currentLayer <= maxLayer; ++currentLayer) {
102 if (currentLayer == lastLayerPerSide) {
103 buildJetAndTracksters(fjInputs, result);
105 const auto &tileOnLayer = input.
tiles[currentLayer];
106 for (
int ieta = 0; ieta <= nEtaBin; ++ieta) {
107 auto offset = ieta * nPhiBin;
111 for (
int iphi = 0; iphi <= nPhiBin; ++iphi) {
114 edm::LogVerbatim(
"PatternRecogntionbyFastJet") <<
"Entries in tileBin: " << tileOnLayer[
offset + iphi].size();
116 for (
auto clusterIdx : tileOnLayer[
offset + iphi]) {
118 if (input.
mask[clusterIdx] == 0.) {
120 edm::LogVerbatim(
"PatternRecogntionbyFastJet") <<
"Skipping masked layerIdx " << clusterIdx;
127 direction = direction.Unit();
128 direction *=
cl.energy();
129 auto fpj = fastjet::PseudoJet(direction.X(), direction.Y(), direction.Z(),
cl.energy());
130 fpj.set_user_index(clusterIdx);
131 fjInputs.push_back(fpj);
138 buildJetAndTracksters(fjInputs, result);
143 rhtools_.getPositionLayer(rhtools_.lastLayerEE(isHFnose), isHFnose).z());
148 for (
auto const &
t : result) {
149 edm::LogVerbatim(
"PatternRecogntionbyFastJet") <<
"Barycenter: " <<
t.barycenter();
150 edm::LogVerbatim(
"PatternRecogntionbyFastJet") <<
"LCs: " <<
t.vertices().size();
152 edm::LogVerbatim(
"PatternRecogntionbyFastJet") <<
"Regressed: " <<
t.regressed_energy();
157 template <
typename TILES>
159 const tensorflow::Session *eidSession,
160 std::vector<Trackster> &tracksters) {
184 std::vector<int> tracksterIndices;
185 for (
int i = 0; i < static_cast<int>(tracksters.size());
i++) {
189 float sumClusterEnergy = 0.;
190 for (
const unsigned int &vertex : tracksters[
i].
vertices()) {
191 sumClusterEnergy +=
static_cast<float>(layerClusters[vertex].energy());
193 if (sumClusterEnergy >= eidMinClusterEnergy_) {
195 tracksters[
i].setRegressedEnergy(0.
f);
196 tracksters[
i].zeroProbabilities();
197 tracksterIndices.push_back(
i);
204 int batchSize =
static_cast<int>(tracksterIndices.size());
205 if (batchSize == 0) {
210 tensorflow::TensorShape shape({batchSize, eidNLayers_, eidNClusters_, eidNFeatures_});
211 tensorflow::Tensor
input(tensorflow::DT_FLOAT, shape);
214 std::vector<tensorflow::Tensor> outputs;
215 std::vector<std::string> outputNames;
216 if (!eidOutputNameEnergy_.empty()) {
217 outputNames.push_back(eidOutputNameEnergy_);
219 if (!eidOutputNameId_.empty()) {
220 outputNames.push_back(eidOutputNameId_);
224 for (
int i = 0;
i < batchSize;
i++) {
225 const Trackster &trackster = tracksters[tracksterIndices[
i]];
230 std::vector<int> clusterIndices(trackster.vertices().size());
231 for (
int k = 0;
k < (int)trackster.vertices().size();
k++) {
232 clusterIndices[
k] =
k;
234 sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](
const int &
a,
const int &
b) {
235 return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(
b)].energy();
239 std::vector<int> seenClusters(eidNLayers_);
242 for (
const int &
k : clusterIndices) {
246 if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
248 float *
features = &input.tensor<float, 4>()(
i, j, seenClusters[j], 0);
251 *(features++) =
float(cluster.
energy() / float(trackster.vertex_multiplicity(
k)));
261 for (
int j = 0;
j < eidNLayers_;
j++) {
262 for (
int k = seenClusters[
j];
k < eidNClusters_;
k++) {
263 float *
features = &input.tensor<float, 4>()(
i,
j,
k, 0);
264 for (
int l = 0;
l < eidNFeatures_;
l++) {
272 tensorflow::run(const_cast<tensorflow::Session *>(eidSession), inputList, outputNames, &outputs);
275 if (!eidOutputNameEnergy_.empty()) {
277 float *
energy = outputs[0].flat<
float>().
data();
279 for (
const int &
i : tracksterIndices) {
280 tracksters[
i].setRegressedEnergy(*(energy++));
285 if (!eidOutputNameId_.empty()) {
287 int probsIdx = eidOutputNameEnergy_.empty() ? 0 : 1;
288 float *probs = outputs[probsIdx].flat<
float>().
data();
290 for (
const int &
i : tracksterIndices) {
291 tracksters[
i].setProbabilities(probs);
292 probs += tracksters[
i].id_probabilities().size();
297 template <
typename TILES>
299 iDesc.
add<
int>(
"algo_verbosity", 0);
300 iDesc.
add<
double>(
"antikt_radius", 0.09)->
setComment(
"Radius to be used while running the Anti-kt clustering");
301 iDesc.
add<
int>(
"minNumLayerCluster", 5)->setComment(
"Not Inclusive");
303 iDesc.
add<
std::string>(
"eid_output_name_energy",
"output/regressed_energy");
304 iDesc.
add<
std::string>(
"eid_output_name_id",
"output/id_probabilities");
305 iDesc.
add<
double>(
"eid_min_cluster_energy", 1.);
306 iDesc.
add<
int>(
"eid_n_layers", 50);
307 iDesc.
add<
int>(
"eid_n_clusters", 10);
Log< level::Info, true > LogVerbatim
void setComment(std::string const &value)
std::vector< NamedTensor > NamedTensorList
void energyRegressionAndID(const std::vector< reco::CaloCluster > &layerClusters, const tensorflow::Session *, std::vector< Trackster > &result)
double eta() const
pseudorapidity of cluster centroid
static std::string const input
bool getData(T &iHolder) const
void assignPCAtoTracksters(std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, const edm::ValueMap< std::pair< float, float >> &, double, bool energyWeight=true)
PatternRecognitionbyFastJet(const edm::ParameterSet &conf, edm::ConsumesCollector)
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
std::vector< float > features(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Abs< T >::type abs(const T &t)
double energy() const
cluster energy
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void makeTracksters(const typename PatternRecognitionAlgoBaseT< TILES >::Inputs &input, std::vector< Trackster > &result, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation) override
void buildJetAndTracksters(std::vector< fastjet::PseudoJet > &, std::vector< ticl::Trackster > &)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
char data[epos_bytes_allocation]
double phi() const
azimuthal angle of cluster centroid