MUQ  0.4.3
muq::Utilities Namespace Reference

Namespaces

 PythonBindings
 
 StringUtilities
 

Classes

class  AnyCast
 Class for easily casting boost::any's in assignment operations. More...
 
class  AnyConstCast
 
struct  VectorLessThan
 
struct  AnyWriter
 
class  Attribute
 
class  AttributeList
 
class  BlockDataset
 
class  H5Object
 
class  HDF5File
 A wrapper around HDF5 (and PHDF5) to handle file I/O. More...
 
struct  MPI_Type
 
struct  MPI_Type< double >
 
struct  MPI_Type< char >
 
struct  MPI_Type< short >
 
struct  MPI_Type< int >
 
struct  MPI_Type< long >
 
struct  MPI_Type< float >
 
struct  MPI_Type< unsigned short >
 
struct  MPI_Type< unsigned int >
 
struct  MPI_Type< unsigned long >
 
struct  HDF5_Type
 
struct  HDF5_Type< double >
 
struct  HDF5_Type< long double >
 
struct  HDF5_Type< int >
 
struct  HDF5_Type< long >
 
struct  HDF5_Type< unsigned long >
 
struct  HDF5_Type< unsigned int >
 
struct  HDF5_Type< float >
 
struct  HDF5_Type< unsigned short >
 
struct  HDF5_Type< short >
 
struct  HDF5_Type< char >
 
struct  HDF5_Type< bool >
 
class  LinearOperatorTypeException
 
struct  LinearOperatorFactory
 
class  LinearOperator
 Generic linear operator base class. More...
 
class  MultiIndexFactory
 A factory class with static methods for generating MultiIndexSets. More...
 
class  MultiIndexLimiter
 An abstract base class for multi index limiters. More...
 
class  TotalOrderLimiter
 Provides a cap on the total-order allowed. More...
 
class  DimensionLimiter
 Provides bounds on what dimensions are allowed to have nonzero values. More...
 
class  AnisotropicLimiter
 Declares multiindices as feasible if their entries for less important dimensions are not too high. More...
 
class  MaxOrderLimiter
 Provides a cap on the maximum value of each component the multiindex. More...
 
class  NoLimiter
 Returns true for an multiindex. More...
 
class  AndLimiter
 Combines two limiters through an AND operation. More...
 
class  OrLimiter
 Combines two limiters through an OR operation. More...
 
class  XorLimiter
 Combines two limiters through an XOR operation. More...
 
class  MultiIndexSet
 A class for holding, sorting, and adapting sets of multiindices. More...
 
class  OTF2TracerBase
 Base interface for OTF2 tracer implemetations. More...
 
class  OTF2TracerDummy
 Fallback dummy implementation not doing anything; Does not require libotf2. More...
 
class  OTF2Tracer
 Tracer implementation writing to OTF2 via libotf2 The result can be viewed by several programs, one example is "vite" developed at INRIA. More...
 
class  RandomGeneratorTemporarySetSeed
 
class  RandomGenerator
 
class  SeedGenerator
 
struct  shared_factory
 
class  VectorSlice
 Enables a subset of a vector to be easily accessed or reversed without copying memory. More...
 
class  WaitBar
 

Typedefs

using OTF2Tracer = OTF2TracerDummy
 

Enumerations

enum  TracerRegions {
  Setup , Finalize , PhonebookBusy , BurnIn ,
  Sampling , CollectorBusy , RetrievingProposal , FetchingProposal
}
 

Functions

std::string demangle (const char *name)
 
template<typename PointerType >
std::string GetTypeName (PointerType const &ptr)
 
H5Object AddChildren (std::shared_ptr< HDF5File > file, std::string const &groupName)
 Recursively add children to create an HDF5 file hierarchy. More...
 
H5Object OpenFile (std::string const &filename)
 Open an HDF5 file and return the root object. More...
 
herr_t H5LTset_attribute (hid_t loc_id, const char *obj_name, const char *attr_name, hid_t mem_type_id, void *buffer, size_t size)
 
std::string H5TypeToString (hid_t type)
 
std::string GetParentPath (std::string const &base)
 
std::pair< std::string, std::string > SplitString (std::string const &path)
 Splits a string on the first forward slash. More...
 
std::shared_ptr< MultiIndexSetoperator+= (std::shared_ptr< MultiIndexSet > x, std::shared_ptr< MultiIndexSet > y)
 
std::shared_ptr< MultiIndexSetoperator+= (std::shared_ptr< MultiIndexSet > x, std::shared_ptr< MultiIndex > y)
 
static OTF2_TimeStamp get_time (void)
 
