Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/trilinos_vector.h>
Public Types | |
typedef ::types::global_dof_index | size_type |
Public Types inherited from TrilinosWrappers::VectorBase | |
typedef TrilinosScalar | value_type |
Public Member Functions | |
Vector () 1 | |
Vector (const size_type n) 1 | |
Vector (const Epetra_Map &partitioning) 1 | |
Vector (const IndexSet &partitioning, const MPI_Comm &communicator=MPI_COMM_WORLD) 1 | |
Vector (const VectorBase &V) 1 | |
template<typename Number > | |
Vector (const ::Vector< Number > &v) 1 | |
void | reinit (const size_type n, const bool omit_zeroing_entries=false) |
void | reinit (const Epetra_Map &input_map, const bool omit_zeroing_entries=false) |
void | reinit (const IndexSet &input_map, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool omit_zeroing_entries=false) |
void | reinit (const VectorBase &V, const bool omit_zeroing_entries=false, const bool allow_different_maps=false) |
Vector & | operator= (const TrilinosScalar s) |
Vector & | operator= (const MPI::Vector &v) |
template<typename Number > | |
Vector & | operator= (const ::Vector< Number > &v) |
Vector & | operator= (const Vector &v) |
void | update_ghost_values () const |
Public Member Functions inherited from TrilinosWrappers::VectorBase | |
VectorBase () | |
VectorBase (const VectorBase &v) | |
virtual | ~VectorBase () |
void | clear () |
void | reinit (const VectorBase &v, const bool omit_zeroing_entries=false) |
void | compress (::VectorOperation::values operation) |
bool | is_compressed () const 1 |
VectorBase & | operator= (const TrilinosScalar s) |
VectorBase & | operator= (const VectorBase &v) |
template<typename Number > | |
VectorBase & | operator= (const ::Vector< Number > &v) |
bool | operator== (const VectorBase &v) const |
bool | operator!= (const VectorBase &v) const |
size_type | size () const |
size_type | local_size () const |
std::pair< size_type, size_type > | local_range () const |
bool | in_local_range (const size_type index) const |
IndexSet | locally_owned_elements () const |
bool | has_ghost_elements () const |
TrilinosScalar | operator* (const VectorBase &vec) const |
real_type | norm_sqr () const |
TrilinosScalar | mean_value () const |
TrilinosScalar | minimal_value () const 1 |
TrilinosScalar | min () const |
TrilinosScalar | max () const |
real_type | l1_norm () const |
real_type | l2_norm () const |
real_type | lp_norm (const TrilinosScalar p) const |
real_type | linfty_norm () const |
TrilinosScalar | add_and_dot (const TrilinosScalar a, const VectorBase &V, const VectorBase &W) |
bool | all_zero () const |
bool | is_non_negative () const |
reference | operator() (const size_type index) |
TrilinosScalar | operator() (const size_type index) const |
reference | operator[] (const size_type index) |
TrilinosScalar | operator[] (const size_type index) const |
void | extract_subvector_to (const std::vector< size_type > &indices, std::vector< TrilinosScalar > &values) const |
template<typename ForwardIterator , typename OutputIterator > | |
void | extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
TrilinosScalar | el (const size_type index) const 1 |
iterator | begin () |
const_iterator | begin () const |
iterator | end () |
const_iterator | end () const |
void | set (const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values) |
void | set (const std::vector< size_type > &indices, const ::Vector< TrilinosScalar > &values) |
void | set (const size_type n_elements, const size_type *indices, const TrilinosScalar *values) |
void | add (const std::vector< size_type > &indices, const std::vector< TrilinosScalar > &values) |
void | add (const std::vector< size_type > &indices, const ::Vector< TrilinosScalar > &values) |
void | add (const size_type n_elements, const size_type *indices, const TrilinosScalar *values) |
VectorBase & | operator*= (const TrilinosScalar factor) |
VectorBase & | operator/= (const TrilinosScalar factor) |
VectorBase & | operator+= (const VectorBase &V) |
VectorBase & | operator-= (const VectorBase &V) |
void | add (const TrilinosScalar s) |
void | add (const VectorBase &V, const bool allow_different_maps=false) |
void | add (const TrilinosScalar a, const VectorBase &V) |
void | add (const TrilinosScalar a, const VectorBase &V, const TrilinosScalar b, const VectorBase &W) |
void | sadd (const TrilinosScalar s, const VectorBase &V) |
void | sadd (const TrilinosScalar s, const TrilinosScalar a, const VectorBase &V) |
void | sadd (const TrilinosScalar s, const TrilinosScalar a, const VectorBase &V, const TrilinosScalar b, const VectorBase &W) 1 |
void | sadd (const TrilinosScalar s, const TrilinosScalar a, const VectorBase &V, const TrilinosScalar b, const VectorBase &W, const TrilinosScalar c, const VectorBase &X) 1 |
void | scale (const VectorBase &scaling_factors) |
void | equ (const TrilinosScalar a, const VectorBase &V) |
void | equ (const TrilinosScalar a, const VectorBase &V, const TrilinosScalar b, const VectorBase &W) 1 |
void | ratio (const VectorBase &a, const VectorBase &b) 1 |
const Epetra_MultiVector & | trilinos_vector () const |
Epetra_FEVector & | trilinos_vector () |
const Epetra_Map & | vector_partitioner () const |
void | print (const char *format=0) const 1 |
void | print (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const |
void | swap (VectorBase &v) |
std::size_t | memory_consumption () const |
const MPI_Comm & | get_mpi_communicator () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Attributes | |
static const bool | supports_distributed_data = false |
Related Functions | |
(Note that these are not member functions.) | |
void | swap (Vector &u, Vector &v) |
Related Functions inherited from TrilinosWrappers::VectorBase | |
void | swap (VectorBase &u, VectorBase &v) |
Additional Inherited Members | |
Static Public Member Functions inherited from TrilinosWrappers::VectorBase | |
static ::ExceptionBase & | ExcDifferentParallelPartitioning () |
static ::ExceptionBase & | ExcTrilinosError (int arg1) |
static ::ExceptionBase & | ExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
This class is a specialization of a Trilinos vector to a localized version. The purpose of this class is to provide a copy interface from the possibly parallel Vector class to a local vector on each processor, in order to be able to access all elements in the vector or to apply certain deal.II functions.
This class is deprecated, use TrilinosWrappers::MPI::Vector instead.
Definition at line 756 of file trilinos_vector.h.
Declare type for container size.
Definition at line 762 of file trilinos_vector.h.
Default constructor that generates an empty (zero size) vector. The function reinit()
will have to give the vector the correct size.
Definition at line 629 of file trilinos_vector.cc.
This constructor takes as input the number of elements in the vector.
Definition at line 638 of file trilinos_vector.cc.
This constructor takes as input the number of elements in the vector. If the map is not localized, i.e., if there are some elements that are not present on all processes, only the global size of the map will be taken and a localized map will be generated internally. In other words, which element of the partitioning
argument are set is in fact ignored, the only thing that matters is the size of the index space described by this argument.
Definition at line 645 of file trilinos_vector.cc.
|
explicit |
This constructor takes as input the number of elements in the vector. If the index set is not localized, i.e., if there are some elements that are not present on all processes, only the global size of the index set will be taken and a localized version will be generated internally. In other words, which element of the partitioning
argument are set is in fact ignored, the only thing that matters is the size of the index space described by this argument.
Definition at line 652 of file trilinos_vector.cc.
|
explicit |
This constructor takes a (possibly parallel) Trilinos Vector and generates a localized version of the whole content on each processor.
Definition at line 660 of file trilinos_vector.cc.
|
explicit |
Copy-constructor from deal.II vectors. Sets the dimension to that of the given vector, and copies all elements.
Reinit function that resizes the vector to the size specified by n
. The flag omit_zeroing_entries
specifies whether the vector entries should be initialized to zero (false). If it is set to true
, the vector entries are in an unspecified state and the user has to set all elements. In the current implementation, this method still sets the entries to zero, but this might change between releases without notification.
Definition at line 680 of file trilinos_vector.cc.
void Vector< Number >::reinit | ( | const Epetra_Map & | input_map, |
const bool | omit_zeroing_entries = false |
||
) |
Initialization with an Epetra_Map. Similar to the call in the other class MPI::Vector, with the difference that now a copy on all processes is generated. This initialization function is appropriate when the data in the localized vector should be imported from a distributed vector that has been initialized with the same communicator. The variable omit_zeroing_entries
determines whether the vector should be filled with zero (false). If it is set to true
, the vector entries are in an unspecified state and the user has to set all elements. In the current implementation, this method still sets the entries to zero, but this might change between releases without notification.
Which element of the input_map
argument are set is in fact ignored, the only thing that matters is the size of the index space described by this argument.
Definition at line 694 of file trilinos_vector.cc.
void Vector< Number >::reinit | ( | const IndexSet & | input_map, |
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const bool | omit_zeroing_entries = false |
||
) |
Initialization with an IndexSet. Similar to the call in the other class MPI::Vector, with the difference that now a copy on all processes is generated. This initialization function is appropriate in case the data in the localized vector should be imported from a distributed vector that has been initialized with the same communicator. The variable omit_zeroing_entries
determines whether the vector should be filled with zero (false). If it is set to true
, the vector entries are in an unspecified state and the user has to set all elements. In the current implementation, this method still sets the entries to zero, but this might change between releases without notification.
Which element of the input_map
argument are set is in fact ignored, the only thing that matters is the size of the index space described by this argument.
Definition at line 709 of file trilinos_vector.cc.
void Vector< Number >::reinit | ( | const VectorBase & | V, |
const bool | omit_zeroing_entries = false , |
||
const bool | allow_different_maps = false |
||
) |
Reinit function. Takes the information of a Vector and copies everything to the calling vector, now also allowing different maps.
If the optional flag omit_zeroing_entries
is not true
, the elements in the vector are initialized with zero. If it is set to true
, the vector entries are in an unspecified state and the user has to set all elements. In the current implementation, this method does not touch the vector entries in case the vector layout is unchanged from before, otherwise entries are set to zero. Note that this behavior might change between releases without notification.
Definition at line 730 of file trilinos_vector.cc.
Vector& TrilinosWrappers::Vector::operator= | ( | const TrilinosScalar | s | ) |
Set all components of the vector to the given number s
. Simply pass this down to the base class, but we still need to declare this function to make the example given in the discussion about making the constructor explicit work.
Vector & Vector< Number >::operator= | ( | const MPI::Vector & | v | ) |
Set the left hand argument to the (parallel) Trilinos Vector. Equivalent to the reinit
function.
Definition at line 807 of file trilinos_vector.cc.
Vector& TrilinosWrappers::Vector::operator= | ( | const ::Vector< Number > & | v | ) |
Set the left hand argument to the deal.II vector.
Copy operator. Copies both the dimension and the content in the right hand argument.
Definition at line 824 of file trilinos_vector.cc.
void TrilinosWrappers::Vector::update_ghost_values | ( | ) | const |
This function does nothing but is there for compatibility with the PETScWrappers::Vector
class.
For the PETSc vector wrapper class, this function updates the ghost values of the PETSc vector. This is necessary after any modification before reading ghost values.
However, for the implementation of this class, it is immaterial and thus an empty function.
Global function swap
which overloads the default implementation of the C++ standard library which uses a temporary object. The function simply exchanges the data of the two vectors.
Definition at line 950 of file trilinos_vector.h.
|
static |
A variable that indicates whether this vector supports distributed data storage. If true, then this vector also needs an appropriate compress() function that allows communicating recent set or add operations to individual elements to be communicated to other processors.
For the current class, the variable equals false, since it does not support parallel data storage. If you do need parallel data storage, use TrilinosWrappers::MPI::Vector.
is_serial_vector< VectorType >::value
Definition at line 777 of file trilinos_vector.h.