libMesh::DenseVector< T > Class Template Reference

#include <analytic_function.h>

Inheritance diagram for libMesh::DenseVector< T >:

Public Member Functions

 DenseVector (const unsigned int n=0)
 
 DenseVector (const unsigned int n, const T &val)
 
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 libmesh_override
 
virtual bool empty () const libmesh_override
 
virtual void zero () libmesh_override
 
const T & operator() (const unsigned int i) const
 
T & operator() (const unsigned int i)
 
virtual T el (const unsigned int i) const libmesh_override
 
virtual T & el (const unsigned int i) libmesh_override
 
template<typename T2 >
DenseVector< T > & operator= (const DenseVector< T2 > &other_vector)
 
void swap (DenseVector< T > &other_vector)
 
void resize (const unsigned int n)
 
template<typename T2 >
void append (const DenseVector< T2 > &other_vector)
 
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, unsigned precision=8) const
 

Private Attributes

std::vector< T > _val
 

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 DenseMatrix class. It has additional capabilities over the std::vector that make it useful for finite elements, particularly for systems of equations. All overridden virtual functions are documented in dense_vector_base.h.

Author
Benjamin S. Kirk
Date
2003

Definition at line 34 of file analytic_function.h.

Constructor & Destructor Documentation

template<typename T >
libMesh::DenseVector< T >::DenseVector ( const unsigned int  n = 0)
inlineexplicit

Constructor. Creates a dense vector of dimension n.

Definition at line 272 of file dense_vector.h.

