[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
#include <NNLS_solver.h>
Public Member Functions | |
NNLSSolver (const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input, const Epetra_CrsMatrix &A, Epetra_MpiComm &Comm, Epetra_Vector &b) | |
Default Constructor. | |
NNLSSolver (const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input, const Epetra_CrsMatrix &A, const bool is_input_A_matrix_transposed, Epetra_MpiComm &Comm, Epetra_Vector &b) | |
Constructor w/ transposed A matrix. | |
NNLSSolver (const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input, const Epetra_CrsMatrix &A, Epetra_MpiComm &Comm, Epetra_Vector &b, bool grad_exit_crit) | |
Constructor w/ Gradient Exit Condition. | |
NNLSSolver (const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input, const Epetra_CrsMatrix &A, const bool is_input_A_matrix_transposed, Epetra_MpiComm &Comm, Epetra_Vector &b, bool grad_exit_crit) | |
Constructor w/ Gradient Exit Condition & transposed A matrix. | |
NNLSSolver (const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input, const Epetra_CrsMatrix &A, Epetra_MpiComm &Comm, Epetra_Vector &b, bool iter_solver, int LS_iter, double LS_tol) | |
Constructor w/ Iterative Linear Solver. | |
NNLSSolver (const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input, const Epetra_CrsMatrix &A, const bool is_input_A_matrix_transposed, Epetra_MpiComm &Comm, Epetra_Vector &b, bool iter_solver, int LS_iter, double LS_tol) | |
Constructor w/ Iterative Linear Solver & transposed A matrix. | |
NNLSSolver (const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input, const Epetra_CrsMatrix &A, const bool is_input_A_matrix_transposed, Epetra_MpiComm &Comm, Epetra_Vector &b, bool grad_exit_crit, bool iter_solver, int LS_iter, double LS_tol) | |
Common Constructor w/ Gradient Exit Condition & Iterative Linear Solver & transposed A matrix. | |
~NNLSSolver () | |
Destructor. | |
bool | solve () |
Call to solve NNLS problem. | |
Epetra_Vector & | get_solution () |
Returns protected approximate solution. | |
void | starting_solution (Epetra_Vector &start) |
Initiliazes the solution vector, must be used before .solve is called. | |
Public Attributes | |
const Parameters::AllParameters *const | all_parameters |
Pointer to all parameters. | |
const dealii::ParameterHandler & | parameter_handler |
Parameter handler for storing the .prm file being ran. | |
bool | iter_solver_ |
Boolean used for iterative solvers. | |
bool | grad_exit_crit_ |
Boolean to use an exit Condition depending on maximum gradent. | |
const bool | is_input_A_matrix_transposed_ |
int | iter_ |
Number of iterations in the NNLS solver. | |
Protected Attributes | |
Epetra_MpiComm | Comm_ |
Epetra Commuicator Object with MPI. | |
const Epetra_CrsMatrix | A_ |
Epetra Matrix A allocated to a single core. | |
Epetra_Vector | b_ |
Epetra Vector b allocated to a single core. | |
Epetra_Vector | x_ |
Epetra Vector x, to be solved. Allocated to a single core. | |
Epetra_Vector | multi_x_ |
Epetra_Vector x, to be solved. Allocated to multiple cores. | |
int | LS_iter_ |
Needed if the an iterative solver is used. | |
double | LS_tol_ |
Needed if the an iterative solver is used. | |
int | numInactive_ |
Number of inactive points. | |
std::vector< bool > | Z |
Vector of booleans representing the columns in the active set. | |
std::vector< bool > | P |
Vector of boolean representing the columns in the inactive set. | |
Eigen::VectorXd | index_set |
Eigen Vector of the index_set. | |
Private Member Functions | |
void | epetra_permutation_matrix (Epetra_CrsMatrix &P_mat) |
Creates square permutation matrix based off the active/inactive set, no longer in use. | |
void | positive_set_matrix (Epetra_CrsMatrix &P_mat) |
Creates a matrix using the columns in A in the set P. | |
void | sub_into_x (Epetra_Vector &temp) |
Replaces the entries with x with the values in temp. | |
void | add_into_x (Epetra_Vector &temp, double alpha) |
Adds the values of temp times alpha into the solution vector x. | |
void | move_to_active_set (int idx) |
Moves the column at idx into the active set (updating the index_set, Z, P, and numInactive_) | |
void | move_to_inactive_set (int idx) |
Moves the column at idx into the inactive set (updating the index_set, Z, P, and numInactive_) | |
Non-Negagive Least Squares Solver for Epetra Structures Based on the NNLS in Eigen unsupported and in MATLAB: https://gitlab.com/libeigen/eigen/-/blob/master/unsupported/Eigen/NNLS https://www.mathworks.com/help/matlab/ref/lsqnonneg.html
Definition at line 77 of file NNLS_solver.h.
const bool PHiLiP::NNLSSolver::is_input_A_matrix_transposed_ |
Boolean to indicate whether A matrix is transposed wrt the dimension of b Note: This is because of the construction of the ECSW training data which when done in parallel requires the matrix A to be transposed
Definition at line 202 of file NNLS_solver.h.