ATLAS Offline Software
Loading...
Searching...
No Matches
DNNCaloSimSvc.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// DNNCaloSimSvc.cxx, (c) ATLAS Detector software //
8
9// class header include
10#include "DNNCaloSimSvc.h"
11
12
13// FastCaloSim includes
20
21//calculators
23
30
31// StoreGate
33
34#include "CaloDetDescr/CaloDetDescrElement.h"
40
42
43#include "lwtnn/parse_json.hh"
44
45#include "TFile.h"
46#include <fstream>
47#include "CLHEP/Random/RandGaussZiggurat.h"
48#include "CLHEP/Random/RandFlat.h"
49#include "CLHEP/Units/SystemOfUnits.h"
50
51using std::abs;
52using std::atan2;
53
55ISF::DNNCaloSimSvc::DNNCaloSimSvc(const std::string& name, ISvcLocator* svc) :
56 BaseSimulationSvc(name, svc),
57 m_graph(nullptr),
58 m_theContainer(nullptr),
59 m_rndGenSvc("AtRndmGenSvc", name),
60 m_randomEngine(nullptr),
61 m_caloDetDescrManager(nullptr),
62 m_caloGeo(nullptr)
63{
64 declareProperty("ParamsInputFilename" , m_paramsFilename="DNNCaloSim/DNNCaloSim_GAN_nn_v0.json"," lwtnn json output describing the trained network");
65 declareProperty("ParamsInputArchitecture" , m_paramsInputArchitecture="GANv0","tag describing additional parameters necessary for the network");
66 declareProperty("CaloCellsOutputName" , m_caloCellsOutputName) ;
67 declareProperty("CaloCellMakerTools_setup" , m_caloCellMakerToolsSetup) ;
68 declareProperty("CaloCellMakerTools_release" , m_caloCellMakerToolsRelease) ;
69 declareProperty("RandomSvc" , m_rndGenSvc );
70 declareProperty("RandomStream" , m_randomEngineName );
71 declareProperty("FastCaloSimCaloExtrapolation" , m_FastCaloSimCaloExtrapolation );
72
73}
74
76= default;
77
80{
81
82 ATH_MSG_INFO(m_screenOutputPrefix << "Initializing ...");
83
84 ATH_CHECK(m_rndGenSvc.retrieve());
86 if (!m_randomEngine)
87 {
88 ATH_MSG_ERROR("Could not get random number engine from RandomNumberService. Abort.");
89 return StatusCode::FAILURE;
90 }
91
94 std::unique_ptr<CaloDetDescrManager> caloMgrPtr = buildCaloDetDescrNoAlign(serviceLocator()
96 ATH_CHECK(detStore()->record(std::move(caloMgrPtr), caloMgrStaticKey));
99 ATH_MSG_ERROR ("Failed to build CaloDetDescrManager.");
100 return StatusCode::FAILURE;
101 }
102 }
103 const FCALDetectorManager * fcalManager=nullptr;
104 ATH_CHECK(detStore()->retrieve(fcalManager));
105
106 const CaloIdManager* caloId_mgr = m_caloDetDescrManager->getCalo_Mgr();
107 m_emID = caloId_mgr->getEM_ID();
108
109 m_caloGeo = std::make_unique<CaloGeometryFromCaloDDM>();
110 m_caloGeo->LoadGeometryFromCaloDDM(m_caloDetDescrManager);
111 if(!m_caloGeo->LoadFCalChannelMapFromFCalDDM(fcalManager) )ATH_MSG_FATAL("Found inconsistency between FCal_Channel map and GEO file. Please, check if they are configured properly.");
112
113
114 // initialize DNN
115 if (initializeNetwork().isFailure())
116 {
117 ATH_MSG_ERROR("Could not initialize network ");
118 return StatusCode::FAILURE;
119
120 }
121
122
123 // Get FastCaloSimCaloExtrapolation
124 if(m_FastCaloSimCaloExtrapolation.retrieve().isFailure())
125 {
126 ATH_MSG_ERROR("FastCaloSimCaloExtrapolation not found ");
127 return StatusCode::FAILURE;
128 }
129
131
132
133
134 return StatusCode::SUCCESS;
135}
136
137// initialize lwtnn network
139{
140
141
142
143 // get neural net JSON file as an std::istream object
144 std::string inputFile=PathResolverFindCalibFile(m_paramsFilename);
145 if (inputFile.empty()){
146 ATH_MSG_ERROR("Could not find json file " << m_paramsFilename );
147 return StatusCode::FAILURE;
148 }
149
150
151
152 ATH_MSG_INFO("Using json file " << inputFile );
153 std::ifstream input(inputFile);
154 // build the graph
155 m_graph=std::make_unique<lwt::LightweightGraph>(lwt::parse_json_graph(input) );
156 if (m_graph==nullptr){
157 ATH_MSG_ERROR("Could not create LightWeightGraph from " << m_paramsFilename );
158 return StatusCode::FAILURE;
159 }
160
161
162
163 // initialize all necessary constants
164 // FIXME eventually all these could be stored in the .json file
165
166 ATH_MSG_INFO("Using ParamsInputArchitecture: " << m_paramsInputArchitecture );
167 if (m_paramsInputArchitecture=="GANv0") // GAN then VAE etc...
168 {
169 m_GANLatentSize = 300;
170 m_logTrueEnergyMean = 9.70406053;
171 m_logTrueEnergyScale = 1.76099569;
172 m_riImpactEtaMean = 3.47603256e-05;
173 m_riImpactEtaScale = 0.00722316;
174 m_riImpactPhiMean = -5.42153684e-05;
175 m_riImpactPhiScale = 0.00708241;
176 }
177
178 if (m_GANLatentSize==0){
179 ATH_MSG_ERROR("m_GANLatentSize uninitialized!.ParamsInputArchitecture= " << m_paramsInputArchitecture );
180 return StatusCode::FAILURE;
181 }
182
183 return StatusCode::SUCCESS;
184}
185
188{
189 ATH_MSG_INFO(m_screenOutputPrefix << "Finalizing ...");
190 return StatusCode::SUCCESS;
191}
192
194{
195 const EventContext& ctx = Gaudi::Hive::currentContext();
196 ATH_MSG_INFO(m_screenOutputPrefix << "setupEvent NEW EVENT! ");
197
199
200 StatusCode sc = evtStore()->record(m_theContainer, m_caloCellsOutputName);
201 if (sc.isFailure())
202 {
203 ATH_MSG_FATAL( m_screenOutputPrefix << "cannot record CaloCellContainer " << m_caloCellsOutputName );
204 return StatusCode::FAILURE;
205 }
206
207 //FIXME why retrieve every event ?
208 CHECK( m_caloCellMakerToolsSetup.retrieve() );
209 ATH_MSG_DEBUG( "Successfully retrieve CaloCellMakerTools: " << m_caloCellMakerToolsSetup );
210 ToolHandleArray<ICaloCellMakerTool>::iterator itrTool = m_caloCellMakerToolsSetup.begin();
211 ToolHandleArray<ICaloCellMakerTool>::iterator endTool = m_caloCellMakerToolsSetup.end();
212 for (; itrTool != endTool; ++itrTool)
213 {
214 std::string chronoName=this->name()+"_"+ itrTool->name();
215 if (m_chrono) m_chrono->chronoStart(chronoName);
216 StatusCode sc = (*itrTool)->process(m_theContainer, ctx);
217 if (m_chrono) {
218 m_chrono->chronoStop(chronoName);
219 ATH_MSG_DEBUG( m_screenOutputPrefix << "Chrono stop : delta " << m_chrono->chronoDelta (chronoName,IChronoStatSvc::USER) * CLHEP::microsecond / CLHEP::second << " second " );
220 }
221
222 if (sc.isFailure())
223 {
224 ATH_MSG_ERROR( m_screenOutputPrefix << "Error executing tool " << itrTool->name() );
225 return StatusCode::FAILURE;
226 }
227 }
228
229
230 // FIXME just for debugging
231 // exercise window building
232 // should be false by default
233 // careful:if enabled change the random number sequence !
234 if (false){
235 TFCSSimulationState testsimulstate(m_randomEngine);
236 const CaloDetDescrElement * testImpactCellDDE;
237 const int ntrial=100;
238 ATH_MSG_INFO ("Trial window building on " << ntrial << " dummy eta phi " );
239 for (int i=0 ; i< ntrial ; i++){
240 const double eta = CLHEP::RandFlat::shoot(testsimulstate.randomEngine(), 0.2, 0.25);
241 const double phi = CLHEP::RandFlat::shoot(testsimulstate.randomEngine(), -CLHEP::pi , CLHEP::pi);
242
243 //randomise eta, phi
244 if (fillWindowCells(eta,phi,testImpactCellDDE).isFailure()){
245 ATH_MSG_WARNING("Could not build trial window cells vector with eta " << eta << " phi " << phi);
246 }
247 }
248 ATH_MSG_INFO ("End of trial window building on " << ntrial << " dummy eta phi " );
249
250 }
251
252
253
254
255 return StatusCode::SUCCESS;
256}
257
259{
260 const EventContext& ctx = Gaudi::Hive::currentContext();
261 ATH_MSG_VERBOSE(m_screenOutputPrefix << "release Event");
262
264
265 //run release tools in a loop
266 ToolHandleArray<ICaloCellMakerTool>::iterator itrTool = m_caloCellMakerToolsRelease.begin();
267 ToolHandleArray<ICaloCellMakerTool>::iterator endTool = m_caloCellMakerToolsRelease.end();
268 for (; itrTool != endTool; ++itrTool)
269 {
270 ATH_MSG_VERBOSE( m_screenOutputPrefix << "Calling tool " << itrTool->name() );
271
272 StatusCode sc = (*itrTool)->process(m_theContainer, ctx);
273
274 if (sc.isFailure())
275 {
276 ATH_MSG_ERROR( m_screenOutputPrefix << "Error executing tool " << itrTool->name() );
277 }
278 }
279
280 return StatusCode::SUCCESS;
281
282}
284{
285 CaloCell_ID::CaloSample aSampling = a->caloDDE()->getSampling();
286 CaloCell_ID::CaloSample bSampling = b->caloDDE()->getSampling();
287
288 float aEtaRaw = a->caloDDE()->eta_raw();
289 float bEtaRaw = b->caloDDE()->eta_raw();
290
291 float aPhiRaw = a->caloDDE()->phi_raw();
292 float bPhiRaw = b->caloDDE()->phi_raw();
293
294 if ((aSampling) < (bSampling))
295 return true;
296 else if ((aSampling) > (bSampling))
297 return false;
298 // reverse sort in eta for left half of detector
299 if ((aEtaRaw) < (bEtaRaw))
300 return false;
301 else if ((aEtaRaw) > (bEtaRaw))
302 return true;
303
304 return CaloPhiRange::diff(aPhiRaw, bPhiRaw) <= 0;
305}
306
308{
309 CaloCell_ID::CaloSample aSampling = a->caloDDE()->getSampling();
310 CaloCell_ID::CaloSample bSampling = b->caloDDE()->getSampling();
311
312 float aEtaRaw = a->caloDDE()->eta_raw();
313 float bEtaRaw = b->caloDDE()->eta_raw();
314
315 float aPhiRaw = a->caloDDE()->phi_raw();
316 float bPhiRaw = b->caloDDE()->phi_raw();
317
318 if ((aSampling) < (bSampling))
319 return true;
320 else if ((aSampling) > (bSampling))
321 return false;
322
323 if ((aEtaRaw) < (bEtaRaw))
324 return true;
325 else if ((aEtaRaw) > (bEtaRaw))
326 return false;
327
328return CaloPhiRange::diff(aPhiRaw, bPhiRaw) <= 0;
329}
330
333{
334
335 ATH_MSG_VERBOSE("NEW PARTICLE! DNNCaloSimSvc called with ISFParticle: " << isfp);
336
337
338 //Don't simulate particles with total energy below 10 MeV
339 if(isfp.ekin() < 10) {
340 ATH_MSG_VERBOSE("Skipping particle with Ekin: " << isfp.ekin() <<" MeV. Below the 10 MeV threshold.");
341 return StatusCode::SUCCESS;
342 }
343
344
345 //Compute all inputs to the network
346 NetworkInputs inputs;
347 double trueEnergy;
348 if (fillNetworkInputs(isfp,inputs,trueEnergy).isFailure()) {
349 ATH_MSG_WARNING("Could not initialize network ");
350 // bail out but do not stop the job
351 return StatusCode::SUCCESS;
352 }
353
354
355 // compute the network output values
356 ATH_MSG_VERBOSE("neural network input = "<<inputs);
357 NetworkOutputs outputs = m_graph->compute(inputs);
358 ATH_MSG_VERBOSE("neural network output = "<<outputs);
359
360
361 // add network output energy in the right place in the calorimeter
362 int itr = 0;
363 for ( auto & windowCell : m_windowCells) {
364 windowCell->addEnergy(trueEnergy * outputs[std::to_string(itr)]);
365 itr++;
366
367 ATH_MSG_VERBOSE(" cell eta_raw " << windowCell->caloDDE()->eta_raw()
368 << " phi_raw " << windowCell->caloDDE()->phi_raw()
369 << " sampling " << windowCell->caloDDE()->getSampling()
370 << " energy " << windowCell->energy());
371 }
372
373
374
375
376
377 return StatusCode::SUCCESS;
378}
379
380
381
382// compute all the necessary inputs to the network
383StatusCode ISF::DNNCaloSimSvc::fillNetworkInputs(const ISF::ISFParticle& isfp, NetworkInputs & inputs, double & trueEnergy)
384{
385 Amg::Vector3D particle_position = isfp.position();
386 Amg::Vector3D particle_direction(isfp.momentum().x(),isfp.momentum().y(),isfp.momentum().z());
387
388
389 TFCSTruthState truth(isfp.momentum().x(),isfp.momentum().y(),isfp.momentum().z(),sqrt(isfp.momentum().mag2()+pow(isfp.mass(),2)),isfp.pdgCode());
390 truth.set_vertex(particle_position[Amg::x], particle_position[Amg::y], particle_position[Amg::z]);
391
392 trueEnergy = isfp.ekin();
393
394
395 TFCSExtrapolationState extrapol;
396 //FIXME this is extrapolating to many layers, when we only need middle layer middle surface
397 //FIXME could have dedicated extrapolation to save time
398 m_FastCaloSimCaloExtrapolation->extrapolate(extrapol,&truth);
399 // extrapol.Print();
400
401
402 if (false)
403 {
404 for (int isam=0; isam< CaloCell_ID_FCS::MaxSample ; isam++){
405 //enum SUBPOS { SUBPOS_MID = 0, SUBPOS_ENT = 1, SUBPOS_EXT = 2}; //MID=middle, ENT=entrance, EXT=exit of cal layer
406
407 for (int isubpos=0; isubpos< 3 ; isubpos++){
408
409 ATH_MSG_VERBOSE("EXTRAPO isam=" << isam <<
410 " isubpos=" << isubpos <<
411 " OK=" << extrapol.OK(isam,isubpos) <<
412 " eta=" << extrapol.eta(isam,isubpos) <<
413 " phi=" << extrapol.phi(isam,isubpos) <<
414 " r=" << extrapol.r(isam,isubpos) );
415
416 }
417 }
418 }
419
420 //FIXME deal with endcap as well
421 int isam=CaloCell_ID_FCS::EMB2;
422 int isubpos=SUBPOS_ENT;
423 double etaExtrap=-999.;
424 double phiExtrap=-999.;
425 if (extrapol.eta(isam,isubpos)) {
426 etaExtrap=extrapol.eta(isam,isubpos);
427 phiExtrap=extrapol.phi(isam,isubpos);
428 }
429
430 ATH_MSG_VERBOSE("Will use isam=" << isam <<
431 " isubpos=" << isubpos <<
432 " eta=" << etaExtrap <<
433 " phi=" << phiExtrap );
434
435
436
437 // fill vector of cells, and DDE of middle cell
438 const CaloDetDescrElement * impactCellDDE;
439 if (fillWindowCells(etaExtrap,phiExtrap,impactCellDDE).isFailure()){
440 ATH_MSG_WARNING("Could not build window cells vector ");
441 return StatusCode::FAILURE;
442 }
443
444
445 // start neural network part
446 // fill a map of input nodes
447 // this is for GAN
448 // most likely it should be specialised as a function of m_ParamsInputArchitecture
449
450 double randGaussz = 0.;
451
452
453 int impactEtaIndex = m_emID->eta(impactCellDDE->identify());
454 int impactPhiIndex = m_emID->phi(impactCellDDE->identify());
455 double etaRawImpactCell=impactCellDDE->eta_raw();
456 double phiRawImpactCell=impactCellDDE->phi_raw();
457
458 ATH_MSG_VERBOSE("impact eta_index " << impactEtaIndex
459 << " phi_index " << impactPhiIndex
460 << " sampling " << m_emID->sampling(impactCellDDE->identify()));
461
462 int pconf = impactPhiIndex % 4 ;
463 int econf = (impactEtaIndex + 1) % 2 ; // ofset corresponds to difference in index calculated for neural net preprocessing
464
465 double riImpactEta = (etaExtrap - etaRawImpactCell);
466 // mirror for left half
467 if (etaExtrap < 0){
468 riImpactEta *= -1.;
469
470 }
471 // scale & centre
472 riImpactEta = (riImpactEta - m_riImpactEtaMean)/m_riImpactEtaScale;
473 double riImpactPhi = (CaloPhiRange::diff(phiExtrap, phiRawImpactCell) - m_riImpactPhiMean)/m_riImpactPhiScale;
474
475
476 // fill randomize latent space
478
479 //FIXME generate in one go
480 for (int i = 0; i< m_GANLatentSize; i ++)
481 {
482 randGaussz = CLHEP::RandGaussZiggurat::shoot(simulstate.randomEngine(), 0., 1.);
483 inputs["Z"].insert ( std::pair<std::string,double>(std::to_string(i), randGaussz) );
484
485 }
486
487 // fill preprocessed true energy
488 inputs["E_true"].insert ( std::pair<std::string,double>("0", (std::log(trueEnergy) - m_logTrueEnergyMean)/m_logTrueEnergyScale) );
489 // fill p,e configurations multi-hot vector
490 for (int i = 0; i< 4; i ++)
491 {
492 if (i == pconf){
493 inputs["pconfig"].insert ( std::pair<std::string,double>(std::to_string(i),1.) );
494 }
495 else{
496 inputs["pconfig"].insert ( std::pair<std::string,double>(std::to_string(i),0.) );
497 }
498 }
499 for (int i = 0; i< 2; i ++){
500 if (i == econf){
501 inputs["econfig"].insert ( std::pair<std::string,double>(std::to_string(i),1.) );
502 }
503 else{
504 inputs["econfig"].insert ( std::pair<std::string,double>(std::to_string(i),0.) );
505 }
506 }
507 // fill position of extrap particle in impact cell
508 inputs["ripos"].insert ( std::pair<std::string,double>("0", riImpactEta) );
509 inputs["ripos"].insert ( std::pair<std::string,double>("1", riImpactPhi ) );
510
511 return StatusCode::SUCCESS;
512}
513
514
515StatusCode ISF::DNNCaloSimSvc::fillWindowCells(const double etaExtrap,const double phiExtrap,const CaloDetDescrElement* & impactCellDDE) {
516
517 //now find the cell it corresponds to
518 //FIXME this is barrel should also look in endcap
519 // (note that is really looking up eta, phi, not raw eta phi
520 impactCellDDE=m_caloDetDescrManager->get_element(CaloCell_ID::EMB2,etaExtrap,phiExtrap);
521 if (impactCellDDE==nullptr){
522 ATH_MSG_WARNING("No cell found for this eta " << etaExtrap << " phi " << phiExtrap);
523 return StatusCode::FAILURE;
524 }
525
526
527
528
529 const int caloHashImpactCell=impactCellDDE->calo_hash();
530 const double etaImpactCell=impactCellDDE->eta();
531 const double phiImpactCell=impactCellDDE->phi();
532 const double etaRawImpactCell=impactCellDDE->eta_raw();
533 const double phiRawImpactCell=impactCellDDE->phi_raw();
534
535
536 ATH_MSG_VERBOSE("impact cell calohash=" << caloHashImpactCell <<
537 " eta=" << etaImpactCell <<
538 " phi=" << phiImpactCell <<
539 " eta raw=" << etaRawImpactCell <<
540 " phi raw=" << phiRawImpactCell );
541
542
543 int nSqCuts = 0;
544
545
546 // select the cells DNN will simulate
547 // note that m_theCellContainer has all the calorimeter cells
548 // this will hold the list of cells to simulate
549 //FIXME this can be sped up certainly
550 m_windowCells.clear();
552 float eta_raw;
553 float phi_raw;
554 for(CaloCell* theCell : * m_theContainer) {
555 sampling = theCell->caloDDE()->getSampling();
556 eta_raw = theCell->caloDDE()->eta_raw();
557 phi_raw = theCell->caloDDE()->phi_raw();
558 if (( eta_raw < etaRawImpactCell + m_etaRawBackCut) && (eta_raw > etaRawImpactCell - m_etaRawBackCut)) {
559 if ((CaloPhiRange::diff(phi_raw, phiRawImpactCell) < m_phiRawStripCut ) && (CaloPhiRange::diff(phi_raw, phiRawImpactCell) > - m_phiRawStripCut) ) {
560
561 }
562 else{
563 continue;
564 }
565 }
566 else{
567 continue;
568 }
569
570 if ((sampling == 0) || (sampling == 1) ){
571 if ((eta_raw < etaRawImpactCell + m_etaRawMiddleCut) && (eta_raw > etaRawImpactCell - m_etaRawMiddleCut)) {
572 nSqCuts ++;
573 // add to vector
574 m_windowCells.push_back(theCell);
575
576 }
577 }
578 else if(sampling == 2) {
579 if ((CaloPhiRange::diff(phi_raw , phiRawImpactCell) < m_phiRawMiddleCut) && (CaloPhiRange::diff(phi_raw, phiRawImpactCell) > - m_phiRawMiddleCut)) {
580 if ((eta_raw < etaRawImpactCell + m_etaRawMiddleCut) && (eta_raw > etaRawImpactCell - m_etaRawMiddleCut)) {
581 nSqCuts ++;
582 m_windowCells.push_back(theCell);
583 }
584 }
585 }
586
587 else if(sampling == 3){
588 if ((CaloPhiRange::diff(phi_raw, phiRawImpactCell) < m_phiRawMiddleCut) && (CaloPhiRange::diff(phi_raw , phiRawImpactCell) > - m_phiRawMiddleCut )) {
589 nSqCuts ++;
590 m_windowCells.push_back(theCell);
591 }
592
593 }
594 }
595
596 if (nSqCuts != m_numberOfCellsForDNN){
597 ATH_MSG_WARNING("Total cells passing DNN selection is " << nSqCuts << " but should be " << m_numberOfCellsForDNN );
598 // bail out, but do not stop the job
599 return StatusCode::FAILURE;
600
601 }
602
603 //sort cells within the cluster like they are fed to DNN
604 if (etaRawImpactCell < 0){
606 }
607 else{
609 }
610
611 return StatusCode::SUCCESS;
612
613}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
std::unique_ptr< CaloDetDescrManager > buildCaloDetDescrNoAlign(ISvcLocator *svcLocator, IMessageSvc *msgSvc)
Definition of CaloDetDescrManager.
static constexpr const char * caloMgrStaticKey
CaloPhiRange class declaration.
#define CHECK(...)
Evaluate an expression and check for errors.
bool compCellsForDNNSortMirror(const CaloCell *a, const CaloCell *b)
bool compCellsForDNNSort(const CaloCell *a, const CaloCell *b)
static Double_t a
static Double_t sc
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
constexpr int pow(int base, int exp) noexcept
Container class for CaloCell.
CaloSampling::CaloSample CaloSample
Definition CaloCell_ID.h:53
Data object for each calorimeter readout cell.
Definition CaloCell.h:57
This class groups all DetDescr information related to a CaloCell.
Identifier identify() const override final
cell identifier
This class provides the client interface for accessing the detector description information common to...
This class initializes the Calo (LAr and Tile) offline identifiers.
const LArEM_ID * getEM_ID(void) const
static double diff(double phi1, double phi2)
simple phi1 - phi2 calculation, but result is fixed to respect range.
A manager class providing access to readout geometry information for the forward calorimeter.
const ServiceHandle< StoreGateSvc > & detStore() const
The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
Gaudi::Property< std::string > m_screenOutputPrefix
Screen output refinement.
const ServiceHandle< StoreGateSvc > & evtStore() const
The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
BaseSimulationSvc(const std::string &name, ISvcLocator *pSvcLocator)
ServiceHandle< IChronoStatSvc > m_chrono
The timing service for general usage.
const CaloDetDescrManager * m_caloDetDescrManager
std::vector< CaloCell * > m_windowCells
virtual StatusCode initialize() override final
Athena algorithm's interface methods.
const double m_phiRawMiddleCut
virtual StatusCode simulate(ISFParticle &isp, McEventCollection *) override final
Simulation Call.
const int m_numberOfCellsForDNN
StatusCode fillWindowCells(const double etaExtrap, const double phiExtrap, const CaloDetDescrElement *&impactCellDDE)
std::string m_caloCellsOutputName
CaloCellContainer * m_theContainer
StatusCode fillNetworkInputs(const ISF::ISFParticle &isfp, NetworkInputs &inputs, double &trueEnergy)
std::string m_paramsInputArchitecture
DNNCaloSimSvc(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters.
std::string m_randomEngineName
const LArEM_ID * m_emID
virtual StatusCode finalize() override final
framework methods
ServiceHandle< IAtRndmGenSvc > m_rndGenSvc
StatusCode initializeNetwork()
helper for initialize
std::unique_ptr< lwt::LightweightGraph > m_graph
virtual StatusCode releaseEvent() override final
Release Event chain - in case of an end-of event action is needed.
const double m_phiRawStripCut
CLHEP::HepRandomEngine * m_randomEngine
ToolHandleArray< ICaloCellMakerTool > m_caloCellMakerToolsSetup
const double m_etaRawBackCut
ToolHandle< IFastCaloSimCaloExtrapolation > m_FastCaloSimCaloExtrapolation
std::unique_ptr< CaloGeometryFromCaloDDM > m_caloGeo
std::map< std::string, std::map< std::string, double > > NetworkInputs
ToolHandleArray< ICaloCellMakerTool > m_caloCellMakerToolsRelease
std::map< std::string, double > NetworkOutputs
virtual ~DNNCaloSimSvc() final
Destructor.
std::string m_paramsFilename
const double m_etaRawMiddleCut
virtual StatusCode setupEvent() override final
Setup Event chain - in case of a begin-of event action is needed.
The generic ISF particle definition,.
Definition ISFParticle.h:42
double ekin() const
Kinetic energy.
const Amg::Vector3D & momentum() const
The current momentum vector of the ISFParticle.
const Amg::Vector3D & position() const
The current position of the ISFParticle.
int pdgCode() const
PDG value.
double mass() const
mass of the particle
This defines the McEventCollection, which is really just an ObjectVector of McEvent objectsFile: Gene...
CLHEP::HepRandomEngine * randomEngine()
void set_vertex(const TLorentzVector &val)
singleton-like access to IMessageSvc via open function and helper
Eigen::Matrix< double, 3, 1 > Vector3D
IMessageSvc * getMessageSvc(bool quiet=false)
@ VIEW_ELEMENTS
this data object is a view, it does not own its elmts
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.