libMesh::DenseVector< T > Class Template Reference

#include <dense_vector.h>

Inheritance diagram for libMesh::DenseVector< T >:

List of all members.

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 void zero ()
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 >
Number dot (const DenseVector< T2 > &vec) const
template<typename T2 >
Number 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)

Detailed Description

template<typename T>
class libMesh::DenseVector< T >

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.

Author:
Benjamin S. Kirk, 2003

Definition at line 51 of file dense_vector.h.


Constructor & Destructor Documentation

template<typename T >
libMesh::DenseVector< T >::DenseVector ( const unsigned int  n = 0  )  [inline, explicit]

Constructor. Creates a dense vector of dimension n.

Definition at line 254 of file dense_vector.h.

00254                                                 :
00255   _val (n, T(0.))
00256 {
00257 }

template<typename T >
template<typename T2 >
libMesh::DenseVector< T >::DenseVector ( const DenseVector< T2 > &  other_vector  )  [inline]

Copy-constructor.

Definition at line 264 of file dense_vector.h.

References libMesh::DenseVector< T >::_val, and libMesh::DenseVector< T >::get_values().

00264                                                                 :
00265   DenseVectorBase<T>()
00266 {
00267   const std::vector<T2> &other_vals = other_vector.get_values();
00268 
00269   _val.clear();
00270   _val.reserve(other_vals.size());
00271 
00272   for (unsigned int i=0; i<other_vals.size(); i++)
00273     _val.push_back(other_vals[i]);
00274 }

template<typename T >
template<typename T2 >
libMesh::DenseVector< T >::DenseVector ( const std::vector< T2 > &  other_vector  )  [inline]

Copy-constructor, from a std::vector.

Definition at line 281 of file dense_vector.h.

00281                                                               :
00282   _val(other_vector)
00283 {
00284 }

template<typename T>
libMesh::DenseVector< T >::~DenseVector (  )  [inline]

Destructor. Does nothing.

Definition at line 76 of file dense_vector.h.

00076 {}


Member Function Documentation

template<typename T >
template<typename T2 , typename T3 >
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 388 of file dense_vector.h.

References libMesh::DenseVector< T >::size().

Referenced by libMesh::EulerSolver::element_residual(), libMesh::EulerSolver::side_residual(), and libMesh::DenseMatrix< T >::vector_mult_add().

00390 {
00391   libmesh_assert_equal_to (this->size(), vec.size());
00392 
00393   for (unsigned int i=0; i<this->size(); i++)
00394     (*this)(i) += factor*vec(i);
00395 }

template<typename T >
template<typename T2 >
Number 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 400 of file dense_vector.h.

References libMesh::libmesh_conj(), and libMesh::DenseVector< T >::size().

00401 {
00402   libmesh_assert_equal_to (this->size(), vec.size());
00403 
00404   Number val = 0.;
00405 
00406   for (unsigned int i=0; i<this->size(); i++)
00407     val += (*this)(i)*libmesh_conj(vec(i));
00408 
00409   return val;
00410 }

template<typename T>
virtual T& libMesh::DenseVector< T >::el ( const unsigned int  i  )  [inline, virtual]
Returns:
the (i) element of the vector as a writeable reference.

Implements libMesh::DenseVectorBase< T >.

Definition at line 108 of file dense_vector.h.

00108 { return (*this)(i); }

template<typename T>
virtual T libMesh::DenseVector< T >::el ( const unsigned int  i  )  const [inline, virtual]
Returns:
the (i) element of the vector.

Implements libMesh::DenseVectorBase< T >.

Definition at line 103 of file dense_vector.h.

00103 { return (*this)(i); }

template<typename T>
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 569 of file dense_vector.h.

References libMesh::DenseVector< T >::resize(), and libMesh::DenseVector< T >::size().

00571 {
00572   libmesh_assert_less_equal ( sub_n, this->size() );
00573 
00574   dest.resize(sub_n);
00575   for(unsigned int i=0; i<sub_n; i++)
00576     dest(i) = (*this)(i);
00577 }

template<typename T>
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 237 of file dense_vector.h.

00237 { return _val; }

template<typename T>
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 230 of file dense_vector.h.

