|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_vector_h
17 #define dealii_vector_h
32 #include <boost/serialization/split_member.hpp>
35 #include <initializer_list>
45 # ifdef DEAL_II_WITH_PETSC
52 # ifdef DEAL_II_WITH_TRILINOS
62 template <
typename number>
110 template <
typename Number>
118 "The Vector class does not support auto-differentiable numbers.");
162 Vector(
const Vector<Number> &v);
168 Vector(Vector<Number> &&v) noexcept =
default;
178 template <
typename OtherNumber>
179 explicit Vector(
const Vector<OtherNumber> &v);
196 template <
typename OtherNumber>
197 explicit Vector(
const std::initializer_list<OtherNumber> &v);
199 #ifdef DEAL_II_WITH_PETSC
214 #ifdef DEAL_II_WITH_TRILINOS
247 template <
typename InputIterator>
248 Vector(
const InputIterator
first,
const InputIterator last);
254 virtual ~Vector()
override =
default;
336 template <
typename Number2>
338 reinit(
const Vector<Number2> &
V,
const bool omit_zeroing_entries =
false);
356 swap(Vector<Number> &v);
385 operator=(Vector<Number> &&v) noexcept =
default;
392 template <
typename Number2>
402 #ifdef DEAL_II_WITH_PETSC
419 #ifdef DEAL_II_WITH_TRILINOS
444 template <
typename Number2>
453 template <
typename Number2>
478 template <
typename Number2>
479 Number
operator*(
const Vector<Number2> &
V)
const;
564 add_and_dot(
const Number a,
const Vector<Number> &
V,
const Vector<Number> &W);
654 template <
typename OtherNumber>
657 std::vector<OtherNumber> &
values)
const;
686 template <
typename ForwardIterator,
typename OutputIterator>
689 const ForwardIterator indices_end,
690 OutputIterator values_begin)
const;
719 template <
typename OtherNumber>
721 add(
const std::vector<size_type> & indices,
722 const std::vector<OtherNumber> &
values);
728 template <
typename OtherNumber>
730 add(
const std::vector<size_type> &indices,
const Vector<OtherNumber> &
values);
737 template <
typename OtherNumber>
741 const OtherNumber *
values);
759 const Vector<Number> &
V,
761 const Vector<Number> &W);
769 add(
const Number a,
const Vector<Number> &
V);
777 sadd(
const Number s,
const Vector<Number> &
V);
785 sadd(
const Number s,
const Number a,
const Vector<Number> &
V);
811 scale(
const Vector<Number> &scaling_factors);
818 template <
typename Number2>
820 scale(
const Vector<Number2> &scaling_factors);
828 equ(
const Number a,
const Vector<Number> &u);
833 template <
typename Number2>
835 equ(
const Number a,
const Vector<Number2> &u);
857 print(std::ostream & out,
858 const unsigned int precision = 3,
859 const bool scientific =
true,
860 const bool across =
true)
const;
888 template <
class Archive>
890 save(Archive &ar,
const unsigned int version)
const;
896 template <
class Archive>
898 load(Archive &ar,
const unsigned int version);
905 template <
class Archive>
907 serialize(Archive &archive,
const unsigned int version);
911 BOOST_SERIALIZATION_SPLIT_MEMBER()
1010 const bool omit_zeroing_entries,
1011 const bool reset_partitioner);
1017 mutable std::shared_ptr<parallel::internal::TBBPartitioner>
1021 template <
typename Number2>
1040 template <
typename Number>
1051 template <
typename Number>
1052 template <
typename OtherNumber>
1059 template <
typename Number>
1060 template <
typename InputIterator>
1071 template <
typename Number>
1082 template <
typename Number>
1090 template <
typename Number>
1099 template <
typename Number>
1108 template <
typename Number>
1117 template <
typename Number>
1126 template <
typename Number>
1135 template <
typename Number>
1144 template <
typename Number>
1153 template <
typename Number>
1163 template <
typename Number>
1173 template <
typename Number>
1181 template <
typename Number>
1189 template <
typename Number>
1190 template <
typename OtherNumber>
1193 std::vector<OtherNumber> &
values)
const
1195 for (
size_type i = 0; i < indices.size(); ++i)
1196 values[i] =
operator()(indices[i]);
1201 template <
typename Number>
1202 template <
typename ForwardIterator,
typename OutputIterator>
1205 const ForwardIterator indices_end,
1206 OutputIterator values_begin)
const
1208 while (indices_begin != indices_end)
1218 template <
typename Number>
1219 inline Vector<Number> &
1231 template <
typename Number>
1232 template <
typename OtherNumber>
1235 const std::vector<OtherNumber> &
values)
1244 template <
typename Number>
1245 template <
typename OtherNumber>
1248 const Vector<OtherNumber> &
values)
1257 template <
typename Number>
1258 template <
typename OtherNumber>
1262 const OtherNumber *
values)
1264 for (
size_type i = 0; i < n_indices; ++i)
1270 "The given value is not finite but either infinite or Not A Number (NaN)"));
1272 this->values[indices[i]] +=
values[i];
1278 template <
typename Number>
1279 template <
typename Number2>
1283 return !(*
this == v);
1288 template <
typename Number>
1293 template <
typename Number>
1300 template <
typename Number>
1310 template <
typename Number>
1320 template <
typename Number>
1321 template <
class Archive>
1332 template <
typename Number>
1333 template <
class Archive>
1359 template <
typename Number>
1361 swap(Vector<Number> &u, Vector<Number> &v)
1374 template <
typename number>
1375 inline std::ostream &
1383 out << v(v.
size() - 1);
1399 template <
typename Number>
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
void add(const std::vector< size_type > &indices, const std::vector< OtherNumber > &values)
void swap(Vector< Number > &u, Vector< Number > &v)
void apply_givens_rotation(const std::array< Number, 3 > &csr, const size_type i, const size_type k)
IndexSet locally_owned_elements() const
const value_type * const_iterator
static ::ExceptionBase & ExcIO()
Vector< Number > & operator+=(const Vector< Number > &V)
bool in_local_range(const size_type global_index) const
const value_type * const_iterator
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
Subscriptor & operator=(const Subscriptor &)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
types::global_dof_index size_type
#define AssertIndexRange(index, range)
void block_write(std::ostream &out) const
void load(Archive &ar, const unsigned int version)
virtual void swap(Vector< Number > &v)
Number operator()(const size_type i) const
real_type norm_sqr() const
__global__ void sadd(const Number s, Number *val, const Number a, const Number *V_val, const size_type N)
void serialize(Archive &ar, const unsigned int version)
typename VectorType::value_type value_type
Vector< Number > & operator/=(const Number factor)
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< OtherNumber > &values) const
void maybe_reset_thread_partitioner()
Number operator[](const size_type i) const
#define AssertIsFinite(number)
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
VectorType::value_type * begin(VectorType &V)
static ::ExceptionBase & ExcEmptyObject()
bool is_finite(const double x)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
virtual ~Vector() override=default
bool has_ghost_elements() const
const value_type & const_reference
#define DEAL_II_NAMESPACE_OPEN
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Number l1_norm(const Tensor< 2, dim, Number > &t)
void compress(::VectorOperation::values operation=::VectorOperation::unknown) const
VectorType::value_type * end(VectorType &V)
void swap(BlockIndices &u, BlockIndices &v)
__global__ void equ(Number *val, const Number a, const Number *V_val, const size_type N)
real_type lp_norm(const real_type p) const
void save(Archive &ar, const unsigned int version) const
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
void swap(AlignedVector< T > &vec)
AlignedVector< Number > values
std::vector< value_type > l2_norm(const typename ::Triangulation< dim, spacedim >::cell_iterator &parent, const value_type parent_value)
bool is_non_negative(const T &t)
types::global_dof_index size_type
#define Assert(cond, exc)
Vector< Number > & operator-=(const Vector< Number > &V)
OutputOperator< VectorType > & operator<<(OutputOperator< VectorType > &out, unsigned int step)
typename numbers::NumberTraits< typename VectorType::value_type >::real_type real_type
bool operator!=(const Vector< Number2 > &v) const
static ::ExceptionBase & ExcZero()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void do_reinit(const size_type new_size, const bool omit_zeroing_entries, const bool reset_partitioner)
std::shared_ptr< parallel::internal::TBBPartitioner > thread_loop_partitioner
const value_type * const_pointer
types::global_dof_index size_type
__global__ void add_and_dot(Number *res, Number *v1, const Number *v2, const Number *v3, const Number a, const size_type N)
Vector< Number > & operator*=(const Number factor)
#define DEAL_II_NAMESPACE_CLOSE
std::string compress(const std::string &input)
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
void grow_or_shrink(const size_type N)
void update_ghost_values() const
#define AssertThrow(cond, exc)
void block_read(std::istream &in)
Number mean_value() const
Number linfty_norm(const Tensor< 2, dim, Number > &t)
void copy(const T *begin, const T *end, U *dest)
std::enable_if< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typename ProductType< std::complex< T >, std::complex< U > >::type >::type operator*(const std::complex< T > &left, const std::complex< U > &right)
const value_type * const_pointer