Reference documentation for deal.II version 8.5.1
trilinos_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__trilinos_vector_h
17 #define dealii__trilinos_vector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
24 # include <deal.II/base/std_cxx11/shared_ptr.h>
25 # include <deal.II/base/subscriptor.h>
26 # include <deal.II/base/index_set.h>
27 # include <deal.II/base/utilities.h>
28 # include <deal.II/lac/exceptions.h>
29 # include <deal.II/lac/vector.h>
30 # include <deal.II/lac/trilinos_vector_base.h>
31 # include <deal.II/lac/vector_type_traits.h>
32 
33 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
34 # include "Epetra_Map.h"
35 # include "Epetra_LocalMap.h"
36 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
41 // forward declaration
42 template <typename> class Vector;
43 
48 namespace TrilinosWrappers
49 {
50  class SparseMatrix;
51 
52  namespace
53  {
54 #ifndef DEAL_II_WITH_64BIT_INDICES
55  // define a helper function that queries the global ID of local ID of
56  // an Epetra_BlockMap object by calling either the 32- or 64-bit
57  // function necessary.
58  inline
59  int gid(const Epetra_BlockMap &map, int i)
60  {
61  return map.GID(i);
62  }
63 #else
64  // define a helper function that queries the global ID of local ID of
65  // an Epetra_BlockMap object by calling either the 32- or 64-bit
66  // function necessary.
67  inline
68  long long int gid(const Epetra_BlockMap &map, int i)
69  {
70  return map.GID64(i);
71  }
72 #endif
73  }
74 
83  namespace MPI
84  {
85  class BlockVector;
86 
254  class Vector : public VectorBase
255  {
256  public:
260  typedef ::types::global_dof_index size_type;
261 
275  static const bool supports_distributed_data DEAL_II_DEPRECATED = true;
276 
286  Vector ();
287 
291  Vector (const Vector &v);
292 
293 #ifdef DEAL_II_WITH_CXX11
294 
301  Vector (Vector &&v);
302 #endif
303 
307  ~Vector ();
308 
332  void reinit (const VectorBase &v,
333  const bool omit_zeroing_entries = false,
334  const bool allow_different_maps = false);
335 
339  void reinit (const BlockVector &v,
340  const bool import_data = false);
341 
348  Vector &operator= (const TrilinosScalar s);
349 
355  Vector &operator= (const Vector &v);
356 
357 #ifdef DEAL_II_WITH_CXX11
358 
365  Vector &operator= (Vector &&v);
366 #endif
367 
375  Vector &operator= (const ::TrilinosWrappers::Vector &v);
376 
385  template <typename Number>
386  Vector &operator= (const ::Vector<Number> &v);
387 
406  (const ::TrilinosWrappers::SparseMatrix &matrix,
407  const Vector &vector);
409 
429  explicit Vector (const Epetra_Map &parallel_partitioning) DEAL_II_DEPRECATED;
430 
447  Vector (const Epetra_Map &parallel_partitioning,
448  const VectorBase &v) DEAL_II_DEPRECATED;
449 
464  template <typename number>
465  void reinit (const Epetra_Map &parallel_partitioner,
466  const ::Vector<number> &v) DEAL_II_DEPRECATED;
467 
482  void reinit (const Epetra_Map &parallel_partitioning,
483  const bool omit_zeroing_entries = false) DEAL_II_DEPRECATED;
484 
499  template <typename Number>
500  Vector (const Epetra_Map &parallel_partitioning,
501  const ::Vector<Number> &v) DEAL_II_DEPRECATED;
503 
525  explicit Vector (const IndexSet &parallel_partitioning,
526  const MPI_Comm &communicator = MPI_COMM_WORLD);
527 
539  Vector (const IndexSet &local,
540  const IndexSet &ghost,
541  const MPI_Comm &communicator = MPI_COMM_WORLD);
542 
557  Vector (const IndexSet &parallel_partitioning,
558  const VectorBase &v,
559  const MPI_Comm &communicator = MPI_COMM_WORLD);
560 
573  template <typename Number>
574  Vector (const IndexSet &parallel_partitioning,
575  const ::Vector<Number> &v,
576  const MPI_Comm &communicator = MPI_COMM_WORLD);
577 
601  void reinit (const IndexSet &parallel_partitioning,
602  const MPI_Comm &communicator = MPI_COMM_WORLD,
603  const bool omit_zeroing_entries = false);
604 
630  void reinit (const IndexSet &locally_owned_entries,
631  const IndexSet &ghost_entries,
632  const MPI_Comm &communicator = MPI_COMM_WORLD,
633  const bool vector_writable = false);
635  };
636 
637 
638 
639 
640 // ------------------- inline and template functions --------------
641 
642 
651  inline
652  void swap (Vector &u, Vector &v)
653  {
654  u.swap (v);
655  }
656 
657 
658 #ifndef DOXYGEN
659 
660  template <typename number>
661  Vector::Vector (const Epetra_Map &input_map,
662  const ::Vector<number> &v)
663  {
664  reinit (input_map, v);
665  }
666 
667 
668 
669  template <typename number>
670  Vector::Vector (const IndexSet &parallel_partitioner,
671  const ::Vector<number> &v,
672  const MPI_Comm &communicator)
673  {
674  *this = Vector(parallel_partitioner.make_trilinos_map (communicator, true),
675  v);
676  owned_elements = parallel_partitioner;
677  }
678 
679 
680 
681 
682  template <typename number>
683  void Vector::reinit (const Epetra_Map &parallel_partitioner,
684  const ::Vector<number> &v)
685  {
686  if (vector.get() == 0 || vector->Map().SameAs(parallel_partitioner) == false)
687  reinit(parallel_partitioner);
688 
689  const int size = parallel_partitioner.NumMyElements();
690 
691  // Need to copy out values, since the deal.II might not use doubles, so
692  // that a direct access is not possible.
693  for (int i=0; i<size; ++i)
694  (*vector)[0][i] = v(gid(parallel_partitioner,i));
695 
696  }
697 
698 
699  inline
700  Vector &
701  Vector::operator= (const TrilinosScalar s)
702  {
704 
705  return *this;
706  }
707 
708 
709  template <typename Number>
710  Vector &
711  Vector::operator= (const ::Vector<Number> &v)
712  {
713  if (size() != v.size())
714  {
715  vector.reset (new Epetra_FEVector(Epetra_Map
716  (static_cast<TrilinosWrappers::types::int_type>(v.size()), 0,
717 #ifdef DEAL_II_WITH_MPI
718  Epetra_MpiComm(MPI_COMM_SELF)
719 #else
720  Epetra_SerialComm()
721 #endif
722  )));
723  }
724 
725  const Epetra_Map &map = vector_partitioner();
726  const int size = map.NumMyElements();
727 
728  // Need to copy out values, since the deal.II might not use doubles, so
729  // that a direct access is not possible.
730  for (int i=0; i<size; ++i)
731  (*vector)[0][i] = v(gid(map,i));
732 
733  return *this;
734  }
735 
736 
737 #endif /* DOXYGEN */
738 
739  } /* end of namespace MPI */
740 
741 
742 
756  class Vector : public VectorBase
757  {
758  public:
762  typedef ::types::global_dof_index size_type;
763 
777  static const bool supports_distributed_data DEAL_II_DEPRECATED = false;
778 
784  Vector () DEAL_II_DEPRECATED;
785 
789  explicit Vector (const size_type n) DEAL_II_DEPRECATED;
790 
800  explicit Vector (const Epetra_Map &partitioning) DEAL_II_DEPRECATED;
801 
811  explicit Vector (const IndexSet &partitioning,
812  const MPI_Comm &communicator = MPI_COMM_WORLD) DEAL_II_DEPRECATED;
813 
818  explicit Vector (const VectorBase &V) DEAL_II_DEPRECATED;
819 
824  template <typename Number>
825  explicit Vector (const ::Vector<Number> &v) DEAL_II_DEPRECATED;
826 
836  void reinit (const size_type n,
837  const bool omit_zeroing_entries = false);
838 
856  void reinit (const Epetra_Map &input_map,
857  const bool omit_zeroing_entries = false);
858 
876  void reinit (const IndexSet &input_map,
877  const MPI_Comm &communicator = MPI_COMM_WORLD,
878  const bool omit_zeroing_entries = false);
879 
893  void reinit (const VectorBase &V,
894  const bool omit_zeroing_entries = false,
895  const bool allow_different_maps = false);
896 
903  Vector &operator= (const TrilinosScalar s);
904 
909  Vector &operator= (const MPI::Vector &v);
910 
914  template <typename Number>
915  Vector &operator= (const ::Vector<Number> &v);
916 
921  Vector &operator= (const Vector &v);
922 
934  void update_ghost_values () const;
935  };
936 
937 
938 // ------------------- inline and template functions --------------
939 
940 
949  inline
950  void swap (Vector &u, Vector &v)
951  {
952  u.swap (v);
953  }
954 
955 
956 #ifndef DOXYGEN
957 
958  template <typename number>
959  Vector::Vector (const ::Vector<number> &v)
960  {
961  Epetra_LocalMap map ((TrilinosWrappers::types::int_type)v.size(), 0, Utilities::Trilinos::comm_self());
962  vector.reset (new Epetra_FEVector(map));
963  *this = v;
964  }
965 
966 
967 
968  inline
969  Vector &
970  Vector::operator= (const TrilinosScalar s)
971  {
973 
974  return *this;
975  }
976 
977 
978 
979  template <typename Number>
980  Vector &
981  Vector::operator= (const ::Vector<Number> &v)
982  {
983  if (size() != v.size())
984  {
985  vector.reset();
986 
987  Epetra_LocalMap map ((TrilinosWrappers::types::int_type)v.size(), 0,
989  vector.reset (new Epetra_FEVector(map));
990  owned_elements = v.locally_owned_elements();
991  }
992 
993  const Epetra_Map &map = vector_partitioner();
994  const TrilinosWrappers::types::int_type size = map.NumMyElements();
995 
996  Assert (map.MaxLID() == size-1,
997  ExcDimensionMismatch(map.MaxLID(), size-1));
998 
999  // Need to copy out values, since the
1000  // deal.II might not use doubles, so
1001  // that a direct access is not possible.
1002  for (TrilinosWrappers::types::int_type i=0; i<size; ++i)
1003  (*vector)[0][i] = v(i);
1004 
1005  return *this;
1006  }
1007 
1008 
1009 
1010  inline
1011  void
1013  {}
1014 
1015 
1016 #endif
1017 
1018 } /* namespace TrilinosWrappers */
1019 
1023 namespace internal
1024 {
1025  namespace LinearOperator
1026  {
1027  template <typename> class ReinitHelper;
1028 
1033  template<>
1034  class ReinitHelper<TrilinosWrappers::MPI::Vector>
1035  {
1036  public:
1037  template <typename Matrix>
1038  static
1039  void reinit_range_vector (const Matrix &matrix,
1041  bool omit_zeroing_entries)
1042  {
1043  v.reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator(), omit_zeroing_entries);
1044  }
1045 
1046  template <typename Matrix>
1047  static
1048  void reinit_domain_vector(const Matrix &matrix,
1050  bool omit_zeroing_entries)
1051  {
1052  v.reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator(), omit_zeroing_entries);
1053  }
1054  };
1055 
1060  template<>
1062  {
1063  public:
1064  template <typename Matrix>
1065  static
1066  void reinit_range_vector (const Matrix &matrix,
1068  bool omit_zeroing_entries)
1069  {
1070  v.reinit(matrix.locally_owned_range_indices(),
1071  matrix.get_mpi_communicator(),
1072  omit_zeroing_entries);
1073  }
1074 
1075  template <typename Matrix>
1076  static
1077  void reinit_domain_vector(const Matrix &matrix,
1079  bool omit_zeroing_entries)
1080  {
1081  v.reinit(matrix.locally_owned_domain_indices(),
1082  matrix.get_mpi_communicator(),
1083  omit_zeroing_entries);
1084  }
1085  };
1086 
1087  } /* namespace LinearOperator */
1088 } /* namespace internal */
1089 
1090 
1096 template <>
1097 struct is_serial_vector< TrilinosWrappers::Vector > : std_cxx11::true_type
1098 {
1099 };
1100 
1101 
1107 template <>
1108 struct is_serial_vector< TrilinosWrappers::MPI::Vector > : std_cxx11::false_type
1109 {
1110 };
1111 
1112 
1113 DEAL_II_NAMESPACE_CLOSE
1114 
1115 #endif // DEAL_II_WITH_TRILINOS
1116 
1117 /*---------------------------- trilinos_vector.h ---------------------------*/
1118 
1119 #endif
1120 /*---------------------------- trilinos_vector.h ---------------------------*/
Vector & operator=(const TrilinosScalar s)
::types::global_dof_index size_type
const Epetra_Comm & comm_self()
Definition: utilities.cc:736
void reinit(const VectorBase &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:535
std_cxx11::shared_ptr< Epetra_FEVector > vector
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const Epetra_Map & vector_partitioner() const
void update_ghost_values() const
void swap(Vector &u, Vector &v)
void import_nonlocal_data_for_fe(const ::TrilinosWrappers::SparseMatrix &matrix, const Vector &vector)
VectorBase & operator=(const TrilinosScalar s)
void reinit(const size_type n, const bool omit_zeroing_entries=false)
::types::global_dof_index size_type
static const bool supports_distributed_data