libmesh解析
本工作只是尝试解析原libmesh的代码,供学习使用
 全部  命名空间 文件 函数 变量 类型定义 枚举 枚举值 友元 
dense_vector.h
浏览该文件的文档.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2023 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_DENSE_VECTOR_H
21 #define LIBMESH_DENSE_VECTOR_H
22 
23 // Local Includes
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"
29 
30 #ifdef LIBMESH_HAVE_EIGEN
31 #include "libmesh/ignore_warnings.h"
32 #include <Eigen/Core>
33 #include "libmesh/restore_warnings.h"
34 #endif
35 
36 #ifdef LIBMESH_HAVE_METAPHYSICL
37 #include "metaphysicl/raw_type.h"
38 #endif
39 
40 // C++ includes
41 #include <vector>
42 #include <initializer_list>
43 
44 namespace libMesh
45 {
46 
55 template<typename T>
56 class DenseVector : public DenseVectorBase<T>
57 {
58 public:
59 
65  explicit
66  DenseVector(const unsigned int n=0);
67 
74  explicit
75  DenseVector(const unsigned int n,
76  const T & val);
77 
83  template <typename T2>
84  DenseVector (const DenseVector<T2> & other_vector);
85 
91  template <typename T2>
92  DenseVector (const std::vector<T2> & other_vector);
93 
99  template <typename T2>
100  DenseVector (std::initializer_list<T2> init_list);
101 
105  DenseVector (DenseVector &&) = default;
106  DenseVector (const DenseVector &) = default;
107  DenseVector & operator= (const DenseVector &) = default;
108  DenseVector & operator= (DenseVector &&) = default;
109  virtual ~DenseVector() = default;
110 
111  virtual unsigned int size() const override final
112  {
113  return cast_int<unsigned int>(_val.size());
114  }
115 
116  virtual bool empty() const override final
117  { return _val.empty(); }
118 
119  virtual void zero() override final;
120 
126  const T & operator() (const unsigned int i) const;
127 
133  T & operator() (const unsigned int i);
134 
135  virtual T el(const unsigned int i) const override final
136  { return (*this)(i); }
137 
138  virtual T & el(const unsigned int i) override final
139  { return (*this)(i); }
140 
148  template <typename T2>
149  DenseVector<T> & operator = (const DenseVector<T2> & other_vector);
150 
156  void swap(DenseVector<T> & other_vector);
157 
163  void resize (const unsigned int n);
164 
170  template <typename T2>
171  void append (const DenseVector<T2> & other_vector);
172 
178  void scale (const T factor);
179 
187  DenseVector<T> & operator*= (const T factor);
188 
198  template <typename T2, typename T3>
199  typename boostcopy::enable_if_c<
200  ScalarTraits<T2>::value, void >::type
201  add (const T2 factor,
202  const DenseVector<T3> & vec);
203 
211  template <typename T2>
212  typename CompareTypes<T, T2>::supertype dot (const DenseVector<T2> & vec) const;
213 
221  template <typename T2>
222  typename CompareTypes<T, T2>::supertype indefinite_dot (const DenseVector<T2> & vec) const;
223 
229  template <typename T2>
230  bool operator== (const DenseVector<T2> & vec) const;
231 
237  template <typename T2>
238  bool operator!= (const DenseVector<T2> & vec) const;
239 
247  template <typename T2>
249 
257  template <typename T2>
259 
263  Real min () const;
264 
268  Real max () const;
269 
273  Real l1_norm () const;
274 
278  Real l2_norm () const;
279 
283  Real linfty_norm () const;
284 
291  void get_principal_subvector (unsigned int sub_n, DenseVector<T> & dest) const;
292 
298  std::vector<T> & get_values() { return _val; }
299 
303  const std::vector<T> & get_values() const { return _val; }
304 
305 private:
306 
310  std::vector<T> _val;
311 };
312 
313 
314 
315 // ------------------------------------------------------------
316 // DenseVector member functions
317 template<typename T>
318 inline
319 DenseVector<T>::DenseVector(const unsigned int n) :
320  _val (n, T{})
321 {
322 }
323 
324 
325 template<typename T>
326 inline
327 DenseVector<T>::DenseVector(const unsigned int n,
328  const T & val) :
329  _val (n, val)
330 {
331 }
332 
333 
334 
335 template<typename T>
336 template<typename T2>
337 inline
339  DenseVectorBase<T>()
340 {
341  const std::vector<T2> & other_vals = other_vector.get_values();
342 
343  _val.clear();
344 
345  const int N = cast_int<int>(other_vals.size());
346  _val.reserve(N);
347 
348  for (int i=0; i<N; i++)
349  _val.push_back(other_vals[i]);
350 }
351 
352 
353 
354 template<typename T>
355 template<typename T2>
356 inline
357 DenseVector<T>::DenseVector (const std::vector<T2> & other_vector) :
358  _val(other_vector)
359 {
360 }
361 
362 
363 template<typename T>
364 template <typename T2>
365 inline
366 DenseVector<T>::DenseVector (std::initializer_list<T2> init_list) :
367  _val(init_list.begin(), init_list.end())
368 {
369 }
370 
371 
372 
373 template<typename T>
374 template<typename T2>
375 inline
377 {
378  const std::vector<T2> & other_vals = other_vector.get_values();
379 
380  _val.clear();
381 
382  const int N = cast_int<int>(other_vals.size());
383  _val.reserve(N);
384 
385  for (int i=0; i<N; i++)
386  _val.push_back(other_vals[i]);
387 
388  return *this;
389 }
390 
391 
392 
393 template<typename T>
394 inline
396 {
397  _val.swap(other_vector._val);
398 }
399 
400 
401 
402 template<typename T>
403 inline
404 void DenseVector<T>::resize(const unsigned int n)
405 {
406  _val.resize(n);
407 
408  zero();
409 }
410 
411 
412 
413 template<typename T>
414 template<typename T2>
415 inline
416 void DenseVector<T>::append (const DenseVector<T2> & other_vector)
417 {
418  const std::vector<T2> & other_vals = other_vector.get_values();
419 
420  _val.reserve(this->size() + other_vals.size());
421  _val.insert(_val.end(), other_vals.begin(), other_vals.end());
422 }
423 
424 
425 
426 template<typename T>
427 inline
429 {
430  std::fill (_val.begin(),
431  _val.end(),
432  T{});
433 }
434 
435 
436 
437 template<typename T>
438 inline
439 const T & DenseVector<T>::operator () (const unsigned int i) const
440 {
441  libmesh_assert_less (i, _val.size());
442 
443  return _val[i];
444 }
445 
446 
447 
448 template<typename T>
449 inline
450 T & DenseVector<T>::operator () (const unsigned int i)
451 {
452  libmesh_assert_less (i, _val.size());
453 
454  return _val[i];
455 }
456 
457 
458 
459 template<typename T>
460 inline
461 void DenseVector<T>::scale (const T factor)
462 {
463  const int N = cast_int<int>(_val.size());
464  for (int i=0; i<N; i++)
465  _val[i] *= factor;
466 }
467 
468 
469 
470 template<typename T>
471 inline
473 {
474  this->scale(factor);
475  return *this;
476 }
477 
478 
479 
480 template<typename T>
481 template<typename T2, typename T3>
482 inline
483 typename boostcopy::enable_if_c<
484  ScalarTraits<T2>::value, void >::type
485 DenseVector<T>::add (const T2 factor,
486  const DenseVector<T3> & vec)
487 {
488  libmesh_assert_equal_to (this->size(), vec.size());
489 
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);
493 }
494 
495 
496 
497 template<typename T>
498 template<typename T2>
499 inline
500 typename CompareTypes<T, T2>::supertype DenseVector<T>::dot (const DenseVector<T2> & vec) const
501 {
502  if (!_val.size())
503  return 0.;
504 
505  libmesh_assert_equal_to (this->size(), vec.size());
506 
507 #ifdef LIBMESH_HAVE_EIGEN
508  // We reverse the order of the arguments to dot() here since
509  // the convention in Eigen is to take the complex conjugate of the
510  // *first* argument, while ours is to take the complex conjugate of
511  // the second.
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()));
514 #else
515  typename CompareTypes<T, T2>::supertype val = 0.;
516 
517  const int N = cast_int<int>(_val.size());
518  // The following pragma tells clang's vectorizer that it is safe to
519  // reorder floating point operations for this loop.
520 #ifdef __clang__
521 #pragma clang loop vectorize(enable)
522 #endif
523  for (int i=0; i<N; i++)
524  val += (*this)(i)*libmesh_conj(vec(i));
525 
526  return val;
527 #endif
528 }
529 
530 template<typename T>
531 template<typename T2>
532 inline
533 typename CompareTypes<T, T2>::supertype DenseVector<T>::indefinite_dot (const DenseVector<T2> & vec) const
534 {
535  libmesh_assert_equal_to (this->size(), vec.size());
536 
537  typename CompareTypes<T, T2>::supertype val = 0.;
538 
539  const int N = cast_int<int>(_val.size());
540  for (int i=0; i<N; i++)
541  val += (*this)(i)*(vec(i));
542 
543  return val;
544 }
545 
546 template<typename T>
547 template<typename T2>
548 inline
550 {
551  libmesh_assert_equal_to (this->size(), vec.size());
552 
553  const int N = cast_int<int>(_val.size());
554  for (int i=0; i<N; i++)
555  if ((*this)(i) != vec(i))
556  return false;
557 
558  return true;
559 }
560 
561 
562 
563 template<typename T>
564 template<typename T2>
565 inline
567 {
568  libmesh_assert_equal_to (this->size(), vec.size());
569 
570  const int N = cast_int<int>(_val.size());
571  for (int i=0; i<N; i++)
572  if ((*this)(i) != vec(i))
573  return true;
574 
575  return false;
576 }
577 
578 
579 
580 template<typename T>
581 template<typename T2>
582 inline
584 {
585  libmesh_assert_equal_to (this->size(), vec.size());
586 
587  const int N = cast_int<int>(_val.size());
588  for (int i=0; i<N; i++)
589  (*this)(i) += vec(i);
590 
591  return *this;
592 }
593 
594 
595 
596 template<typename T>
597 template<typename T2>
598 inline
600 {
601  libmesh_assert_equal_to (this->size(), vec.size());
602 
603  const int N = cast_int<int>(_val.size());
604  for (int i=0; i<N; i++)
605  (*this)(i) -= vec(i);
606 
607  return *this;
608 }
609 
610 
611 
612 template<typename T>
613 inline
615 {
616  libmesh_assert (this->size());
617  Real my_min = libmesh_real((*this)(0));
618 
619  const int N = cast_int<int>(_val.size());
620  for (int i=1; i!=N; i++)
621  {
622  Real current = libmesh_real((*this)(i));
623  my_min = (my_min < current? my_min : current);
624  }
625  return my_min;
626 }
627 
628 
629 
630 template<typename T>
631 inline
633 {
634  libmesh_assert (this->size());
635  Real my_max = libmesh_real((*this)(0));
636 
637  const int N = cast_int<int>(_val.size());
638  for (int i=1; i!=N; i++)
639  {
640  Real current = libmesh_real((*this)(i));
641  my_max = (my_max > current? my_max : current);
642  }
643  return my_max;
644 }
645 
646 
647 
648 template<typename T>
649 inline
651 {
652  if (!_val.size())
653  return 0.;
654 
655 #ifdef LIBMESH_HAVE_EIGEN
656  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).template lpNorm<1>();
657 #else
658  Real my_norm = 0.;
659  const int N = cast_int<int>(_val.size());
660  for (int i=0; i!=N; i++)
661  my_norm += std::abs((*this)(i));
662 
663  return my_norm;
664 #endif
665 }
666 
667 
668 
669 template<typename T>
670 inline
672 {
673  if (!_val.size())
674  return 0.;
675 
676 #ifdef LIBMESH_HAVE_EIGEN
677  return Eigen::Map<const typename Eigen::Matrix<T, Eigen::Dynamic, 1>>(_val.data(), _val.size()).norm();
678 #else
679  Real my_norm = 0.;
680  const int N = cast_int<int>(_val.size());
681  // The following pragma tells clang's vectorizer that it is safe to
682  // reorder floating point operations for this loop.
683 #ifdef __clang__
684 #pragma clang loop vectorize(enable)
685 #endif
686  for (int i=0; i!=N; i++)
687  my_norm += TensorTools::norm_sq((*this)(i));
688 
689  return sqrt(my_norm);
690 #endif
691 }
692 
693 
694 
695 template<typename T>
696 inline
698 {
699  if (!_val.size())
700  return 0.;
701 
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>();
704 #else
705  Real my_norm = TensorTools::norm_sq((*this)(0));
706 
707  const int N = cast_int<int>(_val.size());
708  for (int i=1; i!=N; i++)
709  {
710  Real current = TensorTools::norm_sq((*this)(i));
711  my_norm = (my_norm > current? my_norm : current);
712  }
713  return sqrt(my_norm);
714 #endif
715 }
716 
717 
718 
719 template<typename T>
720 inline
722  DenseVector<T> & dest) const
723 {
724  libmesh_assert_less_equal ( sub_n, this->size() );
725 
726  dest.resize(sub_n);
727  const int N = cast_int<int>(sub_n);
728  for (int i=0; i<N; i++)
729  dest(i) = _val[i];
730 }
731 
732 } // namespace libMesh
733 
734 #ifdef LIBMESH_HAVE_METAPHYSICL
735 namespace MetaPhysicL
736 {
737 template <typename T>
738 struct RawType<libMesh::DenseVector<T>>
739 {
741 
743  {
744  const auto s = in.size();
745  value_type ret(s);
746  for (unsigned int i = 0; i < s; ++i)
747  ret(i) = raw_value(in(i));
748 
749  return ret;
750  }
751 };
752 }
753 #endif
754 
755 #endif // LIBMESH_DENSE_VECTOR_H
T libmesh_real(T a)
DenseVector< T > & operator*=(const T factor)
将向量中的每个元素乘以 factor。
Definition: dense_vector.h:472
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。
Definition: dense_vector.h:485
T libmesh_conj(T a)
void resize(const unsigned int n)
调整向量的大小。将所有元素设置为0。
Definition: dense_vector.h:404
CompareTypes< T, T2 >::supertype dot(const DenseVector< T2 > &vec) const
Definition: dense_vector.h:500
std::vector< T > & get_values()
Definition: dense_vector.h:298
Real min() const
Definition: dense_vector.h:614
static value_type value(const libMesh::DenseVector< T > &in)
Definition: dense_vector.h:742
Real l2_norm() const
Definition: dense_vector.h:671
const T & operator()(const unsigned int i) const
Definition: dense_vector.h:439
virtual bool empty() const overridefinal
Definition: dense_vector.h:116
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的平方根。
Definition: type_vector.h:88
const Number zero
.
Definition: libmesh.h:248
void scale(const T factor)
将向量中的每个元素乘以 factor。
Definition: dense_vector.h:461
void swap(DenseVector< T > &other_vector)
STL 风格的交换方法。
Definition: dense_vector.h:395
virtual T & el(const unsigned int i) overridefinal
Definition: dense_vector.h:138
DenseVector & operator=(const DenseVector &)=default
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的绝对值。
Definition: type_vector.h:112
T norm_sq(std::complex< T > a)
Definition: tensor_tools.h:74
bool operator==(const DenseVector< T2 > &vec) const
Definition: dense_vector.h:549
const std::vector< T > & get_values() const
Definition: dense_vector.h:303
Real max() const
Definition: dense_vector.h:632
libMesh::DenseVector< typename RawType< T >::value_type > value_type
Definition: dense_vector.h:740
Real linfty_norm() const
Definition: dense_vector.h:697
virtual ~DenseVector()=default
Real l1_norm() const
Definition: dense_vector.h:650
std::vector< T > _val
实际数据值,存储为一维数组。
Definition: dense_vector.h:310
virtual T el(const unsigned int i) const overridefinal
Definition: dense_vector.h:135
bool operator!=(const DenseVector< T2 > &vec) const
Definition: dense_vector.h:566
virtual void zero() overridefinal
将向量中的每个元素设置为0。由于派生类中的存储方法可能不同,需要将其声明为纯虚函数。
Definition: dense_vector.h:428
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CompareTypes< T, T2 >::supertype indefinite_dot(const DenseVector< T2 > &vec) const
Definition: dense_vector.h:533
DenseVector(const unsigned int n=0)
构造函数。创建一个维数为 n 的稠密向量。
Definition: dense_vector.h:319
void append(const DenseVector< T2 > &other_vector)
将其他向量的额外条目附加到(调整大小但不改变)向量。
Definition: dense_vector.h:416
定义用于有限元计算的抽象稠密向量基类。 可以从这个类派生出特定的稠密向量,例如 DenseSubVectors。
Definition: dof_map.h:63
定义用于有限元计算的稠密向量类。该类基本上是为了补充 DenseMatrix 类而设计的。 它相对于 std::vector 具有额外的功能,使其在有限元中特别有用,特别是对于方程组。 所有重写的虚拟函...
Definition: dof_map.h:64
ADRealEigenVector< T, D, asd > norm(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的范数。
Definition: type_vector.h:64
DenseVector< T > & operator+=(const DenseVector< T2 > &vec)
将 vec 添加到该向量。
Definition: dense_vector.h:583
void get_principal_subvector(unsigned int sub_n, DenseVector< T > &dest) const
将大小为 sub_n(即前 sub_n 个条目)的主子向量放入 dest。
Definition: dense_vector.h:721
virtual unsigned int size() const overridefinal
Definition: dense_vector.h:111
DenseVector< T > & operator-=(const DenseVector< T2 > &vec)
从该向量中减去 vec。
Definition: dense_vector.h:599