Referenced by libMesh::DenseMatrix< T >::_lu_back_substitute_lapack(), libMesh::DenseMatrix< T >::_matvec_blas(), libMesh::DenseVector< T >::DenseVector(), and libMesh::DenseVector< T >::operator=().

00230 { return _val; }

template<typename T >
template<typename T2 >
Number 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 415 of file dense_vector.h.

References libMesh::DenseVector< T >::size().

00416 {
00417   libmesh_assert_equal_to (this->size(), vec.size());
00418 
00419   Number val = 0.;
00420 
00421   for (unsigned int i=0; i<this->size(); i++)
00422     val += (*this)(i)*(vec(i));
00423 
00424   return val;
00425 }

template<typename T >
Real libMesh::DenseVector< T >::l1_norm (  )  const [inline]
Returns:
the $l_1$-norm of the vector, i.e. the sum of the absolute values.

Definition at line 525 of file dense_vector.h.

References std::abs(), libMesh::Real, and libMesh::DenseVector< T >::size().

00526 {
00527   Real my_norm = 0.;
00528   for (unsigned int i=0; i!=this->size(); i++)
00529     {
00530       my_norm += std::abs((*this)(i));
00531     }
00532   return my_norm;
00533 }

template<typename T >
Real libMesh::DenseVector< T >::l2_norm (  )  const [inline]
Returns:
the $l_2$-norm of the vector, i.e. the square root of the sum of the squares of the elements.

Definition at line 539 of file dense_vector.h.

References libMesh::TensorTools::norm_sq(), libMesh::Real, and libMesh::DenseVector< T >::size().

00540 {
00541   Real my_norm = 0.;
00542   for (unsigned int i=0; i!=this->size(); i++)
00543     {
00544       my_norm += TensorTools::norm_sq((*this)(i));
00545     }
00546   return sqrt(my_norm);
00547 }

template<typename T >
Real libMesh::DenseVector< T >::linfty_norm (  )  const [inline]
Returns:
the maximum absolute value of the elements of this vector, which is the $l_\infty$-norm of a vector.

Definition at line 553 of file dense_vector.h.

References libMesh::TensorTools::norm_sq(), libMesh::Real, and libMesh::DenseVector< T >::size().

00554 {
00555   if (!this->size())
00556     return 0.;
00557   Real my_norm = TensorTools::norm_sq((*this)(0));
00558 
00559   for (unsigned int i=1; i!=this->size(); i++)
00560     {
00561       Real current = TensorTools::norm_sq((*this)(i));
00562       my_norm = (my_norm > current? my_norm : current);
00563     }
00564   return sqrt(my_norm);
00565 }

template<typename T >
Real libMesh::DenseVector< T >::max (  )  const [inline]
Returns:
the maximum element in the vector. In case of complex numbers, this returns the maximum Real part.

Definition at line 508 of file dense_vector.h.

References libMesh::libmesh_real(), libMesh::Real, and libMesh::DenseVector< T >::size().

00509 {
00510   libmesh_assert (this->size());
00511   Real my_max = libmesh_real((*this)(0));
00512 
00513   for (unsigned int i=1; i!=this->size(); i++)
00514     {
00515       Real current = libmesh_real((*this)(i));
00516       my_max = (my_max > current? my_max : current);
00517     }
00518   return my_max;
00519 }

template<typename T >
Real libMesh::DenseVector< T >::min (  )  const [inline]
Returns:
the minimum element in the vector. In case of complex numbers, this returns the minimum Real part.

Definition at line 491 of file dense_vector.h.

References libMesh::libmesh_real(), libMesh::Real, and libMesh::DenseVector< T >::size().

00492 {
00493   libmesh_assert (this->size());
00494   Real my_min = libmesh_real((*this)(0));
00495 
00496   for (unsigned int i=1; i!=this->size(); i++)
00497     {
00498       Real current = libmesh_real((*this)(i));
00499       my_min = (my_min < current? my_min : current);
00500     }
00501   return my_min;
00502 }

template<typename T >
template<typename T2 >
bool libMesh::DenseVector< T >::operator!= ( const DenseVector< T2 > &  vec  )  const [inline]

Tests if vec is not exactly equal to this vector.

Definition at line 446 of file dense_vector.h.

References libMesh::DenseVector< T >::size().

