21 #ifndef LIBMESH_PETSC_VECTOR_H 22 #define LIBMESH_PETSC_VECTOR_H 27 #ifdef LIBMESH_HAVE_PETSC 42 #include <unordered_map> 44 #ifdef LIBMESH_HAVE_CXX11_THREAD 53 template <
typename T>
class SparseMatrix;
100 const std::vector<numeric_index_type> & ghost,
130 virtual void close ()
override;
132 virtual void clear ()
override;
134 virtual void zero ()
override;
136 virtual std::unique_ptr<NumericVector<T>>
zero_clone ()
const override;
138 virtual std::unique_ptr<NumericVector<T>>
clone ()
const override;
142 const bool fast=
false,
146 const bool fast=
false,
151 const std::vector<numeric_index_type> & ghost,
152 const bool fast =
false,
156 const bool fast =
false)
override;
164 virtual Real min ()
const override;
166 virtual Real max ()
const override;
168 virtual T
sum ()
const override;
195 virtual void get(
const std::vector<numeric_index_type> & index,
196 T * values)
const override;
231 const T
value)
override;
234 const T
value)
override;
236 virtual void add (
const T s)
override;
249 const std::vector<numeric_index_type> & dof_indices)
override;
272 virtual void insert (
const T * v,
273 const std::vector<numeric_index_type> & dof_indices)
override;
275 virtual void scale (
const T factor)
override;
279 virtual void abs()
override;
290 virtual void localize (std::vector<T> & v_local)
const override;
295 const std::vector<numeric_index_type> & send_list)
const override;
297 virtual void localize (std::vector<T> & v_local,
298 const std::vector<numeric_index_type> & indices)
const override;
302 const std::vector<numeric_index_type> & send_list)
override;
313 const std::vector<numeric_index_type> & rows)
const override;
340 #ifdef LIBMESH_HAVE_CXX11_THREAD 398 #ifdef LIBMESH_HAVE_CXX11_THREAD 452 template <
typename T>
456 _array_is_present(false),
459 _local_form(nullptr),
461 _global_to_local_map(),
462 _destroy_vec_on_exit(true),
463 _values_manually_retrieved(false),
464 _values_read_only(false)
471 template <
typename T>
477 _array_is_present(false),
478 _local_form(nullptr),
480 _global_to_local_map(),
481 _destroy_vec_on_exit(true),
482 _values_manually_retrieved(false),
483 _values_read_only(false)
485 this->
init(n, n,
false, ptype);
490 template <
typename T>
497 _array_is_present(false),
498 _local_form(nullptr),
500 _global_to_local_map(),
501 _destroy_vec_on_exit(true),
502 _values_manually_retrieved(false),
503 _values_read_only(false)
505 this->
init(n, n_local,
false, ptype);
510 template <
typename T>
515 const std::vector<numeric_index_type> & ghost,
518 _array_is_present(false),
519 _local_form(nullptr),
521 _global_to_local_map(),
522 _destroy_vec_on_exit(true),
523 _values_manually_retrieved(false),
524 _values_read_only(false)
526 this->
init(n, n_local, ghost,
false, ptype);
533 template <
typename T>
538 _array_is_present(false),
539 _local_form(nullptr),
541 _global_to_local_map(),
542 _destroy_vec_on_exit(false),
543 _values_manually_retrieved(false),
544 _values_read_only(false)
552 PetscErrorCode
ierr=0;
553 PetscInt petsc_local_size=0;
554 ierr = VecGetLocalSize(
_vec, &petsc_local_size);
555 LIBMESH_CHKERR(
ierr);
560 #if PETSC_VERSION_LESS_THAN(3,0,0) || !PETSC_VERSION_LESS_THAN(3,4,0) 567 LIBMESH_CHKERR(
ierr);
569 if ((std::strcmp(ptype,VECSHARED) == 0) || (std::strcmp(ptype,VECMPI) == 0))
571 #if PETSC_RELEASE_LESS_THAN(3,1,1) 572 ISLocalToGlobalMapping mapping =
_vec->mapping;
574 ISLocalToGlobalMapping mapping;
575 ierr = VecGetLocalToGlobalMapping(
_vec, &mapping);
576 LIBMESH_CHKERR(
ierr);
584 #if PETSC_RELEASE_LESS_THAN(3,4,0) 588 ierr = ISLocalToGlobalMappingGetSize(mapping, &n);
589 LIBMESH_CHKERR(
ierr);
592 #if PETSC_RELEASE_LESS_THAN(3,1,1) 593 const PetscInt * indices = mapping->indices;
595 const PetscInt * indices;
596 ierr = ISLocalToGlobalMappingGetIndices(mapping,&indices);
597 LIBMESH_CHKERR(
ierr);
602 #if !PETSC_RELEASE_LESS_THAN(3,1,1) 603 ierr = ISLocalToGlobalMappingRestoreIndices(mapping, &indices);
604 LIBMESH_CHKERR(
ierr);
617 template <
typename T>
626 template <
typename T>
633 parallel_object_only();
635 PetscErrorCode
ierr=0;
636 PetscInt petsc_n=
static_cast<PetscInt
>(n);
637 PetscInt petsc_n_local=
static_cast<PetscInt
>(n_local);
654 libmesh_assert ((this->_type==
SERIAL && n==n_local) ||
658 if (this->_type ==
SERIAL)
660 ierr = VecCreateSeq (PETSC_COMM_SELF, petsc_n, &_vec);
661 CHKERRABORT(PETSC_COMM_SELF,
ierr);
663 ierr = VecSetFromOptions (_vec);
664 CHKERRABORT(PETSC_COMM_SELF,
ierr);
669 #ifdef LIBMESH_HAVE_MPI 670 libmesh_assert_less_equal (n_local, n);
671 ierr = VecCreateMPI (this->comm().
get(), petsc_n_local, petsc_n,
673 LIBMESH_CHKERR(
ierr);
675 libmesh_assert_equal_to (n_local, n);
676 ierr = VecCreateSeq (PETSC_COMM_SELF, petsc_n, &_vec);
677 CHKERRABORT(PETSC_COMM_SELF,
ierr);
680 ierr = VecSetFromOptions (_vec);
681 LIBMESH_CHKERR(
ierr);
684 libmesh_error_msg(
"Unsupported type " << this->_type);
687 this->_is_closed =
true;
696 template <
typename T>
702 this->
init(n,n,fast,ptype);
707 template <
typename T>
711 const std::vector<numeric_index_type> & ghost,
715 parallel_object_only();
717 PetscErrorCode
ierr=0;
718 PetscInt petsc_n=
static_cast<PetscInt
>(n);
719 PetscInt petsc_n_local=
static_cast<PetscInt
>(n_local);
720 PetscInt petsc_n_ghost=
static_cast<PetscInt
>(ghost.size());
735 PetscInt * petsc_ghost = ghost.empty() ? PETSC_NULL :
736 const_cast<PetscInt *
>(
reinterpret_cast<const PetscInt *
>(ghost.data()));
748 _global_to_local_map[ghost[i]] = i;
752 ierr = VecCreateGhost (this->comm().
get(), petsc_n_local, petsc_n,
753 petsc_n_ghost, petsc_ghost,
755 LIBMESH_CHKERR(
ierr);
757 ierr = VecSetFromOptions (_vec);
758 LIBMESH_CHKERR(
ierr);
761 this->_is_closed =
true;
769 template <
typename T>
774 parallel_object_only();
780 const PetscVector<T> & v = cast_ref<const PetscVector<T> &>(other);
793 this->_is_closed =
true;
796 this->_type = v.
_type;
799 PetscErrorCode
ierr = VecDuplicate (v.
_vec, &this->_vec);
800 LIBMESH_CHKERR(
ierr);
808 template <
typename T>
812 parallel_object_only();
814 this->_restore_array();
816 PetscErrorCode
ierr=0;
818 ierr = VecAssemblyBegin(_vec);
819 LIBMESH_CHKERR(
ierr);
820 ierr = VecAssemblyEnd(_vec);
821 LIBMESH_CHKERR(
ierr);
825 ierr = VecGhostUpdateBegin(_vec,INSERT_VALUES,SCATTER_FORWARD);
826 LIBMESH_CHKERR(
ierr);
827 ierr = VecGhostUpdateEnd(_vec,INSERT_VALUES,SCATTER_FORWARD);
828 LIBMESH_CHKERR(
ierr);
831 this->_is_closed =
true;
836 template <
typename T>
840 parallel_object_only();
843 this->_restore_array();
845 if ((this->
initialized()) && (this->_destroy_vec_on_exit))
847 PetscErrorCode
ierr=0;
849 ierr = LibMeshVecDestroy(&_vec);
850 LIBMESH_CHKERR(
ierr);
855 _global_to_local_map.clear();
860 template <
typename T>
864 parallel_object_only();
866 libmesh_assert(this->
closed());
868 this->_restore_array();
870 PetscErrorCode
ierr=0;
876 ierr = VecSet (_vec, z);
877 LIBMESH_CHKERR(
ierr);
884 ierr = VecGhostGetLocalForm (_vec,&loc_vec);
885 LIBMESH_CHKERR(
ierr);
887 ierr = VecSet (loc_vec, z);
888 LIBMESH_CHKERR(
ierr);
890 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
891 LIBMESH_CHKERR(
ierr);
897 template <
typename T>
902 cloned_vector->
init(*
this);
903 return std::unique_ptr<NumericVector<T>>(cloned_vector);
908 template <
typename T>
913 cloned_vector->
init(*
this,
true);
914 *cloned_vector = *
this;
915 return std::unique_ptr<NumericVector<T>>(cloned_vector);
920 template <
typename T>
926 PetscErrorCode
ierr=0;
927 PetscInt petsc_size=0;
932 ierr = VecGetSize(_vec, &petsc_size);
933 LIBMESH_CHKERR(
ierr);
940 template <
typename T>
946 PetscErrorCode
ierr=0;
947 PetscInt petsc_size=0;
949 ierr = VecGetLocalSize(_vec, &petsc_size);
950 LIBMESH_CHKERR(
ierr);
957 template <
typename T>
965 if (_array_is_present)
969 PetscErrorCode
ierr=0;
970 PetscInt petsc_first=0, petsc_last=0;
971 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
972 LIBMESH_CHKERR(
ierr);
981 template <
typename T>
989 if (_array_is_present)
993 PetscErrorCode
ierr=0;
994 PetscInt petsc_first=0, petsc_last=0;
995 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
996 LIBMESH_CHKERR(
ierr);
1005 template <
typename T>
1014 if (_array_is_present)
1021 PetscErrorCode
ierr=0;
1022 PetscInt petsc_first=0, petsc_last=0;
1023 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
1024 LIBMESH_CHKERR(
ierr);
1030 if ((i>=first) && (i<last))
1035 GlobalToLocalMap::const_iterator it = _global_to_local_map.find(i);
1037 const GlobalToLocalMap::const_iterator
end = _global_to_local_map.end();
1040 std::ostringstream error_message;
1041 error_message <<
"No index " << i <<
" in ghosted vector.\n" 1042 <<
"Vector contains [" << first <<
',' << last <<
")\n";
1043 GlobalToLocalMap::const_iterator b = _global_to_local_map.begin();
1046 error_message <<
"And empty ghost array.\n";
1050 error_message <<
"And ghost array {" << b->first;
1051 for (++b; b !=
end; ++b)
1052 error_message <<
',' << b->first;
1053 error_message <<
"}\n";
1056 libmesh_error_msg(error_message.str());
1058 libmesh_assert (it != _global_to_local_map.end());
1060 return it->second+last-first;
1065 template <
typename T>
1069 this->_get_array(
true);
1076 libmesh_assert_less (local_index, _local_size);
1080 return static_cast<T
>(_read_only_values[local_index]);
1085 template <
typename T>
1090 this->_get_array(
true);
1092 const std::size_t num = index.size();
1094 for (std::size_t i=0; i<num; i++)
1100 libmesh_assert_less (local_index, _local_size);
1103 values[i] =
static_cast<T
>(_read_only_values[local_index]);
1108 template <
typename T>
1112 _values_manually_retrieved =
true;
1119 template <
typename T>
1123 _values_manually_retrieved =
true;
1126 return _read_only_values;
1129 template <
typename T>
1135 _values_manually_retrieved =
false;
1139 template <
typename T>
1143 parallel_object_only();
1145 this->_restore_array();
1147 PetscErrorCode
ierr=0;
1149 PetscReal returnval=0.;
1151 ierr = VecMin (_vec, &index, &returnval);
1152 LIBMESH_CHKERR(
ierr);
1155 return static_cast<Real>(returnval);
1160 template <
typename T>
1164 parallel_object_only();
1166 this->_restore_array();
1168 PetscErrorCode
ierr=0;
1170 PetscReal returnval=0.;
1172 ierr = VecMax (_vec, &index, &returnval);
1173 LIBMESH_CHKERR(
ierr);
1176 return static_cast<Real>(returnval);
1181 template <
typename T>
1185 parallel_object_only();
1195 #ifdef LIBMESH_HAVE_CXX11_THREAD 1210 #ifdef LIBMESH_HAVE_CXX11 1212 "PETSc and libMesh integer sizes must match!");
1225 #endif // #ifdef LIBMESH_HAVE_PETSC 1226 #endif // LIBMESH_PETSC_VECTOR_H std::string name(const ElemQuality q)
virtual Real l2_norm() const override
PetscVector(const Parallel::Communicator &comm_in, const ParallelType type=AUTOMATIC)
virtual void localize(std::vector< T > &v_local) const override
NumericVector interface to PETSc Vec.
void add_vector_conjugate_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A)
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A) override
virtual void zero() override
bool _destroy_vec_on_exit
virtual numeric_index_type size() const override
virtual Real linfty_norm() const override
numeric_index_type _local_size
virtual numeric_index_type first_local_index() const override
virtual bool initialized() const
void _get_array(bool read_only) const
virtual numeric_index_type local_size() const override
virtual void add(const numeric_index_type i, const T value) override
uint8_t processor_id_type
Threads::spin_mutex _petsc_vector_mutex
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
virtual void get(const std::vector< numeric_index_type > &index, T *values) const override
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
void _restore_array() const
std::atomic< bool > _array_is_present
PetscInt * numeric_petsc_cast(const numeric_index_type *p)
std::mutex _petsc_vector_mutex
virtual void swap(NumericVector< T > &v) override
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices) override
const PetscScalar * _read_only_values
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
virtual Real l1_norm() const override
GlobalToLocalMap _global_to_local_map
virtual void clear() override
virtual void abs() override
dof_id_type numeric_index_type
const PetscScalar * get_array_read() const
virtual NumericVector< T > & operator+=(const NumericVector< T > &v) override
virtual void conjugate() override
void init(triangulateio &t)
virtual T operator()(const numeric_index_type i) const override
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const override
numeric_index_type map_global_to_local_index(const numeric_index_type i) const
virtual void create_subvector(NumericVector< T > &subvector, const std::vector< numeric_index_type > &rows) const override
virtual numeric_index_type last_local_index() const override
virtual T sum() const override
ParallelType type() const
virtual std::unique_ptr< NumericVector< T > > clone() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T indefinite_dot(const NumericVector< T > &v) const
std::unordered_map< numeric_index_type, numeric_index_type > GlobalToLocalMap
virtual void swap(NumericVector< T > &v)
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType type=AUTOMATIC) override
static PetscErrorCode Mat * A
void swap(Iterator &lhs, Iterator &rhs)
virtual Real max() const override
numeric_index_type _first
PetscVector< T > & operator=(const PetscVector< T > &v)
virtual void reciprocal() override
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices) override
bool _values_manually_retrieved
virtual void print_matlab(const std::string &name="") const override
virtual void scale(const T factor) override
PetscScalar * get_array()
virtual std::unique_ptr< NumericVector< T > > zero_clone() const override
virtual Real min() const override
virtual void close() override
virtual NumericVector< T > & operator/=(const NumericVector< T > &v) override
virtual NumericVector< T > & operator-=(const NumericVector< T > &v) override
virtual T dot(const NumericVector< T > &v) const override