Dakota  Version
Explore and Predict with Confidence
Classes | Typedefs | Enumerations | Functions | Variables
dakota::util Namespace Reference

namespace for new Dakota utilities module More...

Classes

class  CholeskySolver
 The CholeskySolver class is used to solve linear systems with a symmetric matrix with a pivoted Cholesky decomposition. More...
 
class  DataScaler
 The DataScaler class computes the scaling coefficients and scales a 2D data matrix with dimensions num_samples by num_features. More...
 
class  LinearSolverBase
 The LinearSolverBase class serves as an API for derived solvers. More...
 
class  LUSolver
 The LUSolver class is used to solve linear systems with the LU decomposition. More...
 
class  NormalizationScaler
 Normalizes the data using max and min feature values. More...
 
class  NoScaler
 Leaves the data unscaled. More...
 
class  QRSolver
 The QRSolver class solves the linear least squares problem with a QR decomposition. More...
 
class  StandardizationScaler
 Standardizes the data so the each feature has zero mean and unit variance. More...
 
class  SVDSolver
 The SVDSolver class is used to solve linear systems with the singular value decomposition. More...
 

Typedefs

using BimapMetrictypeStr = boost::bimap< METRIC_TYPE, std::string >
 alias for Boost Bimap metric type <--> string
 
using SCALER_TYPE = DataScaler::SCALER_TYPE
 alias for DataScaler's SCALER_TYPE
 
using BimapScalertypeStr = boost::bimap< SCALER_TYPE, std::string >
 alias for Boost Bimap scaler type <--> string
 
using SOLVER_TYPE = LinearSolverBase::SOLVER_TYPE
 alias for LinearSolverBase's SOLVER_TYPE
 
using BimapSolvertypeStr = boost::bimap< SOLVER_TYPE, std::string >
 alias for Boost Bimap solver type <--> string
 

Enumerations

enum  METRIC_TYPE {
  SUM_SQUARED, MEAN_SQUARED, ROOT_MEAN_SQUARED, SUM_ABS,
  MEAN_ABS, MAX_ABS, ABS_PERCENTAGE_ERROR, MEAN_ABS_PERCENTAGE_ERROR,
  R_SQUARED
}
 Enumeration for supported metric types.
 

Functions

void error (const std::string &msg)
 Throws a std::runtime_error based on the message argument. More...
 
bool matrix_equals (const MatrixXi &A, const MatrixXi &B)
 Tests whether two Eigen MatrixXi objects are equal. More...
 
bool matrix_equals (const MatrixXd &A, const MatrixXd &B, double tol)
 Tests whether two Eigen MatrixXd objects are equal, within a given tolerance. More...
 
bool matrix_equals (const RealMatrix &A, const RealMatrix &B, double tol)
 Tests whether two Teuchos RealMatrix objects are equal, within a given tolerance. More...
 
bool relative_allclose (const MatrixXd &A, const MatrixXd &B, const double tol)
 Tests whether two Eigen MatrixXd objects relatively equal (element-wise) within a given tolerance. More...
 
double variance (const VectorXd &vec)
 Calculates the variance based on an Eigen VectorXd of double values. More...
 
void populateVectorsFromFile (const std::string &filename, std::vector< VectorXd > &R, int num_datasets, int num_samples)
 Populate a collection of vectors read in a from a text file assuming data layout is one dataset per row. More...
 
void populateMatricesFromFile (const std::string &filename, std::vector< MatrixXd > &S, int num_datasets, int num_vars, int num_samples)
 Populate a collection of matrices read in a from a text file assuming data layout is a "column-major" stack of num_samples by num_vars matrices. More...
 
int n_choose_k (int n, int k)
 Calculate Binomial coefficient n choose k. More...
 
void random_permutation (const int num_pts, const unsigned int seed, VectorXi &permutations)
 Random permutation of int array.
 
void create_cv_folds (const int num_folds, const int num_pts, std::vector< VectorXi > &fold_indices, const int seed=22)
 Generate indices for cross validation folds.
 
MatrixXd create_uniform_random_double_matrix (const int rows, const int cols)
 Generate a real-valued matrix of uniformly distributed random values.
 
MatrixXd create_uniform_random_double_matrix (const int rows, const int cols, const unsigned int seed)
 Generate a real-valued matrix of uniformly distributed random values.
 
MatrixXd create_uniform_random_double_matrix (const int rows, const int cols, const unsigned int seed, bool transform, const double low, const double high)
 Generate a real-valued matrix of uniformly distributed random values.
 
template<typename T >
int num_nonzeros (const T &mat)
 Caclulate and return number of nonzero entries in vector or matrix. More...
 
template<typename T1 , typename T2 >
void nonzero (const T1 &v, T2 &result)
 Create a vector of indices based on nonzero entries in input vector. More...
 
template<typename T1 , typename T2 >
void append_columns (const T1 &new_cols, T2 &target)
 Append columns of input matrix to existing matrix. More...
 
