16 #include <deal.II/lac/petsc_vector_base.h> 18 #ifdef DEAL_II_WITH_PETSC 20 # include <deal.II/base/memory_consumption.h> 21 # include <deal.II/lac/exceptions.h> 22 # include <deal.II/lac/petsc_compatibility.h> 23 # include <deal.II/lac/petsc_parallel_vector.h> 25 # include <deal.II/base/multithread_info.h> 27 DEAL_II_NAMESPACE_OPEN
33 VectorReference::operator PetscScalar ()
const 35 Assert (index < vector.size(),
47 PetscErrorCode ierr = VecGetOwnershipRange (vector.vector, &begin,
51 Vec locally_stored_elements = PETSC_NULL;
52 ierr = VecGhostGetLocalForm(vector.vector, &locally_stored_elements);
56 ierr = VecGetSize(locally_stored_elements, &lsize);
60 ierr = VecGetArray(locally_stored_elements, &ptr);
65 if ( index>=static_cast<size_type>(begin)
66 && index<static_cast<size_type>(end) )
69 value = *(ptr+index-begin);
75 = vector.ghost_indices.index_within_set(index);
78 value = *(ptr+ghostidx+end-begin);
81 ierr = VecRestoreArray(locally_stored_elements, &ptr);
84 ierr = VecGhostRestoreLocalForm(vector.vector, &locally_stored_elements);
96 PetscErrorCode ierr = VecGetOwnershipRange (vector.vector, &begin, &end);
99 AssertThrow ((index >= static_cast<size_type>(begin)) &&
100 (index < static_cast<size_type>(end)),
101 ExcAccessToNonlocalElement (index, begin, end-1));
103 PetscInt idx = index;
105 ierr = VecGetValues(vector.vector, 1, &idx, &value);
117 obtained_ownership(true)
120 ExcMessage(
"PETSc does not support multi-threaded access, set " 121 "the thread limit to 1 in MPI_InitFinalize()."));
130 ghost_indices(v.ghost_indices),
132 obtained_ownership(true)
135 ExcMessage(
"PETSc does not support multi-threaded access, set " 136 "the thread limit to 1 in MPI_InitFinalize()."));
153 obtained_ownership(false)
156 ExcMessage(
"PETSc does not support multi-threaded access, set " 157 "the thread limit to 1 in MPI_InitFinalize()."));
166 const PetscErrorCode ierr = VecDestroy (&
vector);
179 const PetscErrorCode ierr = VecDestroy (&
vector);
198 PetscErrorCode ierr = VecSet (
vector, s);
203 Vec ghost = PETSC_NULL;
204 ierr = VecGhostGetLocalForm(
vector, &ghost);
207 ierr = VecSet (ghost, s);
210 ierr = VecGhostRestoreLocalForm(
vector, &ghost);
226 const PetscErrorCode ierr = VecEqual (
vector, v.
vector, &flag);
229 return (flag == PETSC_TRUE);
241 const PetscErrorCode ierr = VecEqual (
vector, v.
vector, &flag);
244 return (flag == PETSC_FALSE);
249 VectorBase::size_type
253 const PetscErrorCode ierr = VecGetSize (
vector, &sz);
261 VectorBase::size_type
265 const PetscErrorCode ierr = VecGetLocalSize (
vector, &sz);
273 std::pair<VectorBase::size_type, VectorBase::size_type>
277 const PetscErrorCode ierr = VecGetOwnershipRange (static_cast<const Vec &>(
vector),
281 return std::make_pair (begin, end);
288 const std::vector<PetscScalar> &values)
290 Assert (indices.size() == values.size(),
291 ExcMessage (
"Function called with arguments of different sizes"));
299 const std::vector<PetscScalar> &values)
301 Assert (indices.size() == values.size(),
302 ExcMessage (
"Function called with arguments of different sizes"));
310 const ::Vector<PetscScalar> &values)
312 Assert (indices.size() == values.size(),
313 ExcMessage (
"Function called with arguments of different sizes"));
321 const size_type *indices,
322 const PetscScalar *values)
344 const PetscErrorCode ierr = VecDot (vec.
vector,
vector, &result);
368 #ifdef DEAL_II_WITH_MPI 372 int all_int_last_action;
374 const int ierr = MPI_Allreduce
375 (&my_int_last_action, &all_int_last_action, 1, MPI_INT, MPI_BOR,
381 ExcMessage(
"Error: not all processors agree on the last " 382 "VectorOperation before this compress() call."));
389 ExcMessage(
"Missing compress() or calling with wrong VectorOperation argument."));
405 PetscErrorCode ierr = VecAssemblyBegin(
vector);
407 ierr = VecAssemblyEnd(
vector);
418 VectorBase::real_type
432 if (dynamic_cast<const PETScWrappers::MPI::Vector *>(
this) !=
nullptr)
435 const PetscErrorCode ierr = VecSum(
vector, &sum);
437 return sum/
static_cast<PetscReal
>(
size());
442 PetscScalar *start_ptr;
443 PetscErrorCode ierr = VecGetArray (
vector, &start_ptr);
446 PetscScalar mean = 0;
448 PetscScalar sum0 = 0,
456 const PetscScalar *ptr = start_ptr;
457 const PetscScalar *eptr = ptr + (
size()/4)*4;
466 while (ptr != start_ptr+
size())
469 mean = (sum0+sum1+sum2+sum3)/static_cast<PetscReal>(
size());
474 ierr = VecRestoreArray (
vector, &start_ptr);
481 VectorBase::real_type
486 const PetscErrorCode ierr = VecNorm (
vector, NORM_1, &d);
494 VectorBase::real_type
499 const PetscErrorCode ierr = VecNorm (
vector, NORM_2, &d);
507 VectorBase::real_type
512 PetscScalar *start_ptr;
513 PetscErrorCode ierr = VecGetArray (
vector, &start_ptr);
526 const PetscScalar *ptr = start_ptr;
527 const PetscScalar *eptr = ptr + (
size()/4)*4;
536 while (ptr != start_ptr+
size())
539 norm = std::pow(sum0+sum1+sum2+sum3,
545 ierr = VecRestoreArray (
vector, &start_ptr);
553 VectorBase::real_type
558 const PetscErrorCode ierr = VecNorm (
vector, NORM_INFINITY, &d);
566 VectorBase::real_type
572 const PetscErrorCode ierr = VecMin (
vector, &p, &d);
579 VectorBase::real_type
585 const PetscErrorCode ierr = VecMax (
vector, &p, &d);
598 PetscScalar *start_ptr;
599 PetscErrorCode ierr = VecGetArray (
vector, &start_ptr);
602 const PetscScalar *ptr = start_ptr,
617 ierr = VecRestoreArray (
vector, &start_ptr);
626 template <
typename T>
627 bool is_non_negative (
const T &t)
634 template <
typename T>
635 bool is_non_negative (
const std::complex<T> &)
639 "whether it is non-negative."))
651 PetscScalar *start_ptr;
652 PetscErrorCode ierr = VecGetArray (vector, &start_ptr);
655 const PetscScalar *ptr = start_ptr,
656 *eptr = start_ptr + local_size();
660 if (! internal::is_non_negative(*ptr))
670 ierr = VecRestoreArray (vector, &start_ptr);
684 const PetscErrorCode ierr = VecScale (vector, a);
698 const PetscScalar factor = 1./a;
701 const PetscErrorCode ierr = VecScale (vector, factor);
713 const PetscErrorCode ierr = VecAXPY (vector, 1, v);
725 const PetscErrorCode ierr = VecAXPY (vector, -1, v);
739 const PetscErrorCode ierr = VecShift (vector, s);
752 const PetscErrorCode ierr = VecAXPY (vector, a, v);
768 const PetscScalar weights[2] = {a,b};
769 Vec addends[2] = {v.
vector, w.vector};
771 const PetscErrorCode ierr = VecMAXPY (vector, 2, weights, addends);
784 const PetscErrorCode ierr = VecAYPX (vector, s, v);
812 const PetscErrorCode ierr
813 = VecPointwiseMult (vector, factors, vector);
832 const PetscErrorCode ierr = VecCopy (v.
vector, vector);
845 const PetscErrorCode ierr = VecPointwiseDivide (vector, a, b);
857 PetscErrorCode ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
862 ierr = VecView (vector, PETSC_VIEWER_STDOUT_WORLD);
870 const unsigned int precision,
871 const bool scientific,
872 const bool across)
const 879 PetscErrorCode ierr = VecGetArray (vector, &val);
884 const std::ios::fmtflags old_flags = out.flags();
885 const unsigned int old_precision = out.precision (precision);
887 out.precision (precision);
889 out.setf (std::ios::scientific, std::ios::floatfield);
891 out.setf (std::ios::fixed, std::ios::floatfield);
894 for (size_type i=0; i<local_size(); ++i)
895 out << val[i] <<
' ';
897 for (size_type i=0; i<local_size(); ++i)
898 out << val[i] << std::endl;
902 out.flags (old_flags);
903 out.precision(old_precision);
907 ierr = VecRestoreArray (vector, &val);
918 const PetscErrorCode ierr = VecSwap (vector, v.
vector);
924 VectorBase::operator
const Vec &()
const 933 std::size_t mem =
sizeof(Vec)+
sizeof(last_action)
940 mem += local_size()*
sizeof(PetscScalar);
944 mem += ghost_indices.n_elements()*(
sizeof(PetscScalar)+
sizeof(
int));
954 const size_type *indices,
955 const PetscScalar *values,
956 const bool add_values)
961 Assert ((last_action == action)
964 internal::VectorReference::ExcWrongMode (action,
975 const PetscInt *petsc_indices =
reinterpret_cast<const PetscInt *
>(indices);
977 const InsertMode mode = (add_values ? ADD_VALUES : INSERT_VALUES);
978 const PetscErrorCode ierr = VecSetValues (vector, n_elements, petsc_indices,
985 last_action = action;
990 DEAL_II_NAMESPACE_CLOSE
992 #endif // DEAL_II_WITH_PETSC types::global_dof_index size_type
#define AssertNothrow(cond, exc)
VectorBase & operator-=(const VectorBase &V)
static ::ExceptionBase & ExcIO()
void sadd(const PetscScalar s, const VectorBase &V)
real_type l2_norm() const
void ratio(const VectorBase &a, const VectorBase &b)
VectorBase & operator/=(const PetscScalar factor)
#define AssertThrow(cond, exc)
real_type norm_sqr() const
VectorBase & operator=(const PetscScalar s)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void scale(const VectorBase &scaling_factors)
real_type linfty_norm() const
std::pair< size_type, size_type > local_range() const
real_type lp_norm(const real_type p) const
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type local_size() const
void compress(const VectorOperation::values operation)
void do_set_add_operation(const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
VectorOperation::values last_action
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static bool is_running_single_threaded()
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
std::size_t memory_consumption() const
void equ(const PetscScalar a, const VectorBase &V)
VectorBase & operator*=(const PetscScalar factor)
real_type l1_norm() const
static ::ExceptionBase & ExcGhostsPresent()
VectorBase & operator+=(const VectorBase &V)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
#define AssertThrowMPI(error_code)
PetscScalar operator*(const VectorBase &vec) const
bool has_ghost_elements() const
PetscScalar mean_value() const
bool operator==(const VectorBase &v) const
void add(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
virtual const MPI_Comm & get_mpi_communicator() const
#define AssertIsFinite(number)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()
bool operator!=(const VectorBase &v) const
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)