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/base/multithread_info.h> 23 # include <deal.II/lac/exceptions.h> 24 # include <deal.II/lac/petsc_compatibility.h> 25 # include <deal.II/lac/petsc_vector.h> 29 DEAL_II_NAMESPACE_OPEN
35 VectorReference::operator PetscScalar()
const 49 VecGetOwnershipRange(vector.vector, &begin, &end);
52 Vec locally_stored_elements = PETSC_NULL;
53 ierr = VecGhostGetLocalForm(vector.vector, &locally_stored_elements);
57 ierr = VecGetSize(locally_stored_elements, &lsize);
61 ierr = VecGetArray(locally_stored_elements, &ptr);
66 if (index >= static_cast<size_type>(begin) &&
67 index < static_cast<size_type>(end))
70 value = *(ptr + index - begin);
75 const size_type ghostidx =
76 vector.ghost_indices.index_within_set(index);
79 value = *(ptr + ghostidx + end - begin);
82 ierr = VecRestoreArray(locally_stored_elements, &ptr);
86 VecGhostRestoreLocalForm(vector.vector, &locally_stored_elements);
98 PetscErrorCode ierr = VecGetOwnershipRange(vector.vector, &begin, &end);
101 AssertThrow((index >= static_cast<size_type>(begin)) &&
102 (index < static_cast<size_type>(end)),
103 ExcAccessToNonlocalElement(index, begin, end - 1));
105 PetscInt idx = index;
107 ierr = VecGetValues(vector.vector, 1, &idx, &value);
118 , obtained_ownership(true)
121 ExcMessage(
"PETSc does not support multi-threaded access, set " 122 "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()."));
152 , obtained_ownership(false)
155 ExcMessage(
"PETSc does not support multi-threaded access, set " 156 "the thread limit to 1 in MPI_InitFinalize()."));
165 const PetscErrorCode ierr = VecDestroy(&
vector);
178 const PetscErrorCode ierr = VecDestroy(&
vector);
197 PetscErrorCode ierr = VecSet(
vector, s);
202 Vec ghost = PETSC_NULL;
203 ierr = VecGhostGetLocalForm(
vector, &ghost);
206 ierr = VecSet(ghost, s);
209 ierr = VecGhostRestoreLocalForm(
vector, &ghost);
224 const PetscErrorCode ierr = VecEqual(
vector, v.
vector, &flag);
227 return (flag == PETSC_TRUE);
238 const PetscErrorCode ierr = VecEqual(
vector, v.
vector, &flag);
241 return (flag == PETSC_FALSE);
246 VectorBase::size_type
250 const PetscErrorCode ierr = VecGetSize(
vector, &sz);
258 VectorBase::size_type
262 const PetscErrorCode ierr = VecGetLocalSize(
vector, &sz);
270 std::pair<VectorBase::size_type, VectorBase::size_type>
274 const PetscErrorCode ierr =
275 VecGetOwnershipRange(static_cast<const Vec &>(
vector), &begin, &end);
278 return std::make_pair(begin, end);
285 const std::vector<PetscScalar> &values)
287 Assert(indices.size() == values.size(),
288 ExcMessage(
"Function called with arguments of different sizes"));
296 const std::vector<PetscScalar> &values)
298 Assert(indices.size() == values.size(),
299 ExcMessage(
"Function called with arguments of different sizes"));
307 const ::Vector<PetscScalar> &values)
309 Assert(indices.size() == values.size(),
310 ExcMessage(
"Function called with arguments of different sizes"));
318 const size_type * indices,
319 const PetscScalar *values)
339 const PetscErrorCode ierr = VecDot(vec.
vector,
vector, &result);
363 # ifdef DEAL_II_WITH_MPI 367 int all_int_last_action;
369 const int ierr = MPI_Allreduce(&my_int_last_action,
370 &all_int_last_action,
379 ExcMessage(
"Error: not all processors agree on the last " 380 "VectorOperation before this compress() call."));
389 "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, sum1 = 0, sum2 = 0, sum3 = 0;
453 const PetscScalar *ptr = start_ptr;
454 const PetscScalar *eptr = ptr + (
size() / 4) * 4;
463 while (ptr != start_ptr +
size())
466 mean = (sum0 + sum1 + sum2 + sum3) / static_cast<PetscReal>(
size());
471 ierr = VecRestoreArray(
vector, &start_ptr);
478 VectorBase::real_type
483 const PetscErrorCode ierr = VecNorm(
vector, NORM_1, &d);
491 VectorBase::real_type
496 const PetscErrorCode ierr = VecNorm(
vector, NORM_2, &d);
504 VectorBase::real_type
509 PetscScalar * start_ptr;
510 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
515 real_type sum0 = 0, sum1 = 0, sum2 = 0, sum3 = 0;
520 const PetscScalar *ptr = start_ptr;
521 const PetscScalar *eptr = ptr + (
size() / 4) * 4;
530 while (ptr != start_ptr +
size())
533 norm = std::pow(sum0 + sum1 + sum2 + sum3, 1. / p);
538 ierr = VecRestoreArray(
vector, &start_ptr);
546 VectorBase::real_type
551 const PetscErrorCode ierr = VecNorm(
vector, NORM_INFINITY, &d);
559 VectorBase::real_type
565 const PetscErrorCode ierr = VecMin(
vector, &p, &d);
572 VectorBase::real_type
578 const PetscErrorCode ierr = VecMax(
vector, &p, &d);
591 PetscScalar * start_ptr;
592 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
595 const PetscScalar *ptr = start_ptr, *eptr = start_ptr +
local_size();
609 ierr = VecRestoreArray(
vector, &start_ptr);
618 template <
typename T>
620 is_non_negative(
const T &t)
627 template <
typename T>
629 is_non_negative(
const std::complex<T> &)
633 "whether it is non-negative.")) return true;
644 PetscScalar * start_ptr;
645 PetscErrorCode ierr = VecGetArray(vector, &start_ptr);
648 const PetscScalar *ptr = start_ptr, *eptr = start_ptr + local_size();
652 if (!internal::is_non_negative(*ptr))
662 ierr = VecRestoreArray(vector, &start_ptr);
676 const PetscErrorCode ierr = VecScale(vector, a);
690 const PetscScalar factor = 1. / a;
693 const PetscErrorCode ierr = VecScale(vector, factor);
705 const PetscErrorCode ierr = VecAXPY(vector, 1, v);
717 const PetscErrorCode ierr = VecAXPY(vector, -1, v);
731 const PetscErrorCode ierr = VecShift(vector, s);
743 const PetscErrorCode ierr = VecAXPY(vector, a, v);
759 const PetscScalar weights[2] = {a, b};
760 Vec addends[2] = {v.
vector, w.vector};
762 const PetscErrorCode ierr = VecMAXPY(vector, 2, weights, addends);
774 const PetscErrorCode ierr = VecAYPX(vector, s, v);
802 const PetscErrorCode ierr = VecPointwiseMult(vector, factors, vector);
819 const PetscErrorCode ierr = VecCopy(v.
vector, vector);
831 const PetscErrorCode ierr = VecPointwiseDivide(vector, a, b);
843 PetscErrorCode ierr =
844 PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD, format);
848 ierr = VecView(vector, PETSC_VIEWER_STDOUT_WORLD);
856 const unsigned int precision,
857 const bool scientific,
858 const bool across)
const 865 PetscErrorCode ierr = VecGetArray(vector, &val);
870 const std::ios::fmtflags old_flags = out.flags();
871 const unsigned int old_precision = out.precision(precision);
873 out.precision(precision);
875 out.setf(std::ios::scientific, std::ios::floatfield);
877 out.setf(std::ios::fixed, std::ios::floatfield);
880 for (size_type i = 0; i < local_size(); ++i)
881 out << val[i] <<
' ';
883 for (size_type i = 0; i < local_size(); ++i)
884 out << val[i] << std::endl;
888 out.flags(old_flags);
889 out.precision(old_precision);
893 ierr = VecRestoreArray(vector, &val);
904 const PetscErrorCode ierr = VecSwap(vector, v.
vector);
910 VectorBase::operator
const Vec &()
const 919 std::size_t mem =
sizeof(Vec) +
sizeof(last_action) +
926 mem += local_size() *
sizeof(PetscScalar);
930 mem += ghost_indices.n_elements() * (
sizeof(PetscScalar) +
sizeof(
int));
940 const size_type * indices,
941 const PetscScalar *values,
942 const bool add_values)
947 Assert((last_action == action) ||
949 internal::VectorReference::ExcWrongMode(action, last_action));
959 const PetscInt *petsc_indices =
960 reinterpret_cast<const PetscInt *
>(indices);
962 const InsertMode mode = (add_values ? ADD_VALUES : INSERT_VALUES);
963 const PetscErrorCode ierr =
964 VecSetValues(vector, n_elements, petsc_indices, values, mode);
970 last_action = action;
975 DEAL_II_NAMESPACE_CLOSE
977 #endif // DEAL_II_WITH_PETSC
#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)
#define AssertIndexRange(index, range)
VectorBase & operator/=(const PetscScalar factor)
#define AssertThrow(cond, exc)
real_type norm_sqr() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void scale(const VectorBase &scaling_factors)
real_type linfty_norm() const
VectorBase & operator=(const VectorBase &)=delete
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
virtual ~VectorBase() override
#define AssertIsFinite(number)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
bool operator!=(const VectorBase &v) const
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)