Poster of Linux kernelThe best gift for a Linux geek
DenseVector

DenseVector

Section: C Library Functions (3) Updated: Thu Apr 7 2011
Local index Up
 

NAME

DenseVector -  

SYNOPSIS


#include <dense_vector.h>

Inherits DenseVectorBase< T >.  

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 ()

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 > 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 intn = 0) [inline, explicit]Constructor. Creates a dense vector of dimension n.

Definition at line 242 of file dense_vector.h.

                                                :
  _val (n, 0.)
{
}
 

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().

                                                                :
  DenseVectorBase<T>()
{
  const std::vector<T2> &other_vals = other_vector.get_values();

  _val.clear();
  _val.reserve(other_vals.size());

  for (unsigned int i=0; i<other_vals.size(); i++)
    _val.push_back(other_vals[i]);
}
 

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.

                                                              :
  _val(other_vector)
{  
}
 

template<typename T> DenseVector< T >::~DenseVector () [inline]Destructor. Does nothing.

Definition at line 74 of file dense_vector.h.

{}
 

Member Function Documentation

 

template<typename T > template<typename T2 , typename T3 > boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type DenseVector< T >::add (const T2factor, 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().

{
  libmesh_assert (this->size() == vec.size());

  for (unsigned int i=0; i<this->size(); i++)
    (*this)(i) += factor*vec(i);
}
 

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().

{
  libmesh_assert(this->size() == vec.size());

  Number val = 0.;

  for (unsigned int i=0; i<this->size(); i++)
    val += (*this)(i)*vec(i);

  return val;
}
 

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

the (i) element of the vector.

Implements DenseVectorBase< T >.

Definition at line 99 of file dense_vector.h.

{ return (*this)(i); }
 

template<typename T> virtual T& DenseVector< T >::el (const unsigned inti) [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.

{ return (*this)(i); }
 

template<typename T> void DenseVector< T >::get_principal_subvector (unsigned intsub_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().

{
  libmesh_assert( sub_n <= this->size() );

  dest.resize(sub_n);
  for(unsigned int i=0; i<sub_n; i++)
    dest(i) = (*this)(i);
}
 

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 EpetraVector< T >::add_vector(), DenseVector< T >::DenseVector(), and DenseVector< T >::operator=().

{ return _val; }
 

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.

{ 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.

{
  Real my_norm = 0.;
  for (unsigned int i=0; i!=this->size(); i++)
    {
      my_norm += std::abs((*this)(i));
    }
  return my_norm;
}
 

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().

{
  Real my_norm = 0.;
  for (unsigned int i=0; i!=this->size(); i++)
    {
      my_norm += libmesh_norm((*this)(i));
    }
  return sqrt(my_norm);
}
 

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().

{
  if (!this->size())
    return 0.;
  Real my_norm = libmesh_norm((*this)(0));

  for (unsigned int i=1; i!=this->size(); i++)
    {
      Real current = libmesh_norm((*this)(i));
      my_norm = (my_norm > current? my_norm : current);
    }
  return sqrt(my_norm);
}
 

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().

{
  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;
}
 

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().

{
  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;
}
 

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().

{
  libmesh_assert (this->size() == vec.size());

  for (unsigned int i=0; i<this->size(); i++)
    if ((*this)(i) != vec(i))
      return true;

  return false;
}
 

template<typename T > T DenseVector< T >::operator() (const unsigned inti) const [inline]Returns:

the (i) element of the vector.

Definition at line 331 of file dense_vector.h.

{
  libmesh_assert (i < _val.size());
  
  return _val[i]; 
}
 

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

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

Definition at line 342 of file dense_vector.h.

{
  libmesh_assert (i < _val.size());
  
  return _val[i]; 
}
 

template<typename T> DenseVector< T > & DenseVector< T >::operator*= (const Tfactor) [inline]Multiplies every element in the vector by factor.

Definition at line 363 of file dense_vector.h.

References MeshTools::Modification::scale().

{
  this->scale(factor);
  return *this;
}
 

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().

{
  libmesh_assert (this->size() == vec.size());

  for (unsigned int i=0; i<this->size(); i++)
    (*this)(i) += vec(i);

  return *this;
}
 

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().

{
  libmesh_assert (this->size() == vec.size());

  for (unsigned int i=0; i<this->size(); i++)
    (*this)(i) -= vec(i);

  return *this;
}
 

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 >::get_values().

{
  //  _val = other_vector._val;

  const std::vector<T2> &other_vals = other_vector.get_values();

  _val.clear();
  _val.reserve(other_vals.size());

  for (unsigned int i=0; i<other_vals.size(); i++)
    _val.push_back(other_vals[i]);

  return *this;
}
 

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().

{
  libmesh_assert (this->size() == vec.size());

  for (unsigned int i=0; i<this->size(); i++)
    if ((*this)(i) != vec(i))
      return false;

  return true;
}
 

template<typename T > void DenseVectorBase< T >::print (std::ostream &os) const [inherited]Pretty-print the vector to stdout.

Definition at line 61 of file dense_vector_base.C.

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

template<typename T > void 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 29 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
}
 

template<typename T > void DenseVector< T >::resize (const unsigned intn) [inline]Resize the vector. Sets all elements to 0.

Definition at line 309 of file dense_vector.h.

References libMesh::zero.

Referenced by DenseMatrix< T >::_cholesky_back_substitute(), DenseMatrix< T >::_lu_back_substitute(), HPCoarsenTest::add_projection(), FEBase::coarsened_dof_values(), FEBase::compute_periodic_constraints(), FEBase::compute_proj_constraints(), DofMap::constrain_element_dyad_matrix(), DofMap::constrain_element_matrix_and_vector(), DofMap::constrain_element_vector(), JumpErrorEstimator::estimate_error(), DenseVector< T >::get_principal_subvector(), System::ProjectVector::operator()(), PatchRecoveryErrorEstimator::EstimateError::operator()(), MeshFunction::operator()(), FEMContext::reinit(), and HPCoarsenTest::select_refinement().

{
  _val.resize(n);

  zero();
}
 

template<typename T> void DenseVector< T >::scale (const Tfactor) [inline]Multiplies every element in the vector by factor.

Definition at line 353 of file dense_vector.h.

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

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

the size of the vector.

Implements DenseVectorBase< T >.

Definition at line 79 of file dense_vector.h.

Referenced by DenseMatrix< T >::_lu_back_substitute(), DenseVector< T >::add(), HPCoarsenTest::add_projection(), EpetraVector< T >::add_vector(), PetscVector< T >::add_vector(), LaspackVector< T >::add_vector(), DistributedVector< T >::add_vector(), KellyErrorEstimator::boundary_side_integration(), DiscontinuityMeasure::boundary_side_integration(), DofMap::constrain_element_dyad_matrix(), DofMap::constrain_element_matrix_and_vector(), DofMap::constrain_element_vector(), DenseVector< T >::dot(), EulerSolver::element_residual(), Euler2Solver::element_residual(), ExactErrorEstimator::find_squared_element_error(), EpetraVector< T >::insert(), PetscVector< T >::insert(), LaspackVector< T >::insert(), DistributedVector< T >::insert(), KellyErrorEstimator::internal_side_integration(), LaplacianErrorEstimator::internal_side_integration(), DiscontinuityMeasure::internal_side_integration(), DenseVector< T >::operator!=(), DenseVector< T >::operator+=(), DenseVector< T >::operator-=(), DenseVector< T >::operator==(), HPCoarsenTest::select_refinement(), EulerSolver::side_residual(), Euler2Solver::side_residual(), DenseMatrix< T >::vector_mult(), and DenseMatrix< T >::vector_mult_add().

{ return _val.size(); }
 

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 EulerSolver::element_residual(), Euler2Solver::element_residual(), EulerSolver::side_residual(), and Euler2Solver::side_residual().

{
  _val.swap(other_vector._val);
}
 

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.

Referenced by HPCoarsenTest::add_projection(), FEBase::coarsened_dof_values(), DofMap::constrain_element_matrix_and_vector(), DofMap::constrain_element_vector(), EulerSolver::element_residual(), Euler2Solver::element_residual(), FEMSystem::numerical_jacobian(), System::ProjectVector::operator()(), HPCoarsenTest::select_refinement(), EulerSolver::side_residual(), Euler2Solver::side_residual(), and DenseMatrix< T >::vector_mult().

{
  std::fill (_val.begin(),
             _val.end(),
             0.);
}
 

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.

  {
    v.print(os);
    return os;
  }
 

Member Data Documentation

 

template<typename T> std::vector<T> DenseVector< T >::_val [private]The actual data values, stored as a 1D array.

Definition at line 232 of file dense_vector.h.

Referenced by DenseVector< T >::DenseVector(), DenseVector< Number >::get_values(), DenseVector< Number >::size(), and DenseVector< T >::swap().

 

Author

Generated automatically by Doxygen for libMesh from the source code.


 

Index

NAME
SYNOPSIS
Public Member Functions
Private Attributes
Friends
Detailed Description
template<typename T> class DenseVector< T >
Constructor & Destructor Documentation
template<typename T > DenseVector< T >::DenseVector (const unsigned intn = 0) [inline, explicit]Constructor. Creates a dense vector of dimension n.
template<typename T > template<typename T2 > DenseVector< T >::DenseVector (const DenseVector< T2 > &other_vector) [inline]Copy-constructor.
template<typename T > template<typename T2 > DenseVector< T >::DenseVector (const std::vector< T2 > &other_vector) [inline]Copy-constructor, from a std::vector.
template<typename T> DenseVector< T >::~DenseVector () [inline]Destructor. Does nothing.
Member Function Documentation
template<typename T > template<typename T2 , typename T3 > boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type DenseVector< T >::add (const T2factor, 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
template<typename T > template<typename T2 > Number DenseVector< T >::dot (const DenseVector< T2 > &vec) const [inline]Evaluate dot product with vec.
template<typename T> virtual T DenseVector< T >::el (const unsigned inti) const [inline, virtual]Returns:
template<typename T> virtual T& DenseVector< T >::el (const unsigned inti) [inline, virtual]Returns:
template<typename T> void DenseVector< T >::get_principal_subvector (unsigned intsub_n, DenseVector< T > &dest) const [inline]Puts the principal subvector of size sub_n (i.e. first sub_n entries) into dest.
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.
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.
template<typename T > Real DenseVector< T >::l1_norm () const [inline]Returns:
template<typename T > Real DenseVector< T >::l2_norm () const [inline]Returns:
template<typename T > Real DenseVector< T >::linfty_norm () const [inline]Returns:
template<typename T > Real DenseVector< T >::max () const [inline]Returns:
template<typename T > Real DenseVector< T >::min () const [inline]Returns:
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.
template<typename T > T DenseVector< T >::operator() (const unsigned inti) const [inline]Returns:
template<typename T > T & DenseVector< T >::operator() (const unsigned inti) [inline]Returns:
template<typename T> DenseVector< T > & DenseVector< T >::operator*= (const Tfactor) [inline]Multiplies every element in the vector by factor.
template<typename T > template<typename T2 > DenseVector< T > & DenseVector< T >::operator+= (const DenseVector< T2 > &vec) [inline]Adds vec to this vector.
template<typename T > template<typename T2 > DenseVector< T > & DenseVector< T >::operator-= (const DenseVector< T2 > &vec) [inline]Subtracts vec from this vector.
template<typename T > template<typename T2 > DenseVector< T > & DenseVector< T >::operator= (const DenseVector< T2 > &other_vector) [inline]Assignment operator.
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.
template<typename T > void DenseVectorBase< T >::print (std::ostream &os) const [inherited]Pretty-print the vector to stdout.
template<typename T > void DenseVectorBase< T >::print_scientific (std::ostream &os) const [inherited]Prints the entries of the vector with additional decimal places in scientific notation.
template<typename T > void DenseVector< T >::resize (const unsigned intn) [inline]Resize the vector. Sets all elements to 0.
template<typename T> void DenseVector< T >::scale (const Tfactor) [inline]Multiplies every element in the vector by factor.
template<typename T> virtual unsigned int DenseVector< T >::size () const [inline, virtual]Returns:
template<typename T> void DenseVector< T >::swap (DenseVector< T > &other_vector) [inline]STL-like swap method
template<typename T > void DenseVector< T >::zero () [inline, virtual]Set every element in the vector to 0.
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.
Member Data Documentation
template<typename T> std::vector<T> DenseVector< T >::_val [private]The actual data values, stored as a 1D array.
Author

This document was created by man2html, using the manual pages.
Time: 21:44:52 GMT, April 16, 2011