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/petsc_vector.h> 22 # include <deal.II/lac/petsc_parallel_vector.h> 24 # include <deal.II/base/multithread_info.h> 26 DEAL_II_NAMESPACE_OPEN
32 VectorReference::operator PetscScalar ()
const 34 Assert (index < vector.size(),
35 ExcIndexRange (index, 0, vector.size()));
40 if (dynamic_cast<const PETScWrappers::Vector *>(&vector) != 0)
44 int ierr = VecGetValues(vector.vector, 1, &idx, &value);
50 else if (dynamic_cast<const PETScWrappers::MPI::Vector *>(&vector) != 0)
74 ierr = VecGetOwnershipRange (vector.vector, &begin, &end);
77 Vec locally_stored_elements = PETSC_NULL;
78 ierr = VecGhostGetLocalForm(vector.vector, &locally_stored_elements);
82 ierr = VecGetSize(locally_stored_elements, &lsize);
86 ierr = VecGetArray(locally_stored_elements, &ptr);
91 if ( index>=static_cast<size_type>(begin)
92 && index<static_cast<size_type>(end) )
95 value = *(ptr+index-begin);
101 = vector.ghost_indices.index_within_set(index);
103 Assert(ghostidx+end-begin<(size_type)lsize, ExcInternalError());
104 value = *(ptr+ghostidx+end-begin);
107 ierr = VecRestoreArray(locally_stored_elements, &ptr);
110 ierr = VecGhostRestoreLocalForm(vector.vector, &locally_stored_elements);
122 ierr = VecGetOwnershipRange (vector.vector, &begin, &end);
127 AssertThrow ((index >= static_cast<size_type>(begin)) &&
128 (index < static_cast<size_type>(end)),
129 ExcAccessToNonlocalElement (index, begin, end-1));
133 PetscInt idx = index;
135 ierr = VecGetValues(vector.vector, 1, &idx, &value);
143 Assert (
false, ExcInternalError());
153 attained_ownership(true)
156 ExcMessage(
"PETSc does not support multi-threaded access, set " 157 "the thread limit to 1 in MPI_InitFinalize()."));
171 ExcMessage(
"PETSc does not support multi-threaded access, set " 172 "the thread limit to 1 in MPI_InitFinalize()."));
192 ExcMessage(
"PETSc does not support multi-threaded access, set " 193 "the thread limit to 1 in MPI_InitFinalize()."));
202 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 203 const int ierr = VecDestroy (
vector);
205 const int ierr = VecDestroy (&
vector);
218 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 219 const int ierr = VecDestroy (
vector);
221 const int ierr = VecDestroy (&
vector);
241 int ierr = VecSet (
vector, s);
246 Vec ghost = PETSC_NULL;
247 ierr = VecGhostGetLocalForm(
vector, &ghost);
250 ierr = VecSet (ghost, s);
253 ierr = VecGhostRestoreLocalForm(
vector, &ghost);
266 ExcDimensionMismatch(
size(), v.
size()));
268 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 278 return (flag == PETSC_TRUE);
287 ExcDimensionMismatch(
size(), v.
size()));
289 #if DEAL_II_PETSC_VERSION_LT(3,2,0) 299 return (flag == PETSC_FALSE);
304 VectorBase::size_type
308 const int ierr = VecGetSize (
vector, &sz);
316 VectorBase::size_type
320 const int ierr = VecGetLocalSize (
vector, &sz);
328 std::pair<VectorBase::size_type, VectorBase::size_type>
332 const int ierr = VecGetOwnershipRange (static_cast<const Vec &>(
vector),
336 return std::make_pair (begin, end);
343 const std::vector<PetscScalar> &values)
345 Assert (indices.size() == values.size(),
346 ExcMessage (
"Function called with arguments of different sizes"));
354 const std::vector<PetscScalar> &values)
356 Assert (indices.size() == values.size(),
357 ExcMessage (
"Function called with arguments of different sizes"));
365 const ::Vector<PetscScalar> &values)
367 Assert (indices.size() == values.size(),
368 ExcMessage (
"Function called with arguments of different sizes"));
376 const size_type *indices,
377 const PetscScalar *values)
388 ExcDimensionMismatch(
size(), vec.
size()));
419 #ifdef DEAL_II_WITH_MPI 423 int all_int_last_action;
425 MPI_Allreduce(&my_int_last_action, &all_int_last_action, 1, MPI_INT,
428 AssertThrow(all_int_last_action != (::VectorOperation::add | ::VectorOperation::insert),
429 ExcMessage(
"Error: not all processors agree on the last VectorOperation before this compress() call."));
435 ExcMessage(
"Missing compress() or calling with wrong VectorOperation argument."));
452 ierr = VecAssemblyBegin(
vector);
454 ierr = VecAssemblyEnd(
vector);
465 VectorBase::real_type
481 if (dynamic_cast<const PETScWrappers::MPI::Vector *>(
this) != 0)
484 ierr = VecSum(
vector, &sum);
486 return sum/
static_cast<PetscReal
>(
size());
491 PetscScalar *start_ptr;
492 ierr = VecGetArray (
vector, &start_ptr);
495 PetscScalar mean = 0;
497 PetscScalar sum0 = 0,
505 const PetscScalar *ptr = start_ptr;
506 const PetscScalar *eptr = ptr + (
size()/4)*4;
515 while (ptr != start_ptr+
size())
518 mean = (sum0+sum1+sum2+sum3)/static_cast<PetscReal>(
size());
523 ierr = VecRestoreArray (
vector, &start_ptr);
530 VectorBase::real_type
535 const int ierr = VecNorm (
vector, NORM_1, &d);
543 VectorBase::real_type
548 const int ierr = VecNorm (
vector, NORM_2, &d);
556 VectorBase::real_type
561 PetscScalar *start_ptr;
562 int ierr = VecGetArray (
vector, &start_ptr);
575 const PetscScalar *ptr = start_ptr;
576 const PetscScalar *eptr = ptr + (
size()/4)*4;
585 while (ptr != start_ptr+
size())
588 norm = std::pow(sum0+sum1+sum2+sum3,
594 ierr = VecRestoreArray (
vector, &start_ptr);
602 VectorBase::real_type
607 const int ierr = VecNorm (
vector, NORM_INFINITY, &d);
615 VectorBase::real_type
620 const int ierr = VecNormalize (
vector, &d);
627 VectorBase::real_type
633 const int ierr = VecMin (
vector, &p, &d);
640 VectorBase::real_type
646 const int ierr = VecMax (
vector, &p, &d);
658 const int ierr = VecAbs (
vector);
671 const int ierr = VecConjugate (
vector);
707 const int ierr = VecPointwiseMult (
vector,u,v);
718 PetscScalar *start_ptr;
719 int ierr = VecGetArray (
vector, &start_ptr);
722 const PetscScalar *ptr = start_ptr,
737 ierr = VecRestoreArray (
vector, &start_ptr);
746 template <
typename T>
754 template <
typename T>
759 "whether it is non-negative."))
771 PetscScalar *start_ptr;
772 int ierr = VecGetArray (
vector, &start_ptr);
775 const PetscScalar *ptr = start_ptr,
780 if (! internal::is_non_negative(*ptr))
790 ierr = VecRestoreArray (
vector, &start_ptr);
804 const int ierr = VecScale (
vector, a);
818 const PetscScalar factor = 1./a;
821 const int ierr = VecScale (
vector, factor);
833 const int ierr = VecAXPY (
vector, 1, v);
845 const int ierr = VecAXPY (
vector, -1, v);
859 const int ierr = VecShift (
vector, s);
880 const int ierr = VecAXPY (
vector, a, v);
896 const PetscScalar weights[2] = {a,b};
899 const int ierr = VecMAXPY (
vector, 2, weights, addends);
912 const int ierr = VecAYPX (
vector, s, v);
952 const PetscScalar weights[2] = {a,b};
955 const int ierr = VecMAXPY (
vector, 2, weights, addends);
980 const PetscScalar weights[3] = {a,b,c};
983 const int ierr = VecMAXPY (
vector, 3, weights, addends);
1008 ExcDimensionMismatch (
size(), v.
size()));
1024 const PetscScalar b,
1032 ExcDimensionMismatch (
size(), v.
size()));
1050 const int ierr = VecPointwiseDivide (
vector, a, b);
1062 PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
1066 VecView (
vector, PETSC_VIEWER_STDOUT_WORLD);
1073 const unsigned int precision,
1074 const bool scientific,
1075 const bool across)
const 1082 int ierr = VecGetArray (
vector, &val);
1087 const std::ios::fmtflags old_flags = out.flags();
1088 const unsigned int old_precision = out.precision (precision);
1090 out.precision (precision);
1092 out.setf (std::ios::scientific, std::ios::floatfield);
1094 out.setf (std::ios::fixed, std::ios::floatfield);
1098 out << val[i] <<
' ';
1101 out << val[i] << std::endl;
1105 out.flags (old_flags);
1106 out.precision(old_precision);
1110 ierr = VecRestoreArray (
vector, &val);
1127 VectorBase::operator
const Vec &()
const 1157 const size_type *indices,
1158 const PetscScalar *values,
1159 const bool add_values)
1161 ::VectorOperation::values action = (add_values ?
1162 ::VectorOperation::add :
1163 ::VectorOperation::insert);
1167 internal::VectorReference::ExcWrongMode (action,
1176 if (n_elements != 0)
1178 #ifdef PETSC_USE_64BIT_INDICES 1179 std::vector<PetscInt> petsc_ind (n_elements);
1180 for (size_type i=0; i<n_elements; ++i)
1181 petsc_ind[i] = indices[i];
1182 const PetscInt *petsc_indices = &petsc_ind[0];
1184 const int *petsc_indices = (
const int *)indices;
1187 const InsertMode mode = (add_values ? ADD_VALUES : INSERT_VALUES);
1189 = VecSetValues (
vector, n_elements, petsc_indices, values,
1201 DEAL_II_NAMESPACE_CLOSE
1203 #endif // DEAL_II_WITH_PETSC types::global_dof_index size_type
VectorBase & operator-=(const VectorBase &V)
void sadd(const PetscScalar s, const VectorBase &V)
real_type l2_norm() const
::ExceptionBase & ExcMessage(std::string arg1)
VectorBase & mult() DEAL_II_DEPRECATED
bool is_non_negative() const
VectorBase & operator/=(const PetscScalar factor)
#define AssertThrow(cond, exc)
real_type norm_sqr() const
VectorBase & operator=(const PetscScalar s)
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
VectorBase & abs() DEAL_II_DEPRECATED
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 bool is_running_single_threaded()
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
real_type normalize() const DEAL_II_DEPRECATED
std::size_t memory_consumption() const
void equ(const PetscScalar a, const VectorBase &V)
VectorBase & operator*=(const PetscScalar factor)
real_type l1_norm() const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
VectorBase & operator+=(const VectorBase &V)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
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)
void ratio(const VectorBase &a, const VectorBase &b) DEAL_II_DEPRECATED
virtual const MPI_Comm & get_mpi_communicator() const
size_type n_elements() const
#define AssertIsFinite(number)
VectorBase & conjugate() DEAL_II_DEPRECATED
bool operator!=(const VectorBase &v) const
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)