$extrastylesheet
#include <dense_vector.h>

Public Member Functions | |
| DenseVector (const unsigned int n=0) | |
| template<typename T2 > | |
| DenseVector (const DenseVector< T2 > &other_vector) | |
| template<typename T2 > | |
| DenseVector (const std::vector< T2 > &other_vector) | |
| ~DenseVector () | |
| virtual unsigned int | size () const |
| virtual bool | empty () const |
| virtual void | zero () |
| const T & | operator() (const unsigned int i) const |
| T & | operator() (const unsigned int i) |
| virtual T | el (const unsigned int i) const |
| virtual T & | el (const unsigned int i) |
| template<typename T2 > | |
| DenseVector< T > & | operator= (const DenseVector< T2 > &other_vector) |
| void | swap (DenseVector< T > &other_vector) |
| void | resize (const unsigned int n) |
| void | scale (const T factor) |
| DenseVector< T > & | operator*= (const T factor) |
| template<typename T2 , typename T3 > | |
| boostcopy::enable_if_c < ScalarTraits< T2 >::value, void >::type | add (const T2 factor, const DenseVector< T3 > &vec) |
| template<typename T2 > | |
| CompareTypes< T, T2 >::supertype | dot (const DenseVector< T2 > &vec) const |
| template<typename T2 > | |
| CompareTypes< T, T2 >::supertype | indefinite_dot (const DenseVector< T2 > &vec) const |
| template<typename T2 > | |
| bool | operator== (const DenseVector< T2 > &vec) const |
| template<typename T2 > | |
| bool | operator!= (const DenseVector< T2 > &vec) const |
| template<typename T2 > | |
| DenseVector< T > & | operator+= (const DenseVector< T2 > &vec) |
| template<typename T2 > | |
| DenseVector< T > & | operator-= (const DenseVector< T2 > &vec) |
| Real | min () const |
| Real | max () const |
| Real | l1_norm () const |
| Real | l2_norm () const |
| Real | linfty_norm () const |
| void | get_principal_subvector (unsigned int sub_n, DenseVector< T > &dest) const |
| std::vector< T > & | get_values () |
| const std::vector< T > & | get_values () const |
| void | print (std::ostream &os) const |
| void | print_scientific (std::ostream &os) const |
Private Attributes | |
| std::vector< T > | _val |
Friends | |
| std::ostream & | operator<< (std::ostream &os, const DenseVectorBase< T > &v) |
Defines a dense vector for use in Finite Element-type computations. This class is to basically compliment the DenseMatix class. It has additional capabilities over the std::vector that make it useful for finite elements, particulary for systems of equations.
Definition at line 51 of file dense_vector.h.
| libMesh::DenseVector< T >::DenseVector | ( | const unsigned int | n = 0 | ) | [inline, explicit] |
Constructor. Creates a dense vector of dimension n.
Definition at line 259 of file dense_vector.h.
: _val (n, T(0.)) { }
| libMesh::DenseVector< T >::DenseVector | ( | const DenseVector< T2 > & | other_vector | ) | [inline] |
Copy-constructor.
Definition at line 269 of file dense_vector.h.
References libMesh::DenseVector< T >::_val, and libMesh::DenseVector< T >::get_values().
| libMesh::DenseVector< T >::DenseVector | ( | const std::vector< T2 > & | other_vector | ) | [inline] |
Copy-constructor, from a std::vector.
Definition at line 286 of file dense_vector.h.
: _val(other_vector) { }
| libMesh::DenseVector< T >::~DenseVector | ( | ) | [inline] |
| boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type libMesh::DenseVector< T >::add | ( | const T2 | factor, |
| const DenseVector< T3 > & | vec | ||
| ) | [inline] |
Adds factor times vec to this vector. This should only work if T += T2 * T3 is valid C++ and if T2 is scalar. Return type is void
Definition at line 393 of file dense_vector.h.
References libMesh::DenseVector< T >::size().
Referenced by libMesh::EulerSolver::_general_residual(), and libMesh::DenseMatrix< T >::vector_mult_add().
| CompareTypes< T, T2 >::supertype libMesh::DenseVector< T >::dot | ( | const DenseVector< T2 > & | vec | ) | const [inline] |
Evaluate dot product with vec. In the complex-valued case, use the complex conjugate of vec.
Definition at line 405 of file dense_vector.h.
References libMesh::libmesh_conj(), and libMesh::DenseVector< T >::size().
{
libmesh_assert_equal_to (this->size(), vec.size());
typename CompareTypes<T, T2>::supertype val = 0.;
for (unsigned int i=0; i<this->size(); i++)
val += (*this)(i)*libmesh_conj(vec(i));
return val;
}
| virtual T libMesh::DenseVector< T >::el | ( | const unsigned int | i | ) | const [inline, virtual] |
(i) element of the vector. Implements libMesh::DenseVectorBase< T >.
Definition at line 108 of file dense_vector.h.
{ return (*this)(i); }
| virtual T& libMesh::DenseVector< T >::el | ( | const unsigned int | i | ) | [inline, virtual] |
(i) element of the vector as a writeable reference. Implements libMesh::DenseVectorBase< T >.
Definition at line 113 of file dense_vector.h.
{ return (*this)(i); }
| virtual bool libMesh::DenseVector< T >::empty | ( | ) | const [inline, virtual] |
Reimplemented from libMesh::DenseVectorBase< T >.
Definition at line 88 of file dense_vector.h.
Referenced by libMesh::NumericVector< T >::add_vector(), and libMesh::NumericVector< T >::insert().
{ return _val.empty(); }
| void libMesh::DenseVector< T >::get_principal_subvector | ( | unsigned int | sub_n, |
| DenseVector< T > & | dest | ||
| ) | const [inline] |
Puts the principal subvector of size sub_n (i.e. first sub_n entries) into dest.
Definition at line 574 of file dense_vector.h.
References libMesh::DenseVector< T >::resize().
{
libmesh_assert_less_equal ( sub_n, this->size() );
dest.resize(sub_n);
for(unsigned int i=0; i<sub_n; i++)
dest(i) = (*this)(i);
}
| std::vector<T>& libMesh::DenseVector< T >::get_values | ( | ) | [inline] |
Access to the values array. This should be used with caution but can be used to speed up code compilation significantly.
Definition at line 235 of file dense_vector.h.
Referenced by libMesh::DenseMatrix< T >::_evd_lapack(), libMesh::DenseMatrix< T >::_lu_back_substitute_lapack(), libMesh::DenseMatrix< T >::_matvec_blas(), libMesh::DenseVector< T >::DenseVector(), libMesh::DenseVector< T >::operator=(), and libMesh::FEMContext::pre_fe_reinit().
{ return _val; }
| const std::vector<T>& libMesh::DenseVector< T >::get_values | ( | ) | const [inline] |
Access to the values array. This should be used with caution but can be used to speed up code compilation significantly.
Definition at line 242 of file dense_vector.h.
{ return _val; }
| CompareTypes< T, T2 >::supertype libMesh::DenseVector< T >::indefinite_dot | ( | const DenseVector< T2 > & | vec | ) | const [inline] |
Evaluate dot product with vec. In the complex-valued case, do not use the complex conjugate of vec.
Definition at line 420 of file dense_vector.h.
References libMesh::DenseVector< T >::size().
| Real libMesh::DenseVector< T >::l1_norm | ( | ) | const [inline] |
-norm of the vector, i.e. the sum of the absolute values. Definition at line 530 of file dense_vector.h.
References std::abs(), and libMesh::Real.
| Real libMesh::DenseVector< T >::l2_norm | ( | ) | const [inline] |
-norm of the vector, i.e. the square root of the sum of the squares of the elements. Definition at line 544 of file dense_vector.h.
References libMesh::TensorTools::norm_sq(), and libMesh::Real.
{
Real my_norm = 0.;
for (unsigned int i=0; i!=this->size(); i++)
{
my_norm += TensorTools::norm_sq((*this)(i));
}
return sqrt(my_norm);
}
| Real libMesh::DenseVector< T >::linfty_norm | ( | ) | const [inline] |
-norm of a vector. Definition at line 558 of file dense_vector.h.
References libMesh::TensorTools::norm_sq(), and libMesh::Real.
{
if (!this->size())
return 0.;
Real my_norm = TensorTools::norm_sq((*this)(0));
for (unsigned int i=1; i!=this->size(); i++)
{
Real current = TensorTools::norm_sq((*this)(i));
my_norm = (my_norm > current? my_norm : current);
}
return sqrt(my_norm);
}
| Real libMesh::DenseVector< T >::max | ( | ) | const [inline] |
Definition at line 513 of file dense_vector.h.
References libMesh::libmesh_assert(), libMesh::libmesh_real(), and libMesh::Real.
{
libmesh_assert (this->size());
Real my_max = libmesh_real((*this)(0));
for (unsigned int i=1; i!=this->size(); i++)
{
Real current = libmesh_real((*this)(i));
my_max = (my_max > current? my_max : current);
}
return my_max;
}
| Real libMesh::DenseVector< T >::min | ( | ) | const [inline] |
Definition at line 496 of file dense_vector.h.
References libMesh::libmesh_assert(), libMesh::libmesh_real(), and libMesh::Real.
{
libmesh_assert (this->size());
Real my_min = libmesh_real((*this)(0));
for (unsigned int i=1; i!=this->size(); i++)
{
Real current = libmesh_real((*this)(i));
my_min = (my_min < current? my_min : current);
}
return my_min;
}
| bool libMesh::DenseVector< T >::operator!= | ( | const DenseVector< T2 > & | vec | ) | const [inline] |
Tests if vec is not exactly equal to this vector.
Definition at line 451 of file dense_vector.h.
References libMesh::DenseVector< T >::size().
| const T & libMesh::DenseVector< T >::operator() | ( | const unsigned int | i | ) | const [inline] |
(i) element of the vector as a const reference. Definition at line 348 of file dense_vector.h.
| T & libMesh::DenseVector< T >::operator() | ( | const unsigned int | i | ) | [inline] |
(i,j) element of the vector as a writeable reference. Definition at line 359 of file dense_vector.h.
| DenseVector< T > & libMesh::DenseVector< T >::operator*= | ( | const T | factor | ) | [inline] |
Multiplies every element in the vector by factor.
Definition at line 380 of file dense_vector.h.
References libMesh::MeshTools::Modification::scale().
{
this->scale(factor);
return *this;
}
| DenseVector< T > & libMesh::DenseVector< T >::operator+= | ( | const DenseVector< T2 > & | vec | ) | [inline] |
Adds vec to this vector.
Definition at line 467 of file dense_vector.h.
References libMesh::DenseVector< T >::size().
| DenseVector< T > & libMesh::DenseVector< T >::operator-= | ( | const DenseVector< T2 > & | vec | ) | [inline] |
Subtracts vec from this vector.
Definition at line 482 of file dense_vector.h.
References libMesh::DenseVector< T >::size().
| DenseVector< T > & libMesh::DenseVector< T >::operator= | ( | const DenseVector< T2 > & | other_vector | ) | [inline] |
Assignment operator.
Definition at line 298 of file dense_vector.h.
References libMesh::DenseVector< T >::get_values().
| bool libMesh::DenseVector< T >::operator== | ( | const DenseVector< T2 > & | vec | ) | const [inline] |
Tests if vec is exactly equal to this vector.
Definition at line 435 of file dense_vector.h.
References libMesh::DenseVector< T >::size().
| void libMesh::DenseVectorBase< T >::print | ( | std::ostream & | os | ) | const [inherited] |
Pretty-print the vector to stdout.
Definition at line 62 of file dense_vector_base.C.
| void libMesh::DenseVectorBase< T >::print_scientific | ( | std::ostream & | os | ) | const [inherited] |
Prints the entries of the vector with additional decimal places in scientific notation.
Definition at line 30 of file dense_vector_base.C.
{
#ifndef LIBMESH_BROKEN_IOSTREAM
// save the initial format flags
std::ios_base::fmtflags os_flags = os.flags();
// Print the vector entries.
for (unsigned int i=0; i<this->size(); i++)
os << std::setw(10)
<< std::scientific
<< std::setprecision(8)
<< this->el(i)
<< std::endl;
// reset the original format flags
os.flags(os_flags);
#else
// Print the matrix entries.
for (unsigned int i=0; i<this->size(); i++)
os << std::setprecision(8)
<< this->el(i)
<< std::endl;
#endif
}
| void libMesh::DenseVector< T >::resize | ( | const unsigned int | n | ) | [inline] |
Resize the vector. Sets all elements to 0.
Definition at line 326 of file dense_vector.h.
References libMesh::zero.
Referenced by libMesh::DenseMatrix< T >::_cholesky_back_substitute(), libMesh::DenseMatrix< T >::_evd_lapack(), libMesh::DenseMatrix< T >::_lu_back_substitute(), libMesh::DenseMatrix< T >::_svd_lapack(), libMesh::HPCoarsenTest::add_projection(), libMesh::DofMap::build_constraint_matrix_and_vector(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::FEGenericBase< OutputType >::compute_periodic_constraints(), libMesh::FEGenericBase< OutputType >::compute_proj_constraints(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::DenseVector< T >::get_principal_subvector(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), libMesh::CompositeFunction< Output >::operator()(), libMesh::CompositeFEMFunction< Output >::operator()(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::ProjectSolution::operator()(), libMesh::MeshFunction::operator()(), libMesh::ProjectFEMSolution::operator()(), libMesh::BoundaryProjectSolution::operator()(), libMesh::FEMContext::pre_fe_reinit(), libMesh::HPCoarsenTest::select_refinement(), libMesh::DenseMatrix< T >::vector_mult(), libMesh::DenseMatrix< T >::vector_mult_add(), and libMesh::DenseMatrix< T >::vector_mult_transpose().
| void libMesh::DenseVector< T >::scale | ( | const T | factor | ) | [inline] |
Multiplies every element in the vector by factor.
Definition at line 370 of file dense_vector.h.
| virtual unsigned int libMesh::DenseVector< T >::size | ( | ) | const [inline, virtual] |
Implements libMesh::DenseVectorBase< T >.
Definition at line 81 of file dense_vector.h.
Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::DenseMatrix< T >::_lu_back_substitute(), libMesh::DenseMatrix< T >::_matvec_blas(), libMesh::DenseVector< T >::add(), libMesh::HPCoarsenTest::add_projection(), libMesh::NumericVector< T >::add_vector(), libMesh::FEMSystem::assembly(), libMesh::DiscontinuityMeasure::boundary_side_integration(), libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::DofMap::constrain_element_dyad_matrix(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::DofMap::constrain_element_vector(), libMesh::DenseVector< T >::dot(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector(), libMesh::DofMap::heterogenously_constrain_element_vector(), libMesh::DenseVector< T >::indefinite_dot(), libMesh::NumericVector< T >::insert(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), libMesh::KellyErrorEstimator::internal_side_integration(), libMesh::DenseVector< T >::operator!=(), libMesh::ConstFunction< Output >::operator()(), libMesh::ConstFEMFunction< Output >::operator()(), libMesh::WrappedFunction< Output >::operator()(), libMesh::CompositeFunction< Output >::operator()(), libMesh::CompositeFEMFunction< Output >::operator()(), libMesh::ParsedFunction< Output, OutputGradient >::operator()(), libMesh::ParsedFEMFunction< Output >::operator()(), libMesh::DenseVector< T >::operator+=(), libMesh::DenseVector< T >::operator-=(), libMesh::DenseVector< T >::operator==(), libMesh::HPCoarsenTest::select_refinement(), libMesh::DenseMatrix< T >::vector_mult(), libMesh::DenseMatrix< T >::vector_mult_add(), and libMesh::DenseMatrix< T >::vector_mult_transpose().
{
return cast_int<unsigned int>(_val.size());
}
| void libMesh::DenseVector< T >::swap | ( | DenseVector< T > & | other_vector | ) | [inline] |
STL-like swap method
Definition at line 317 of file dense_vector.h.
References libMesh::DenseVector< T >::_val.
Referenced by libMesh::EulerSolver::_general_residual(), and libMesh::Euler2Solver::_general_residual().
{
_val.swap(other_vector._val);
}
| void libMesh::DenseVector< T >::zero | ( | ) | [inline, virtual] |
Set every element in the vector to 0.
Implements libMesh::DenseVectorBase< T >.
Definition at line 337 of file dense_vector.h.
Referenced by libMesh::HPCoarsenTest::add_projection(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::FEMSystem::numerical_jacobian(), libMesh::CompositeFunction< Output >::operator()(), libMesh::CompositeFEMFunction< Output >::operator()(), libMesh::ProjectSolution::operator()(), libMesh::ProjectFEMSolution::operator()(), libMesh::BoundaryProjectSolution::operator()(), and libMesh::HPCoarsenTest::select_refinement().
| std::ostream& operator<< | ( | std::ostream & | os, |
| const DenseVectorBase< T > & | v | ||
| ) | [friend, inherited] |
Same as above, but allows you to print using the usual stream syntax.
Definition at line 96 of file dense_vector_base.h.
{
v.print(os);
return os;
}
std::vector<T> libMesh::DenseVector< T >::_val [private] |
The actual data values, stored as a 1D array.
Definition at line 249 of file dense_vector.h.
Referenced by libMesh::DenseVector< T >::DenseVector(), libMesh::DenseVector< Number >::empty(), libMesh::DenseVector< Number >::get_values(), libMesh::DenseVector< Number >::size(), and libMesh::DenseVector< T >::swap().