DenseVector< T > Class Template Reference

#include <dense_vector.h>

Inheritance diagram for 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 >
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 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 49 of file dense_vector.h.


Constructor & Destructor Documentation

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

Constructor. Creates a dense vector of dimension n.

Definition at line 242 of file dense_vector.h.

00242                                                 :
00243   _val (n, 0.)
00244 {
00245 }

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

Copy-constructor.

Definition at line 252 of file dense_vector.h.

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

00252                                                                 :
00253   DenseVectorBase<T>()
00254 {
00255   const std::vector<T2> &other_vals = other_vector.get_values();
00256 
00257   _val.clear();
00258   _val.reserve(other_vals.size());
00259 
00260   for (unsigned int i=0; i<other_vals.size(); i++)
00261     _val.push_back(other_vals[i]);
00262 }

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

Copy-constructor, from a std::vector.

Definition at line 269 of file dense_vector.h.

00269                                                               :
00270   _val(other_vector)
00271 {  
00272 }

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

Destructor. Does nothing.

Definition at line 74 of file dense_vector.h.

00074 {}


Member Function Documentation

template<typename T >
template<typename T2 , typename T3 >
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type 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 376 of file dense_vector.h.

References DenseVector< T >::size().

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

00378 {
00379   libmesh_assert (this->size() == vec.size());
00380 
00381   for (unsigned int i=0; i<this->size(); i++)
00382     (*this)(i) += factor*vec(i);
00383 }

template<typename T >
template<typename T2 >
Number DenseVector< T >::dot ( const DenseVector< T2 > &  vec  )  const [inline]

Evaluate dot product with vec.

Definition at line 388 of file dense_vector.h.

References DenseVector< T >::size().

00389 {
00390   libmesh_assert(this->size() == vec.size());
00391 
00392   Number val = 0.;
00393 
00394   for (unsigned int i=0; i<this->size(); i++)
00395     val += (*this)(i)*vec(i);
00396 
00397   return val;
00398 }

template<typename T >
virtual T& DenseVector< T >::el ( const unsigned int  i  )  [inline, virtual]

Returns:
the (i) element of the vector as a writeable reference.

Implements DenseVectorBase< T >.

Definition at line 104 of file dense_vector.h.

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

template<typename T >
virtual T DenseVector< T >::el ( const unsigned int  i  )  const [inline, virtual]

Returns:
the (i) element of the vector.

Implements DenseVectorBase< T >.

Definition at line 99 of file dense_vector.h.

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

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

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

00544 {
00545   libmesh_assert( sub_n <= this->size() );
00546 
00547   dest.resize(sub_n);
00548   for(unsigned int i=0; i<sub_n; i++)
00549     dest(i) = (*this)(i);
00550 }

template<typename T >
const std::vector<T>& 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 225 of file dense_vector.h.

00225 { return _val; }

template<typename T >
std::vector<T>& 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 218 of file dense_vector.h.

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

00218 { return _val; }

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

References DenseVector< T >::size().

00499 {
00500   Real my_norm = 0.;
00501   for (unsigned int i=0; i!=this->size(); i++)
00502     {
00503       my_norm += std::abs((*this)(i));
00504     }
00505   return my_norm;
00506 }

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

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

00513 {
00514   Real my_norm = 0.;
00515   for (unsigned int i=0; i!=this->size(); i++)
00516     {
00517       my_norm += libmesh_norm((*this)(i));
00518     }
00519   return sqrt(my_norm);
00520 }

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

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

00527 {
00528   if (!this->size())
00529     return 0.;
00530   Real my_norm = libmesh_norm((*this)(0));
00531 
00532   for (unsigned int i=1; i!=this->size(); i++)
00533     {
00534       Real current = libmesh_norm((*this)(i));
00535       my_norm = (my_norm > current? my_norm : current);
00536     }
00537   return sqrt(my_norm);
00538 }

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

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

00482 {
00483   libmesh_assert (this->size());
00484   Real my_max = libmesh_real((*this)(0));
00485 
00486   for (unsigned int i=1; i!=this->size(); i++)
00487     {
00488       Real current = libmesh_real((*this)(i));
00489       my_max = (my_max > current? my_max : current);
00490     }
00491   return my_max;
00492 }

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

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

00465 {
00466   libmesh_assert (this->size());
00467   Real my_min = libmesh_real((*this)(0));
00468 
00469   for (unsigned int i=1; i!=this->size(); i++)
00470     {
00471       Real current = libmesh_real((*this)(i));
00472       my_min = (my_min < current? my_min : current);
00473     }
00474   return my_min;
00475 }

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

