libmesh解析
本工作只是尝试解析原libmesh的代码,供学习使用
 全部  命名空间 文件 函数 变量 类型定义 枚举 枚举值 友元 
laspack_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 
21 #ifndef LIBMESH_LASPACK_VECTOR_H
22 #define LIBMESH_LASPACK_VECTOR_H
23 
24 
25 
26 #include "libmesh/libmesh_common.h"
27 
28 #ifdef LIBMESH_HAVE_LASPACK
29 
30 // Local includes
31 #include "libmesh/numeric_vector.h"
32 
33 // Laspack includes
34 #include <operats.h>
35 #include <qvector.h>
36 
37 // C++ includes
38 #include <limits>
39 #include <mutex>
40 
41 namespace libMesh
42 {
43 
44 // Forward declarations
45 template <typename T> class LaspackLinearSolver;
46 template <typename T> class SparseMatrix;
47 
55 template <typename T>
56 class LaspackVector final : public NumericVector<T>
57 {
58 public:
59 
67  explicit
68  LaspackVector (const Parallel::Communicator & comm,
69  const ParallelType = AUTOMATIC);
70 
79  explicit
80  LaspackVector (const Parallel::Communicator & comm,
81  const numeric_index_type n,
82  const ParallelType = AUTOMATIC);
83 
94  LaspackVector (const Parallel::Communicator & comm,
95  const numeric_index_type n,
96  const numeric_index_type n_local,
97  const ParallelType = AUTOMATIC);
98 
111  LaspackVector (const Parallel::Communicator & comm,
112  const numeric_index_type N,
113  const numeric_index_type n_local,
114  const std::vector<numeric_index_type> & ghost,
115  const ParallelType = AUTOMATIC);
116 
125  LaspackVector<T> & operator= (const LaspackVector<T> & v);
126 
133  LaspackVector (LaspackVector &&) = delete;
134  LaspackVector (const LaspackVector &) = delete;
135  LaspackVector & operator= (LaspackVector &&) = delete;
136  virtual ~LaspackVector ();
137 
142  virtual void close () override;
143 
148  virtual void clear () override;
149 
154  virtual void zero () override;
155 
156 
161  virtual std::unique_ptr<NumericVector<T>> zero_clone () const override;
162 
167  virtual std::unique_ptr<NumericVector<T>> clone () const override;
168 
178  virtual void init (const numeric_index_type N,
179  const numeric_index_type n_local,
180  const bool fast=false,
181  const ParallelType ptype=AUTOMATIC) override;
182 
191  virtual void init (const numeric_index_type N,
192  const bool fast=false,
193  const ParallelType ptype=AUTOMATIC) override;
194 
205  virtual void init (const numeric_index_type N,
206  const numeric_index_type n_local,
207  const std::vector<numeric_index_type> & ghost,
208  const bool fast = false,
209  const ParallelType = AUTOMATIC) override;
210 
218  virtual void init (const NumericVector<T> & other,
219  const bool fast = false) override;
220 
228  virtual NumericVector<T> & operator= (const T s) override;
229 
237  virtual NumericVector<T> & operator= (const NumericVector<T> & v) override;
238 
246  virtual NumericVector<T> & operator= (const std::vector<T> & v) override;
247 
254  virtual Real min () const override;
255 
261  virtual Real max () const override;
262 
269  virtual T sum () const override;
270 
277  virtual Real l1_norm () const override;
278 
285  virtual Real l2_norm () const override;
286 
293  virtual Real linfty_norm () const override;
294 
300  virtual numeric_index_type size () const override;
301 
307  virtual numeric_index_type local_size() const override;
308 
314  virtual numeric_index_type first_local_index() const override;
315 
321  virtual numeric_index_type last_local_index() const override;
322 
329  virtual T operator() (const numeric_index_type i) const override;
330 
337  virtual NumericVector<T> & operator += (const NumericVector<T> & v) override;
338 
345  virtual NumericVector<T> & operator -= (const NumericVector<T> & v) override;
346 
353  virtual NumericVector<T> & operator *= (const NumericVector<T> & v) override;
354 
361  virtual NumericVector<T> & operator /= (const NumericVector<T> & v) override;
362 
366  virtual void reciprocal() override;
367 
371  virtual void conjugate() override;
372 
379  virtual void set (const numeric_index_type i, const T value) override;
380 
387  virtual void add (const numeric_index_type i, const T value) override;
388 
394  virtual void add (const T s) override;
395 
401  virtual void add (const NumericVector<T> & v) override;
402 
409  virtual void add (const T a, const NumericVector<T> & v) override;
410 
415  using NumericVector<T>::add_vector;
416 
423  virtual void add_vector (const NumericVector<T> & v, const SparseMatrix<T> & A) override;
424 
431  virtual void add_vector_transpose (const NumericVector<T> & v, const SparseMatrix<T> & A) override;
432 
438  virtual void scale (const T factor) override;
439 
443  virtual void abs() override;
444 
451  virtual T dot(const NumericVector<T> & v) const override;
452 
458  virtual void localize (std::vector<T> & v_local) const override;
459 
465  virtual void localize (NumericVector<T> & v_local) const override;
466 
473  virtual void localize (NumericVector<T> & v_local, const std::vector<numeric_index_type> & send_list) const override;
474 
481  virtual void localize (std::vector<T> & v_local, const std::vector<numeric_index_type> & indices) const override;
482 
490  virtual void localize (const numeric_index_type first_local_idx, const numeric_index_type last_local_idx, const std::vector<numeric_index_type> & send_list) override;
491 
498  virtual void localize_to_one (std::vector<T> & v_local, const processor_id_type proc_id=0) const override;
499 
506  virtual void pointwise_mult (const NumericVector<T> & vec1, const NumericVector<T> & vec2) override;
507 
514  virtual void pointwise_divide (const NumericVector<T> & vec1, const NumericVector<T> & vec2) override;
515 
521  virtual void swap (NumericVector<T> & v) override;
522 
528  virtual std::size_t max_allowed_id() const override;
529 
530 private:
531 
535  QVector _vec;
536 
540  friend class LaspackLinearSolver<T>;
541 };
542 
543 
544 
545 //----------------------------------------------------------
546 // LaspackVector inline methods
550 template <typename T>
551 inline
552 LaspackVector<T>::LaspackVector (const Parallel::Communicator & comm,
553  const ParallelType ptype)
554  : NumericVector<T>(comm, ptype)
555 {
556  this->_type = ptype;
557 }
558 
559 
563 template <typename T>
564 inline
565 LaspackVector<T>::LaspackVector (const Parallel::Communicator & comm,
566  const numeric_index_type n,
567  const ParallelType ptype)
568  : NumericVector<T>(comm, ptype)
569 {
570  this->init(n, n, false, ptype);
571 }
572 
573 
577 template <typename T>
578 inline
579 LaspackVector<T>::LaspackVector (const Parallel::Communicator & comm,
580  const numeric_index_type n,
581  const numeric_index_type n_local,
582  const ParallelType ptype)
583  : NumericVector<T>(comm, ptype)
584 {
585  this->init(n, n_local, false, ptype);
586 }
587 
588 
592 template <typename T>
593 inline
594 LaspackVector<T>::LaspackVector (const Parallel::Communicator & comm,
595  const numeric_index_type N,
596  const numeric_index_type n_local,
597  const std::vector<numeric_index_type> & ghost,
598  const ParallelType ptype)
599  : NumericVector<T>(comm, ptype)
600 {
601  this->init(N, n_local, ghost, false, ptype);
602 }
603 
604 
608 template <typename T>
609 inline
611 {
612  this->clear ();
613 }
614 
615 
626 template <typename T>
627 inline
629  const numeric_index_type libmesh_dbg_var(n_local),
630  const bool fast,
631  const ParallelType ptype)
632 {
633  //Laspack向量仅用于串行情况,但可以在一个处理器上提供“并行”向量。
634  libmesh_assert_equal_to (n, n_local);
635 
636  this->_type = SERIAL;
637 
638  // 清除已初始化的向量
639  if (this->initialized())
640  this->clear();
641 
642  // 创建一个顺序向量
643 
644  static int cnt = 0;
645  std::string foo{"Vec-"};
646  foo += std::to_string(cnt++);
647 
648  // 在 Laspack 中构造一个向量对象
649  V_Constr(&_vec, const_cast<char *>(foo.c_str()), n, Normal, _LPTrue);
650 
651  // 标记对象已初始化
652  this->_is_initialized = true;
653 #ifndef NDEBUG
654  this->_is_closed = true;
655 #endif
656 
657  // 可选地将所有组件清零
658  if (fast == false)
659  this->zero ();
660 
661  return;
662 }
663 
664 
665 template <typename T>
666 inline
668  const bool fast,
669  const ParallelType ptype)
670 {
671  this->init(n,n,fast,ptype);
672 }
673 
674 
675 template <typename T>
676 inline
678  const numeric_index_type n_local,
679  const std::vector<numeric_index_type> & libmesh_dbg_var(ghost),
680  const bool fast,
681  const ParallelType ptype)
682 {
683  libmesh_assert(ghost.empty());
684  this->init(n,n_local,fast,ptype);
685 }
686 
687 
688 
689 /* Default implementation for solver packages for which ghosted
690  vectors are not yet implemented. */
691 template <class T>
693  const bool fast)
694 {
695  this->init(other.size(),other.local_size(),fast,other.type());
696 }
697 
698 
699 
700 template <typename T>
701 inline
703 {
704  libmesh_assert (this->initialized());
705 
706 #ifndef NDEBUG
707  this->_is_closed = true;
708 #endif
709 }
710 
711 
712 
713 template <typename T>
714 inline
716 {
717  if (this->initialized())
718  {
719  V_Destr (&_vec);
720  }
721 
722  this->_is_initialized = false;
723 #ifndef NDEBUG
724  this->_is_closed = false;
725 #endif
726 }
727 
728 
729 
730 template <typename T> inline
732 {
733  libmesh_assert (this->initialized());
734  libmesh_assert (this->closed());
735 
736  V_SetAllCmp (&_vec, 0.);
737 }
738 
739 
740 
741 template <typename T>
742 inline
743 std::unique_ptr<NumericVector<T>> LaspackVector<T>::zero_clone () const
744 {
745  NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
746 
747  cloned_vector->init(*this);
748 
749  return std::unique_ptr<NumericVector<T>>(cloned_vector);
750 }
751 
752 
753 
754 template <typename T>
755 inline
756 std::unique_ptr<NumericVector<T>> LaspackVector<T>::clone () const
757 {
758  NumericVector<T> * cloned_vector = new LaspackVector<T>(this->comm());
759 
760  cloned_vector->init(*this, true);
761 
762  *cloned_vector = *this;
763 
764  return std::unique_ptr<NumericVector<T>>(cloned_vector);
765 }
766 
767 
768 
769 template <typename T>
770 inline
772 {
773  libmesh_assert (this->initialized());
774 
775  return static_cast<numeric_index_type>(V_GetDim(const_cast<QVector*>(&_vec)));
776 }
777 
778 
779 
780 template <typename T>
781 inline
783 {
784  libmesh_assert (this->initialized());
785 
786  return this->size();
787 }
788 
789 
790 
791 template <typename T>
792 inline
794 {
795  libmesh_assert (this->initialized());
796 
797  return 0;
798 }
799 
800 
801 
802 template <typename T>
803 inline
805 {
806  libmesh_assert (this->initialized());
807 
808  return this->size();
809 }
810 
811 
812 
813 template <typename T>
814 inline
815 void LaspackVector<T>::set (const numeric_index_type i, const T value)
816 {
817  libmesh_assert (this->initialized());
818  libmesh_assert_less (i, this->size());
819 
820  std::scoped_lock lock(this->_numeric_vector_mutex);
821  V_SetCmp (&_vec, i+1, value);
822 
823 #ifndef NDEBUG
824  this->_is_closed = false;
825 #endif
826 }
827 
828 
829 
830 template <typename T>
831 inline
832 void LaspackVector<T>::add (const numeric_index_type i, const T value)
833 {
834  libmesh_assert (this->initialized());
835  libmesh_assert_less (i, this->size());
836 
837  std::scoped_lock lock(this->_numeric_vector_mutex);
838  V_AddCmp (&_vec, i+1, value);
839 
840 #ifndef NDEBUG
841  this->_is_closed = false;
842 #endif
843 }
844 
845 
846 
847 template <typename T>
848 inline
850 {
851  libmesh_assert (this->initialized());
852  libmesh_assert ( ((i >= this->first_local_index()) &&
853  (i < this->last_local_index())) );
854 
855 
856  return static_cast<T>(V_GetCmp(const_cast<QVector*>(&_vec), i+1));
857 }
858 
859 
860 
861 template <typename T>
862 inline
864 {
865  LaspackVector<T> & v = cast_ref<LaspackVector<T> &>(other);
866 
867  // This is all grossly dependent on Laspack version...
868 
869  std::swap(_vec.Name, v._vec.Name);
870  std::swap(_vec.Dim, v._vec.Dim);
871  std::swap(_vec.Instance, v._vec.Instance);
872  std::swap(_vec.LockLevel, v._vec.LockLevel);
873  std::swap(_vec.Multipl, v._vec.Multipl);
874  std::swap(_vec.OwnData, v._vec.OwnData);
875 
876  // This should still be O(1), since _vec.Cmp is just a pointer to
877  // data on the heap
878 
879  std::swap(_vec.Cmp, v._vec.Cmp);
880 }
881 
882 
883 
884 template <typename T>
885 inline
887 {
888  // The QVector type declares a "size_t Dim;"
889  return std::numeric_limits<std::size_t>::max();
890 }
891 
892 
893 } // namespace libMesh
894 
895 
896 #endif // #ifdef LIBMESH_HAVE_LASPACK
897 #endif // LIBMESH_LASPACK_VECTOR_H
virtual void add_vector(const NumericVector< T > &v, const SparseMatrix< T > &A) override
将一个向量添加到该向量并将结果存储在其中。
bool closed()
Checks that the library has been closed.
Definition: libmesh.C:268
virtual void swap(NumericVector< T > &v) override
交换两个向量的内容。
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
使用给定的全局大小和局部大小初始化向量。
virtual void close() override
Close the vector.
virtual numeric_index_type size() const override
返回向量的大小(全局维度)。
virtual NumericVector< T > & operator-=(const NumericVector< T > &v) override
将另一个向量的元素按元素从该向量中减去。
virtual void abs() override
计算向量的每个元素的绝对值。
virtual void scale(const T factor) override
缩放向量的所有元素,将它们乘以指定的标量因子。
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A) override
将一个向量的转置添加到该向量并将结果存储在其中。
virtual T dot(const NumericVector< T > &v) const override
计算向量与另一个向量的点积。
virtual void set(const numeric_index_type i, const T value) override
设置指定索引处的元素值。
virtual numeric_index_type size() const =0
获取向量的大小。
LaspackVector< T > & operator=(const LaspackVector< T > &v)
Copy assignment operator.
LaspackVector(const Parallel::Communicator &comm, const ParallelType=AUTOMATIC)
Dummy-Constructor.
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
对向量的每个元素执行逐元素乘法。
提供了不同线性代数库的向量存储方案的统一接口。
Definition: dof_map.h:67
virtual T operator()(const numeric_index_type i) const override
返回指定索引处向量的值。
const Number zero
.
Definition: libmesh.h:248
virtual NumericVector< T > & operator+=(const NumericVector< T > &v) override
将另一个向量的元素按元素加到该向量。
virtual NumericVector< T > & operator*=(const NumericVector< T > &v) override
将该向量的元素按元素与另一个向量的元素相乘。
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
更改向量的维度为 n 。如果可能的话 ,该向量的保留内存保持不变。 如果 n==0 ,所有内存都将被释放。因此,如果要调整向量的大小并释放不需要的内存, 必须首先调用 init(0) ,然后调用 ini...
uint8_t processor_id_type
Definition: id_types.h:104
virtual std::size_t max_allowed_id() const override
获取允许的最大ID(索引)。
virtual std::unique_ptr< NumericVector< T > > zero_clone() const override
创建并返回一个指向该对象的零克隆的 unique_ptr。
virtual std::unique_ptr< NumericVector< T > > clone() const override
创建并返回一个指向该对象的克隆的 unique_ptr。
virtual Real max() const override
返回向量的最大元素。
virtual Real l2_norm() const override
返回向量的 L2 范数。
dof_id_type numeric_index_type
Definition: id_types.h:99
bool _is_initialized
Flag that tells if init() has been called.
Definition: libmesh.C:242
virtual void add(const numeric_index_type i, const T value) override
将指定值添加到指定索引处的元素。
virtual Real min() const override
返回向量的最小元素。
virtual void localize(std::vector< T > &v_local) const override
将向量的元素本地化,以便在不同处理器之间交换。
virtual numeric_index_type last_local_index() const override
返回向量的最后一个局部索引。
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const override
将向量本地化到指定处理器。
ParallelType _type
向量的类型。
virtual numeric_index_type first_local_index() const override
返回向量的第一个局部索引。
这个类为基于laspackc的串行向量数据结构提供了一个很好的接口。 所有被覆盖的虚函数都记录在numeric_vector.h中。
virtual void zero() override
Set all elements of the vector to zero.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void conjugate() override
计算向量中所有元素的复共轭。
virtual Real linfty_norm() const override
返回向量的 Linfinity 范数。
virtual numeric_index_type local_size() const =0
获取向量的本地大小,即 index_stop - index_start。
virtual numeric_index_type local_size() const override
返回向量的局部大小(局部维度)。
virtual void clear() override
Clear the vector.
ParallelType type() const
获取向量的类型。
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:261
virtual NumericVector< T > & operator/=(const NumericVector< T > &v) override
将该向量的元素按元素除以另一个向量的元素。
virtual T sum() const override
返回向量的所有元素之和。
QVector _vec
用于保存向量条目的实际 Laspack 向量数据类型。
virtual void reciprocal() override
计算向量中所有元素的倒数(1/x)。
virtual Real l1_norm() const override
返回向量的 L1 范数。
virtual ~LaspackVector()
LaspackVector 对象的析构函数。
virtual void pointwise_divide(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
对向量的每个元素执行逐元素除法。