00447 {
00448   libmesh_assert_equal_to (this->size(), vec.size());
00449 
00450   for (unsigned int i=0; i<this->size(); i++)
00451     if ((*this)(i) != vec(i))
00452       return true;
00453 
00454   return false;
00455 }

template<typename T >
T & libMesh::DenseVector< T >::operator() ( const unsigned int  i  )  [inline]
Returns:
the (i,j) element of the vector as a writeable reference.

Definition at line 354 of file dense_vector.h.

References libMesh::DenseVector< T >::_val.

00355 {
00356   libmesh_assert_less (i, _val.size());
00357 
00358   return _val[i];
00359 }

template<typename T >
T libMesh::DenseVector< T >::operator() ( const unsigned int  i  )  const [inline]
Returns:
the (i) element of the vector.

Definition at line 343 of file dense_vector.h.

References libMesh::DenseVector< T >::_val.

00344 {
00345   libmesh_assert_less (i, _val.size());
00346 
00347   return _val[i];
00348 }

template<typename T>
DenseVector< T > & libMesh::DenseVector< T >::operator*= ( const T  factor  )  [inline]

Multiplies every element in the vector by factor.

Definition at line 375 of file dense_vector.h.

References libMesh::DenseVector< T >::scale().

00376 {
00377   this->scale(factor);
00378   return *this;
00379 }

template<typename T >
template<typename T2 >
DenseVector< T > & libMesh::DenseVector< T >::operator+= ( const DenseVector< T2 > &  vec  )  [inline]

Adds vec to this vector.

Definition at line 462 of file dense_vector.h.

References libMesh::DenseVector< T >::size().

00463 {
00464   libmesh_assert_equal_to (this->size(), vec.size());
00465 
00466   for (unsigned int i=0; i<this->size(); i++)
00467     (*this)(i) += vec(i);
00468 
00469   return *this;
00470 }

template<typename T >
template<typename T2 >
DenseVector< T > & libMesh::DenseVector< T >::operator-= ( const DenseVector< T2 > &  vec  )  [inline]

Subtracts vec from this vector.

Definition at line 477 of file dense_vector.h.

References libMesh::DenseVector< T >::size().

00478 {
00479   libmesh_assert_equal_to (this->size(), vec.size());
00480 
00481   for (unsigned int i=0; i<this->size(); i++)
00482     (*this)(i) -= vec(i);
00483 
00484   return *this;
00485 }

template<typename T >
template<typename T2 >
DenseVector< T > & libMesh::DenseVector< T >::operator= ( const DenseVector< T2 > &  other_vector  )  [inline]

Assignment operator.

Definition at line 293 of file dense_vector.h.

References libMesh::DenseVector< T >::_val, and libMesh::DenseVector< T >::get_values().

00294 {
00295   //  _val = other_vector._val;
00296 
00297   const std::vector<T2> &other_vals = other_vector.get_values();
00298 
00299   _val.clear();
00300   _val.reserve(other_vals.size());
00301 
00302   for (unsigned int i=0; i<other_vals.size(); i++)
00303     _val.push_back(other_vals[i]);
00304 
00305   return *this;
00306 }

template<typename T >
template<typename T2 >
bool libMesh::DenseVector< T >::operator== ( const DenseVector< T2 > &  vec  )  const [inline]

Tests if vec is exactly equal to this vector.

Definition at line 430 of file dense_vector.h.

References libMesh::DenseVector< T >::size().

00431 {
00432   libmesh_assert_equal_to (this->size(), vec.size());
00433 
00434   for (unsigned int i=0; i<this->size(); i++)
00435     if ((*this)(i) != vec(i))
00436       return false;
00437 
00438   return true;
00439 }

template<typename T >
void libMesh::DenseVectorBase< T >::print ( std::ostream &  os  )  const [inline, inherited]

Pretty-print the vector to stdout.

Definition at line 62 of file dense_vector_base.C.

References libMesh::DenseVectorBase< T >::el(), and libMesh::DenseVectorBase< T >::size().

00063 {
00064   for (unsigned int i=0; i<this->size(); i++)
00065     os << std::setw(8)
00066        << this->el(i)
00067        << std::endl;
00068 }

template<typename T >
void libMesh::DenseVectorBase< T >::print_scientific ( std::ostream &  os  )  const [inline, inherited]

