20 #ifndef LIBMESH_DENSE_VECTOR_H
21 #define LIBMESH_DENSE_VECTOR_H
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/dense_vector_base.h"
26 #include "libmesh/compare_types.h"
27 #include "libmesh/int_range.h"
28 #include "libmesh/tensor_tools.h"
30 #ifdef LIBMESH_HAVE_EIGEN
31 #include "libmesh/ignore_warnings.h"
33 #include "libmesh/restore_warnings.h"
36 #ifdef LIBMESH_HAVE_METAPHYSICL
37 #include "metaphysicl/raw_type.h"
42 #include <initializer_list>
56 class DenseVector :
public DenseVectorBase<T>
83 template <
typename T2>
91 template <
typename T2>
99 template <
typename T2>
111 virtual unsigned int size() const override final
113 return cast_int<unsigned int>(
_val.size());
116 virtual bool empty() const override final
117 {
return _val.empty(); }
119 virtual void zero() override final;
126 const T & operator() (const
unsigned int i) const;
133 T & operator() (const
unsigned int i);
135 virtual T
el(const
unsigned int i) const override final
136 {
return (*
this)(i); }
138 virtual T &
el(
const unsigned int i)
override final
139 {
return (*
this)(i); }
148 template <
typename T2>
163 void resize (
const unsigned int n);
170 template <
typename T2>
178 void scale (
const T factor);
198 template <
typename T2,
typename T3>
199 typename boostcopy::enable_if_c<
200 ScalarTraits<T2>::value,
void >::type
201 add (
const T2 factor,
211 template <
typename T2>
221 template <
typename T2>
229 template <
typename T2>
237 template <
typename T2>
247 template <
typename T2>
257 template <
typename T2>
336 template<
typename T2>
341 const std::vector<T2> & other_vals = other_vector.get_values();
345 const int N = cast_int<int>(other_vals.size());
348 for (
int i=0; i<N; i++)
349 _val.push_back(other_vals[i]);
355 template<
typename T2>
364 template <
typename T2>
367 _val(init_list.begin(), init_list.end())
374 template<
typename T2>
378 const std::vector<T2> & other_vals = other_vector.get_values();
382 const int N = cast_int<int>(other_vals.size());
385 for (
int i=0; i<N; i++)
386 _val.push_back(other_vals[i]);
397 _val.swap(other_vector._val);
414 template<
typename T2>
418 const std::vector<T2> & other_vals = other_vector.get_values();
420 _val.reserve(this->size() + other_vals.size());
421 _val.insert(_val.end(), other_vals.begin(), other_vals.end());
430 std::fill (_val.begin(),
441 libmesh_assert_less (i, _val.size());
452 libmesh_assert_less (i, _val.size());
463 const int N = cast_int<int>(_val.size());
464 for (
int i=0; i<N; i++)
481 template<
typename T2,
typename T3>
483 typename boostcopy::enable_if_c<
484 ScalarTraits<T2>::value,
void >::type
488 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 (*
this)(i) += static_cast<T>(factor)*vec(i);
498 template<
typename T2>
505 libmesh_assert_equal_to (this->size(), vec.size());
507 #ifdef LIBMESH_HAVE_EIGEN
512 return Eigen::Map<const typename Eigen::Matrix<T2, Eigen::Dynamic, 1>>(vec.get_values().data(), vec.size())
513 .dot(Eigen::Map<
const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()));
515 typename CompareTypes<T, T2>::supertype val = 0.;
517 const int N = cast_int<int>(_val.size());
521 #pragma clang loop vectorize(enable)
523 for (
int i=0; i<N; i++)
531 template<
typename T2>
535 libmesh_assert_equal_to (this->size(), vec.size());
537 typename CompareTypes<T, T2>::supertype val = 0.;
539 const int N = cast_int<int>(_val.size());
540 for (
int i=0; i<N; i++)
541 val += (*
this)(i)*(vec(i));
547 template<
typename T2>
551 libmesh_assert_equal_to (this->size(), vec.size());
553 const int N = cast_int<int>(_val.size());
554 for (
int i=0; i<N; i++)
555 if ((*
this)(i) != vec(i))
564 template<
typename T2>
568 libmesh_assert_equal_to (this->size(), vec.size());
570 const int N = cast_int<int>(_val.size());
571 for (
int i=0; i<N; i++)
572 if ((*
this)(i) != vec(i))
581 template<
typename T2>
585 libmesh_assert_equal_to (this->size(), vec.size());
587 const int N = cast_int<int>(_val.size());
588 for (
int i=0; i<N; i++)
589 (*
this)(i) += vec(i);
597 template<
typename T2>
601 libmesh_assert_equal_to (this->size(), vec.size());
603 const int N = cast_int<int>(_val.size());
604 for (
int i=0; i<N; i++)
605 (*
this)(i) -= vec(i);
616 libmesh_assert (this->size());
619 const int N = cast_int<int>(_val.size());
620 for (
int i=1; i!=N; i++)
623 my_min = (my_min < current? my_min : current);
634 libmesh_assert (this->size());
637 const int N = cast_int<int>(_val.size());
638 for (
int i=1; i!=N; i++)
641 my_max = (my_max > current? my_max : current);
655 #ifdef LIBMESH_HAVE_EIGEN
656 return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).
template lpNorm<1>();
659 const int N = cast_int<int>(_val.size());
660 for (
int i=0; i!=N; i++)
676 #ifdef LIBMESH_HAVE_EIGEN
677 return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).
norm();
680 const int N = cast_int<int>(_val.size());
684 #pragma clang loop vectorize(enable)
686 for (
int i=0; i!=N; i++)
689 return sqrt(my_norm);
702 #ifdef LIBMESH_HAVE_EIGEN
703 return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).
template lpNorm<Eigen::Infinity>();
707 const int N = cast_int<int>(_val.size());
708 for (
int i=1; i!=N; i++)
711 my_norm = (my_norm > current? my_norm : current);
713 return sqrt(my_norm);
724 libmesh_assert_less_equal ( sub_n, this->size() );
727 const int N = cast_int<int>(sub_n);
728 for (
int i=0; i<N; i++)
734 #ifdef LIBMESH_HAVE_METAPHYSICL
735 namespace MetaPhysicL
737 template <
typename T>
738 struct RawType<libMesh::DenseVector<T>>
744 const auto s = in.
size();
746 for (
unsigned int i = 0; i < s; ++i)
747 ret(i) = raw_value(in(i));
755 #endif // LIBMESH_DENSE_VECTOR_H
DenseVector< T > & operator*=(const T factor)
将向量中的每个元素乘以 factor。
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseVector< T3 > &vec)
将 factor 乘以 vec 添加到该向量。这应该仅在 T += T2 * T3 是有效的 C++ 且 T2 是标量的情况下起作用。 返回类型是 void。
void resize(const unsigned int n)
调整向量的大小。将所有元素设置为0。
CompareTypes< T, T2 >::supertype dot(const DenseVector< T2 > &vec) const
std::vector< T > & get_values()
const T & operator()(const unsigned int i) const
virtual bool empty() const overridefinal
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的平方根。
void scale(const T factor)
将向量中的每个元素乘以 factor。
void swap(DenseVector< T > &other_vector)
STL 风格的交换方法。
virtual T & el(const unsigned int i) overridefinal
DenseVector & operator=(const DenseVector &)=default
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的绝对值。
bool operator==(const DenseVector< T2 > &vec) const
const std::vector< T > & get_values() const
virtual ~DenseVector()=default
std::vector< T > _val
实际数据值,存储为一维数组。
virtual T el(const unsigned int i) const overridefinal
bool operator!=(const DenseVector< T2 > &vec) const
virtual void zero() overridefinal
将向量中的每个元素设置为0。由于派生类中的存储方法可能不同,需要将其声明为纯虚函数。
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CompareTypes< T, T2 >::supertype indefinite_dot(const DenseVector< T2 > &vec) const
DenseVector(const unsigned int n=0)
构造函数。创建一个维数为 n 的稠密向量。
void append(const DenseVector< T2 > &other_vector)
将其他向量的额外条目附加到(调整大小但不改变)向量。
定义用于有限元计算的抽象稠密向量基类。 可以从这个类派生出特定的稠密向量,例如 DenseSubVectors。
定义用于有限元计算的稠密向量类。该类基本上是为了补充 DenseMatrix 类而设计的。 它相对于 std::vector 具有额外的功能,使其在有限元中特别有用,特别是对于方程组。 所有重写的虚拟函...
ADRealEigenVector< T, D, asd > norm(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的范数。
DenseVector< T > & operator+=(const DenseVector< T2 > &vec)
将 vec 添加到该向量。
void get_principal_subvector(unsigned int sub_n, DenseVector< T > &dest) const
将大小为 sub_n(即前 sub_n 个条目)的主子向量放入 dest。
virtual unsigned int size() const overridefinal
DenseVector< T > & operator-=(const DenseVector< T2 > &vec)
从该向量中减去 vec。