272  :
273  _val (n, T(0.))
274 {
275 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
libMesh::DenseVector< T >::DenseVector ( const unsigned int  n,
const T &  val 
)
inlineexplicit

Constructor. Creates a dense vector of dimension n where all entries have value val.

Definition at line 280 of file dense_vector.h.

281  :
282  _val (n, val)
283 {
284 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
template<typename T2 >
libMesh::DenseVector< T >::DenseVector ( const DenseVector< T2 > &  other_vector)
inline

Copy-constructor.

Definition at line 291 of file dense_vector.h.

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

291  :
292  DenseVectorBase<T>()
293 {
294  const std::vector<T2> & other_vals = other_vector.get_values();
295 
296  _val.clear();
297 
298  const int N = cast_int<int>(other_vals.size());
299  _val.reserve(N);
300 
301  for (int i=0; i<N; i++)
302  _val.push_back(other_vals[i]);
303 }
std::vector< T > _val
Definition: dense_vector.h:263
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 310 of file dense_vector.h.

310  :
311  _val(other_vector)
312 {
313 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
libMesh::DenseVector< T >::~DenseVector ( )
inline

Destructor. Does nothing.

Definition at line 85 of file dense_vector.h.

85 {}

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

Returns
A reference to *this.

Definition at line 431 of file dense_vector.h.

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

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::DenseVector< Number >::el(), and libMesh::DenseMatrix< T >::vector_mult_add().

433 {
434  libmesh_assert_equal_to (this->size(), vec.size());
435 
436  const int N = cast_int<int>(_val.size());
437  for (int i=0; i<N; i++)
438  (*this)(i) += static_cast<T>(factor)*vec(i);
439 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
template<typename T2 >
void libMesh::DenseVector< T >::append ( const DenseVector< T2 > &  other_vector)
inline

Append additional entries to (resizing, but unchanging) the vector.

Definition at line 362 of file dense_vector.h.

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

Referenced by libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), and libMesh::DenseVector< Number >::el().

363 {
364  const std::vector<T2> & other_vals = other_vector.get_values();
365 
366  _val.reserve(this->size() + other_vals.size());
367  _val.insert(_val.end(), other_vals.begin(), other_vals.end());
368 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
template<typename T2 >
CompareTypes< T, T2 >::supertype libMesh::DenseVector< T >::dot ( const DenseVector< T2 > &  vec) const
inline
Returns
The dot product of *this with vec.

In the complex-valued case, uses the complex conjugate of vec.

Definition at line 446 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

447 {
448  if (!_val.size())
449  return 0.;
450 
451  libmesh_assert_equal_to (this->size(), vec.size());
452 
453 #ifdef LIBMESH_HAVE_EIGEN
454  // We reverse the order of the arguments to dot() here since
455  // the convention in Eigen is to take the complex conjugate of the
456  // *first* argument, while ours is to take the complex conjugate of
457  // the second.
458  return Eigen::Map<const typename Eigen::Matrix<T2, Eigen::Dynamic, 1> >(&(vec.get_values()[0]), vec.size())
459  .dot(Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()));
460 #else
461  typename CompareTypes<T, T2>::supertype val = 0.;
462 
463  const int N = cast_int<int>(_val.size());
464  // The following pragma tells clang's vectorizer that it is safe to
465  // reorder floating point operations for this loop.
466 #ifdef __clang__
467 #pragma clang loop vectorize(enable)
468 #endif
469  for (int i=0; i<N; i++)
470  val += (*this)(i)*libmesh_conj(vec(i));
471 
472  return val;
473 #endif
474 }
T libmesh_conj(T a)
CompareTypes< T, T2 >::supertype dot(const DenseVector< T2 > &vec) const
Definition: dense_vector.h:446
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
virtual T libMesh::DenseVector< T >::el ( const unsigned int  i) const
inlinevirtual
Returns
The (i) element of the vector.

Implements libMesh::DenseVectorBase< T >.

Definition at line 107 of file dense_vector.h.

108  { return (*this)(i); }
template<typename T>
virtual T& libMesh::DenseVector< T >::el ( const unsigned int  i)
inlinevirtual
Returns
The (i) element of the vector as a writable reference.

Implements libMesh::DenseVectorBase< T >.

Definition at line 110 of file dense_vector.h.

111  { return (*this)(i); }
template<typename T>
virtual bool libMesh::DenseVector< T >::empty ( ) const
inlinevirtual
Returns
true iff size() is 0.

Reimplemented from libMesh::DenseVectorBase< T >.

Definition at line 92 of file dense_vector.h.

Referenced by libMesh::NumericVector< T >::add_vector(), and libMesh::NumericVector< T >::insert().

93  { return _val.empty(); }
std::vector< T > _val
Definition: dense_vector.h:263
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 667 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

669 {
670  libmesh_assert_less_equal ( sub_n, this->size() );
671 
672  dest.resize(sub_n);
673  const int N = cast_int<int>(sub_n);
674  for (int i=0; i<N; i++)
675  dest(i) = _val[i];
676 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
std::vector<T>& libMesh::DenseVector< T >::get_values ( )
inline
Returns
A reference to the underlying data storage vector.

This should be used with caution (i.e. one should not change the size of the vector, etc.) but is useful for interoperating with low level BLAS routines which expect a simple array.

Definition at line 251 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::DenseMatrix< T >::_svd_solve_lapack(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), and libMesh::FEMContext::pre_fe_reinit().

251 { return _val; }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
const std::vector<T>& libMesh::DenseVector< T >::get_values ( ) const
inline
Returns
A constant reference to the underlying data storage vector.

Definition at line 256 of file dense_vector.h.

256 { return _val; }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
template<typename T2 >
CompareTypes< T, T2 >::supertype libMesh::DenseVector< T >::indefinite_dot ( const DenseVector< T2 > &  vec) const
inline
Returns
The dot product of *this with vec.

In the complex-valued case, does not use the complex conjugate of vec.

Definition at line 479 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

480 {
481  libmesh_assert_equal_to (this->size(), vec.size());
482 
483  typename CompareTypes<T, T2>::supertype val = 0.;
484 
485  const int N = cast_int<int>(_val.size());
486  for (int i=0; i<N; i++)
487  val += (*this)(i)*(vec(i));
488 
489  return val;
490 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
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 of the entries.

Definition at line 596 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

597 {
598  if (!_val.size())
599  return 0.;
600 
601 #ifdef LIBMESH_HAVE_EIGEN
602  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()).template lpNorm<1>();
603 #else
604  Real my_norm = 0.;
605  const int N = cast_int<int>(_val.size());
606  for (int i=0; i!=N; i++)
607  my_norm += std::abs((*this)(i));
608 
609  return my_norm;
610 #endif
611 }
double abs(double a)
std::vector< T > _val
Definition: dense_vector.h:263
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 entries.

Definition at line 617 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

618 {
619  if (!_val.size())
620  return 0.;
621 
622 #ifdef LIBMESH_HAVE_EIGEN
623  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()).norm();
624 #else
625  Real my_norm = 0.;
626  const int N = cast_int<int>(_val.size());
627  // The following pragma tells clang's vectorizer that it is safe to
628  // reorder floating point operations for this loop.
629 #ifdef __clang__
630 #pragma clang loop vectorize(enable)
631 #endif
632  for (int i=0; i!=N; i++)
633  my_norm += TensorTools::norm_sq((*this)(i));
634 
635  return sqrt(my_norm);
636 #endif
637 }
std::vector< T > _val
Definition: dense_vector.h:263
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template<typename T >
Real libMesh::DenseVector< T >::linfty_norm ( ) const
inline
Returns
The $l_\infty$-norm of the vector, i.e. the maximum absolute value of the entries.

Definition at line 643 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

644 {
645  if (!_val.size())
646  return 0.;
647 
648 #ifdef LIBMESH_HAVE_EIGEN
649  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()).template lpNorm<Eigen::Infinity>();
650 #else
651  Real my_norm = TensorTools::norm_sq((*this)(0));
652 
653  const int N = cast_int<int>(_val.size());
654  for (int i=1; i!=N; i++)
655  {
656  Real current = TensorTools::norm_sq((*this)(i));
657  my_norm = (my_norm > current? my_norm : current);
658  }
659  return sqrt(my_norm);
660 #endif
661 }
std::vector< T > _val
Definition: dense_vector.h:263
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template<typename T >
Real libMesh::DenseVector< T >::max ( ) const
inline
Returns
The maximum entry of the vector, or the maximum real part in the case of complex numbers.

Definition at line 578 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

579 {
580  libmesh_assert (this->size());
581  Real my_max = libmesh_real((*this)(0));
582 
583  const int N = cast_int<int>(_val.size());
584  for (int i=1; i!=N; i++)
585  {
586  Real current = libmesh_real((*this)(i));
587  my_max = (my_max > current? my_max : current);
588  }
589  return my_max;
590 }
T libmesh_real(T a)
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
libmesh_assert(j)
std::vector< T > _val
Definition: dense_vector.h:263
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template<typename T >
Real libMesh::DenseVector< T >::min ( ) const
inline
Returns
The minimum entry of the vector, or the minimum real part in the case of complex numbers.

Definition at line 560 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

561 {
562  libmesh_assert (this->size());
563  Real my_min = libmesh_real((*this)(0));
564 
565  const int N = cast_int<int>(_val.size());
566  for (int i=1; i!=N; i++)
567  {
568  Real current = libmesh_real((*this)(i));
569  my_min = (my_min < current? my_min : current);
570  }
571  return my_min;
572 }
T libmesh_real(T a)
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
libmesh_assert(j)
std::vector< T > _val
Definition: dense_vector.h:263
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
template<typename T >
template<typename T2 >
bool libMesh::DenseVector< T >::operator!= ( const DenseVector< T2 > &  vec) const
inline
Returns
true if vec is not exactly equal to this vector, false otherwise.

Definition at line 512 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

513 {
514  libmesh_assert_equal_to (this->size(), vec.size());
515 
516  const int N = cast_int<int>(_val.size());
517  for (int i=0; i<N; i++)
518  if ((*this)(i) != vec(i))
519  return true;
520 
521  return false;
522 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
const T & libMesh::DenseVector< T >::operator() ( const unsigned int  i) const
inline
Returns
Entry i of the vector as a const reference.

Definition at line 385 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::empty().

386 {
387  libmesh_assert_less (i, _val.size());
388 
389  return _val[i];
390 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
T & libMesh::DenseVector< T >::operator() ( const unsigned int  i)
inline
Returns
Entry i of the vector as a writable reference.

Definition at line 396 of file dense_vector.h.

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

397 {
398  libmesh_assert_less (i, _val.size());
399 
400  return _val[i];
401 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
DenseVector< T > & libMesh::DenseVector< T >::operator*= ( const T  factor)
inline

Multiplies every element in the vector by factor.

Returns
A reference to *this.

Definition at line 418 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

419 {
420  this->scale(factor);
421  return *this;
422 }
void scale(const T factor)
Definition: dense_vector.h:407
template<typename T >
template<typename T2 >
DenseVector< T > & libMesh::DenseVector< T >::operator+= ( const DenseVector< T2 > &  vec)
inline

Adds vec to this vector.

Returns
A reference to *this.

Definition at line 529 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

530 {
531  libmesh_assert_equal_to (this->size(), vec.size());
532 
533  const int N = cast_int<int>(_val.size());
534  for (int i=0; i<N; i++)
535  (*this)(i) += vec(i);
536 
537  return *this;
538 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
template<typename T2 >
DenseVector< T > & libMesh::DenseVector< T >::operator-= ( const DenseVector< T2 > &  vec)
inline

Subtracts vec from this vector.

Returns
A reference to *this.

Definition at line 545 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

546 {
547  libmesh_assert_equal_to (this->size(), vec.size());
548 
549  const int N = cast_int<int>(_val.size());
550  for (int i=0; i<N; i++)
551  (*this)(i) -= vec(i);
552 
553  return *this;
554 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
template<typename T2 >
DenseVector< T > & libMesh::DenseVector< T >::operator= ( const DenseVector< T2 > &  other_vector)
inline

Assignment operator.

Returns
A reference to *this.

Definition at line 322 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

323 {
324  const std::vector<T2> & other_vals = other_vector.get_values();
325 
326  _val.clear();
327 
328  const int N = cast_int<int>(other_vals.size());
329  _val.reserve(N);
330 
331  for (int i=0; i<N; i++)
332  _val.push_back(other_vals[i]);
333 
334  return *this;
335 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
template<typename T2 >
bool libMesh::DenseVector< T >::operator== ( const DenseVector< T2 > &  vec) const
inline
Returns
true if vec is exactly equal to this vector, false otherwise.

Definition at line 495 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el().

496 {
497  libmesh_assert_equal_to (this->size(), vec.size());
498 
499  const int N = cast_int<int>(_val.size());
500  for (int i=0; i<N; i++)
501  if ((*this)(i) != vec(i))
502  return false;
503 
504  return true;
505 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
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.

Referenced by libMesh::DenseVectorBase< T >::empty().

63 {
64  for (unsigned int i=0; i<this->size(); i++)
65  os << std::setw(8)
66  << this->el(i)
67  << std::endl;
68 }
virtual T el(const unsigned int i) const =0
virtual unsigned int size() const =0
template<typename T >
void libMesh::DenseVectorBase< T >::print_scientific ( std::ostream &  os,
unsigned  precision = 8 
) 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.

31 {
32 #ifndef LIBMESH_BROKEN_IOSTREAM
33 
34  // save the initial format flags
35  std::ios_base::fmtflags os_flags = os.flags();
36 
37  // Print the vector entries.
38  for (unsigned int i=0; i<this->size(); i++)
39  os << std::setw(10)
40  << std::scientific
41  << std::setprecision(precision)
42  << this->el(i)
43  << std::endl;
44 
45  // reset the original format flags
46  os.flags(os_flags);
47 
48 #else
49 
50  // Print the matrix entries.
51  for (unsigned int i=0; i<this->size(); i++)
52  os << std::setprecision(precision)
53  << this->el(i)
54  << std::endl;
55 
56 #endif
57 }
virtual T el(const unsigned int i) const =0
virtual unsigned int size() const =0
template<typename T >
void libMesh::DenseVector< T >::resize ( const unsigned int  n)
inline

Resize the vector. Sets all elements to 0.

Definition at line 350 of file dense_vector.h.

References libMesh::DenseVector< T >::_val, and libMesh::DenseVector< T >::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::DenseMatrix< T >::_svd_solve_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::DenseVector< Number >::el(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), libMesh::CompositeFEMFunction< Output >::operator()(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::CompositeFunction< Output >::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::MeshFunction::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().

351 {
352  _val.resize(n);
353 
354  zero();
355 }
virtual void zero() libmesh_override
Definition: dense_vector.h:374
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
void libMesh::DenseVector< T >::scale ( const T  factor)
inline

Multiplies every element in the vector by factor.

Definition at line 407 of file dense_vector.h.

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

Referenced by libMesh::DenseVector< Number >::el(), and libMesh::DenseVector< T >::operator*=().

408 {
409  const int N = cast_int<int>(_val.size());
410  for (int i=0; i<N; i++)
411  _val[i] *= factor;
412 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
virtual unsigned int libMesh::DenseVector< T >::size ( ) const
inlinevirtual
Returns
The size of the vector.

Implements libMesh::DenseVectorBase< T >.

Definition at line 87 of file dense_vector.h.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::DenseMatrix< T >::_lu_back_substitute(), libMesh::DenseMatrix< T >::_matvec_blas(), libMesh::DenseMatrix< T >::_svd_lapack(), libMesh::DenseMatrix< T >::_svd_solve_lapack(), libMesh::DenseVector< T >::add(), libMesh::HPCoarsenTest::add_projection(), libMesh::NumericVector< T >::add_vector(), libMesh::DenseVector< T >::append(), libMesh::FEMSystem::assembly(), libMesh::DofMap::constrain_element_dyad_matrix(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::DofMap::constrain_element_vector(), libMesh::MeshFunction::discontinuous_gradient(), libMesh::MeshFunction::discontinuous_value(), libMesh::DenseVector< T >::dot(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::DenseVector< T >::get_principal_subvector(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector(), libMesh::DofMap::heterogenously_constrain_element_vector(), libMesh::DenseVector< T >::indefinite_dot(), libMesh::VectorSetAction< Val >::insert(), libMesh::NumericVector< T >::insert(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DenseVector< T >::max(), libMesh::DenseVector< T >::min(), libMesh::DenseVector< T >::operator!=(), libMesh::ConstFEMFunction< Output >::operator()(), libMesh::WrappedFunction< Output >::operator()(), libMesh::ParsedFEMFunction< Output >::operator()(), libMesh::CompositeFEMFunction< Output >::operator()(), libMesh::CompositeFunction< Output >::operator()(), libMesh::MeshFunction::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().

88  {
89  return cast_int<unsigned int>(_val.size());
90  }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T>
void libMesh::DenseVector< T >::swap ( DenseVector< T > &  other_vector)
inline

STL-like swap method

Definition at line 341 of file dense_vector.h.

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

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), and libMesh::DenseVector< Number >::el().

342 {
343  _val.swap(other_vector._val);
344 }
std::vector< T > _val
Definition: dense_vector.h:263
template<typename T >
void libMesh::DenseVector< T >::zero ( )
inlinevirtual

Member Data Documentation


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