18#ifdef DEAL_II_WITH_PETSC
36 VectorReference::operator PetscScalar()
const
50 VecGetOwnershipRange(vector.vector, &
begin, &
end);
53 Vec locally_stored_elements = PETSC_NULL;
54 ierr = VecGhostGetLocalForm(vector.vector, &locally_stored_elements);
58 ierr = VecGetSize(locally_stored_elements, &lsize);
62 ierr = VecGetArray(locally_stored_elements, &ptr);
71 value = *(ptr + index -
begin);
77 vector.ghost_indices.index_within_set(index);
80 value = *(ptr + ghostidx +
end -
begin);
83 ierr = VecRestoreArray(locally_stored_elements, &ptr);
87 VecGhostRestoreLocalForm(vector.vector, &locally_stored_elements);
99 PetscErrorCode ierr = VecGetOwnershipRange(vector.vector, &
begin, &
end);
104 ExcAccessToNonlocalElement(index,
begin,
end - 1));
106 PetscInt idx = index;
108 ierr = VecGetValues(vector.vector, 1, &idx, &value);
120 , obtained_ownership(true)
123 ExcMessage(
"PETSc does not support multi-threaded access, set "
124 "the thread limit to 1 in MPI_InitFinalize()."));
132 , ghost_indices(v.ghost_indices)
134 , obtained_ownership(true)
137 ExcMessage(
"PETSc does not support multi-threaded access, set "
138 "the thread limit to 1 in MPI_InitFinalize()."));
154 , obtained_ownership(false)
157 ExcMessage(
"PETSc does not support multi-threaded access, set "
158 "the thread limit to 1 in MPI_InitFinalize()."));
167 const PetscErrorCode ierr = VecDestroy(&
vector);
180 const PetscErrorCode ierr = VecDestroy(&
vector);
199 PetscErrorCode ierr = VecSet(
vector, s);
204 Vec ghost = PETSC_NULL;
205 ierr = VecGhostGetLocalForm(
vector, &ghost);
208 ierr = VecSet(ghost, s);
211 ierr = VecGhostRestoreLocalForm(
vector, &ghost);
226 const PetscErrorCode ierr = VecEqual(
vector, v.
vector, &flag);
229 return (flag == PETSC_TRUE);
240 const PetscErrorCode ierr = VecEqual(
vector, v.
vector, &flag);
243 return (flag == PETSC_FALSE);
252 const PetscErrorCode ierr = VecGetSize(
vector, &sz);
264 const PetscErrorCode ierr = VecGetLocalSize(
vector, &sz);
276 const PetscErrorCode ierr = VecGetLocalSize(
vector, &sz);
284 std::pair<VectorBase::size_type, VectorBase::size_type>
288 const PetscErrorCode ierr =
289 VecGetOwnershipRange(
static_cast<const Vec &
>(
vector), &
begin, &
end);
299 const std::vector<PetscScalar> &
values)
302 ExcMessage(
"Function called with arguments of different sizes"));
310 const std::vector<PetscScalar> &
values)
313 ExcMessage(
"Function called with arguments of different sizes"));
321 const ::Vector<PetscScalar> &
values)
324 ExcMessage(
"Function called with arguments of different sizes"));
333 const PetscScalar *
values)
353 const PetscErrorCode ierr = VecDot(vec.
vector,
vector, &result);
377# ifdef DEAL_II_WITH_MPI
381 int all_int_last_action;
383 const int ierr = MPI_Allreduce(&my_int_last_action,
384 &all_int_last_action,
393 ExcMessage(
"Error: not all processors agree on the last "
394 "VectorOperation before this compress() call."));
403 "Missing compress() or calling with wrong VectorOperation argument."));
419 PetscErrorCode ierr = VecAssemblyBegin(
vector);
421 ierr = VecAssemblyEnd(
vector);
449 const PetscErrorCode ierr = VecSum(
vector, &
sum);
451 return sum /
static_cast<PetscReal
>(
size());
456 PetscScalar * start_ptr;
457 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
460 PetscScalar
mean = 0;
462 PetscScalar sum0 = 0, sum1 = 0, sum2 = 0, sum3 = 0;
467 const PetscScalar *ptr = start_ptr;
468 const PetscScalar *eptr = ptr + (
size() / 4) * 4;
477 while (ptr != start_ptr +
size())
480 mean = (sum0 + sum1 + sum2 + sum3) /
static_cast<PetscReal
>(
size());
485 ierr = VecRestoreArray(
vector, &start_ptr);
497 const PetscErrorCode ierr = VecNorm(
vector, NORM_1, &
d);
510 const PetscErrorCode ierr = VecNorm(
vector, NORM_2, &
d);
523 PetscScalar * start_ptr;
524 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
529 real_type sum0 = 0, sum1 = 0, sum2 = 0, sum3 = 0;
534 const PetscScalar *ptr = start_ptr;
535 const PetscScalar *eptr = ptr + (
size() / 4) * 4;
544 while (ptr != start_ptr +
size())
552 ierr = VecRestoreArray(
vector, &start_ptr);
565 const PetscErrorCode ierr = VecNorm(
vector, NORM_INFINITY, &
d);
579 const PetscErrorCode ierr = VecMin(
vector, &p, &
d);
592 const PetscErrorCode ierr = VecMax(
vector, &p, &
d);
605 PetscScalar * start_ptr;
606 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
609 const PetscScalar *ptr = start_ptr,
624 ierr = VecRestoreArray(
vector, &start_ptr);
633 template <
typename T>
642 template <
typename T>
648 "whether it is non-negative."))
return true;
659 PetscScalar * start_ptr;
660 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
663 const PetscScalar *ptr = start_ptr,
678 ierr = VecRestoreArray(
vector, &start_ptr);
692 const PetscErrorCode ierr = VecScale(
vector, a);
706 const PetscScalar factor = 1. / a;
709 const PetscErrorCode ierr = VecScale(
vector, factor);
721 const PetscErrorCode ierr = VecAXPY(
vector, 1, v);
733 const PetscErrorCode ierr = VecAXPY(
vector, -1, v);
747 const PetscErrorCode ierr = VecShift(
vector, s);
759 const PetscErrorCode ierr = VecAXPY(
vector, a, v);
775 const PetscScalar weights[2] = {a,
b};
776 Vec addends[2] = {v.
vector,
w.vector};
778 const PetscErrorCode ierr = VecMAXPY(
vector, 2, weights, addends);
790 const PetscErrorCode ierr = VecAYPX(
vector, s, v);
818 const PetscErrorCode ierr = VecPointwiseMult(
vector, factors,
vector);
849 PetscErrorCode ierr =
850 PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD, format);
854 ierr = VecView(
vector, PETSC_VIEWER_STDOUT_WORLD);
862 const unsigned int precision,
863 const bool scientific,
864 const bool across)
const
871 PetscErrorCode ierr = VecGetArray(
vector, &val);
876 const std::ios::fmtflags old_flags = out.flags();
877 const unsigned int old_precision = out.precision(precision);
879 out.precision(precision);
881 out.setf(std::ios::scientific, std::ios::floatfield);
883 out.setf(std::ios::fixed, std::ios::floatfield);
887 out << val[i] <<
' ';
890 out << val[i] << std::endl;
894 out.flags(old_flags);
895 out.precision(old_precision);
899 ierr = VecRestoreArray(
vector, &val);
916 VectorBase::operator
const Vec &()
const
925 std::size_t mem =
sizeof(Vec) +
sizeof(
last_action) +
947 const PetscScalar *
values,
948 const bool add_values)
955 internal::VectorReference::ExcWrongMode(action,
last_action));
965 const PetscInt *petsc_indices =
966 reinterpret_cast<const PetscInt *
>(indices);
968 const InsertMode mode = (add_values ? ADD_VALUES : INSERT_VALUES);
969 const PetscErrorCode ierr =
970 VecSetValues(
vector, n_elements, petsc_indices,
values, mode);
size_type n_elements() const
static bool is_running_single_threaded()
real_type lp_norm(const real_type p) const
real_type l1_norm() const
VectorBase & operator+=(const VectorBase &V)
PetscScalar mean_value() const
size_type local_size() const
VectorOperation::values last_action
PetscScalar operator*(const VectorBase &vec) const
bool operator==(const VectorBase &v) const
VectorBase & operator*=(const PetscScalar factor)
VectorBase & operator-=(const VectorBase &V)
bool operator!=(const VectorBase &v) const
std::size_t memory_consumption() const
VectorBase & operator/=(const PetscScalar factor)
bool is_non_negative() const
virtual const MPI_Comm & get_mpi_communicator() const
void scale(const VectorBase &scaling_factors)
std::pair< size_type, size_type > local_range() const
void sadd(const PetscScalar s, const VectorBase &V)
real_type norm_sqr() const
void compress(const VectorOperation::values operation)
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
virtual ~VectorBase() override
real_type linfty_norm() const
void add(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
bool has_ghost_elements() const
size_type locally_owned_size() const
void equ(const PetscScalar a, const VectorBase &V)
void do_set_add_operation(const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
real_type l2_norm() const
VectorBase & operator=(const VectorBase &)=delete
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcGhostsPresent()
#define Assert(cond, exc)
#define AssertIsFinite(number)
#define AssertThrowMPI(error_code)
#define AssertNothrow(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
types::global_dof_index size_type
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
bool is_non_negative(const T &t)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
T sum(const T &t, const MPI_Comm &mpi_communicator)
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)