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)
 
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 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
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 269 of file dense_vector.h.

269  :
270  _val (n, T(0.))
271 {
272 }
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T >
template<typename T2 >
libMesh::DenseVector< T >::DenseVector ( const DenseVector< T2 > &  other_vector)
inline

Copy-constructor.

Definition at line 279 of file dense_vector.h.

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

279  :
280  DenseVectorBase<T>()
281 {
282  const std::vector<T2> & other_vals = other_vector.get_values();
283 
284  _val.clear();
285 
286  const int N = cast_int<int>(other_vals.size());
287  _val.reserve(N);
288 
289  for (int i=0; i<N; i++)
290  _val.push_back(other_vals[i]);
291 }
std::vector< T > _val
Definition: dense_vector.h:260
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 298 of file dense_vector.h.

298  :
299  _val(other_vector)
300 {
301 }
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T>
libMesh::DenseVector< T >::~DenseVector ( )
inline

Destructor. Does nothing.

Definition at line 76 of file dense_vector.h.

76 {}

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

421 {
422  libmesh_assert_equal_to (this->size(), vec.size());
423 
424  const int N = cast_int<int>(_val.size());
425  for (int i=0; i<N; i++)
426  (*this)(i) += static_cast<T>(factor)*vec(i);
427 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
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 350 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().

351 {
352  const std::vector<T2> & other_vals = other_vector.get_values();
353 
354  _val.reserve(this->size() + other_vals.size());
355  _val.insert(_val.end(), other_vals.begin(), other_vals.end());
356 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T >
template<typename T2 >
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 434 of file dense_vector.h.

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

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

435 {
436  if (!_val.size())
437  return 0.;
438 
439  libmesh_assert_equal_to (this->size(), vec.size());
440 
441 #ifdef LIBMESH_HAVE_EIGEN
442  // Note: we reverse the order of the arguments to dot() here since
443  // the convention in Eigen is to take the complex conjugate of the
444  // *first* argument, while ours is to take the complex conjugate of
445  // the second.
446  return Eigen::Map<const typename Eigen::Matrix<T2, Eigen::Dynamic, 1> >(&(vec.get_values()[0]), vec.size())
447  .dot(Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()));
448 #else
449  typename CompareTypes<T, T2>::supertype val = 0.;
450 
451  const int N = cast_int<int>(_val.size());
452  // The following pragma tells clang's vectorizer that it is safe to
453  // reorder floating point operations for this loop.
454 #pragma clang loop vectorize(enable)
455  for (int i=0; i<N; i++)
456  val += (*this)(i)*libmesh_conj(vec(i));
457 
458  return val;
459 #endif
460 }
T libmesh_conj(T a)
CompareTypes< T, T2 >::supertype dot(const DenseVector< T2 > &vec) const
Definition: dense_vector.h:434
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
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 110 of file dense_vector.h.

111  { 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 writeable reference.

Implements libMesh::DenseVectorBase< T >.

Definition at line 116 of file dense_vector.h.

117  { 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 89 of file dense_vector.h.

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

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

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

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

653 {
654  libmesh_assert_less_equal ( sub_n, this->size() );
655 
656  dest.resize(sub_n);
657  const int N = cast_int<int>(sub_n);
658  for (int i=0; i<N; i++)
659  dest(i) = _val[i];
660 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T>
std::vector<T>& libMesh::DenseVector< T >::get_values ( )
inline
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 253 of file dense_vector.h.

253 { return _val; }
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T >
template<typename T2 >
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 465 of file dense_vector.h.

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

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

466 {
467  libmesh_assert_equal_to (this->size(), vec.size());
468 
469  typename CompareTypes<T, T2>::supertype val = 0.;
470 
471  const int N = cast_int<int>(_val.size());
472  for (int i=0; i<N; i++)
473  val += (*this)(i)*(vec(i));
474 
475  return val;
476 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
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 582 of file dense_vector.h.

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

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

583 {
584  if (!_val.size())
585  return 0.;
586 
587 #ifdef LIBMESH_HAVE_EIGEN
588  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()).template lpNorm<1>();
589 #else
590  Real my_norm = 0.;
591  const int N = cast_int<int>(_val.size());
592  for (int i=0; i!=N; i++)
593  my_norm += std::abs((*this)(i));
594 
595  return my_norm;
596 #endif
597 }
double abs(double a)
std::vector< T > _val
Definition: dense_vector.h:260
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 elements.

Definition at line 603 of file dense_vector.h.

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

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

604 {
605  if (!_val.size())
606  return 0.;
607 
608 #ifdef LIBMESH_HAVE_EIGEN
609  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()).norm();
610 #else
611  Real my_norm = 0.;
612  const int N = cast_int<int>(_val.size());
613  // The following pragma tells clang's vectorizer that it is safe to
614  // reorder floating point operations for this loop.
615 #pragma clang loop vectorize(enable)
616  for (int i=0; i!=N; i++)
617  my_norm += TensorTools::norm_sq((*this)(i));
618 
619  return sqrt(my_norm);
620 #endif
621 }
std::vector< T > _val
Definition: dense_vector.h:260
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 627 of file dense_vector.h.

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

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

628 {
629  if (!_val.size())
630  return 0.;
631 
632 #ifdef LIBMESH_HAVE_EIGEN
633  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1> >(&_val[0], _val.size()).template lpNorm<Eigen::Infinity>();
634 #else
635  Real my_norm = TensorTools::norm_sq((*this)(0));
636 
637  const int N = cast_int<int>(_val.size());
638  for (int i=1; i!=N; i++)
639  {
640  Real current = TensorTools::norm_sq((*this)(i));
641  my_norm = (my_norm > current? my_norm : current);
642  }
643  return sqrt(my_norm);
644 #endif
645 }
std::vector< T > _val
Definition: dense_vector.h:260
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 564 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().

565 {
566  libmesh_assert (this->size());
567  Real my_max = libmesh_real((*this)(0));
568 
569  const int N = cast_int<int>(_val.size());
570  for (int i=1; i!=N; i++)
571  {
572  Real current = libmesh_real((*this)(i));
573  my_max = (my_max > current? my_max : current);
574  }
575  return my_max;
576 }
T libmesh_real(T a)
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
libmesh_assert(j)
std::vector< T > _val
Definition: dense_vector.h:260
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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 546 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().

547 {
548  libmesh_assert (this->size());
549  Real my_min = libmesh_real((*this)(0));
550 
551  const int N = cast_int<int>(_val.size());
552  for (int i=1; i!=N; i++)
553  {
554  Real current = libmesh_real((*this)(i));
555  my_min = (my_min < current? my_min : current);
556  }
557  return my_min;
558 }
T libmesh_real(T a)
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
libmesh_assert(j)
std::vector< T > _val
Definition: dense_vector.h:260
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

Tests if vec is not exactly equal to this vector.

Definition at line 498 of file dense_vector.h.

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

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

499 {
500  libmesh_assert_equal_to (this->size(), vec.size());
501 
502  const int N = cast_int<int>(_val.size());
503  for (int i=0; i<N; i++)
504  if ((*this)(i) != vec(i))
505  return true;
506 
507  return false;
508 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T >
const T & libMesh::DenseVector< T >::operator() ( const unsigned int  i) const
inline
Returns
the (i) element of the vector as a const reference.

Definition at line 373 of file dense_vector.h.

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

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

374 {
375  libmesh_assert_less (i, _val.size());
376 
377  return _val[i];
378 }
std::vector< T > _val
Definition: dense_vector.h:260
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 384 of file dense_vector.h.

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

385 {
386  libmesh_assert_less (i, _val.size());
387 
388  return _val[i];
389 }
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T>
DenseVector< T > & libMesh::DenseVector< T >::operator*= ( const T  factor)
inline

Multiplies every element in the vector by factor.

Definition at line 406 of file dense_vector.h.

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

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

407 {
408  this->scale(factor);
409  return *this;
410 }
void scale(const T factor)
Definition: dense_vector.h:395
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 515 of file dense_vector.h.

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

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

516 {
517  libmesh_assert_equal_to (this->size(), vec.size());
518 
519  const int N = cast_int<int>(_val.size());
520  for (int i=0; i<N; i++)
521  (*this)(i) += vec(i);
522 
523  return *this;
524 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
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 531 of file dense_vector.h.

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

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

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

Assignment operator.

Definition at line 310 of file dense_vector.h.

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

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

311 {
312  const std::vector<T2> & other_vals = other_vector.get_values();
313 
314  _val.clear();
315 
316  const int N = cast_int<int>(other_vals.size());
317  _val.reserve(N);
318 
319  for (int i=0; i<N; i++)
320  _val.push_back(other_vals[i]);
321 
322  return *this;
323 }
std::vector< T > _val
Definition: dense_vector.h:260
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 481 of file dense_vector.h.

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

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

482 {
483  libmesh_assert_equal_to (this->size(), vec.size());
484 
485  const int N = cast_int<int>(_val.size());
486  for (int i=0; i<N; i++)
487  if ((*this)(i) != vec(i))
488  return false;
489 
490  return true;
491 }
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:81
std::vector< T > _val
Definition: dense_vector.h:260
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 338 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::CompositeFEMFunction< Output >::operator()(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::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().

339 {
340  _val.resize(n);
341 
342  zero();
343 }
virtual void zero() libmesh_override
Definition: dense_vector.h:362
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T>
void libMesh::DenseVector< T >::scale ( const T  factor)
inline

Multiplies every element in the vector by factor.

Definition at line 395 of file dense_vector.h.

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

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

396 {
397  const int N = cast_int<int>(_val.size());
398  for (int i=0; i<N; i++)
399  _val[i] *= factor;
400 }
std::vector< T > _val
Definition: dense_vector.h:260
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 81 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::DenseVector< T >::dot(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::DenseVector< T >::get_principal_subvector(), 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::CompositeFEMFunction< Output >::operator()(), libMesh::ParsedFEMFunction< Output >::operator()(), libMesh::CompositeFunction< 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().

82  {
83  return cast_int<unsigned int>(_val.size());
84  }
std::vector< T > _val
Definition: dense_vector.h:260
template<typename T>
void libMesh::DenseVector< T >::swap ( DenseVector< T > &  other_vector)
inline

STL-like swap method

Definition at line 329 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().

330 {
331  _val.swap(other_vector._val);
332 }
std::vector< T > _val
Definition: dense_vector.h:260

Member Data Documentation


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