static OTF2_FlushType pre_flush (void *userData, OTF2_FileType fileType, OTF2_LocationRef location, void *callerData, bool final)
 
static OTF2_TimeStamp post_flush (void *userData, OTF2_FileType fileType, OTF2_LocationRef location)
 
void AddDictToPtree (pybind11::dict dict, std::string basePath, boost::property_tree::ptree &pt)
 
boost::property_tree::ptree ConvertDictToPtree (pybind11::dict dict)
 
template<typename ScalarType >
VectorSlice< std::vector< ScalarType >, ScalarType > GetSlice (std::vector< ScalarType > &dataIn, int startIndIn, int endIndIn, int skipIn=1)
 
template<typename ScalarType >
VectorSlice< Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 >, ScalarType > GetSlice (Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 > &dataIn, int startIndIn, int endIndIn, int skipIn=1)
 
template<typename VectorType , typename ScalarType >
VectorSlice< VectorType, ScalarType > GetSlice (VectorSlice< VectorType, ScalarType > &dataIn, int startIndIn, int endIndIn, int skipIn=1)
 

Variables

static OTF2_FlushCallbacks flush_callbacks
 

Typedef Documentation

◆ OTF2Tracer

Definition at line 287 of file OTF2Tracer.h.

Enumeration Type Documentation

◆ TracerRegions

Enumerator
Setup 
Finalize 
PhonebookBusy 
BurnIn 
Sampling 
CollectorBusy 
RetrievingProposal 
FetchingProposal 

Definition at line 12 of file OTF2Tracer.h.

Function Documentation

◆ AddChildren()

H5Object muq::Utilities::AddChildren ( std::shared_ptr< HDF5File file,
std::string const &  groupName 
)

Recursively add children to create an HDF5 file hierarchy.

Definition at line 290 of file H5Object.cpp.

References muq::Utilities::H5Object::children.

Referenced by OpenFile().

◆ AddDictToPtree()

void muq::Utilities::AddDictToPtree ( pybind11::dict  dict,
std::string  basePath,
boost::property_tree::ptree &  pt 
)
inline

Definition at line 16 of file PyDictConversion.h.

Referenced by ConvertDictToPtree().

◆ ConvertDictToPtree()

boost::property_tree::ptree muq::Utilities::ConvertDictToPtree ( pybind11::dict  dict)
inline

This function is useful for definining python wrappers of functions with ptree arguments. It can be used in the python bindings to convert a python dictionary object into a ptree before calling the c++ function. For example, it may be used like:

struct PtreePrinter{
void PrintPtree(boost::property_tree::ptree pt){
boost::property_tree::write_json(std::cout, pt, true);
}
};
// ... other stuff and module initialization
py::class_<PtreePrinter, std::shared_ptr<PtreePrinter>> pp(m, "PtreePrinter");
pp
.def(py::init<>())
.def("PrintPtree", [](PtreePrinter* p, py::dict d){p->PrintPtree(muq::Utilities::ConvertDictToPtree(d));});
boost::property_tree::ptree ConvertDictToPtree(pybind11::dict dict)

Definition at line 67 of file PyDictConversion.h.

References AddDictToPtree().

◆ demangle()

◆ get_time()

static OTF2_TimeStamp muq::Utilities::get_time ( void  )
static

◆ GetParentPath()

std::string muq::Utilities::GetParentPath ( std::string const &  base)
   @brief Get the parent folder of a path
   @details Given a path of the form /a/b/c, return the parent path /a/b
std::string base = "/some/path";
std::string parent = GetParent(base);

Parent now holds "/some"

Definition at line 6 of file PathTools.cpp.

Referenced by muq::Utilities::HDF5File::CreateDataset(), muq::Utilities::HDF5File::CreateGroup(), muq::Utilities::HDF5File::DoesDataSetExist(), muq::Utilities::HDF5File::DoesGroupExist(), and muq::Utilities::HDF5File::WriteMatrix().

◆ GetSlice() [1/3]

template<typename ScalarType >
VectorSlice<Eigen::Matrix<ScalarType,Eigen::Dynamic,1>, ScalarType> muq::Utilities::GetSlice ( Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 > &  dataIn,
int  startIndIn,
int  endIndIn,
int  skipIn = 1 
)

Definition at line 118 of file VectorSlice.h.

◆ GetSlice() [2/3]

template<typename ScalarType >
VectorSlice<std::vector<ScalarType>, ScalarType> muq::Utilities::GetSlice ( std::vector< ScalarType > &  dataIn,
int  startIndIn,
int  endIndIn,
int  skipIn = 1 
)

Definition at line 109 of file VectorSlice.h.

