C++ API Reference¶
The C++ API reference is generated from Doxygen comments in the source code using Breathe.
Namespace¶
All OpenImpala classes live in the OpenImpala namespace.
Key classes¶
I/O Readers¶
-
class TiffReader¶
Public Functions
-
int bitsPerSample() const¶
-
amrex::Box box() const¶
-
int depth() const¶
-
inline uint16_t getFillOrder() const¶
-
int height() const¶
-
inline bool isRead() const¶
-
bool readFile(const std::string &filename)¶
-
bool readFileSequence(const std::string &base_pattern, int num_files, int start_index = 0, int digits = 1, const std::string &suffix = ".tif")¶
-
int sampleFormat() const¶
-
int samplesPerPixel() const¶
-
void threshold(double raw_threshold, amrex::iMultiFab &mf) const¶
-
void threshold(double raw_threshold, int value_if_true, int value_if_false, amrex::iMultiFab &mf) const¶
-
TiffReader()¶
-
TiffReader(const std::string &base_pattern, int num_files, int start_index = 0, int digits = 1, const std::string &suffix = ".tif")¶
-
explicit TiffReader(const std::string &filename)¶
-
int width() const¶
-
virtual ~TiffReader() = default¶
-
int bitsPerSample() const¶
-
class HDF5Reader¶
Public Functions
-
amrex::Box box() const¶
-
int depth() const¶
-
std::map<std::string, std::string> getAllAttributes() const¶
-
std::string getAttribute(const std::string &attr_name) const¶
-
inline H5::DataType getNativeDataType() const¶
-
HDF5Reader()¶
-
explicit HDF5Reader(const std::string &filename, const std::string &hdf5dataset)¶
-
int height() const¶
-
inline bool isRead() const¶
-
bool readFile(const std::string &filename, const std::string &hdf5dataset)¶
-
void threshold(double raw_threshold, amrex::iMultiFab &mf) const¶
-
void threshold(double raw_threshold, int value_if_true, int value_if_false, amrex::iMultiFab &mf) const¶
-
int width() const¶
-
virtual ~HDF5Reader() = default¶
-
amrex::Box box() const¶
-
class RawReader¶
Public Types
-
using ByteType = unsigned char¶
Public Functions
-
amrex::Box box() const¶
-
int depth() const¶
-
RawDataType getDataType() const¶
-
double getValue(int i, int j, int k) const¶
-
int height() const¶
-
bool isRead() const¶
-
RawReader()¶
-
explicit RawReader(const std::string &filename, int width, int height, int depth, RawDataType data_type)¶
-
bool readFile(const std::string &filename, int width, int height, int depth, RawDataType data_type)¶
-
void threshold(double threshold_value, amrex::iMultiFab &mf) const¶
-
void threshold(double threshold_value, int value_if_true, int value_if_false, amrex::iMultiFab &mf) const¶
-
int width() const¶
-
virtual ~RawReader() = default¶
-
using ByteType = unsigned char¶
Transport Solvers¶
-
class TortuosityHypre : public OpenImpala::Tortuosity, public OpenImpala::HypreStructSolver¶
Public Types
-
using SolverType = OpenImpala::SolverType¶
Public Functions
-
bool checkMatrixProperties()¶
-
inline BCType getInletOutletBCType() const¶
-
inline const std::map<int, amrex::Real> &getPhaseCoeffMap() const¶
-
inline const std::vector<amrex::Real> &getPlaneFluxes() const¶
-
inline amrex::Real getPlaneFluxMaxDeviation() const¶
-
inline BCType getSidesBCType() const¶
-
inline bool isMultiPhase() const¶
-
TortuosityHypre(const amrex::Geometry &geom, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::iMultiFab &mf_phase_input, const amrex::Real vf, const int phase, const OpenImpala::Direction dir, const SolverType solvertype, const std::string &resultspath, const amrex::Real vlo = 0.0, const amrex::Real vhi = 1.0, int verbose = 0, bool write_plotfile = false, const PrecondType precond_type = PrecondType::SMG)¶
-
virtual amrex::Real value(const bool refresh = false) override¶
-
virtual ~TortuosityHypre() override = default¶
-
using SolverType = OpenImpala::SolverType¶
-
class TortuosityMLMG : public OpenImpala::TortuositySolverBase¶
Public Functions
-
virtual bool solve() override¶
-
TortuosityMLMG(const amrex::Geometry &geom, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::iMultiFab &mf_phase_input, const amrex::Real vf, const int phase, const OpenImpala::Direction dir, const std::string &resultspath, const amrex::Real vlo = 0.0, const amrex::Real vhi = 1.0, int verbose = 0, bool write_plotfile = false, amrex::Real eps = 1.0e-9, int maxiter = 200, int max_coarsening_level = 30)¶
-
virtual bool solve() override¶
-
class EffectiveDiffusivityHypre : public OpenImpala::HypreStructSolver¶
Public Types
-
using SolverType = OpenImpala::SolverType¶
Public Functions
-
EffectiveDiffusivityHypre(const amrex::Geometry &geom, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::iMultiFab &mf_phase_input, const int phase_id, const OpenImpala::Direction dir_of_chi_k, const SolverType solver_type, const std::string &resultspath, int verbose_level = 1, bool write_plotfile_flag = false)¶
-
void getChiSolution(amrex::MultiFab &chi_field)¶
-
bool solve()¶
-
~EffectiveDiffusivityHypre() override = default¶
-
using SolverType = OpenImpala::SolverType¶
Utilities¶
-
class VolumeFraction¶
Public Functions
-
VolumeFraction &operator=(const VolumeFraction&) = delete¶
-
VolumeFraction &operator=(VolumeFraction&&) = delete¶
-
void value(long long &phase_count, long long &total_count, bool local = false) const¶
-
explicit VolumeFraction(const amrex::iMultiFab &fm, const int phase = 0, int comp = 0)¶
-
VolumeFraction(const VolumeFraction&) = delete¶
-
VolumeFraction(VolumeFraction&&) = delete¶
-
virtual ~VolumeFraction() = default¶
-
VolumeFraction &operator=(const VolumeFraction&) = delete¶
-
class PercolationCheck¶
Public Functions
-
inline amrex::Real activeVolumeFraction() const¶
-
inline bool percolates() const¶
-
PercolationCheck(const amrex::Geometry &geom, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::iMultiFab &mf_phase, int phase_id, OpenImpala::Direction dir, int verbose = 0)¶
Public Static Functions
-
static std::string directionString(OpenImpala::Direction dir)¶
-
inline amrex::Real activeVolumeFraction() const¶
-
class TortuositySolverBase : public OpenImpala::Tortuosity¶
Subclassed by OpenImpala::TortuosityMLMG
Public Functions
-
virtual bool solve() = 0¶
-
TortuositySolverBase(const amrex::Geometry &geom, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, const amrex::iMultiFab &mf_phase_input, const amrex::Real vf, const int phase, const OpenImpala::Direction dir, const std::string &resultspath, const amrex::Real vlo, const amrex::Real vhi, int verbose, bool write_plotfile)¶
-
virtual amrex::Real value(const bool refresh = false) override¶
-
virtual bool solve() = 0¶
Configuration¶
-
struct PhysicsConfig¶
Full Doxygen output¶
For the complete class hierarchy, include dependency graphs, and file-level documentation, see the Doxygen pages.