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 (DenseVector &&)=default
 
 DenseVector (const DenseVector &)=default
 
DenseVectoroperator= (const DenseVector &)=default
 
DenseVectoroperator= (DenseVector &&)=default
 
virtual ~DenseVector ()=default
 
virtual unsigned int size () const override
 
virtual bool empty () const override
 
virtual void zero () override
 
const T & operator() (const unsigned int i) const
 
T & operator() (const unsigned int i)
 
virtual T el (const unsigned int i) const override
 
virtual T & el (const unsigned int i) 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

◆ DenseVector() [1/6]

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

Constructor. Creates a dense vector of dimension n.

Definition at line 277 of file dense_vector.h.

277  :
278  _val (n, T(0.))
279 {
280 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ DenseVector() [2/6]

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

286  :
287  _val (n, val)
288 {
289 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ DenseVector() [3/6]

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

Copy-constructor.

Definition at line 296 of file dense_vector.h.

296  :
297  DenseVectorBase<T>()
298 {
299  const std::vector<T2> & other_vals = other_vector.get_values();
300 
301  _val.clear();
302 
303  const int N = cast_int<int>(other_vals.size());
304  _val.reserve(N);
305 
306  for (int i=0; i<N; i++)
307  _val.push_back(other_vals[i]);
308 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ DenseVector() [4/6]

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

315  :
316  _val(other_vector)
317 {
318 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ DenseVector() [5/6]

template<typename T>
libMesh::DenseVector< T >::DenseVector ( DenseVector< T > &&  )
default

The 5 special functions can be defaulted for this class, as it does not manage any memory itself.

◆ DenseVector() [6/6]

template<typename T>
libMesh::DenseVector< T >::DenseVector ( const DenseVector< T > &  )
default

◆ ~DenseVector()

template<typename T>
virtual libMesh::DenseVector< T >::~DenseVector ( )
virtualdefault

Member Function Documentation

◆ add()

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

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

438 {
439  libmesh_assert_equal_to (this->size(), vec.size());
440 
441  const int N = cast_int<int>(_val.size());
442  for (int i=0; i<N; i++)
443  (*this)(i) += static_cast<T>(factor)*vec(i);
444 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ append()

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

Referenced by libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values().

368 {
369  const std::vector<T2> & other_vals = other_vector.get_values();
370 
371  _val.reserve(this->size() + other_vals.size());
372  _val.insert(_val.end(), other_vals.begin(), other_vals.end());
373 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ dot()

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

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

◆ el() [1/2]

template<typename T>
virtual T libMesh::DenseVector< T >::el ( const unsigned int  i) const
inlineoverridevirtual
Returns
The (i) element of the vector.

Implements libMesh::DenseVectorBase< T >.

Definition at line 112 of file dense_vector.h.

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

◆ el() [2/2]

template<typename T>
virtual T& libMesh::DenseVector< T >::el ( const unsigned int  i)
inlineoverridevirtual
Returns
The (i) element of the vector as a writable reference.

Implements libMesh::DenseVectorBase< T >.

Definition at line 115 of file dense_vector.h.

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

◆ empty()

template<typename T>
virtual bool libMesh::DenseVector< T >::empty ( ) const
inlineoverridevirtual
Returns
true iff size() is 0.

Reimplemented from libMesh::DenseVectorBase< T >.

Definition at line 97 of file dense_vector.h.

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

98  { return _val.empty(); }
std::vector< T > _val
Definition: dense_vector.h:268

◆ get_principal_subvector()

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

674 {
675  libmesh_assert_less_equal ( sub_n, this->size() );
676 
677  dest.resize(sub_n);
678  const int N = cast_int<int>(sub_n);
679  for (int i=0; i<N; i++)
680  dest(i) = _val[i];
681 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ get_values() [1/2]

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

Referenced by libMesh::DenseMatrix< Number >::_evd_lapack(), libMesh::DenseMatrix< Number >::_lu_back_substitute_lapack(), libMesh::DenseMatrix< Number >::_matvec_blas(), libMesh::DenseMatrix< Number >::_svd_solve_lapack(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), and libMesh::FEMContext::pre_fe_reinit().

256 { return _val; }
std::vector< T > _val
Definition: dense_vector.h:268

◆ get_values() [2/2]

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

261 { return _val; }
std::vector< T > _val
Definition: dense_vector.h:268

◆ indefinite_dot()

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

485 {
486  libmesh_assert_equal_to (this->size(), vec.size());
487 
488  typename CompareTypes<T, T2>::supertype val = 0.;
489 
490  const int N = cast_int<int>(_val.size());
491  for (int i=0; i<N; i++)
492  val += (*this)(i)*(vec(i));
493 
494  return val;
495 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ l1_norm()

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

602 {
603  if (!_val.size())
604  return 0.;
605 
606 #ifdef LIBMESH_HAVE_EIGEN
607  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).template lpNorm<1>();
608 #else
609  Real my_norm = 0.;
610  const int N = cast_int<int>(_val.size());
611  for (int i=0; i!=N; i++)
612  my_norm += std::abs((*this)(i));
613 
614  return my_norm;
615 #endif
616 }
double abs(double a)
std::vector< T > _val
Definition: dense_vector.h:268
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ l2_norm()

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

623 {
624  if (!_val.size())
625  return 0.;
626 
627 #ifdef LIBMESH_HAVE_EIGEN
628  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).norm();
629 #else
630  Real my_norm = 0.;
631  const int N = cast_int<int>(_val.size());
632  // The following pragma tells clang's vectorizer that it is safe to
633  // reorder floating point operations for this loop.
634 #ifdef __clang__
635 #pragma clang loop vectorize(enable)
636 #endif
637  for (int i=0; i!=N; i++)
638  my_norm += TensorTools::norm_sq((*this)(i));
639 
640  return sqrt(my_norm);
641 #endif
642 }
std::vector< T > _val
Definition: dense_vector.h:268
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ linfty_norm()

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

649 {
650  if (!_val.size())
651  return 0.;
652 
653 #ifdef LIBMESH_HAVE_EIGEN
654  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).template lpNorm<Eigen::Infinity>();
655 #else
656  Real my_norm = TensorTools::norm_sq((*this)(0));
657 
658  const int N = cast_int<int>(_val.size());
659  for (int i=1; i!=N; i++)
660  {
661  Real current = TensorTools::norm_sq((*this)(i));
662  my_norm = (my_norm > current? my_norm : current);
663  }
664  return sqrt(my_norm);
665 #endif
666 }
std::vector< T > _val
Definition: dense_vector.h:268
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ max()

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

584 {
585  libmesh_assert (this->size());
586  Real my_max = libmesh_real((*this)(0));
587 
588  const int N = cast_int<int>(_val.size());
589  for (int i=1; i!=N; i++)
590  {
591  Real current = libmesh_real((*this)(i));
592  my_max = (my_max > current? my_max : current);
593  }
594  return my_max;
595 }
T libmesh_real(T a)
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ min()

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

566 {
567  libmesh_assert (this->size());
568  Real my_min = libmesh_real((*this)(0));
569 
570  const int N = cast_int<int>(_val.size());
571  for (int i=1; i!=N; i++)
572  {
573  Real current = libmesh_real((*this)(i));
574  my_min = (my_min < current? my_min : current);
575  }
576  return my_min;
577 }
T libmesh_real(T a)
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ operator!=()

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

518 {
519  libmesh_assert_equal_to (this->size(), vec.size());
520 
521  const int N = cast_int<int>(_val.size());
522  for (int i=0; i<N; i++)
523  if ((*this)(i) != vec(i))
524  return true;
525 
526  return false;
527 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ operator()() [1/2]

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

391 {
392  libmesh_assert_less (i, _val.size());
393 
394  return _val[i];
395 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ operator()() [2/2]

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

402 {
403  libmesh_assert_less (i, _val.size());
404 
405  return _val[i];
406 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ operator*=()

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

424 {
425  this->scale(factor);
426  return *this;
427 }
void scale(const T factor)
Definition: dense_vector.h:412

◆ operator+=()

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

535 {
536  libmesh_assert_equal_to (this->size(), vec.size());
537 
538  const int N = cast_int<int>(_val.size());
539  for (int i=0; i<N; i++)
540  (*this)(i) += vec(i);
541 
542  return *this;
543 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ operator-=()

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

551 {
552  libmesh_assert_equal_to (this->size(), vec.size());
553 
554  const int N = cast_int<int>(_val.size());
555  for (int i=0; i<N; i++)
556  (*this)(i) -= vec(i);
557 
558  return *this;
559 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ operator=() [1/3]

template<typename T>
DenseVector& libMesh::DenseVector< T >::operator= ( const DenseVector< T > &  )
default

◆ operator=() [2/3]

template<typename T>
DenseVector& libMesh::DenseVector< T >::operator= ( DenseVector< T > &&  )
default

◆ operator=() [3/3]

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

328 {
329  const std::vector<T2> & other_vals = other_vector.get_values();
330 
331  _val.clear();
332 
333  const int N = cast_int<int>(other_vals.size());
334  _val.reserve(N);
335 
336  for (int i=0; i<N; i++)
337  _val.push_back(other_vals[i]);
338 
339  return *this;
340 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ operator==()

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

501 {
502  libmesh_assert_equal_to (this->size(), vec.size());
503 
504  const int N = cast_int<int>(_val.size());
505  for (int i=0; i<N; i++)
506  if ((*this)(i) != vec(i))
507  return false;
508 
509  return true;
510 }
virtual unsigned int size() const override
Definition: dense_vector.h:92
std::vector< T > _val
Definition: dense_vector.h:268

◆ print()

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

Pretty-print the vector to stdout.

Definition at line 51 of file dense_vector_base.C.

52 {
53  for (auto i : IntRange<int>(0, this->size()))
54  os << std::setw(8)
55  << this->el(i)
56  << std::endl;
57 }
virtual unsigned int size() const =0
virtual T el(const unsigned int i) const =0

◆ print_scientific()

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 31 of file dense_vector_base.C.

32 {
33  // save the initial format flags
34  std::ios_base::fmtflags os_flags = os.flags();
35 
36  // Print the vector entries.
37  for (auto i : IntRange<int>(0, this->size()))
38  os << std::setw(10)
39  << std::scientific
40  << std::setprecision(precision)
41  << this->el(i)
42  << std::endl;
43 
44  // reset the original format flags
45  os.flags(os_flags);
46 }
virtual unsigned int size() const =0
virtual T el(const unsigned int i) const =0

◆ resize()

template<typename T >
void libMesh::DenseVector< T >::resize ( const unsigned int  n)
inline

Resize the vector. Sets all elements to 0.

Definition at line 355 of file dense_vector.h.

Referenced by libMesh::DenseMatrix< Number >::_cholesky_back_substitute(), libMesh::DenseMatrix< Number >::_evd_lapack(), libMesh::DenseMatrix< Number >::_lu_back_substitute(), libMesh::DenseMatrix< Number >::_svd_lapack(), libMesh::DenseMatrix< Number >::_svd_solve_lapack(), libMesh::HPCoarsenTest::add_projection(), libMesh::DofMap::build_constraint_matrix_and_vector(), libMesh::FEGenericBase< FEOutputType< T >::type >::coarsened_dof_values(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_periodic_constraints(), libMesh::FEGenericBase< FEOutputType< T >::type >::compute_proj_constraints(), libMesh::DGFEMContext::neighbor_side_fe_reinit(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::operator()(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::CompositeFEMFunction< Output >::operator()(), libMesh::CompositeFunction< Output >::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::MeshFunction::operator()(), libMesh::BoundaryProjectSolution::operator()(), libMesh::MeshTools::Generation::Private::GaussLobattoRedistributionFunction::operator()(), libMesh::FEMContext::pre_fe_reinit(), libMesh::DenseVector< Output >::resize(), libMesh::HPCoarsenTest::select_refinement(), libMesh::DenseMatrix< Number >::vector_mult(), libMesh::DenseMatrix< Number >::vector_mult_add(), and libMesh::DenseMatrix< Number >::vector_mult_transpose().

356 {
357  _val.resize(n);
358 
359  zero();
360 }
std::vector< T > _val
Definition: dense_vector.h:268
virtual void zero() override
Definition: dense_vector.h:379

◆ scale()

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

Multiplies every element in the vector by factor.

Definition at line 412 of file dense_vector.h.

413 {
414  const int N = cast_int<int>(_val.size());
415  for (int i=0; i<N; i++)
416  _val[i] *= factor;
417 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ size()

template<typename T>
virtual unsigned int libMesh::DenseVector< T >::size ( ) const
inlineoverridevirtual
Returns
The size of the vector.

Implements libMesh::DenseVectorBase< T >.

Definition at line 92 of file dense_vector.h.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::DenseMatrix< Number >::_lu_back_substitute(), libMesh::DenseMatrix< Number >::_matvec_blas(), libMesh::DenseMatrix< Number >::_svd_lapack(), libMesh::DenseMatrix< Number >::_svd_solve_lapack(), libMesh::HPCoarsenTest::add_projection(), libMesh::NumericVector< Number >::add_vector(), 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::JumpErrorEstimator::estimate_error(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::MeshFunction::gradient(), libMesh::MeshFunction::hessian(), libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector(), libMesh::DofMap::heterogenously_constrain_element_vector(), libMesh::VectorSetAction< Val >::insert(), libMesh::NumericVector< Number >::insert(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::ConstFunction< Output >::operator()(), libMesh::ConstFEMFunction< Output >::operator()(), libMesh::WrappedFunction< Output >::operator()(), libMesh::ParsedFEMFunction< T >::operator()(), libMesh::CompositeFEMFunction< Output >::operator()(), libMesh::CompositeFunction< Output >::operator()(), libMesh::MeshFunction::operator()(), libMesh::DenseMatrix< Number >::outer_product(), libMesh::HPCoarsenTest::select_refinement(), libMesh::DenseMatrix< Number >::vector_mult(), libMesh::DenseMatrix< Number >::vector_mult_add(), and libMesh::DenseMatrix< Number >::vector_mult_transpose().

93  {
94  return cast_int<unsigned int>(_val.size());
95  }
std::vector< T > _val
Definition: dense_vector.h:268

◆ swap()

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

STL-like swap method

Definition at line 346 of file dense_vector.h.

Referenced by libMesh::EulerSolver::_general_residual(), and libMesh::Euler2Solver::_general_residual().

347 {
348  _val.swap(other_vector._val);
349 }
std::vector< T > _val
Definition: dense_vector.h:268

◆ zero()

template<typename T >
void libMesh::DenseVector< T >::zero ( )
inlineoverridevirtual

Member Data Documentation

◆ _val

template<typename T>
std::vector<T> libMesh::DenseVector< T >::_val
private

The actual data values, stored as a 1D array.

Definition at line 268 of file dense_vector.h.

Referenced by libMesh::DenseVector< Output >::empty(), libMesh::DenseVector< Output >::get_values(), and libMesh::DenseVector< Output >::size().


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