Skip to content
Snippets Groups Projects
Commit 78e4dc4b authored by Maik Nijhuis's avatar Maik Nijhuis
Browse files

Remove singletons

parent ab5ae22f
No related branches found
Tags v0.4.0
1 merge request!258Remove singletons
Pipeline #37525 passed with warnings
// Copyright (C) 2020 ASTRON (Netherlands Institute for Radio Astronomy)
// SPDX-License-Identifier: GPL-3.0-or-later
#ifndef EVERYBEAM_SINGLETON_H
#define EVERYBEAM_SINGLETON_H
namespace everybeam {
namespace common {
template <typename T>
struct Singleton {
// Factory function to obtain the one-and-only instance of the Singleton
static std::shared_ptr<T> GetInstance() {
// Static variable, initialized on first call to GetInstance()
static std::shared_ptr<T> instance = std::make_shared<T>();
return instance;
}
private:
// Make the constructor private, to prevent direct instantiation
// outside of the factory function GetInstance()
Singleton() {}
// Forbid to make copies of the Singleton, by deleting the
// copy and assignment constructors
Singleton(Singleton const&) = delete;
void operator=(Singleton const&) = delete;
};
} // namespace common
} // namespace everybeam
#endif
\ No newline at end of file
......@@ -13,6 +13,8 @@
#include <algorithm>
#include "config.h"
#include "common/mathutils.h"
#include "elementresponsefixeddirection.h"
......@@ -80,4 +82,10 @@ std::shared_ptr<ElementResponse> ElementResponse::FixateDirection(
shared_from_this(), thetaphi[0], thetaphi[1]);
}
std::filesystem::path ElementResponse::GetPath(
const std::filesystem::path& relative_path) {
// TODO (AST-803) Allow specifying a different data directory at run time.
return std::filesystem::path(EVERYBEAM_DATA_DIR) / relative_path;
}
} // namespace everybeam
......@@ -5,6 +5,7 @@
#define EVERYBEAM_ELEMENTRESPONSE_H
#include <complex>
#include <filesystem>
#include <memory>
#include <ostream>
#include <vector>
......@@ -97,9 +98,18 @@ class ElementResponse : public std::enable_shared_from_this<ElementResponse> {
return Response(freq, theta, phi);
}
static std::shared_ptr<ElementResponse> GetInstance(
static std::shared_ptr<const ElementResponse> GetInstance(
ElementResponseModel model, const std::string& name,
const Options& options);
protected:
/**
* Get the path to an EveryBeam data file or directory.
* @param relative_path A path relative to the EveryBeam data directory.
* @return The full path of the data file or data directory.
*/
static std::filesystem::path GetPath(
const std::filesystem::path& relative_path);
};
} // namespace everybeam
#endif
......@@ -15,18 +15,18 @@
namespace everybeam {
std::shared_ptr<ElementResponse> ElementResponse::GetInstance(
std::shared_ptr<const ElementResponse> ElementResponse::GetInstance(
ElementResponseModel model, const std::string& name,
const Options& options) {
switch (model) {
case ElementResponseModel::kHamaker:
return HamakerElementResponse::GetInstance(name);
return std::make_shared<HamakerElementResponse>(name);
case ElementResponseModel::kHamakerLba:
return HamakerElementResponse::GetLbaInstance();
return std::make_shared<HamakerElementResponse>("LBA");
case ElementResponseModel::kOSKARDipole:
return OSKARElementResponseDipole::GetInstance();
return std::make_shared<OSKARElementResponseDipole>();
case ElementResponseModel::kOSKARSphericalWave:
return OSKARElementResponseSphericalWave::GetInstance();
return std::make_shared<OSKARElementResponseSphericalWave>();
case ElementResponseModel::kLOBES:
try {
return LOBESElementResponse::GetInstance(name, options);
......@@ -35,7 +35,7 @@ std::shared_ptr<ElementResponse> ElementResponse::GetInstance(
<< " failed because: " << std::endl;
std::cout << e.what() << std::endl;
std::cout << "Switching to HamakerElementResponse instead" << std::endl;
return GetInstance(ElementResponseModel::kHamaker, name, options);
return std::make_shared<HamakerElementResponse>(name);
}
default:
std::stringstream message;
......
......@@ -12,12 +12,12 @@ H5::CompType GetComplexDoubleType() {
}
size_t HamakerCoefficients::GetIndex(unsigned int n, unsigned int t,
unsigned int f) {
unsigned int f) const {
return (n * nPowerTheta_ + t) * nPowerFreq_ * nInner_ + f * nInner_;
}
// Constructor for reading coeff from file
HamakerCoefficients::HamakerCoefficients(std::string& filename) {
HamakerCoefficients::HamakerCoefficients(const std::string& filename) {
ReadCoefficients(filename);
}
......@@ -58,13 +58,13 @@ void HamakerCoefficients::SetCoefficients(
void HamakerCoefficients::GetCoefficient(
unsigned int n, unsigned int t, unsigned int f,
std::pair<std::complex<double>, std::complex<double>>& value) {
std::pair<std::complex<double>, std::complex<double>>& value) const {
size_t index = GetIndex(n, t, f);
value.first = coeff_[index];
value.second = coeff_[index + 1];
}
void HamakerCoefficients::ReadCoefficients(std::string& filename) {
void HamakerCoefficients::ReadCoefficients(const std::string& filename) {
// Open file
H5::H5File file(filename, H5F_ACC_RDONLY);
......@@ -94,7 +94,7 @@ void HamakerCoefficients::ReadCoefficients(std::string& filename) {
dataset.read(coeff_.data(), data_type, dataspace);
}
void HamakerCoefficients::WriteCoefficients(std::string& filename) {
void HamakerCoefficients::WriteCoefficients(const std::string& filename) {
// Open file
H5::H5File file(filename, H5F_ACC_TRUNC);
......
......@@ -21,7 +21,7 @@ class HamakerCoefficients {
HamakerCoefficients();
//! Constructor for reading coeff from file
HamakerCoefficients(std::string& filename);
HamakerCoefficients(const std::string& filename);
//! Constructor for writing coeff to file
HamakerCoefficients(const double freq_center, const double freq_range,
......@@ -60,19 +60,19 @@ class HamakerCoefficients {
void GetCoefficient(
unsigned int n, unsigned int t, unsigned int f,
std::pair<std::complex<double>, std::complex<double>>& value);
std::pair<std::complex<double>, std::complex<double>>& value) const;
// HDF5 I/O
void ReadCoefficients(std::string& filename);
void ReadCoefficients(const std::string& filename);
void WriteCoefficients(std::string& filename);
void WriteCoefficients(const std::string& filename);
// Debugging
void PrintCoefficients();
private:
// Methods
size_t GetIndex(unsigned int n, unsigned int t, unsigned int f);
size_t GetIndex(unsigned int n, unsigned int t, unsigned int f) const;
// Parameters
double freq_center_;
......
......@@ -10,35 +10,30 @@
#include "config.h"
#include "hamakerelementresponse.h"
#include "../common/singleton.h"
#include <aocommon/matrix2x2.h>
namespace everybeam {
std::shared_ptr<HamakerElementResponse>
HamakerElementResponse::GetLbaInstance() {
return common::Singleton<HamakerElementResponseLBA>::GetInstance();
}
std::shared_ptr<HamakerElementResponse> HamakerElementResponse::GetInstance(
const std::string& name) {
HamakerElementResponse::HamakerElementResponse(const std::string& name) {
if (name.find("LBA") != std::string::npos) {
return common::Singleton<HamakerElementResponseLBA>::GetInstance();
}
if (name.find("HBA") != std::string::npos) {
return common::Singleton<HamakerElementResponseHBA>::GetInstance();
coefficients_ = cached_lba_coefficients_.lock();
if (!coefficients_) {
const std::string path = GetPath("HamakerLBACoeff.h5");
coefficients_ = std::make_shared<HamakerCoefficients>(path);
cached_lba_coefficients_ = coefficients_;
}
} else if (name.find("HBA") != std::string::npos) {
coefficients_ = cached_hba_coefficients_.lock();
if (!coefficients_) {
const std::string path = GetPath("HamakerHBACoeff.h5");
coefficients_ = std::make_shared<HamakerCoefficients>(path);
cached_hba_coefficients_ = coefficients_;
}
} else {
throw std::invalid_argument(
"HamakerElementResponse: name should end in either 'LBA' or 'HBA'");
}
throw std::invalid_argument(
"HamakerElementResponse::GetInstance: name should end in either 'LBA' or "
"'HBA'");
}
std::string HamakerElementResponse::GetPath(const char* filename) const {
std::stringstream ss;
ss << EVERYBEAM_DATA_DIR << "/";
ss << filename;
return ss.str();
}
aocommon::MC2x2 HamakerElementResponse::Response(double freq, double theta,
......@@ -50,11 +45,11 @@ aocommon::MC2x2 HamakerElementResponse::Response(double freq, double theta,
return response;
}
const double freq_center = coeffs_->GetFreqCenter();
const double freq_range = coeffs_->GetFreqRange();
const unsigned int nHarmonics = coeffs_->Get_nHarmonics();
const unsigned int nPowerTheta = coeffs_->Get_nPowerTheta();
const unsigned int nPowerFreq = coeffs_->Get_nPowerFreq();
const double freq_center = coefficients_->GetFreqCenter();
const double freq_range = coefficients_->GetFreqRange();
const unsigned int nHarmonics = coefficients_->Get_nHarmonics();
const unsigned int nPowerTheta = coefficients_->Get_nPowerTheta();
const unsigned int nPowerFreq = coefficients_->Get_nPowerFreq();
// The model is parameterized in terms of a normalized frequency in the
// range [-1, 1]. The appropriate conversion is taken care of below.
......@@ -78,19 +73,20 @@ aocommon::MC2x2 HamakerElementResponse::Response(double freq, double theta,
// start indexing the block of coefficients at the last element
// Evaluate the highest order term.
coeffs_->GetCoefficient(k, nPowerTheta - 1, nPowerFreq - 1, P);
coefficients_->GetCoefficient(k, nPowerTheta - 1, nPowerFreq - 1, P);
for (unsigned int i = 0; i < nPowerFreq - 1; ++i) {
coeffs_->GetCoefficient(k, nPowerTheta - 1, nPowerFreq - i - 2, Pk);
coefficients_->GetCoefficient(k, nPowerTheta - 1, nPowerFreq - i - 2, Pk);
P.first = P.first * freq + Pk.first;
P.second = P.second * freq + Pk.second;
}
// Evaluate the remaining terms.
for (unsigned int j = 0; j < nPowerTheta - 1; ++j) {
coeffs_->GetCoefficient(k, nPowerTheta - j - 2, nPowerFreq - 1, Pj);
coefficients_->GetCoefficient(k, nPowerTheta - j - 2, nPowerFreq - 1, Pj);
for (unsigned int i = 0; i < nPowerFreq - 1; ++i) {
coeffs_->GetCoefficient(k, nPowerTheta - j - 2, nPowerFreq - i - 2, Pk);
coefficients_->GetCoefficient(k, nPowerTheta - j - 2,
nPowerFreq - i - 2, Pk);
Pj.first = Pj.first * freq + Pk.first;
Pj.second = Pj.second * freq + Pk.second;
}
......@@ -115,13 +111,9 @@ aocommon::MC2x2 HamakerElementResponse::Response(double freq, double theta,
return response;
}
HamakerElementResponseHBA::HamakerElementResponseHBA() {
std::string path = GetPath("HamakerHBACoeff.h5");
coeffs_ = std::make_unique<HamakerCoefficients>(path);
}
std::weak_ptr<const HamakerCoefficients>
HamakerElementResponse::cached_lba_coefficients_;
std::weak_ptr<const HamakerCoefficients>
HamakerElementResponse::cached_hba_coefficients_;
HamakerElementResponseLBA::HamakerElementResponseLBA() {
std::string path = GetPath("HamakerLBACoeff.h5");
coeffs_ = std::make_unique<HamakerCoefficients>(path);
}
} // namespace everybeam
......@@ -14,6 +14,8 @@ namespace everybeam {
//! Implementation of the Hamaker response model
class HamakerElementResponse : public ElementResponse {
public:
explicit HamakerElementResponse(const std::string& name);
ElementResponseModel GetModel() const final override {
return ElementResponseModel::kHamaker;
}
......@@ -26,29 +28,11 @@ class HamakerElementResponse : public ElementResponse {
*/
static std::shared_ptr<HamakerElementResponse> GetLbaInstance();
/**
* @brief Get instance of HamakerElementResponse, infer type (LBA/HBA)
* from station name.
*
* @param name Station Name
*/
static std::shared_ptr<HamakerElementResponse> GetInstance(
const std::string& name);
protected:
std::string GetPath(const char*) const;
std::unique_ptr<HamakerCoefficients> coeffs_;
};
class HamakerElementResponseHBA : public HamakerElementResponse {
public:
HamakerElementResponseHBA();
};
class HamakerElementResponseLBA : public HamakerElementResponse {
public:
HamakerElementResponseLBA();
private:
// Since coefficients are equal for all LBA and HBA instances, share them.
static std::weak_ptr<const HamakerCoefficients> cached_lba_coefficients_;
static std::weak_ptr<const HamakerCoefficients> cached_hba_coefficients_;
std::shared_ptr<const HamakerCoefficients> coefficients_;
};
} // namespace everybeam
......
......@@ -71,27 +71,6 @@ static std::optional<AartfaacStation> GetAartfaacStation(
namespace everybeam {
namespace {
/**
* @brief Search for LOBES h5 coefficient file
* on the suggested path \param search_path. Returns
* an empty string if the file cannot be found.
*
* @param search_path Search path
* @param station_name Station name, as read from MS
* @return std::string Path to file or empty string if file cannot be found
*/
std::filesystem::path FindCoeffFile(const std::string& search_path,
std::string_view station_name) {
const std::string station_file = "LOBES_" + std::string{station_name} + ".h5";
return search_path.empty()
? std::filesystem::path(std::string{EVERYBEAM_DATA_DIR} +
std::string{"/lobes"}) /
station_file
: std::filesystem::path(search_path) / station_file;
}
} // namespace
static const H5::CompType kH5Dcomplex = [] {
const std::string REAL("r");
const std::string IMAG("i");
......@@ -163,8 +142,13 @@ LOBESElementResponse::LOBESElementResponse(const std::string& name,
const std::optional<AartfaacStation> aartfaac_station =
GetAartfaacStation(name, AartfaacElements::kInner);
std::filesystem::path coeff_file_path = FindCoeffFile(
options.coeff_path, aartfaac_station ? aartfaac_station->station : name);
const std::filesystem::path search_path =
options.coeff_path.empty() ? GetPath("lobes")
: std::filesystem::path(options.coeff_path);
const std::string_view station_name =
aartfaac_station ? aartfaac_station->station : name;
const std::string station_file = "LOBES_" + std::string(station_name) + ".h5";
const std::filesystem::path coeff_file_path = search_path / station_file;
H5::H5File h5file;
if (!std::filesystem::exists(coeff_file_path)) {
......@@ -278,17 +262,28 @@ aocommon::MC2x2 LOBESElementResponse::Response(
return response;
}
std::shared_ptr<LOBESElementResponse> LOBESElementResponse::GetInstance(
std::shared_ptr<const LOBESElementResponse> LOBESElementResponse::GetInstance(
const std::string& name, const Options& options) {
static std::map<std::string, std::shared_ptr<LOBESElementResponse>>
// Using a single LOBESElementResponse object for each name reduces memory
// usage since the coefficients are only loaded once.
// Using weak pointers in this map ensures that LOBESElementResponse objects,
// are deleted when they are no longer used, which saves memory.
static std::map<std::string, std::weak_ptr<const LOBESElementResponse>>
name_response_map;
std::shared_ptr<const LOBESElementResponse> instance;
auto entry = name_response_map.find(name);
if (entry == name_response_map.end()) {
entry = name_response_map.insert(
entry, {name, std::make_shared<LOBESElementResponse>(name, options)});
instance = std::make_shared<const LOBESElementResponse>(name, options);
name_response_map.insert({name, instance});
} else {
instance = entry->second.lock();
if (!instance) {
instance = std::make_shared<const LOBESElementResponse>(name, options);
entry->second = instance;
}
}
return entry->second;
return instance;
}
} // namespace everybeam
......@@ -68,9 +68,8 @@ class LOBESElementResponse : public ElementResponse {
* @brief Create LOBESElementResponse
*
* @param name Station name, e.g. CS302LBA
* @return std::shared_ptr<LOBESElementResponse>
*/
static std::shared_ptr<LOBESElementResponse> GetInstance(
static std::shared_ptr<const LOBESElementResponse> GetInstance(
const std::string& name, const Options& options);
/**
......
......@@ -28,15 +28,17 @@ aocommon::MC2x2 OSKARElementResponseDipole::Response(double freq, double theta,
return response;
}
OSKARElementResponseSphericalWave::OSKARElementResponseSphericalWave() {
std::string path = GetPath("oskar.h5");
datafile_.reset(new Datafile(path));
OSKARElementResponseSphericalWave::OSKARElementResponseSphericalWave()
: datafile_(cached_datafile_.lock()) {
if (!datafile_) {
datafile_ = std::make_shared<Datafile>(GetPath("oskar.h5"));
cached_datafile_ = datafile_;
}
}
OSKARElementResponseSphericalWave::OSKARElementResponseSphericalWave(
const std::string& path) {
datafile_.reset(new Datafile(path));
}
const std::string& filename)
: datafile_(std::make_shared<Datafile>(filename)) {}
aocommon::MC2x2 OSKARElementResponseSphericalWave::Response(
[[maybe_unused]] double freq, [[maybe_unused]] double theta,
......@@ -76,12 +78,6 @@ aocommon::MC2x2 OSKARElementResponseSphericalWave::Response(int element_id,
return response;
}
std::string OSKARElementResponseSphericalWave::GetPath(
const char* filename) const {
std::stringstream ss;
ss << EVERYBEAM_DATA_DIR << "/";
ss << filename;
return ss.str();
}
std::weak_ptr<Datafile> OSKARElementResponseSphericalWave::cached_datafile_;
} // namespace everybeam
......@@ -5,7 +5,6 @@
#define OSKAR_ELEMENTRESPONSE_H
#include "../elementresponse.h"
#include "../common/singleton.h"
#include "oskardatafile.h"
......@@ -16,10 +15,6 @@ namespace everybeam {
//! Implementation of the OSKAR dipole response model
class OSKARElementResponseDipole : public ElementResponse {
public:
static std::shared_ptr<OSKARElementResponseDipole> GetInstance() {
return common::Singleton<OSKARElementResponseDipole>::GetInstance();
}
ElementResponseModel GetModel() const final override {
return ElementResponseModel::kOSKARDipole;
}
......@@ -31,24 +26,14 @@ class OSKARElementResponseDipole : public ElementResponse {
//! Implementation of the OSKAR spherical wave response model
class OSKARElementResponseSphericalWave : public ElementResponse {
public:
/**
* A constructor-like static method to instantiate the class
*
* returns a globally shared instance of the class that is instantiated
* in the first call
*/
static std::shared_ptr<OSKARElementResponseSphericalWave> GetInstance() {
return common::Singleton<OSKARElementResponseSphericalWave>::GetInstance();
}
/** Constructor loading the default coefficients file */
OSKARElementResponseSphericalWave();
/** Constructor loading a custom coefficients file
*
* @param path Path to the coefficients file to load
/**
* Constructor loading a custom coefficients file
* @param filename Filename of HDF5 file with coefficients.
*/
OSKARElementResponseSphericalWave(const std::string& path);
OSKARElementResponseSphericalWave(const std::string& filename);
ElementResponseModel GetModel() const final override {
return ElementResponseModel::kOSKARSphericalWave;
......@@ -60,10 +45,10 @@ class OSKARElementResponseSphericalWave : public ElementResponse {
aocommon::MC2x2 Response(int element_id, double freq, double theta,
double phi) const final override;
protected:
std::string GetPath(const char*) const;
std::unique_ptr<Datafile> datafile_;
private:
// This weak pointer allows reusing the default coefficients.
static std::weak_ptr<Datafile> cached_datafile_;
std::shared_ptr<Datafile> datafile_;
};
} // namespace everybeam
......
......@@ -90,7 +90,7 @@ for em_idx in range(2):
)
subprocess.check_call(
["oskar_csv_to_hdf5.py", "telescope.tm", "oskar.h5"]
["oskar_csv_to_hdf5.py", "telescope.tm", "oskar-comparison.h5"]
)
subprocess.check_call(["make_element_response_image", str(npixels)])
......
......@@ -10,7 +10,8 @@
#include "../../external/npy.hpp" // to save arrays in numpy format
int main(int argc, char** argv) {
everybeam::OSKARElementResponseSphericalWave element_response("oskar.h5");
everybeam::OSKARElementResponseSphericalWave element_response(
"oskar-comparison.h5");
double freq = 50e6;
int N;
......
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