Referenced by muq::SamplingAlgorithms::DRKernel::Alpha().

◆ GetSlice() [3/3]

template<typename VectorType , typename ScalarType >
VectorSlice<VectorType, ScalarType> muq::Utilities::GetSlice ( VectorSlice< VectorType, ScalarType > &  dataIn,
int  startIndIn,
int  endIndIn,
int  skipIn = 1 
)

◆ GetTypeName()

template<typename PointerType >
std::string muq::Utilities::GetTypeName ( PointerType const &  ptr)

◆ H5LTset_attribute()

herr_t muq::Utilities::H5LTset_attribute ( hid_t  loc_id,
const char *  obj_name,
const char *  attr_name,
hid_t  mem_type_id,
void *  buffer,
size_t  size 
)

Definition at line 19 of file HDF5Types.cpp.

Referenced by muq::Utilities::HDF5File::WriteVectorAttribute().

◆ H5TypeToString()

std::string muq::Utilities::H5TypeToString ( hid_t  type)

Definition at line 45 of file HDF5Types.cpp.

◆ OpenFile()

H5Object muq::Utilities::OpenFile ( std::string const &  filename)

Open an HDF5 file and return the root object.

Primary way to open an HDF5file. It returns the root group.

// Open the HDF5 file for reading and writing
auto f = muq::Utilities::OpenFile("Data.h5");
// Create a new group
f.CreateGroup("/NewGroup");
// Add an attribute to the group
f["/NewGroup"].attrs["Some Metadata"] = "Created with MH5!";
// Store a vector in the group
f["/NewGroup/Ones"] = Eigen::VectorXd::Ones(10);
// Store another vector in a different group. The group is automatically generated
f["/AnotherGroup/Zeros"] = Eigen::MatrixXd::Zero(5,5);
H5Object OpenFile(std::string const &filename)
Open an HDF5 file and return the root object.
Definition: H5Object.cpp:321

Definition at line 321 of file H5Object.cpp.

References AddChildren().

Referenced by muq::Utilities::MultiIndexSet::FromHDF5(), muq::Approximation::BasisExpansion::FromHDF5(), main(), ReadData(), muq::Utilities::MultiIndexSet::ToHDF5(), muq::Approximation::BasisExpansion::ToHDF5(), and WriteResults().

◆ operator+=() [1/2]

std::shared_ptr< MultiIndexSet > muq::Utilities::operator+= ( std::shared_ptr< MultiIndexSet x,
std::shared_ptr< MultiIndex >  y 
)

Definition at line 16 of file MultiIndexSet.cpp.

◆ operator+=() [2/2]

std::shared_ptr< MultiIndexSet > muq::Utilities::operator+= ( std::shared_ptr< MultiIndexSet x,
std::shared_ptr< MultiIndexSet y 
)

◆ post_flush()

static OTF2_TimeStamp muq::Utilities::post_flush ( void *  userData,
OTF2_FileType  fileType,
OTF2_LocationRef  location 
)
static

Definition at line 75 of file OTF2Tracer.h.

References get_time().

◆ pre_flush()

static OTF2_FlushType muq::Utilities::pre_flush ( void *  userData,
OTF2_FileType  fileType,
OTF2_LocationRef  location,
void *  callerData,
bool  final 
)
static

Definition at line 66 of file OTF2Tracer.h.

◆ SplitString()

std::pair< std::string, std::string > muq::Utilities::SplitString ( std::string const &  path)

Splits a string on the first forward slash.

Given a path, like "/a/b/c", this function returns two strings, "/a", and "/b/c". @params[in] path The original plan to split.

Returns
A pair of std::string holind the two components of the path.

Definition at line 19 of file PathTools.cpp.

Referenced by muq::Utilities::H5Object::CreateGroup(), muq::Utilities::H5Object::CreatePlaceholder(), and muq::Utilities::H5Object::operator[]().

Variable Documentation

◆ flush_callbacks

OTF2_FlushCallbacks muq::Utilities::flush_callbacks
static
Initial value:
=
{
.otf2_pre_flush = pre_flush,
.otf2_post_flush = post_flush
}
static OTF2_TimeStamp post_flush(void *userData, OTF2_FileType fileType, OTF2_LocationRef location)
Definition: OTF2Tracer.h:75
static OTF2_FlushType pre_flush(void *userData, OTF2_FileType fileType, OTF2_LocationRef location, void *callerData, bool final)
Definition: OTF2Tracer.h:66

Definition at line 81 of file OTF2Tracer.h.

Referenced by muq::Utilities::OTF2Tracer::OTF2Tracer().