Tests if vec is not exactly equal to this vector.

Definition at line 419 of file dense_vector.h.

References DenseVector< T >::size().

00420 {
00421   libmesh_assert (this->size() == vec.size());
00422 
00423   for (unsigned int i=0; i<this->size(); i++)
00424     if ((*this)(i) != vec(i))
00425       return true;
00426 
00427   return false;
00428 }

template<typename T >
T & DenseVector< T >::operator() ( const unsigned int  i  )  [inline]

Returns:
the (i,j) element of the vector as a writeable reference.

Definition at line 342 of file dense_vector.h.

References DenseVector< T >::_val.

00343 {
00344   libmesh_assert (i < _val.size());
00345   
00346   return _val[i]; 
00347 }

template<typename T >
T DenseVector< T >::operator() ( const unsigned int  i  )  const [inline]

Returns:
the (i) element of the vector.

Definition at line 331 of file dense_vector.h.

References DenseVector< T >::_val.

00332 {
00333   libmesh_assert (i < _val.size());
00334   
00335   return _val[i]; 
00336 }

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

Multiplies every element in the vector by factor.

Definition at line 363 of file dense_vector.h.

References DenseVector< T >::scale().

00364 {
00365   this->scale(factor);
00366   return *this;
00367 }

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

Adds vec to this vector.

Definition at line 435 of file dense_vector.h.

References DenseVector< T >::size().

00436 {
00437   libmesh_assert (this->size() == vec.size());
00438 
00439   for (unsigned int i=0; i<this->size(); i++)
00440     (*this)(i) += vec(i);
00441 
00442   return *this;
00443 }

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

Subtracts vec from this vector.

Definition at line 450 of file dense_vector.h.

References DenseVector< T >::size().

00451 {
00452   libmesh_assert (this->size() == vec.size());
00453 
00454   for (unsigned int i=0; i<this->size(); i++)
00455     (*this)(i) -= vec(i);
00456 
00457   return *this;
00458 }

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

Assignment operator.

Definition at line 281 of file dense_vector.h.

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

00282 {
00283   //  _val = other_vector._val;
00284 
00285   const std::vector<T2> &other_vals = other_vector.get_values();
00286 
00287   _val.clear();
00288   _val.reserve(other_vals.size());
00289 
00290   for (unsigned int i=0; i<other_vals.size(); i++)
00291     _val.push_back(other_vals[i]);
00292 
00293   return *this;
00294 }

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

Tests if vec is exactly equal to this vector.

Definition at line 403 of file dense_vector.h.

References DenseVector< T >::size().

00404 {
00405   libmesh_assert (this->size() == vec.size());
00406 
00407   for (unsigned int i=0; i<this->size(); i++)
00408     if ((*this)(i) != vec(i))
00409       return false;
00410 
00411   return true;
00412 }

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

Pretty-print the vector to stdout.

Definition at line 61 of file dense_vector_base.C.

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

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

template<typename T >
void 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 29 of file dense_vector_base.C.

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

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

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

Multiplies every element in the vector by factor.

Definition at line 353 of file dense_vector.h.

References DenseVector< T >::_val.

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

00354 {
00355   for (unsigned int i=0; i<_val.size(); i++)
00356     _val[i] *= factor;
00357 }

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

STL-like swap method

Definition at line 300 of file dense_vector.h.

References DenseVector< T >::_val.

Referenced by DenseMatrix< T >::_lu_back_substitute_lapack().

00301 {
00302   _val.swap(other_vector._val);
00303 }

template<typename T >
void DenseVector< T >::zero (  )  [inline, virtual]

Set every element in the vector to 0.

Implements DenseVectorBase< T >.

Definition at line 320 of file dense_vector.h.

References DenseVector< T >::_val.

Referenced by FEBase::coarsened_dof_values(), DofMap::constrain_element_matrix_and_vector(), DofMap::constrain_element_vector(), System::ProjectVector::operator()(), and DenseVector< T >::resize().

00321 {
00322   std::fill (_val.begin(),
00323              _val.end(),
00324              0.);
00325 }


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 88 of file dense_vector_base.h.

00089   {
00090     v.print(os);
00091     return os;
00092   }


Member Data Documentation


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

Site Created By: libMesh Developers
Last modified: November 25 2009 03:44:03.

Hosted By:
SourceForge.net Logo