21 #ifndef LIBMESH_PETSC_VECTOR_H
22 #define LIBMESH_PETSC_VECTOR_H
25 #include "libmesh/libmesh_config.h"
27 #ifdef LIBMESH_HAVE_PETSC
30 #include "libmesh/numeric_vector.h"
31 #include "libmesh/petsc_macro.h"
32 #include "libmesh/int_range.h"
33 #include "libmesh/libmesh_common.h"
34 #include "libmesh/petsc_solver_exception.h"
35 #include "libmesh/parallel_only.h"
39 # define LIBMESH_SAW_I
43 # undef I // Avoid complex.h contamination
50 #include <unordered_map>
52 #ifdef LIBMESH_HAVE_CXX11_THREAD
61 template <
typename T>
class SparseMatrix;
80 PetscVector (
const Parallel::Communicator & comm_in,
87 PetscVector (
const Parallel::Communicator & comm_in,
95 PetscVector (
const Parallel::Communicator & comm_in,
105 PetscVector (
const Parallel::Communicator & comm_in,
108 const std::vector<numeric_index_type> & ghost,
119 const Parallel::Communicator & comm_in);
138 virtual void close ()
override;
140 virtual void clear ()
override;
142 virtual void zero ()
override;
144 virtual std::unique_ptr<NumericVector<T>>
zero_clone ()
const override;
146 virtual std::unique_ptr<NumericVector<T>>
clone ()
const override;
150 const bool fast=
false,
154 const bool fast=
false,
159 const std::vector<numeric_index_type> & ghost,
160 const bool fast =
false,
164 const bool fast =
false)
override;
172 virtual Real min ()
const override;
174 virtual Real max ()
const override;
176 virtual T
sum ()
const override;
203 virtual void get(
const std::vector<numeric_index_type> & index,
204 T * values)
const override;
239 const T
value)
override;
242 const T
value)
override;
244 virtual void add (
const T s)
override;
257 const std::vector<numeric_index_type> & dof_indices)
override;
280 virtual void insert (
const T * v,
281 const std::vector<numeric_index_type> & dof_indices)
override;
283 virtual void scale (
const T factor)
override;
289 virtual void abs()
override;
300 virtual void localize (std::vector<T> & v_local)
const override;
305 const std::vector<numeric_index_type> & send_list)
const override;
307 virtual void localize (std::vector<T> & v_local,
308 const std::vector<numeric_index_type> & indices)
const override;
312 const std::vector<numeric_index_type> & send_list)
override;
323 const std::vector<numeric_index_type> & rows)
const override;
352 #ifdef LIBMESH_HAVE_CXX11_THREAD
410 #ifdef LIBMESH_HAVE_CXX11_THREAD
464 template <
typename T>
468 _array_is_present(false),
471 _local_form(nullptr),
473 _global_to_local_map(),
474 _destroy_vec_on_exit(true),
475 _values_manually_retrieved(false),
476 _values_read_only(false)
483 template <
typename T>
489 _array_is_present(false),
490 _local_form(nullptr),
492 _global_to_local_map(),
493 _destroy_vec_on_exit(true),
494 _values_manually_retrieved(false),
495 _values_read_only(false)
497 this->
init(n, n,
false, ptype);
502 template <
typename T>
509 _array_is_present(false),
510 _local_form(nullptr),
512 _global_to_local_map(),
513 _destroy_vec_on_exit(true),
514 _values_manually_retrieved(false),
515 _values_read_only(false)
517 this->
init(n, n_local,
false, ptype);
522 template <
typename T>
527 const std::vector<numeric_index_type> & ghost,
530 _array_is_present(false),
531 _local_form(nullptr),
533 _global_to_local_map(),
534 _destroy_vec_on_exit(true),
535 _values_manually_retrieved(false),
536 _values_read_only(false)
538 this->
init(n, n_local, ghost,
false, ptype);
545 template <
typename T>
548 const Parallel::Communicator & comm_in) :
550 _array_is_present(false),
551 _local_form(nullptr),
553 _global_to_local_map(),
554 _destroy_vec_on_exit(false),
555 _values_manually_retrieved(false),
556 _values_read_only(false)
559 this->_is_closed =
true;
564 PetscErrorCode
ierr=0;
565 PetscInt petsc_local_size=0;
566 ierr = VecGetLocalSize(_vec, &petsc_local_size);
567 LIBMESH_CHKERR(
ierr);
572 #if PETSC_VERSION_LESS_THAN(3,0,0) || !PETSC_VERSION_LESS_THAN(3,4,0)
578 ierr = VecGetType(_vec, &ptype);
579 LIBMESH_CHKERR(
ierr);
581 if ((std::strcmp(ptype,VECSHARED) == 0) || (std::strcmp(ptype,VECMPI) == 0))
583 #if PETSC_RELEASE_LESS_THAN(3,1,1)
584 ISLocalToGlobalMapping mapping = _vec->mapping;
586 ISLocalToGlobalMapping mapping;
587 ierr = VecGetLocalToGlobalMapping(_vec, &mapping);
588 LIBMESH_CHKERR(
ierr);
594 const numeric_index_type my_local_size = static_cast<numeric_index_type>(petsc_local_size);
595 const numeric_index_type ghost_begin = static_cast<numeric_index_type>(petsc_local_size);
596 #if PETSC_RELEASE_LESS_THAN(3,4,0)
600 ierr = ISLocalToGlobalMappingGetSize(mapping, &n);
601 LIBMESH_CHKERR(
ierr);
604 #if PETSC_RELEASE_LESS_THAN(3,1,1)
605 const PetscInt * indices = mapping->indices;
607 const PetscInt * indices;
608 ierr = ISLocalToGlobalMappingGetIndices(mapping,&indices);
609 LIBMESH_CHKERR(
ierr);
612 _global_to_local_map[indices[i]] = i-my_local_size;
614 #if !PETSC_RELEASE_LESS_THAN(3,1,1)
615 ierr = ISLocalToGlobalMappingRestoreIndices(mapping, &indices);
616 LIBMESH_CHKERR(
ierr);
629 template <
typename T>
638 template <
typename T>
645 parallel_object_only();
647 PetscErrorCode
ierr=0;
648 PetscInt petsc_n=static_cast<PetscInt>(n);
668 if (this->_type ==
SERIAL)
670 ierr = VecCreateSeq (PETSC_COMM_SELF, petsc_n, &_vec);
671 CHKERRABORT(PETSC_COMM_SELF,
ierr);
673 ierr = VecSetFromOptions (_vec);
674 CHKERRABORT(PETSC_COMM_SELF,
ierr);
679 #ifdef LIBMESH_HAVE_MPI
680 PetscInt petsc_n_local=cast_int<PetscInt>(n_local);
681 libmesh_assert_less_equal (n_local, n);
682 ierr = VecCreateMPI (this->comm().
get(), petsc_n_local, petsc_n,
684 LIBMESH_CHKERR(
ierr);
686 libmesh_assert_equal_to (n_local, n);
687 ierr = VecCreateSeq (PETSC_COMM_SELF, petsc_n, &_vec);
688 CHKERRABORT(PETSC_COMM_SELF,
ierr);
691 ierr = VecSetFromOptions (_vec);
692 LIBMESH_CHKERR(
ierr);
695 libmesh_error_msg(
"Unsupported type " << this->_type);
698 this->_is_closed =
true;
707 template <
typename T>
713 this->
init(n,n,fast,ptype);
718 template <
typename T>
722 const std::vector<numeric_index_type> & ghost,
726 parallel_object_only();
728 PetscErrorCode
ierr=0;
729 PetscInt petsc_n=static_cast<PetscInt>(n);
730 PetscInt petsc_n_local=static_cast<PetscInt>(n_local);
731 PetscInt petsc_n_ghost=static_cast<PetscInt>(ghost.size());
746 PetscInt * petsc_ghost = ghost.empty() ? PETSC_NULL :
747 const_cast<PetscInt *>(reinterpret_cast<const PetscInt *>(ghost.data()));
759 _global_to_local_map[ghost[i]] = i;
763 ierr = VecCreateGhost (this->comm().
get(), petsc_n_local, petsc_n,
764 petsc_n_ghost, petsc_ghost,
766 LIBMESH_CHKERR(
ierr);
768 ierr = VecSetFromOptions (_vec);
769 LIBMESH_CHKERR(
ierr);
772 this->_is_closed =
true;
780 template <
typename T>
785 parallel_object_only();
791 const PetscVector<T> & v = cast_ref<const PetscVector<T> &>(other);
804 this->_is_closed =
true;
807 this->_type = v.
_type;
810 PetscErrorCode
ierr = VecDuplicate (v.
_vec, &this->_vec);
811 LIBMESH_CHKERR(
ierr);
819 template <
typename T>
823 parallel_object_only();
825 this->_restore_array();
827 PetscErrorCode
ierr=0;
829 ierr = VecAssemblyBegin(_vec);
830 LIBMESH_CHKERR(
ierr);
831 ierr = VecAssemblyEnd(_vec);
832 LIBMESH_CHKERR(
ierr);
836 ierr = VecGhostUpdateBegin(_vec,INSERT_VALUES,SCATTER_FORWARD);
837 LIBMESH_CHKERR(
ierr);
838 ierr = VecGhostUpdateEnd(_vec,INSERT_VALUES,SCATTER_FORWARD);
839 LIBMESH_CHKERR(
ierr);
842 this->_is_closed =
true;
847 template <
typename T>
851 parallel_object_only();
854 this->_restore_array();
856 if ((this->
initialized()) && (this->_destroy_vec_on_exit))
858 PetscErrorCode
ierr=0;
860 ierr = LibMeshVecDestroy(&_vec);
861 LIBMESH_CHKERR(
ierr);
866 _global_to_local_map.clear();
871 template <
typename T>
875 parallel_object_only();
879 this->_restore_array();
881 PetscErrorCode
ierr=0;
887 ierr = VecSet (_vec, z);
888 LIBMESH_CHKERR(
ierr);
895 ierr = VecGhostGetLocalForm (_vec,&loc_vec);
896 LIBMESH_CHKERR(
ierr);
898 ierr = VecSet (loc_vec, z);
899 LIBMESH_CHKERR(
ierr);
901 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
902 LIBMESH_CHKERR(
ierr);
908 template <
typename T>
913 cloned_vector->
init(*
this);
914 return std::unique_ptr<NumericVector<T>>(cloned_vector);
919 template <
typename T>
924 cloned_vector->
init(*
this,
true);
925 *cloned_vector = *
this;
926 return std::unique_ptr<NumericVector<T>>(cloned_vector);
931 template <
typename T>
937 PetscErrorCode
ierr=0;
938 PetscInt petsc_size=0;
943 ierr = VecGetSize(_vec, &petsc_size);
944 LIBMESH_CHKERR(
ierr);
946 return static_cast<numeric_index_type>(petsc_size);
951 template <
typename T>
957 PetscErrorCode
ierr=0;
958 PetscInt petsc_size=0;
960 ierr = VecGetLocalSize(_vec, &petsc_size);
961 LIBMESH_CHKERR(
ierr);
963 return static_cast<numeric_index_type>(petsc_size);
968 template <
typename T>
976 if (_array_is_present)
980 PetscErrorCode
ierr=0;
981 PetscInt petsc_first=0, petsc_last=0;
982 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
983 LIBMESH_CHKERR(
ierr);
984 first = static_cast<numeric_index_type>(petsc_first);
992 template <
typename T>
1000 if (_array_is_present)
1004 PetscErrorCode
ierr=0;
1005 PetscInt petsc_first=0, petsc_last=0;
1006 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
1007 LIBMESH_CHKERR(
ierr);
1008 last = static_cast<numeric_index_type>(petsc_last);
1016 template <
typename T>
1025 if (_array_is_present)
1032 PetscErrorCode
ierr=0;
1033 PetscInt petsc_first=0, petsc_last=0;
1034 ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
1035 LIBMESH_CHKERR(
ierr);
1036 first = static_cast<numeric_index_type>(petsc_first);
1037 last = static_cast<numeric_index_type>(petsc_last);
1041 if ((i>=first) && (i<last))
1046 GlobalToLocalMap::const_iterator it = _global_to_local_map.find(i);
1048 const GlobalToLocalMap::const_iterator
end = _global_to_local_map.end();
1051 std::ostringstream error_message;
1052 error_message <<
"No index " << i <<
" in ghosted vector.\n"
1053 <<
"Vector contains [" << first <<
',' << last <<
")\n";
1054 GlobalToLocalMap::const_iterator b = _global_to_local_map.begin();
1057 error_message <<
"And empty ghost array.\n";
1061 error_message <<
"And ghost array {" << b->first;
1062 for (++b; b !=
end; ++b)
1063 error_message <<
',' << b->first;
1064 error_message <<
"}\n";
1067 libmesh_error_msg(error_message.str());
1071 return it->second+last-first;
1076 template <
typename T>
1080 this->_get_array(
true);
1087 libmesh_assert_less (local_index, _local_size);
1091 return static_cast<T>(_read_only_values[local_index]);
1096 template <
typename T>
1101 this->_get_array(
true);
1103 const std::size_t num = index.size();
1105 for (std::size_t i=0; i<num; i++)
1111 libmesh_assert_less (local_index, _local_size);
1114 values[i] = static_cast<T>(_read_only_values[local_index]);
1119 template <
typename T>
1123 _values_manually_retrieved =
true;
1130 template <
typename T>
1134 _values_manually_retrieved =
true;
1137 return _read_only_values;
1140 template <
typename T>
1146 _values_manually_retrieved =
false;
1150 template <
typename T>
1154 parallel_object_only();
1156 this->_restore_array();
1158 PetscErrorCode
ierr=0;
1160 PetscReal returnval=0.;
1162 ierr = VecMin (_vec, &index, &returnval);
1163 LIBMESH_CHKERR(
ierr);
1166 return static_cast<Real>(returnval);
1171 template <
typename T>
1175 parallel_object_only();
1177 this->_restore_array();
1179 PetscErrorCode
ierr=0;
1181 PetscReal returnval=0.;
1183 ierr = VecMax (_vec, &index, &returnval);
1184 LIBMESH_CHKERR(
ierr);
1187 return static_cast<Real>(returnval);
1192 template <
typename T>
1196 parallel_object_only();
1206 #ifdef LIBMESH_HAVE_CXX11_THREAD
1221 #ifdef LIBMESH_HAVE_CXX11
1223 "PETSc and libMesh integer sizes must match!");
1230 return reinterpret_cast<PetscInt *>(const_cast<numeric_index_type *>(p));
1236 #endif // #ifdef LIBMESH_HAVE_PETSC
1237 #endif // LIBMESH_PETSC_VECTOR_H