Skip to content
Snippets Groups Projects
Commit 9215c5bb authored by Brian Petersen's avatar Brian Petersen
Browse files

Merge branch 'L1CaloFEXSim_reentrant' into '23.0'

L1CaloFEXSim: make jTowerMakerFromJfexTowers algorithm reentrant

See merge request !64836
parents 04526c29 a8a7b748
No related branches found
No related tags found
2 merge requests!648892023-08-08: daily sweep of 23.0 into main,!64836L1CaloFEXSim: make jTowerMakerFromJfexTowers algorithm reentrant
......@@ -34,8 +34,8 @@ class jSuperCellTowerMapper: public AthAlgTool, virtual public IjSuperCellTowerM
/** standard Athena-Algorithm method */
virtual StatusCode initialize() override;
virtual StatusCode AssignSuperCellsToTowers(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) override;
virtual StatusCode AssignTriggerTowerMapper(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) override;
virtual StatusCode AssignSuperCellsToTowers(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) const override;
virtual StatusCode AssignTriggerTowerMapper(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) const override;
virtual void reset() override;
......@@ -47,12 +47,12 @@ class jSuperCellTowerMapper: public AthAlgTool, virtual public IjSuperCellTowerM
Gaudi::Property<bool> m_apply_masking {this, "SCellMasking", false, "Applies masking. Only use for data"};
virtual int FindAndConnectTower(/*jTowerContainer**/std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,CaloSampling::CaloSample sample,const int region, int layer, const int pos_neg, const int eta_index, const int phi_index, Identifier ID, float et, int prov, bool doPrint,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0) override;
virtual void ConnectSuperCellToTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw, int iETower, Identifier ID, int iCell, float et, int layer) override;
virtual int FindTowerIDForSuperCell(int towereta, int towerphi) override;
virtual void PrintCellSpec(const CaloSampling::CaloSample sample, int layer, const int region, const int eta_index, const int phi_index, const int pos_neg, int iETower, int iCell, int prov, Identifier ID, bool doenergysplit,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0,bool cellValid=true) override;
virtual int FindAndConnectTower(/*jTowerContainer**/std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,CaloSampling::CaloSample sample,const int region, int layer, const int pos_neg, const int eta_index, const int phi_index, Identifier ID, float et, int prov, bool doPrint,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0) const override;
virtual void ConnectSuperCellToTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw, int iETower, Identifier ID, int iCell, float et, int layer) const override;
virtual int FindTowerIDForSuperCell(int towereta, int towerphi) const override;
virtual void PrintCellSpec(const CaloSampling::CaloSample sample, int layer, const int region, const int eta_index, const int phi_index, const int pos_neg, int iETower, int iCell, int prov, Identifier ID, bool doenergysplit,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0,bool cellValid=true) const override;
std::string DectectorName(const CaloSampling::CaloSample sample);
std::string DectectorName(const CaloSampling::CaloSample sample) const;
};
......
......@@ -32,9 +32,9 @@ class jTowerBuilder: public AthAlgTool, virtual public IjTowerBuilder {
virtual ~jTowerBuilder() = default;
virtual StatusCode initialize() override;
virtual void init(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) override ;
virtual void execute(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) override ;
virtual void reset() override ;
virtual void init(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const override ;
virtual void execute(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const override ;
virtual void reset() const override ;
......@@ -56,14 +56,14 @@ class jTowerBuilder: public AthAlgTool, virtual public IjTowerBuilder {
//property for jFEX mapping
Gaudi::Property<std::string> m_PileupWeigthFile {this, "PileupWeigthFile", "Run3L1CaloSimulation/Noise/jTowerCorrection.20210308.r12406.root", "Root file for the pileup weight"};
Gaudi::Property<std::string> m_PileupHelperFile {this, "PileupHelperFile", "Run3L1CaloSimulation/Calibrations/jFEX_MatchedMapping.2022Mar10.r12406.root", "Root file to set the jTower coordinated (float eta/phi)"};
Gaudi::Property<std::string> m_PileupHelperFile {this, "PileupHelperFile", "Run3L1CaloSimulation/Calibrations/jFEX_MatchedMapping.2022Mar10.r12406.root", "Root file to set the jTower coordinates (float eta/phi)"};
//histograms need to set coordinates and noise subtraction
TH1F* m_jTowerArea_hist;
TH1I* m_Firmware2BitwiseID;
TH1I* m_BinLayer;
TH1F* m_EtaCoords;
TH1F* m_PhiCoords;
TH1F* m_jTowerArea_hist = nullptr;
TH1I* m_Firmware2BitwiseID = nullptr;
TH1I* m_BinLayer = nullptr;
TH1F* m_EtaCoords = nullptr;
TH1F* m_PhiCoords = nullptr;
};
......
......@@ -58,7 +58,7 @@ StatusCode jSuperCellTowerMapper::initialize()
}
StatusCode jSuperCellTowerMapper::AssignTriggerTowerMapper(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) {
StatusCode jSuperCellTowerMapper::AssignTriggerTowerMapper(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) const {
static constexpr float pi_over_32 = M_PI/32;
......@@ -104,7 +104,7 @@ void jSuperCellTowerMapper::reset(){
}
// works for real supercells from MC
StatusCode jSuperCellTowerMapper::AssignSuperCellsToTowers(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw)
StatusCode jSuperCellTowerMapper::AssignSuperCellsToTowers(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) const
{
bool doPrint = false;
......@@ -306,7 +306,7 @@ void jSuperCellTowerMapper::reset(){
}
void jSuperCellTowerMapper::ConnectSuperCellToTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,int iJTower, Identifier ID, int iCell, float et, int layer){
void jSuperCellTowerMapper::ConnectSuperCellToTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,int iJTower, Identifier ID, int iCell, float et, int layer) const {
LVL1::jTower * tmpTower = my_jTowerContainerRaw->findTower(iJTower);
if(tmpTower){
......@@ -315,7 +315,7 @@ void jSuperCellTowerMapper::ConnectSuperCellToTower(std::unique_ptr<jTowerContai
}
int jSuperCellTowerMapper::FindAndConnectTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,CaloSampling::CaloSample sample,const int region, int layer, const int pos_neg, const int eta_index, const int phi_index, Identifier ID, float et, int prov,bool doPrint, float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0)
int jSuperCellTowerMapper::FindAndConnectTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,CaloSampling::CaloSample sample,const int region, int layer, const int pos_neg, const int eta_index, const int phi_index, Identifier ID, float et, int prov,bool doPrint, float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0) const
{
// bool for the special case of 1.8 < eta < 2.0 only in the front layer
......@@ -1108,21 +1108,16 @@ int jSuperCellTowerMapper::FindAndConnectTower(std::unique_ptr<jTowerContainer>
// END ITERATING OVER SUPER CELLS+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ +++++++++++++++++++++++++++++++++++++++++++++
return 1;
}
int jSuperCellTowerMapper::FindTowerIDForSuperCell(int towereta, int towerphi)
int jSuperCellTowerMapper::FindTowerIDForSuperCell(int towereta, int towerphi) const
{
return (towerphi + (64 * towereta));
}
void jSuperCellTowerMapper::PrintCellSpec(const CaloSampling::CaloSample sample, int layer, const int region, const int eta_index, const int phi_index, const int pos_neg, int iJTower, int iCell, int prov, Identifier ID ,bool doenergysplit, float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0,bool cellValid)
void jSuperCellTowerMapper::PrintCellSpec(const CaloSampling::CaloSample sample, int layer, const int region, const int eta_index, const int phi_index, const int pos_neg, int iJTower, int iCell, int prov, Identifier ID ,bool doenergysplit, float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0,bool cellValid) const
{
std::string sampleName = "";
......@@ -1196,7 +1191,7 @@ int jSuperCellTowerMapper::FindTowerIDForSuperCell(int towereta, int towerphi)
}
std::string jSuperCellTowerMapper::DectectorName(const CaloSampling::CaloSample sample){
std::string jSuperCellTowerMapper::DectectorName(const CaloSampling::CaloSample sample) const {
std::string sampleName ="";
switch (sample) {
case CaloSampling::PreSamplerB: { sampleName = "PreSamplerB"; break; }
......
......@@ -63,7 +63,7 @@ StatusCode jTowerBuilder::initialize()
void jTowerBuilder::init(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) {
void jTowerBuilder::init(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
execute(jTowerContainerRaw);
......@@ -74,12 +74,12 @@ void jTowerBuilder::init(std::unique_ptr<jTowerContainer> & jTowerContainerRaw)
}
void jTowerBuilder::reset() {
void jTowerBuilder::reset() const {
}
void jTowerBuilder::execute(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) {
void jTowerBuilder::execute(std::unique_ptr<jTowerContainer> & jTowerContainerRaw) const {
BuildAllTowers(jTowerContainerRaw);
}
......@@ -386,7 +386,6 @@ void jTowerBuilder::BuildAllTowers(std::unique_ptr<jTowerContainer> & jTowerCont
BuildEMEjTowers(jTowerContainerRaw);
BuildEMIEjTowers(jTowerContainerRaw);
BuildFCALjTowers(jTowerContainerRaw);
}
} // end of LVL1 namespace
......
......@@ -34,7 +34,7 @@
namespace LVL1 {
jTowerMakerFromJfexTowers::jTowerMakerFromJfexTowers(const std::string& name, ISvcLocator* pSvcLocator)
: AthAlgorithm(name, pSvcLocator)
: AthReentrantAlgorithm(name, pSvcLocator)
{}
......@@ -53,7 +53,7 @@ StatusCode jTowerMakerFromJfexTowers::initialize()
}
StatusCode jTowerMakerFromJfexTowers::execute()
StatusCode jTowerMakerFromJfexTowers::execute(const EventContext& ctx) const
{
ATH_MSG_DEBUG("Executing " << name() << ", processing event number " );
......@@ -61,7 +61,7 @@ StatusCode jTowerMakerFromJfexTowers::execute()
SG::ReadHandle<xAOD::jFexTowerContainer> jDataTowerContainer;
bool jDataTowerFilled = false;
if(!m_isMC){
jDataTowerContainer = SG::ReadHandle<xAOD::jFexTowerContainer>(m_DataTowerKey);
jDataTowerContainer = SG::ReadHandle<xAOD::jFexTowerContainer>(m_DataTowerKey, ctx);
if(!jDataTowerContainer.isValid()) {
ATH_MSG_FATAL("Could not retrieve collection " << jDataTowerContainer.key() );
return StatusCode::FAILURE;
......@@ -73,7 +73,7 @@ StatusCode jTowerMakerFromJfexTowers::execute()
SG::ReadHandle<xAOD::jFexTowerContainer> jEmulatedTowerContainer;
if(m_UseEmulated){
jEmulatedTowerContainer = SG::ReadHandle<xAOD::jFexTowerContainer>(m_EmulTowerKey);
jEmulatedTowerContainer = SG::ReadHandle<xAOD::jFexTowerContainer>(m_EmulTowerKey, ctx);
if(!jEmulatedTowerContainer.isValid()) {
ATH_MSG_FATAL("Could not retrieve collection " << jEmulatedTowerContainer.key() );
return StatusCode::FAILURE;
......@@ -139,7 +139,7 @@ StatusCode jTowerMakerFromJfexTowers::execute()
}
// STEP 3 - Write the completed jTowerContainer into StoreGate (move the local copy in memory)
SG::WriteHandle<LVL1::jTowerContainer> jTowerContainerSG(m_jTowerContainerSGKey);
SG::WriteHandle<LVL1::jTowerContainer> jTowerContainerSG(m_jTowerContainerSGKey, ctx);
ATH_CHECK(jTowerContainerSG.record(std::move( local_jTowerContainerRaw ) ) );
// STEP 4 - Close and clean the event
......
......@@ -9,7 +9,7 @@
#include <string>
// Athena/Gaudi
#include "AthenaBaseComps/AthAlgorithm.h"
#include "AthenaBaseComps/AthReentrantAlgorithm.h"
#include "L1CaloFEXSim/jTowerContainer.h"
#include "xAODTrigL1Calo/jFexTowerContainer.h"
#include "L1CaloFEXSim/jSuperCellTowerMapper.h"
......@@ -19,7 +19,7 @@ class CaloIdManager;
namespace LVL1 {
class jTowerMakerFromJfexTowers : public AthAlgorithm
class jTowerMakerFromJfexTowers : public AthReentrantAlgorithm
{
public:
......@@ -27,7 +27,8 @@ class jTowerMakerFromJfexTowers : public AthAlgorithm
~jTowerMakerFromJfexTowers() = default;
virtual StatusCode initialize() override;
virtual StatusCode execute() override;
virtual StatusCode execute(const EventContext& ctx) const override;
private:
......
/*
Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
*/
//***************************************************************************
......@@ -36,16 +36,16 @@ Interface definition for jSuperCellTowerMapper
public:
static const InterfaceID& interfaceID( ) ;
virtual StatusCode AssignSuperCellsToTowers(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) = 0;
virtual StatusCode AssignTriggerTowerMapper(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) = 0;
virtual StatusCode AssignSuperCellsToTowers(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) const = 0;
virtual StatusCode AssignTriggerTowerMapper(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw) const = 0;
virtual void reset() = 0;
virtual int FindAndConnectTower(/*jTowerContainer**/std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,CaloSampling::CaloSample sample,const int region, int layer, const int pos_neg, const int eta_index, const int phi_index, Identifier ID, float et, int prov, bool doPrint,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0) = 0;
virtual void ConnectSuperCellToTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw, int iETower, Identifier ID, int iCell, float et, int layer) = 0;
virtual int FindTowerIDForSuperCell(int towereta, int towerphi) = 0;
virtual void PrintCellSpec(const CaloSampling::CaloSample sample, int layer, const int region, const int eta_index, const int phi_index, const int pos_neg, int iETower, int iCell, int prov, Identifier ID, bool doenergysplit,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0, bool cellValid = false) = 0;
virtual int FindAndConnectTower(/*jTowerContainer**/std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw,CaloSampling::CaloSample sample,const int region, int layer, const int pos_neg, const int eta_index, const int phi_index, Identifier ID, float et, int prov, bool doPrint,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0) const = 0;
virtual void ConnectSuperCellToTower(std::unique_ptr<jTowerContainer> & my_jTowerContainerRaw, int iETower, Identifier ID, int iCell, float et, int layer) const = 0;
virtual int FindTowerIDForSuperCell(int towereta, int towerphi) const = 0;
virtual void PrintCellSpec(const CaloSampling::CaloSample sample, int layer, const int region, const int eta_index, const int phi_index, const int pos_neg, int iETower, int iCell, int prov, Identifier ID, bool doenergysplit,float eta_min, float eta_max, float eta0, float phi_min, float phi_max, float phi0, bool cellValid = false) const = 0;
private:
......
......@@ -31,9 +31,9 @@ Interface definition for jTowerBuilder
virtual StatusCode initialize() = 0;
virtual void init(std::unique_ptr<jTowerContainer> & jTowerContainer) = 0;
virtual void execute(std::unique_ptr<jTowerContainer> & jTowerContainer) = 0;
virtual void reset() = 0;
virtual void init(std::unique_ptr<jTowerContainer> & jTowerContainer) const = 0;
virtual void execute(std::unique_ptr<jTowerContainer> & jTowerContainer) const = 0;
virtual void reset() const = 0;
private:
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment