16 #ifndef dealii__read_write_vector_h 17 #define dealii__read_write_vector_h 19 #include <deal.II/base/config.h> 20 #include <deal.II/base/index_set.h> 21 #include <deal.II/base/mpi.h> 22 #include <deal.II/base/template_constraints.h> 23 #include <deal.II/base/types.h> 24 #include <deal.II/base/utilities.h> 25 #include <deal.II/base/memory_consumption.h> 26 #include <deal.II/base/thread_management.h> 27 #include <deal.II/lac/vector_view.h> 32 #ifdef DEAL_II_WITH_TRILINOS 33 #include <deal.II/lac/trilinos_epetra_communication_pattern.h> 34 #include "Epetra_MultiVector.h" 37 DEAL_II_NAMESPACE_OPEN
41 class CommunicationPatternBase;
44 template <
typename>
class Vector;
48 #ifdef DEAL_II_WITH_PETSC 58 #ifdef DEAL_II_WITH_TRILINOS 69 namespace EpetraWrappers
76 #ifdef DEAL_II_WITH_CUDA 79 namespace CUDAWrappers
81 template <
typename>
class Vector;
127 template <
typename Number>
186 const bool omit_zeroing_entries =
false);
196 template <
typename Number2>
198 const bool omit_zeroing_entries =
false);
210 const bool omit_zeroing_entries =
false);
212 #ifdef DEAL_II_WITH_CXX11 226 template <
typename Functor>
227 void apply(
const Functor &func);
253 template <
typename Number2>
273 std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
274 std_cxx11::shared_ptr<const CommunicationPatternBase> ());
276 #ifdef DEAL_II_WITH_PETSC 287 std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
288 std_cxx11::shared_ptr<const CommunicationPatternBase> ());
291 #ifdef DEAL_II_WITH_TRILINOS 302 std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
303 std_cxx11::shared_ptr<const CommunicationPatternBase> ());
315 std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
316 std_cxx11::shared_ptr<const CommunicationPatternBase> ());
319 #ifdef DEAL_II_WITH_CUDA 326 void import(
const CUDAWrappers::Vector<Number> &cuda_vec,
328 std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
329 std_cxx11::shared_ptr<const CommunicationPatternBase> ());
340 size_type
size()
const;
365 const_iterator
begin ()
const;
377 const_iterator
end ()
const;
391 Number
operator () (
const size_type global_index)
const;
398 Number &
operator () (
const size_type global_index);
407 Number
operator [] (
const size_type global_index)
const;
416 Number &
operator [] (
const size_type global_index);
424 template <
typename Number2>
426 std::vector<Number2> &values)
const;
432 template <
typename ForwardIterator,
typename OutputIterator>
434 const ForwardIterator indices_end,
435 OutputIterator values_begin)
const;
472 template <
typename Number2>
473 void add (
const std::vector<size_type> &indices,
474 const std::vector<Number2> &values);
480 template <
typename Number2>
481 void add (
const std::vector<size_type> &indices,
489 template <
typename Number2>
491 const size_type *indices,
492 const Number2 *values);
497 void print (std::ostream &out,
498 const unsigned int precision = 3,
499 const bool scientific =
true)
const;
508 #ifdef DEAL_II_WITH_TRILINOS 514 void import(
const Epetra_MultiVector &multivector,
515 const IndexSet &locally_owned_elements,
517 const MPI_Comm &mpi_comm,
518 std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern);
535 void resize_val (
const size_type new_allocated_size);
537 #if defined(DEAL_II_WITH_TRILINOS) && defined(DEAL_II_WITH_MPI) 544 const MPI_Comm &mpi_comm);
579 #ifdef DEAL_II_WITH_CXX11 586 template <
typename Functor>
600 const size_type
end);
623 template <
typename Number>
634 template <
typename Number>
646 template <
typename Number>
657 template <
typename Number>
663 reinit (locally_stored_indices);
668 template <
typename Number>
677 template <
typename Number>
679 typename ReadWriteVector<Number>::size_type
687 template <
typename Number>
689 typename ReadWriteVector<Number>::size_type
697 template <
typename Number>
707 template <
typename Number>
709 typename ReadWriteVector<Number>::iterator
717 template <
typename Number>
719 typename ReadWriteVector<Number>::const_iterator
727 template <
typename Number>
729 typename ReadWriteVector<Number>::iterator
737 template <
typename Number>
739 typename ReadWriteVector<Number>::const_iterator
747 template <
typename Number>
757 template <
typename Number>
767 template <
typename Number>
777 template <
typename Number>
787 template <
typename Number>
788 template <
typename Number2>
791 std::vector<Number2> &values)
const 793 for (size_type i = 0; i < indices.size(); ++i)
794 values[i] =
operator()(indices[i]);
799 template <
typename Number>
800 template <
typename ForwardIterator,
typename OutputIterator>
803 const ForwardIterator indices_end,
804 OutputIterator values_begin)
const 806 while (indices_begin != indices_end)
816 template <
typename Number>
823 return val[local_index];
828 template <
typename Number>
835 return val[local_index];
840 template <
typename Number>
841 template <
typename Number2>
845 const std::vector<Number2> &values)
848 add (indices.size(), &indices[0], &values[0]);
853 template <
typename Number>
854 template <
typename Number2>
858 const ReadWriteVector<Number2> &values)
860 const size_type
size = indices.size();
861 for (size_type i=0; i<
size; ++i)
864 ExcMessage(
"The given value is not finite but either infinite or Not A Number (NaN)"));
865 this->
operator()(indices[i]) += values[indices[i]];
871 template <
typename Number>
872 template <
typename Number2>
876 const size_type *indices,
877 const Number2 *values)
879 for (size_type i=0; i<n_indices; ++i)
882 ExcMessage(
"The given value is not finite but either infinite or Not A Number (NaN)"));
889 #ifdef DEAL_II_WITH_CXX11 890 template <
typename Number>
891 template <
typename Functor>
894 ReadWriteVector<Number> &parent,
895 const Functor &functor)
903 template <
typename Number>
904 template <
typename Functor>
910 functor(parent.val[i]);
914 #endif // ifndef DOXYGEN 928 template <
typename Number>
937 DEAL_II_NAMESPACE_CLOSE
IndexSet source_stored_elements
#define AssertDimension(dim1, dim2)
Number operator[](const size_type global_index) const
void swap(ReadWriteVector< Number > &v)
Number local_element(const size_type local_index) const
FunctorTemplate(ReadWriteVector< Number > &parent, const Functor &functor)
void resize_val(const size_type new_allocated_size)
#define AssertIndexRange(index, range)
std::size_t memory_consumption() const
void swap(LinearAlgebra::ReadWriteVector< Number > &u, LinearAlgebra::ReadWriteVector< Number > &v)
bool is_finite(const double x)
static ::ExceptionBase & ExcMessage(std::string arg1)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
unsigned int global_dof_index
#define Assert(cond, exc)
size_type index_within_set(const size_type global_index) const
Number operator()(const size_type global_index) const
ReadWriteVector< Number > & operator=(const ReadWriteVector< Number > &in_vector)
EpetraWrappers::CommunicationPattern create_epetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm &mpi_comm)
std_cxx11::shared_ptr< ::parallel::internal::TBBPartitioner > thread_loop_partitioner
void add(const std::vector< size_type > &indices, const std::vector< Number2 > &values)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true) const
virtual void operator()(const size_type begin, const size_type end)
void extract_subvector_to(const std::vector< size_type > &indices, std::vector< Number2 > &values) const
void apply(const Functor &func)
unsigned int global_to_local(const types::global_dof_index global_index) const
size_type n_elements() const
size_type n_elements() const
const IndexSet & get_stored_elements() const
std_cxx11::shared_ptr< CommunicationPatternBase > comm_pattern