template<typename T >
double p_norm (const T &v, double p)
 Caclulate and return p-norm of a vector. More...
 
METRIC_TYPE metric_type (const std::string &metric_name)
 Convert the metric from string to enum. More...
 
double compute_metric (const VectorXd &p, const VectorXd &d, const std::string &metric_name)
 Computes the difference between prediction and data vectors. More...
 
std::shared_ptr< DataScalerscaler_factory (DataScaler::SCALER_TYPE scaler_type, const MatrixXd &unscaled_matrix)
 Free function to construct DataScaler. More...
 
std::shared_ptr< LinearSolverBasesolver_factory (LinearSolverBase::SOLVER_TYPE type)
 Free function to construct LinearSolverBase. More...
 

Variables

static BimapMetrictypeStr type_name_bimap
 Bimap between metric types and names. More...
 
static BimapScalertypeStr type_name_bimap
 Bimap between scaler types and names. More...
 
static BimapSolvertypeStr type_name_bimap
 Bimap between solver types and names. More...
 

Detailed Description

namespace for new Dakota utilities module

Function Documentation

◆ error()

void error ( const std::string &  msg)

Throws a std::runtime_error based on the message argument.

Parameters
[in]msgThe error message to throw

Referenced by compute_metric(), create_cv_folds(), and matrix_equals().

◆ matrix_equals() [1/3]

bool matrix_equals ( const MatrixXi A,
const MatrixXi B 
)

Tests whether two Eigen MatrixXi objects are equal.

Parameters
[in]AThe first matrix to test
[in]BThe second matrix to test
Returns
Whether the matrices are equal

◆ matrix_equals() [2/3]

bool matrix_equals ( const MatrixXd A,
const MatrixXd B,
double  tol 
)

Tests whether two Eigen MatrixXd objects are equal, within a given tolerance.

Parameters
[in]AThe first matrix to test
[in]BThe second matrix to test
[in]tolThe tolerance to use when comparing double values
Returns
Whether the matrices are equal to within tolerance

References error().

◆ matrix_equals() [3/3]

bool matrix_equals ( const RealMatrix A,
const RealMatrix B,
double  tol 
)

Tests whether two Teuchos RealMatrix objects are equal, within a given tolerance.

Parameters
[in]AThe first matrix to test
[in]BThe second matrix to test
[in]tolThe tolerance to use when comparing double values
Returns
Whether the matrices are equal to within tolerance

References error().

◆ relative_allclose()

bool relative_allclose ( const MatrixXd A,
const MatrixXd B,
const double  tol 
)

Tests whether two Eigen MatrixXd objects relatively equal (element-wise) within a given tolerance.

Parameters
[in]AThe first matrix to test
[in]BThe second matrix to test
[in]tolThe relative tolerance to use when comparing double values
Returns
Whether the matrices are relatively equal (within the tolerance)

◆ variance()

double variance ( const VectorXd vec)

Calculates the variance based on an Eigen VectorXd of double values.

Parameters
[in]vecThe vector
Returns
The calculated variance

◆ populateVectorsFromFile()

void populateVectorsFromFile ( const std::string &  filename,
std::vector< VectorXd > &  R,
int  num_datasets,
int  num_samples 
)

Populate a collection of vectors read in a from a text file assuming data layout is one dataset per row.

Parameters
[in]filenameThe file that contains the data
[out]RThe collection of vectors to be populated
[in]num_datasetsThe number of datasets to read in
[in]num_samplesThe number of data points (e.g. function values, build points) per dataset

◆ populateMatricesFromFile()

void populateMatricesFromFile ( const std::string &  filename,
std::vector< MatrixXd > &  S,
int  num_datasets,
int  num_vars,
int  num_samples 
)

Populate a collection of matrices read in a from a text file assuming data layout is a "column-major" stack of num_samples by num_vars matrices.

Parameters
[in]filenameThe file that contains the data
[out]SThe collection of vectors to be populated
[in]num_datasetsThe number of datasets to read in
[in]num_samplesThe number of data points (e.g. function values, build points) per dataset (row dim)
[in]num_varsThe number of variables per dataset (column dim)

◆ n_choose_k()

int n_choose_k ( int  n,
int  k 
)

Calculate Binomial coefficient n choose k.

Parameters
[in]nNumber of elements in set
[in]kNumber of elements in subset k where n >= k >= 0
Returns
Number of ways to choose an (unordered) subset of k elements from a fixed set of n elements

Referenced by dakota::surrogates::size_level_index_vector().

◆ num_nonzeros()

int dakota::util::num_nonzeros ( const T &  mat)

Caclulate and return number of nonzero entries in vector or matrix.

Parameters
[in]matIncoming vector or matrix
Returns
Number of nonzeros

Referenced by dakota::surrogates::compute_hyperbolic_level_indices(), dakota::surrogates::compute_hyperbolic_subdim_level_indices(), and nonzero().

◆ nonzero()

void dakota::util::nonzero ( const T1 &  v,
T2 &  result 
)

Create a vector of indices based on nonzero entries in input vector.