Prints the entries of the vector with additional decimal places in scientific notation.

Definition at line 30 of file dense_vector_base.C.

References libMesh::DenseVectorBase< T >::el(), and libMesh::DenseVectorBase< T >::size().

00031 {
00032 #ifndef LIBMESH_BROKEN_IOSTREAM
00033 
00034   // save the initial format flags
00035   std::ios_base::fmtflags os_flags = os.flags();
00036 
00037   // Print the vector entries.
00038   for (unsigned int i=0; i<this->size(); i++)
00039     os << std::setw(10)
00040        << std::scientific
00041        << std::setprecision(8)
00042        << this->el(i)
00043        << std::endl;
00044 
00045   // reset the original format flags
00046   os.flags(os_flags);
00047 
00048 #else
00049 
00050   // Print the matrix entries.
00051   for (unsigned int i=0; i<this->size(); i++)
00052     os << std::setprecision(8)
00053        << this->el(i)
00054        << std::endl;
00055 
00056 #endif
00057 }

template<typename T>
void libMesh::DenseVector< T >::scale ( const T  factor  )  [inline]

Multiplies every element in the vector by factor.

Definition at line 365 of file dense_vector.h.

References libMesh::DenseVector< T >::_val.

Referenced by libMesh::DenseVector< T >::operator*=().

00366 {
00367   for (unsigned int i=0; i<_val.size(); i++)
00368     _val[i] *= factor;
00369 }

template<typename T>
virtual unsigned int libMesh::DenseVector< T >::size (  )  const [inline, virtual]
Returns:
the size of the vector.

Implements libMesh::DenseVectorBase< T >.

Definition at line 81 of file dense_vector.h.

Referenced by libMesh::DenseMatrix< T >::_lu_back_substitute(), libMesh::DenseMatrix< T >::_matvec_blas(), libMesh::DenseVector< T >::add(), libMesh::HPCoarsenTest::add_projection(), libMesh::PetscVector< T >::add_vector(), libMesh::LaspackVector< T >::add_vector(), libMesh::DistributedVector< T >::add_vector(), libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::DiscontinuityMeasure::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::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::DenseVector< T >::get_principal_subvector(), libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector(), libMesh::DenseVector< T >::indefinite_dot(), libMesh::PetscVector< T >::insert(), libMesh::LaspackVector< T >::insert(), libMesh::DistributedVector< T >::insert(), libMesh::KellyErrorEstimator::internal_side_integration(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), libMesh::DenseVector< T >::l1_norm(), libMesh::DenseVector< T >::l2_norm(), libMesh::DenseVector< T >::linfty_norm(), libMesh::DenseVector< T >::max(), libMesh::DenseVector< T >::min(), libMesh::DenseVector< T >::operator!=(), libMesh::WrappedFunction< Output >::operator()(), libMesh::ParsedFunction< Output >::operator()(), libMesh::ConstFunction< Output >::operator()(), libMesh::ConstFEMFunction< Output >::operator()(), libMesh::DenseVector< T >::operator+=(), libMesh::DenseVector< T >::operator-=(), libMesh::DenseVector< T >::operator==(), libMesh::HPCoarsenTest::select_refinement(), libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), libMesh::DenseMatrix< T >::vector_mult(), libMesh::DenseMatrix< T >::vector_mult_add(), and libMesh::DenseMatrix< T >::vector_mult_transpose().

00081                                     {
00082     return libmesh_cast_int<unsigned int>(_val.size());
00083   }

template<typename T>
void libMesh::DenseVector< T >::swap ( DenseVector< T > &  other_vector  )  [inline]

STL-like swap method

Definition at line 312 of file dense_vector.h.

References libMesh::DenseVector< T >::_val.

Referenced by libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::EulerSolver::side_residual(), and libMesh::Euler2Solver::side_residual().

00313 {
00314   _val.swap(other_vector._val);
00315 }


Friends And Related Function Documentation

template<typename T>
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 91 of file dense_vector_base.h.

00092   {
00093     v.print(os);
00094     return os;
00095   }


Member Data Documentation


The documentation for this class was generated from the following file:

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:13 UTC

Hosted By:
SourceForge.net Logo