20 #ifndef LIBMESH_DENSE_VECTOR_H 21 #define LIBMESH_DENSE_VECTOR_H 29 #ifdef LIBMESH_HAVE_EIGEN 52 class DenseVector :
public DenseVectorBase<T>
73 template <
typename T2>
79 template <
typename T2>
92 virtual unsigned int size()
const override 94 return cast_int<unsigned int>(
_val.size());
97 virtual bool empty()
const override 98 {
return _val.empty(); }
100 virtual void zero()
override;
105 const T &
operator() (
const unsigned int i)
const;
112 virtual T
el(
const unsigned int i)
const override 113 {
return (*
this)(i); }
115 virtual T &
el(
const unsigned int i)
override 116 {
return (*
this)(i); }
123 template <
typename T2>
134 void resize (
const unsigned int n);
140 template <
typename T2>
146 void scale (
const T factor);
162 template <
typename T2,
typename T3>
165 add (
const T2 factor,
173 template <
typename T2>
182 template <
typename T2>
188 template <
typename T2>
194 template <
typename T2>
202 template <
typename T2>
210 template <
typename T2>
294 template<
typename T2>
299 const std::vector<T2> & other_vals = other_vector.get_values();
303 const int N = cast_int<int>(other_vals.size());
306 for (
int i=0; i<N; i++)
307 _val.push_back(other_vals[i]);
313 template<
typename T2>
325 template<
typename T2>
329 const std::vector<T2> & other_vals = other_vector.get_values();
333 const int N = cast_int<int>(other_vals.size());
336 for (
int i=0; i<N; i++)
337 _val.push_back(other_vals[i]);
348 _val.swap(other_vector._val);
365 template<
typename T2>
369 const std::vector<T2> & other_vals = other_vector.get_values();
371 _val.reserve(this->size() + other_vals.size());
372 _val.insert(_val.end(), other_vals.begin(), other_vals.end());
381 std::fill (_val.begin(),
392 libmesh_assert_less (i, _val.size());
403 libmesh_assert_less (i, _val.size());
414 const int N = cast_int<int>(_val.size());
415 for (
int i=0; i<N; i++)
432 template<
typename T2,
typename T3>
439 libmesh_assert_equal_to (this->size(), vec.size());
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);
449 template<
typename T2>
456 libmesh_assert_equal_to (this->size(), vec.size());
458 #ifdef LIBMESH_HAVE_EIGEN 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()));
468 const int N = cast_int<int>(_val.size());
472 #pragma clang loop vectorize(enable) 474 for (
int i=0; i<N; i++)
482 template<
typename T2>
486 libmesh_assert_equal_to (this->size(), vec.size());
490 const int N = cast_int<int>(_val.size());
491 for (
int i=0; i<N; i++)
492 val += (*
this)(i)*(vec(i));
498 template<
typename T2>
502 libmesh_assert_equal_to (this->size(), vec.size());
504 const int N = cast_int<int>(_val.size());
505 for (
int i=0; i<N; i++)
506 if ((*
this)(i) != vec(i))
515 template<
typename T2>
519 libmesh_assert_equal_to (this->size(), vec.size());
521 const int N = cast_int<int>(_val.size());
522 for (
int i=0; i<N; i++)
523 if ((*
this)(i) != vec(i))
532 template<
typename T2>
536 libmesh_assert_equal_to (this->size(), vec.size());
538 const int N = cast_int<int>(_val.size());
539 for (
int i=0; i<N; i++)
540 (*
this)(i) += vec(i);
548 template<
typename T2>
552 libmesh_assert_equal_to (this->size(), vec.size());
554 const int N = cast_int<int>(_val.size());
555 for (
int i=0; i<N; i++)
556 (*
this)(i) -= vec(i);
567 libmesh_assert (this->size());
570 const int N = cast_int<int>(_val.size());
571 for (
int i=1; i!=N; i++)
574 my_min = (my_min < current? my_min : current);
585 libmesh_assert (this->size());
588 const int N = cast_int<int>(_val.size());
589 for (
int i=1; i!=N; i++)
592 my_max = (my_max > current? my_max : current);
606 #ifdef LIBMESH_HAVE_EIGEN 607 return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).
template lpNorm<1>();
610 const int N = cast_int<int>(_val.size());
611 for (
int i=0; i!=N; i++)
627 #ifdef LIBMESH_HAVE_EIGEN 628 return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).norm();
631 const int N = cast_int<int>(_val.size());
635 #pragma clang loop vectorize(enable) 637 for (
int i=0; i!=N; i++)
640 return sqrt(my_norm);
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>();
658 const int N = cast_int<int>(_val.size());
659 for (
int i=1; i!=N; i++)
662 my_norm = (my_norm > current? my_norm : current);
664 return sqrt(my_norm);
675 libmesh_assert_less_equal ( sub_n, this->size() );
678 const int N = cast_int<int>(sub_n);
679 for (
int i=0; i<N; i++)
685 #endif // LIBMESH_DENSE_VECTOR_H
virtual unsigned int size() const override
DenseVector< T > & operator*=(const T factor)
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseVector< T3 > &vec)
void resize(const unsigned int n)
std::vector< T > & get_values()
bool operator==(const DenseVector< T2 > &vec) const
virtual T el(const unsigned int i) const override
void scale(const T factor)
long double max(long double a, double b)
void get_principal_subvector(unsigned int sub_n, DenseVector< T > &dest) const
const std::vector< T > & get_values() const
virtual bool empty() const override
void swap(DenseVector< T > &other_vector)
DenseVector & operator=(const DenseVector &)=default
bool operator!=(const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &x, const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &y)
virtual ~DenseVector()=default
virtual T & el(const unsigned int i) override
bool operator!=(const DenseVector< T2 > &vec) const
const T & operator()(const unsigned int i) const
CompareTypes< T, T2 >::supertype indefinite_dot(const DenseVector< T2 > &vec) const
CompareTypes< T, T2 >::supertype dot(const DenseVector< T2 > &vec) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void swap(Iterator &lhs, Iterator &rhs)
DenseVector(const unsigned int n=0)
bool operator==(const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &x, const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &y)
void append(const DenseVector< T2 > &other_vector)
Iterator & operator=(const Iterator &rhs)
DenseVector< T > & operator+=(const DenseVector< T2 > &vec)
long double min(long double a, double b)
DenseVector< T > & operator-=(const DenseVector< T2 > &vec)
virtual void zero() override