ATLAS Offline Software
Loading...
Searching...
No Matches
TFCSPredictExtrapWeights.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
11
12#include "TFile.h"
13#include "TClass.h"
14
15#include "HepPDT/ParticleData.hh"
16#include "HepPDT/ParticleDataTable.hh"
17
18#include "CLHEP/Random/RandGauss.h"
19#include "CLHEP/Random/RandFlat.h"
20
21#if defined(__FastCaloSimStandAlone__)
22#include "CLHEP/Random/TRandomEngine.h"
23#else
24#include <CLHEP/Random/RanluxEngine.h>
25#endif
26
27#include <iostream>
28#include <fstream>
29
30// LWTNN
31#include "lwtnn/LightweightGraph.hh"
32#include "lwtnn/LightweightNeuralNetwork.hh"
33#include "lwtnn/parse_json.hh"
34
35//=============================================
36//======= TFCSPredictExtrapWeights =========
37//=============================================
38
45
46// Destructor
48 if (m_input != nullptr) {
49 delete m_input;
50 }
51 if (m_relevantLayers != nullptr) {
52 delete m_relevantLayers;
53 }
54 if (m_normLayers != nullptr) {
55 delete m_normLayers;
56 }
57 if (m_normMeans != nullptr) {
58 delete m_normMeans;
59 }
60 if (m_normStdDevs != nullptr) {
61 delete m_normStdDevs;
62 }
63 if (m_nn != nullptr) {
64 delete m_nn;
65 }
66}
67
69 const TFCSParametrizationBase &ref) const {
70 if (IsA() != ref.IsA()) {
71 ATH_MSG_DEBUG("operator==: different class types "
72 << IsA()->GetName() << " != " << ref.IsA()->GetName());
73 return false;
74 }
75 const TFCSPredictExtrapWeights &ref_typed =
76 static_cast<const TFCSPredictExtrapWeights &>(ref);
77
79 return true;
81 return false;
83 return false;
84
85 return (m_input->compare(*ref_typed.m_input) == 0);
86}
87
88// getNormInputs()
89// Get values needed to normalize inputs
91 const std::string &etaBin, const std::string &FastCaloTXTInputFolderName) {
92 ATH_MSG_DEBUG(" Getting normalization inputs... ");
93
94 // Open corresponding TXT file and extract mean/std dev for each variable
95 if (m_normLayers != nullptr) {
96 m_normLayers->clear();
97 } else {
98 m_normLayers = new std::vector<int>();
99 }
100 if (m_normMeans != nullptr) {
101 m_normMeans->clear();
102 } else {
103 m_normMeans = new std::vector<float>();
104 }
105 if (m_normStdDevs != nullptr) {
106 m_normStdDevs->clear();
107 } else {
108 m_normStdDevs = new std::vector<float>();
109 }
110 std::string inputFileName = FastCaloTXTInputFolderName +
111 "MeanStdDevEnergyFractions_eta_" + etaBin +
112 ".txt";
113 ATH_MSG_DEBUG(" Opening " << inputFileName);
114 std::ifstream inputTXT(inputFileName);
115 if (inputTXT.is_open()) {
116 std::string line;
117 while (getline(inputTXT, line)) {
118 std::stringstream ss(line);
119 unsigned int counter = 0;
120 while (ss.good()) {
121 std::string substr;
122 getline(ss, substr, ' ');
123 if (counter == 0) { // Get index (#layer or -1 if var == etrue)
124 if (substr != "etrue") {
125 int index = std::stoi(substr.substr(substr.find('_') + 1));
126 m_normLayers->push_back(index);
127 } else { // etrue
128 m_normLayers->push_back(-1);
129 }
130 } else if (counter == 1) {
131 m_normMeans->push_back(std::stof(substr));
132 } else if (counter == 2) {
133 m_normStdDevs->push_back(std::stof(substr));
134 }
135 counter++;
136 }
137 }
138 inputTXT.close();
139 } else {
140 ATH_MSG_ERROR(" Unable to open file " << inputFileName);
141 return false;
142 }
143
144 return true;
145}
146
147// prepareInputs()
148// Prepare input variables to the Neural Network
149std::map<std::string, double>
151 const float truthE) const {
152 std::map<std::string, double> inputVariables;
153 for (int ilayer = 0; ilayer < CaloCell_ID_FCS::MaxSample; ++ilayer) {
154 if (std::find(m_relevantLayers->cbegin(), m_relevantLayers->cend(),
155 ilayer) != m_relevantLayers->cend()) {
156 const std::string layer = std::to_string(ilayer);
157 // Find index
158 auto itr =
159 std::find(m_normLayers->cbegin(), m_normLayers->cend(), ilayer);
160 if (itr != m_normLayers->cend()) {
161 const int index = std::distance(m_normLayers->cbegin(), itr);
162 inputVariables["ef_" + layer] =
163 (simulstate.Efrac(ilayer) - std::as_const(m_normMeans)->at(index)) /
164 std::as_const(m_normStdDevs)->at(index);
165 } else {
166 ATH_MSG_ERROR("Normalization information not found for layer "
167 << ilayer);
168 }
169 }
170 }
171 // Find index for truth energy
172 auto itr = std::find(m_normLayers->cbegin(), m_normLayers->cend(), -1);
173 int index = std::distance(m_normLayers->cbegin(), itr);
174 inputVariables["etrue"] = (truthE - std::as_const(m_normMeans)->at(index)) /
175 std::as_const(m_normStdDevs)->at(index);
176 if (is_match_pdgid(22)) {
177 inputVariables["pdgId"] = 1; // one hot enconding
178 } else if (is_match_pdgid(11) || is_match_pdgid(-11)) {
179 inputVariables["pdgId"] = 0; // one hot enconding
180 }
181 return inputVariables;
182}
183
184// simulate()
185// get predicted extrapolation weights and save them as AuxInfo in simulstate
187 TFCSSimulationState &simulstate, const TFCSTruthState *truth,
188 const TFCSExtrapolationState * /*extrapol*/) const {
189
190 // Get inputs to Neural Network
191 std::map<std::string, double> inputVariables =
192 prepareInputs(simulstate, truth->E() * 0.001);
193
194 // Get predicted extrapolation weights
195 auto outputs = m_nn->compute(inputVariables);
196 for (int ilayer = 0; ilayer < CaloCell_ID_FCS::MaxSample; ++ilayer) {
197 if (std::find(m_relevantLayers->cbegin(), m_relevantLayers->cend(),
198 ilayer) != m_relevantLayers->cend()) {
199 ATH_MSG_DEBUG("TFCSPredictExtrapWeights::simulate: layer: "
200 << ilayer << " weight: "
201 << outputs["extrapWeight_" + std::to_string(ilayer)]);
202 float weight = outputs["extrapWeight_" + std::to_string(ilayer)];
203 // Protections
204 if (weight < 0) {
205 weight = 0;
206 } else if (weight > 1) {
207 weight = 1;
208 }
209 simulstate.setAuxInfo<float>(ilayer, weight);
210 } else { // use weight=0.5 for non-relevant layers
212 "Setting weight=0.5 for layer = " << std::to_string(ilayer));
213 simulstate.setAuxInfo<float>(ilayer, float(0.5));
214 }
215 }
216 return FCSSuccess;
217}
218
219// simulate_hit()
221 Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState * /*truth*/,
222 const TFCSExtrapolationState *extrapol) {
223
224 const int cs = calosample();
225
226 // Get corresponding predicted extrapolation weight from simulstate
227 float extrapWeight;
228 if (simulstate.hasAuxInfo(cs)) {
229 extrapWeight = simulstate.getAuxInfo<float>(cs);
230 } else { // missing AuxInfo
232 "Simulstate is not decorated with extrapolation weights for cs = "
233 << std::to_string(cs));
234 return FCSFatal;
235 }
236
237 double eta = (1. - extrapWeight) * extrapol->eta(cs, SUBPOS_ENT) +
238 extrapWeight * extrapol->eta(cs, SUBPOS_EXT);
239 double phi = (1. - extrapWeight) * extrapol->phi(cs, SUBPOS_ENT) +
240 extrapWeight * extrapol->phi(cs, SUBPOS_EXT);
241 float extrapWeight_for_r_z = extrapWeight;
242 if (UseHardcodedWeight()) {
243 extrapWeight_for_r_z = 0.5;
245 "Will use extrapWeight=0.5 for r and z when constructing a hit");
246 } else {
247 ATH_MSG_DEBUG("Will use predicted extrapWeight also for r and z when "
248 "constructing a hit");
249 }
250 double r = (1. - extrapWeight_for_r_z) * extrapol->r(cs, SUBPOS_ENT) +
251 extrapWeight_for_r_z * extrapol->r(cs, SUBPOS_EXT);
252 double z = (1. - extrapWeight_for_r_z) * extrapol->z(cs, SUBPOS_ENT) +
253 extrapWeight_for_r_z * extrapol->z(cs, SUBPOS_EXT);
254
255 if (!std::isfinite(r) || !std::isfinite(z) || !std::isfinite(eta) ||
256 !std::isfinite(phi)) {
257 ATH_MSG_WARNING("Extrapolator contains NaN or infinite number.\nSetting "
258 "center position to calo boundary.");
259 ATH_MSG_WARNING("Before fix: center_r: "
260 << r << " center_z: " << z << " center_phi: " << phi
261 << " center_eta: " << eta << " weight: " << extrapWeight
262 << " cs: " << cs);
263 // If extrapolator fails we can set position to calo boundary
264 r = extrapol->IDCaloBoundary_r();
265 z = extrapol->IDCaloBoundary_z();
266 eta = extrapol->IDCaloBoundary_eta();
267 phi = extrapol->IDCaloBoundary_phi();
268 ATH_MSG_WARNING("After fix: center_r: "
269 << r << " center_z: " << z << " center_phi: " << phi
270 << " center_eta: " << eta << " weight: " << extrapWeight
271 << " cs: " << cs);
272 }
273
274 hit.setCenter_r(r);
275 hit.setCenter_z(z);
276 hit.setCenter_eta(eta);
277 hit.setCenter_phi(phi);
278
279 ATH_MSG_DEBUG("TFCSPredictExtrapWeights: center_r: "
280 << hit.center_r() << " center_z: " << hit.center_z()
281 << " center_phi: " << hit.center_phi()
282 << " center_eta: " << hit.center_eta()
283 << " weight: " << extrapWeight << " cs: " << cs);
284
285 return FCSSuccess;
286}
287
288// initializeNetwork()
289// Initialize lwtnn network
291 int pid, const std::string &etaBin,
292 const std::string &FastCaloNNInputFolderName) {
293
295 "Using FastCaloNNInputFolderName: " << FastCaloNNInputFolderName);
296 set_pdgid(pid);
297
298 std::string inputFileName =
299 FastCaloNNInputFolderName + "NN_" + etaBin + ".json";
300 ATH_MSG_DEBUG("Will read JSON file: " << inputFileName);
301 if (inputFileName.empty()) {
302 ATH_MSG_ERROR("Could not find json file " << inputFileName);
303 return false;
304 } else {
305 ATH_MSG_INFO("For pid: " << pid << " and etaBin" << etaBin
306 << ", loading json file " << inputFileName);
307 std::ifstream input(inputFileName);
308 std::stringstream sin;
309 sin << input.rdbuf();
310 input.close();
311 auto config = lwt::parse_json(sin);
312 m_nn = new lwt::LightweightNeuralNetwork(config.inputs, config.layers,
313 config.outputs);
314 if (m_nn == nullptr) {
315 ATH_MSG_ERROR("Could not create LightWeightNeuralNetwork from "
316 << inputFileName);
317 return false;
318 }
319 if (m_input != nullptr) {
320 delete m_input;
321 }
322 m_input = new std::string(sin.str());
323 // Extract relevant layers from the outputs
324 m_relevantLayers = new std::vector<int>();
325 for (auto name : config.outputs) {
326 int layer = std::stoi(
327 name.erase(0, 13)); // remove "extrapWeight_" and convert to int
328 m_relevantLayers->push_back(layer);
329 }
330 }
331 return true;
332}
333
334// Streamer()
335void TFCSPredictExtrapWeights::Streamer(TBuffer &R__b) {
336 // Stream an object of class TFCSPredictExtrapWeights
337
338 if (R__b.IsReading()) {
339 R__b.ReadClassBuffer(TFCSPredictExtrapWeights::Class(), this);
340 if (m_nn != nullptr) {
341 delete m_nn;
342 m_nn = nullptr;
343 }
344 if (m_input && !m_input->empty()) {
345 std::stringstream sin;
346 sin.str(*m_input);
347 auto config = lwt::parse_json(sin);
348 m_nn = new lwt::LightweightNeuralNetwork(config.inputs, config.layers,
349 config.outputs);
350 }
351#ifndef __FastCaloSimStandAlone__
352 // When running inside Athena, delete input/config/normInputs to free the
353 // memory
354 if (freemem()) {
355 delete m_input;
356 m_input = nullptr;
357 }
358#endif
359 } else {
360 R__b.WriteClassBuffer(TFCSPredictExtrapWeights::Class(), this);
361 }
362}
363
364// unit_test()
365// Function for testing
367 TFCSSimulationState *simulstate, const TFCSTruthState *truth,
368 const TFCSExtrapolationState *extrapol) {
369 const std::string this_file = __FILE__;
370 const std::string parent_dir = this_file.substr(0, this_file.find("/src/"));
371 const std::string norm_path = parent_dir + "/share/NormPredExtrapSample/";
372 std::string net_path = "/cvmfs/atlas-nightlies.cern.ch/repo/data/data-art/"
373 "FastCaloSim/LWTNNPredExtrapSample/";
374 test_path(net_path, norm_path, simulstate, truth, extrapol);
375 //net_path = "/cvmfs/atlas-nightlies.cern.ch/repo/data/data-art/FastCaloSim/"
376 // "ONNXPredExtrapSample/";
377 //test_path(net_path, norm_path, simulstate, truth, extrapol);
378}
379
380// test_path()
381// Function for testing
383 std::string &net_path, std::string const &norm_path,
384 TFCSSimulationState *simulstate, const TFCSTruthState *truth,
385 const TFCSExtrapolationState *extrapol) {
387 ATH_MSG_NOCLASS(logger, "Testing net path ..."
388 << net_path.substr(net_path.length() - 20)
389 << " and norm path ..."
390 << norm_path.substr(norm_path.length() - 20));
391 if (!simulstate) {
392 simulstate = new TFCSSimulationState();
393#if defined(__FastCaloSimStandAlone__)
394 simulstate->setRandomEngine(new CLHEP::TRandomEngine());
395#else
396 simulstate->setRandomEngine(new CLHEP::RanluxEngine());
397#endif
398 }
399 if (!truth) {
401 t->SetPtEtaPhiM(524288000, 0, 0, 130); // 524288 GeV
402 t->set_pdgid(22); // photon
403 truth = t;
404 }
405 if (!extrapol) {
407 e->set_IDCaloBoundary_eta(truth->Eta());
408 for (int i = 0; i < 24; ++i) {
409 e->set_eta(i, TFCSExtrapolationState::SUBPOS_ENT, truth->Eta());
410 e->set_eta(i, TFCSExtrapolationState::SUBPOS_EXT, truth->Eta());
411 e->set_eta(i, TFCSExtrapolationState::SUBPOS_MID, truth->Eta());
412 e->set_phi(i, TFCSExtrapolationState::SUBPOS_ENT, 0);
413 e->set_phi(i, TFCSExtrapolationState::SUBPOS_EXT, 0);
414 e->set_phi(i, TFCSExtrapolationState::SUBPOS_MID, 0);
415 e->set_r(i, TFCSExtrapolationState::SUBPOS_ENT, 1500 + i * 10);
416 e->set_r(i, TFCSExtrapolationState::SUBPOS_EXT, 1510 + i * 10);
417 e->set_r(i, TFCSExtrapolationState::SUBPOS_MID, 1505 + i * 10);
418 e->set_z(i, TFCSExtrapolationState::SUBPOS_ENT, 3500 + i * 10);
419 e->set_z(i, TFCSExtrapolationState::SUBPOS_EXT, 3510 + i * 10);
420 e->set_z(i, TFCSExtrapolationState::SUBPOS_MID, 3505 + i * 10);
421 }
422 extrapol = e;
423 }
424
425 // Set energy in layers which then will be retrieved in simulate_hit()
426 simulstate->set_E(0, 1028.77124023);
427 simulstate->set_E(1, 68199.0625);
428 simulstate->set_E(2, 438270.78125);
429 simulstate->set_E(3, 3024.02929688);
430 simulstate->set_E(12, 1330.10131836);
431 simulstate->set_E(1028.77124023 + 68199.0625 + 438270.78125 + 3024.02929688 +
432 1330.10131836);
433 simulstate->set_Efrac(0, simulstate->E(0) / simulstate->E());
434 simulstate->set_Efrac(1, simulstate->E(1) / simulstate->E());
435 simulstate->set_Efrac(2, simulstate->E(2) / simulstate->E());
436 simulstate->set_Efrac(3, simulstate->E(3) / simulstate->E());
437 simulstate->set_Efrac(12, simulstate->E(12) / simulstate->E());
438
439 const int pdgId = truth->pdgid();
440 const float Ekin = truth->Ekin();
441 const float eta = truth->Eta();
442
443 ATH_MSG_NOCLASS(logger, "True energy " << Ekin << " pdgId " << pdgId
444 << " eta " << eta);
445
446 // Find eta bin
447 const int Eta = eta * 10;
448 std::string etaBin = "";
449 for (int i = 0; i <= 25; ++i) {
450 int etaTmp = i * 5;
451 if (Eta >= etaTmp && Eta < (etaTmp + 5)) {
452 etaBin = std::to_string(i * 5) + "_" + std::to_string((i + 1) * 5);
453 }
454 }
455
456 ATH_MSG_NOCLASS(logger, "etaBin = " << etaBin);
457
458 TFCSPredictExtrapWeights NN("NN", "NN");
459 NN.setLevel(MSG::INFO);
460 const int pid = truth->pdgid();
461 NN.initializeNetwork(pid, etaBin, net_path);
462 NN.getNormInputs(etaBin, norm_path);
463
464 // Get extrapWeights and save them as AuxInfo in simulstate
465
466 // Get inputs to Neural Network
467 std::map<std::string, double> inputVariables =
468 NN.prepareInputs(*simulstate, truth->E() * 0.001);
469
470 // Get predicted extrapolation weights
471 ATH_MSG_NOCLASS(logger, "computing with m_nn");
472 auto outputs = NN.m_nn->compute(inputVariables);
473 const std::vector<int> layers = {0, 1, 2, 3, 12};
474 for (int ilayer : layers) {
475 simulstate->setAuxInfo<float>(
476 ilayer, outputs["extrapWeight_" + std::to_string(ilayer)]);
477 }
478
479 // Simulate
480 const int layer = 0;
481 NN.set_calosample(layer);
483 NN.simulate_hit(hit, *simulstate, truth, extrapol);
484
485 // Write
486 TFile *fNN = new TFile("FCSNNtest.root", "RECREATE");
487 NN.Write();
488 fNN->ls();
489 fNN->Close();
490 delete fNN;
491
492 // Open
493 fNN = TFile::Open("FCSNNtest.root");
494 TFCSPredictExtrapWeights *NN2 = (TFCSPredictExtrapWeights *)(fNN->Get("NN"));
495
496 NN2->setLevel(MSG::INFO);
497 NN2->simulate_hit(hit, *simulstate, truth, extrapol);
498 //simulstate->Print();
499
500 return;
501}
502
503void TFCSPredictExtrapWeights::Print(Option_t *option) const {
504 TString opt(option);
505 bool shortprint = opt.Index("short") >= 0;
506 bool longprint = msgLvl(MSG::DEBUG) || (msgLvl(MSG::INFO) && !shortprint);
507 TString optprint = opt;
508 optprint.ReplaceAll("short", "");
510
511 if (longprint)
512 ATH_MSG_INFO(optprint << " m_input (TFCSPredictExtrapWeights): "
514 if (longprint)
515 ATH_MSG_INFO(optprint << " Address of m_nn: " << (void *)m_nn);
516}
const boost::regex ref(r_ef)
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static Double_t ss
#define ATH_MSG_NOCLASS(logger_name, x)
Definition MLogging.h:52
FCSReturnCode
Base class for all FastCaloSim parametrizations Functionality in derivde classes is provided through ...
#define z
Helper for getting a const version of a pointer.
Cut down AthMessaging.
Definition MLogging.h:176
bool msgLvl(const MSG::Level lvl) const
Check whether the logging system is active at the provided verbosity level.
Definition MLogging.h:222
virtual void setLevel(MSG::Level lvl)
Update outputlevel.
Definition MLogging.cxx:105
double phi(int layer, int subpos) const
double z(int layer, int subpos) const
double r(int layer, int subpos) const
double eta(int layer, int subpos) const
TFCSLateralShapeParametrizationHitBase(const char *name=nullptr, const char *title=nullptr)
void Print(Option_t *option="") const override
bool compare(const TFCSParametrizationBase &ref) const
bool compare(const TFCSParametrizationBase &ref) const
Do not persistify!
virtual bool is_match_pdgid(int id) const override
virtual void set_pdgid(int id)
bool compare(const TFCSParametrizationBase &ref) const
virtual FCSReturnCode simulate(TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) const override
Method in all derived classes to do some simulation.
bool getNormInputs(const std::string &etaBin, const std::string &FastCaloTXTInputFolderName)
std::vector< float > * m_normMeans
lwt::LightweightNeuralNetwork * m_nn
virtual bool operator==(const TFCSParametrizationBase &ref) const override
The == operator compares the content of instances.
TFCSPredictExtrapWeights(const char *name=nullptr, const char *title=nullptr)
std::vector< float > * m_normStdDevs
std::vector< int > * m_relevantLayers
void Print(Option_t *option="") const override
std::vector< int > * m_normLayers
Do not persistify.
std::map< std::string, double > prepareInputs(TFCSSimulationState &simulstate, const float truthE) const
virtual FCSReturnCode simulate_hit(Hit &hit, TFCSSimulationState &simulstate, const TFCSTruthState *truth, const TFCSExtrapolationState *extrapol) override
simulated one hit position with some energy.
static void test_path(std::string &net_path, std::string const &norm_path, TFCSSimulationState *simulstate=nullptr, const TFCSTruthState *truth=nullptr, const TFCSExtrapolationState *extrapol=nullptr)
static void unit_test(TFCSSimulationState *simulstate=nullptr, const TFCSTruthState *truth=nullptr, const TFCSExtrapolationState *extrapol=nullptr)
bool initializeNetwork(int pid, const std::string &etaBin, const std::string &FastCaloNNInputFolderName)
double Efrac(int sample) const
void set_E(int sample, double Esample)
bool hasAuxInfo(std::uint32_t index) const
const T getAuxInfo(std::uint32_t index) const
void set_Efrac(int sample, double Efracsample)
void setAuxInfo(std::uint32_t index, const T &val)
void setRandomEngine(CLHEP::HepRandomEngine *engine)
int pdgid() const
double Ekin() const
int r
Definition globals.cxx:22
static Root::TMsgLogger logger("iLumiCalc")
const T * as_const_ptr(const T *p)
Helper for getting a const version of a pointer.
Definition index.py:1
#define IsA
Declare the TObject style functions.