Parameters
[in]vIncoming vector
[out]resultVector having values at nonzero locations of incoming vector and value equal to ordinal of occurrence

References num_nonzeros().

Referenced by dakota::surrogates::compute_hyperbolic_level_indices().

◆ append_columns()

void dakota::util::append_columns ( const T1 &  new_cols,
T2 &  target 
)

Append columns of input matrix to existing matrix.

Parameters
[in]new_colsIncoming matrix of column vectors to append
[out]targetMatrix to augment with appended columns

Referenced by dakota::surrogates::compute_hyperbolic_indices(), dakota::surrogates::compute_hyperbolic_level_indices(), and dakota::surrogates::compute_reduced_indices().

◆ p_norm()

double dakota::util::p_norm ( const T &  v,
double  p 
)

Caclulate and return p-norm of a vector.

Parameters
[in]vIncoming vector
[in]pOrder or norm to compute
Returns
p-norm of incoming vector

Referenced by dakota::surrogates::compute_hyperbolic_subdim_level_indices().

◆ metric_type()

METRIC_TYPE metric_type ( const std::string &  metric_name)

Convert the metric from string to enum.

Parameters
[in]metric_namemetric
Returns
converted metric

References type_name_bimap.

Referenced by compute_metric().

◆ compute_metric()

double compute_metric ( const VectorXd p,
const VectorXd d,
const std::string &  metric_name 
)

Computes the difference between prediction and data vectors.

Parameters
[in]pprediction vector.
[in]ddata vector.
[in]metric_namemetric to compute.
Returns
the value of the computed metric.

References error(), and metric_type().

Referenced by Dakota::compute_regression_coeffs(), and Surrogate::evaluate_metrics().

◆ scaler_factory()

std::shared_ptr< DataScaler > scaler_factory ( DataScaler::SCALER_TYPE  scaler_type,
const MatrixXd unscaled_matrix 
)

Free function to construct DataScaler.

Parameters
[in]scaler_typeWhich scaler to construct
[in]unscaled_matrixUnscaled data matrix - (num_samples by num_features)
Returns
Shared pointer to a DataScaler

Referenced by GaussianProcess::build(), and PolynomialRegression::build().

◆ solver_factory()

std::shared_ptr< LinearSolverBase > solver_factory ( LinearSolverBase::SOLVER_TYPE  type)

Free function to construct LinearSolverBase.

Parameters
[in]typeWhich solver to construct
Returns
Shared pointer to a LinearSolverBase

Referenced by PolynomialRegression::build().

Variable Documentation

◆ type_name_bimap [1/3]

BimapMetrictypeStr type_name_bimap
static
Initial value:
=
boost::assign::list_of<BimapMetrictypeStr::relation>
(METRIC_TYPE::SUM_SQUARED, "sum_squared")
(METRIC_TYPE::MEAN_SQUARED, "mean_squared")
(METRIC_TYPE::ROOT_MEAN_SQUARED, "root_mean_squared")
(METRIC_TYPE::SUM_ABS, "sum_abs")
(METRIC_TYPE::MEAN_ABS, "mean_abs")
(METRIC_TYPE::MAX_ABS, "max_abs")
(METRIC_TYPE::ABS_PERCENTAGE_ERROR, "ape")
(METRIC_TYPE::MEAN_ABS_PERCENTAGE_ERROR, "mape")
(METRIC_TYPE::R_SQUARED, "rsquared")

Bimap between metric types and names.

Referenced by metric_type(), DataScaler::scaler_type(), and LinearSolverBase::solver_type().

◆ type_name_bimap [2/3]

BimapScalertypeStr type_name_bimap
static
Initial value:
=
boost::assign::list_of<BimapScalertypeStr::relation>
(SCALER_TYPE::NONE, "none")
(SCALER_TYPE::STANDARDIZATION, "standardization")
(SCALER_TYPE::MEAN_NORMALIZATION, "mean normalization")
(SCALER_TYPE::MINMAX_NORMALIZATION, "min-max normalization")

Bimap between scaler types and names.

◆ type_name_bimap [3/3]

BimapSolvertypeStr type_name_bimap
static
Initial value:
=
boost::assign::list_of<BimapSolvertypeStr::relation>
(SOLVER_TYPE::CHOLESKY, "cholesky")
(SOLVER_TYPE::EQ_CONS_LEAST_SQ_REGRESSION, "equality-constrained lsq regression")
(SOLVER_TYPE::LASSO_REGRESSION, "lasso regression")
(SOLVER_TYPE::LEAST_ANGLE_REGRESSION, "least angle regression")
(SOLVER_TYPE::LU, "LU")
(SOLVER_TYPE::ORTHOG_MATCH_PURSUIT, "orthogonal matching pursuit")
(SOLVER_TYPE::QR_LEAST_SQ_REGRESSION, "QR lsq regression")
(SOLVER_TYPE::SVD_LEAST_SQ_REGRESSION, "SVD")

Bimap between solver types and names.