Reference documentation for deal.II version 8.5.1
Public Types | Static Public Member Functions | Private Attributes | Friends | Related Functions | List of all members
TrilinosWrappers::VectorBase Class Reference

#include <deal.II/lac/trilinos_vector_base.h>

Inheritance diagram for TrilinosWrappers::VectorBase:
[legend]

Public Types

typedef TrilinosScalar value_type
 

Public Member Functions

1: Basic Object-handling
 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
 
VectorBaseoperator= (const TrilinosScalar s)
 
VectorBaseoperator= (const VectorBase &v)
 
template<typename Number >
VectorBaseoperator= (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
 
2: Data-Access
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
 
3: Modification of vectors
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)
 
VectorBaseoperator*= (const TrilinosScalar factor)
 
VectorBaseoperator/= (const TrilinosScalar factor)
 
VectorBaseoperator+= (const VectorBase &V)
 
VectorBaseoperator-= (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
 
4: Mixed stuff
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 ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (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 Member Functions

static ::ExceptionBaseExcDifferentParallelPartitioning ()
 
static ::ExceptionBaseExcTrilinosError (int arg1)
 
static ::ExceptionBaseExcAccessToNonLocalElement (size_type arg1, size_type arg2, size_type arg3, size_type arg4)
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *arg2)
 

Private Attributes

Epetra_CombineMode last_action
 
bool compressed
 
bool has_ghosts
 
std_cxx11::shared_ptr< Epetra_FEVector > vector
 
std_cxx11::shared_ptr< Epetra_MultiVector > nonlocal_vector
 
IndexSet owned_elements
 

Friends

class internal::VectorReference
 
class Vector
 

Related Functions

(Note that these are not member functions.)

void swap (VectorBase &u, VectorBase &v)
 

Detailed Description

Base class for the two types of Trilinos vectors, the distributed memory vector MPI::Vector and a localized vector Vector. The latter is designed for use in either serial implementations or as a localized copy on each processor. The implementation of this class is based on the Trilinos vector class Epetra_FEVector, the (parallel) partitioning of which is governed by an Epetra_Map. This means that the vector type is generic and can be done in this base class, while the definition of the partition map (and hence, the constructor and reinit function) will have to be done in the derived classes. The Epetra_FEVector is precisely the kind of vector we deal with all the time - we probably get it from some assembly process, where also entries not locally owned might need to written and hence need to be forwarded to the owner. The only requirement for this class to work is that Trilinos is installed with the same compiler as is used for compilation of deal.II.

The interface of this class is modeled after the existing Vector class in deal.II. It has almost the same member functions, and is often exchangeable. However, since Trilinos only supports a single scalar type (double), it is not templated, and only works with that type.

Note that Trilinos only guarantees that operations do what you expect if the function GlobalAssemble has been called after vector assembly in order to distribute the data. Therefore, you need to call Vector::compress() before you actually use the vectors.

Author
Martin Kronbichler, 2008

Definition at line 222 of file trilinos_vector_base.h.

Member Typedef Documentation

◆ value_type

Declare some of the standard types used in all containers. These types parallel those in the C standard libraries vector<...> class.

Definition at line 230 of file trilinos_vector_base.h.

Constructor & Destructor Documentation

◆ VectorBase() [1/2]

TrilinosWrappers::VectorBase::VectorBase ( )

Default constructor that generates an empty (zero size) vector. The function reinit() will have to give the vector the correct size and distribution among processes in case of an MPI run.

Definition at line 83 of file trilinos_vector_base.cc.

◆ VectorBase() [2/2]

TrilinosWrappers::VectorBase::VectorBase ( const VectorBase v)

Copy constructor. Sets the dimension to that of the given vector, and copies all the elements.

Definition at line 99 of file trilinos_vector_base.cc.

◆ ~VectorBase()

TrilinosWrappers::VectorBase::~VectorBase ( )
virtual

Destructor

Definition at line 111 of file trilinos_vector_base.cc.

Member Function Documentation

◆ clear()

void TrilinosWrappers::VectorBase::clear ( )

Release all memory and return to a state just like after having called the default constructor.

Definition at line 117 of file trilinos_vector_base.cc.

◆ reinit()

void TrilinosWrappers::VectorBase::reinit ( const VectorBase v,
const bool  omit_zeroing_entries = false 
)

Reinit functionality, sets the dimension and possibly the parallel partitioning (Epetra_Map) of the calling vector to the settings of the input vector.

◆ compress()

void TrilinosWrappers::VectorBase::compress ( ::VectorOperation::values  operation)

Compress the underlying representation of the Trilinos object, i.e. flush the buffers of the vector object if it has any. This function is necessary after writing into a vector element-by-element and before anything else can be done on it.

The (defaulted) argument can be used to specify the compress mode (Add or Insert) in case the vector has not been written to since the last time this function was called. The argument is ignored if the vector has been added or written to since the last time compress() was called.

See Compressing distributed objects for more information.

Definition at line 194 of file trilinos_vector_base.cc.

◆ is_compressed()

bool TrilinosWrappers::VectorBase::is_compressed ( ) const

Return the state of the vector, i.e., whether compress() has already been called after an operation requiring data exchange.

This function is deprecated.

◆ operator=() [1/3]

VectorBase& TrilinosWrappers::VectorBase::operator= ( const TrilinosScalar  s)

Set all components of the vector to the given number s. Simply pass this down to the Trilinos Epetra object, but we still need to declare this function to make the example given in the discussion about making the constructor explicit work.

Since the semantics of assigning a scalar to a vector are not immediately clear, this operator should really only be used if you want to set the entire vector to zero. This allows the intuitive notation v=0. Assigning other values is deprecated and may be disallowed in the future.

◆ operator=() [2/3]

VectorBase & TrilinosWrappers::VectorBase::operator= ( const VectorBase v)

Copy function. This function takes a VectorBase vector and copies all the elements. The target vector will have the same parallel distribution as the calling vector.

Definition at line 135 of file trilinos_vector_base.cc.

◆ operator=() [3/3]

template<typename Number >
VectorBase& TrilinosWrappers::VectorBase::operator= ( const ::Vector< Number > &  v)

Another copy function. This one takes a deal.II vector and copies it into a TrilinosWrapper vector. Note that since we do not provide any Epetra_map that tells about the partitioning of the vector among the MPI processes, the size of the TrilinosWrapper vector has to be the same as the size of the input vector. In order to change the map, use the reinit(const Epetra_Map &input_map) function.

◆ operator==()

bool TrilinosWrappers::VectorBase::operator== ( const VectorBase v) const

Test for equality. This function assumes that the present vector and the one to compare with have the same size already, since comparing vectors of different sizes makes not much sense anyway.

Definition at line 324 of file trilinos_vector_base.cc.

◆ operator!=()

bool TrilinosWrappers::VectorBase::operator!= ( const VectorBase v) const

Test for inequality. This function assumes that the present vector and the one to compare with have the same size already, since comparing vectors of different sizes makes not much sense anyway.

Definition at line 341 of file trilinos_vector_base.cc.

◆ size()

size_type TrilinosWrappers::VectorBase::size ( ) const

Return the global dimension of the vector.

◆ local_size()

size_type TrilinosWrappers::VectorBase::local_size ( ) const

Return the local dimension of the vector, i.e. the number of elements stored on the present MPI process. For sequential vectors, this number is the same as size(), but for parallel vectors it may be smaller.

To figure out which elements exactly are stored locally, use local_range().

If the vector contains ghost elements, they are included in this number.

◆ local_range()

std::pair<size_type, size_type> TrilinosWrappers::VectorBase::local_range ( ) const

Return a pair of indices indicating which elements of this vector are stored locally. The first number is the index of the first element stored, the second the index of the one past the last one that is stored locally. If this is a sequential vector, then the result will be the pair (0,N), otherwise it will be a pair (i,i+n), where n=local_size() and i is the first element of the vector stored on this processor, corresponding to the half open interval \([i,i+n)\)

Note
The description above is true most of the time, but not always. In particular, Trilinos vectors need not store contiguous ranges of elements such as \([i,i+n)\). Rather, it can store vectors where the elements are distributed in an arbitrary way across all processors and each processor simply stores a particular subset, not necessarily contiguous. In this case, this function clearly makes no sense since it could, at best, return a range that includes all elements that are stored locally. Thus, the function only succeeds if the locally stored range is indeed contiguous. It will trigger an assertion if the local portion of the vector is not contiguous.

◆ in_local_range()

bool TrilinosWrappers::VectorBase::in_local_range ( const size_type  index) const

Return whether index is in the local range or not, see also local_range().

Note
The same limitation for the applicability of this function applies as listed in the documentation of local_range().

◆ locally_owned_elements()

IndexSet TrilinosWrappers::VectorBase::locally_owned_elements ( ) const

Return an index set that describes which elements of this vector are owned by the current processor. Note that this index set does not include elements this vector may store locally as ghost elements but that are in fact owned by another processor. As a consequence, the index sets returned on different processors if this is a distributed vector will form disjoint sets that add up to the complete index set. Obviously, if a vector is created on only one processor, then the result would satisfy

vec.locally_owned_elements() == complete_index_set (vec.size())

◆ has_ghost_elements()

bool TrilinosWrappers::VectorBase::has_ghost_elements ( ) const

Return if the vector contains ghost elements. This answer is true if there are ghost elements on at least one process.

See also
vectors with ghost elements

◆ operator*()

TrilinosScalar TrilinosWrappers::VectorBase::operator* ( const VectorBase vec) const

Return the scalar (inner) product of two vectors. The vectors must have the same size.

◆ norm_sqr()

real_type TrilinosWrappers::VectorBase::norm_sqr ( ) const

Return square of the \(l_2\)-norm.

◆ mean_value()

TrilinosScalar TrilinosWrappers::VectorBase::mean_value ( ) const

Mean value of the elements of this vector.

◆ minimal_value()

TrilinosScalar TrilinosWrappers::VectorBase::minimal_value ( ) const

Compute the minimal value of the elements of this vector.

This function is deprecated use min() instead.

◆ min()

TrilinosScalar TrilinosWrappers::VectorBase::min ( ) const

Compute the minimal value of the elements of this vector.

◆ max()

TrilinosScalar TrilinosWrappers::VectorBase::max ( ) const

Compute the maximal value of the elements of this vector.

◆ l1_norm()

real_type TrilinosWrappers::VectorBase::l1_norm ( ) const

\(l_1\)-norm of the vector. The sum of the absolute values.

◆ l2_norm()

real_type TrilinosWrappers::VectorBase::l2_norm ( ) const

\(l_2\)-norm of the vector. The square root of the sum of the squares of the elements.

◆ lp_norm()

real_type TrilinosWrappers::VectorBase::lp_norm ( const TrilinosScalar  p) const

\(l_p\)-norm of the vector. The pth root of the sum of the pth powers of the absolute values of the elements.

◆ linfty_norm()

real_type TrilinosWrappers::VectorBase::linfty_norm ( ) const

Maximum absolute value of the elements.

◆ add_and_dot()

TrilinosScalar TrilinosWrappers::VectorBase::add_and_dot ( const TrilinosScalar  a,
const VectorBase V,
const VectorBase W 
)

Performs a combined operation of a vector addition and a subsequent inner product, returning the value of the inner product. In other words, the result of this function is the same as if the user called

this->add(a, V);
return_value = *this * W;

The reason this function exists is for compatibility with deal.II's own vector classes which can implement this functionality with less memory transfer. However, for Trilinos vectors such a combined operation is not natively supported and thus the cost is completely equivalent as calling the two methods separately.

◆ all_zero()

bool TrilinosWrappers::VectorBase::all_zero ( ) const

Return whether the vector contains only elements with value zero. This is a collective operation. This function is expensive, because potentially all elements have to be checked.

Definition at line 352 of file trilinos_vector_base.cc.

◆ is_non_negative()

bool TrilinosWrappers::VectorBase::is_non_negative ( ) const

Return true if the vector has no negative entries, i.e. all entries are zero or positive. This function is used, for example, to check whether refinement indicators are really all positive (or zero).

Definition at line 387 of file trilinos_vector_base.cc.

◆ operator()() [1/2]

reference TrilinosWrappers::VectorBase::operator() ( const size_type  index)

Provide access to a given element, both read and write.

When using a vector distributed with MPI, this operation only makes sense for elements that are actually present on the calling processor. Otherwise, an exception is thrown. This is different from the el() function below that always succeeds (but returns zero on non-local elements).

◆ operator()() [2/2]

TrilinosScalar TrilinosWrappers::VectorBase::operator() ( const size_type  index) const

Provide read-only access to an element.

When using a vector distributed with MPI, this operation only makes sense for elements that are actually present on the calling processor. Otherwise, an exception is thrown. This is different from the el() function below that always succeeds (but returns zero on non-local elements).

Definition at line 265 of file trilinos_vector_base.cc.

◆ operator[]() [1/2]

reference TrilinosWrappers::VectorBase::operator[] ( const size_type  index)

Provide access to a given element, both read and write.

Exactly the same as operator().

◆ operator[]() [2/2]

TrilinosScalar TrilinosWrappers::VectorBase::operator[] ( const size_type  index) const

Provide read-only access to an element.

Exactly the same as operator().

◆ extract_subvector_to() [1/2]

void TrilinosWrappers::VectorBase::extract_subvector_to ( const std::vector< size_type > &  indices,
std::vector< TrilinosScalar > &  values 
) const

A collective get operation: instead of getting individual elements of a vector, this function allows to get a whole set of elements at once. The indices of the elements to be read are stated in the first argument, the corresponding values are returned in the second.

◆ extract_subvector_to() [2/2]

template<typename ForwardIterator , typename OutputIterator >
void TrilinosWrappers::VectorBase::extract_subvector_to ( ForwardIterator  indices_begin,
const ForwardIterator  indices_end,
OutputIterator  values_begin 
) const

Just as the above, but with pointers. Useful in minimizing copying of data around.

◆ el()

TrilinosScalar TrilinosWrappers::VectorBase::el ( const size_type  index) const

Return the value of the vector entry i. Note that this function does only work properly when we request a data stored on the local processor. In case the elements sits on another process, this function returns 0 which might or might not be appropriate in a given situation. If you rely on consistent results, use the access functions () or [] that throw an assertion in case a non-local element is used.

This function is deprecated.

Definition at line 248 of file trilinos_vector_base.cc.

◆ begin() [1/2]

iterator TrilinosWrappers::VectorBase::begin ( )

Make the Vector class a bit like the vector<> class of the C++ standard library by returning iterators to the start and end of the locally owned elements of this vector. The ordering of local elements corresponds to the one given by the global indices in case the vector is constructed from an IndexSet or other methods in deal.II (note that an Epetra_Map can contain elements in arbitrary orders, though).

It holds that end() - begin() == local_size().

◆ begin() [2/2]

const_iterator TrilinosWrappers::VectorBase::begin ( ) const

Return constant iterator to the start of the locally owned elements of the vector.

◆ end() [1/2]

iterator TrilinosWrappers::VectorBase::end ( )

Return an iterator pointing to the element past the end of the array of locally owned entries.

◆ end() [2/2]

const_iterator TrilinosWrappers::VectorBase::end ( ) const

Return a constant iterator pointing to the element past the end of the array of the locally owned entries.

◆ set() [1/3]

void TrilinosWrappers::VectorBase::set ( const std::vector< size_type > &  indices,
const std::vector< TrilinosScalar > &  values 
)

A collective set operation: instead of setting individual elements of a vector, this function allows to set a whole set of elements at once. The indices of the elements to be set are stated in the first argument, the corresponding values in the second.

◆ set() [2/3]

void TrilinosWrappers::VectorBase::set ( const std::vector< size_type > &  indices,
const ::Vector< TrilinosScalar > &  values 
)

This is a second collective set operation. As a difference, this function takes a deal.II vector of values.

◆ set() [3/3]

void TrilinosWrappers::VectorBase::set ( const size_type  n_elements,
const size_type *  indices,
const TrilinosScalar *  values 
)

This collective set operation is of lower level and can handle anything else — the only thing you have to provide is an address where all the indices are stored and the number of elements to be set.

◆ add() [1/7]

void TrilinosWrappers::VectorBase::add ( const std::vector< size_type > &  indices,
const std::vector< TrilinosScalar > &  values 
)

A collective add operation: This function adds a whole set of values stored in values to the vector components specified by indices.

◆ add() [2/7]

void TrilinosWrappers::VectorBase::add ( const std::vector< size_type > &  indices,
const ::Vector< TrilinosScalar > &  values 
)

This is a second collective add operation. As a difference, this function takes a deal.II vector of values.

◆ add() [3/7]

void TrilinosWrappers::VectorBase::add ( const size_type  n_elements,
const size_type *  indices,
const TrilinosScalar *  values 
)

Take an address where n_elements are stored contiguously and add them into the vector. Handles all cases which are not covered by the other two add() functions above.

◆ operator*=()

VectorBase& TrilinosWrappers::VectorBase::operator*= ( const TrilinosScalar  factor)

Multiply the entire vector by a fixed factor.

◆ operator/=()

VectorBase& TrilinosWrappers::VectorBase::operator/= ( const TrilinosScalar  factor)

Divide the entire vector by a fixed factor.

◆ operator+=()

VectorBase& TrilinosWrappers::VectorBase::operator+= ( const VectorBase V)

Add the given vector to the present one.

◆ operator-=()

VectorBase& TrilinosWrappers::VectorBase::operator-= ( const VectorBase V)

Subtract the given vector from the present one.

◆ add() [4/7]

void TrilinosWrappers::VectorBase::add ( const TrilinosScalar  s)

Addition of s to all components. Note that s is a scalar and not a vector.

◆ add() [5/7]

void TrilinosWrappers::VectorBase::add ( const VectorBase V,
const bool  allow_different_maps = false 
)

Simple vector addition, equal to the operator +=.

Though, if the second argument allow_different_maps is set, then it is possible to add data from a vector that uses a different map, i.e., a vector whose elements are split across processors differently. This may include vectors with ghost elements, for example. In general, however, adding vectors with a different element-to- processor map requires communicating data among processors and, consequently, is a slower operation than when using vectors using the same map.

Definition at line 289 of file trilinos_vector_base.cc.

◆ add() [6/7]

void TrilinosWrappers::VectorBase::add ( const TrilinosScalar  a,
const VectorBase V 
)

Simple addition of a multiple of a vector, i.e. *this += a*V.

◆ add() [7/7]

void TrilinosWrappers::VectorBase::add ( const TrilinosScalar  a,
const VectorBase V,
const TrilinosScalar  b,
const VectorBase W 
)

Multiple addition of scaled vectors, i.e. *this += a*V + b*W.

◆ sadd() [1/4]

void TrilinosWrappers::VectorBase::sadd ( const TrilinosScalar  s,
const VectorBase V 
)

Scaling and simple vector addition, i.e. *this = s*(*this) + V.

◆ sadd() [2/4]

void TrilinosWrappers::VectorBase::sadd ( const TrilinosScalar  s,
const TrilinosScalar  a,
const VectorBase V 
)

Scaling and simple addition, i.e. *this = s*(*this) + a*V.

◆ sadd() [3/4]

void TrilinosWrappers::VectorBase::sadd ( const TrilinosScalar  s,
const TrilinosScalar  a,
const VectorBase V,
const TrilinosScalar  b,
const VectorBase W 
)

Scaling and multiple addition.

This function is deprecated.

◆ sadd() [4/4]

void TrilinosWrappers::VectorBase::sadd ( const TrilinosScalar  s,
const TrilinosScalar  a,
const VectorBase V,
const TrilinosScalar  b,
const VectorBase W,
const TrilinosScalar  c,
const VectorBase X 
)

Scaling and multiple addition. *this = s*(*this) + a*V + b*W + c*X.

This function is deprecated.

◆ scale()

void TrilinosWrappers::VectorBase::scale ( const VectorBase scaling_factors)

Scale each element of this vector by the corresponding element in the argument. This function is mostly meant to simulate multiplication (and immediate re-assignment) by a diagonal scaling matrix.

◆ equ() [1/2]

void TrilinosWrappers::VectorBase::equ ( const TrilinosScalar  a,
const VectorBase V 
)

Assignment *this = a*V.

◆ equ() [2/2]

void TrilinosWrappers::VectorBase::equ ( const TrilinosScalar  a,
const VectorBase V,
const TrilinosScalar  b,
const VectorBase W 
)

Assignment *this = a*V + b*W.

This function is deprecated.

Definition at line 427 of file trilinos_vector_base.cc.

◆ ratio()

void TrilinosWrappers::VectorBase::ratio ( const VectorBase a,
const VectorBase b 
)

Compute the elementwise ratio of the two given vectors, that is let this[i] = a[i]/b[i]. This is useful for example if you want to compute the cellwise ratio of true to estimated error.

This vector is appropriately scaled to hold the result.

If any of the b[i] is zero, the result is undefined. No attempt is made to catch such situations.

◆ trilinos_vector() [1/2]

const Epetra_MultiVector& TrilinosWrappers::VectorBase::trilinos_vector ( ) const

Return a const reference to the underlying Trilinos Epetra_MultiVector class.

◆ trilinos_vector() [2/2]

Epetra_FEVector& TrilinosWrappers::VectorBase::trilinos_vector ( )

Return a (modifyable) reference to the underlying Trilinos Epetra_FEVector class.

◆ vector_partitioner()

const Epetra_Map& TrilinosWrappers::VectorBase::vector_partitioner ( ) const

Return a const reference to the underlying Trilinos Epetra_Map that sets the parallel partitioning of the vector.

◆ print() [1/2]

void TrilinosWrappers::VectorBase::print ( const char *  format = 0) const

Output of vector in user-defined format in analogy to the Vector class.

This function is deprecated.

Definition at line 469 of file trilinos_vector_base.cc.

◆ print() [2/2]

void TrilinosWrappers::VectorBase::print ( std::ostream &  out,
const unsigned int  precision = 3,
const bool  scientific = true,
const bool  across = true 
) const

Print to a stream. precision denotes the desired precision with which values shall be printed, scientific whether scientific notation shall be used. If across is true then the vector is printed in a line, while if false then the elements are printed on a separate line each.

Definition at line 489 of file trilinos_vector_base.cc.

◆ swap()

void TrilinosWrappers::VectorBase::swap ( VectorBase v)

Swap the contents of this vector and the other vector v. One could do this operation with a temporary variable and copying over the data elements, but this function is significantly more efficient since it only swaps the pointers to the data of the two vectors and therefore does not need to allocate temporary storage and move data around. Note that the vectors need to be of the same size and base on the same map.

This function is analogous to the the swap function of all C++ standard containers. Also, there is a global function swap(u,v) that simply calls u.swap(v), again in analogy to standard functions.

Definition at line 530 of file trilinos_vector_base.cc.

◆ memory_consumption()

std::size_t TrilinosWrappers::VectorBase::memory_consumption ( ) const

Estimate for the memory consumption in bytes.

Definition at line 540 of file trilinos_vector_base.cc.

◆ get_mpi_communicator()

const MPI_Comm& TrilinosWrappers::VectorBase::get_mpi_communicator ( ) const

Return a reference to the MPI communicator object in use with this object.

Friends And Related Function Documentation

◆ internal::VectorReference

friend class internal::VectorReference
friend

Make the reference class a friend.

Definition at line 954 of file trilinos_vector_base.h.

◆ Vector

friend class Vector
friend

Make all other ReadWriteVector types friends.

Typedef for the vector type used.

Definition at line 955 of file trilinos_vector_base.h.

◆ swap()

void swap ( VectorBase u,
VectorBase v 
)
related

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.

Author
Martin Kronbichler, Wolfgang Bangerth, 2008

Definition at line 973 of file trilinos_vector_base.h.

Member Data Documentation

◆ last_action

Epetra_CombineMode TrilinosWrappers::VectorBase::last_action
private

Trilinos doesn't allow to mix additions to matrix entries and overwriting them (to make synchronisation of parallel computations simpler). The way we do it is to, for each access operation, store whether it is an insertion or an addition. If the previous one was of different type, then we first have to flush the Trilinos buffers; otherwise, we can simply go on. Luckily, Trilinos has an object for this which does already all the parallel communications in such a case, so we simply use their model, which stores whether the last operation was an addition or an insertion.

Definition at line 918 of file trilinos_vector_base.h.

◆ compressed

bool TrilinosWrappers::VectorBase::compressed
private

A boolean variable to hold information on whether the vector is compressed or not.

Definition at line 924 of file trilinos_vector_base.h.

◆ has_ghosts

bool TrilinosWrappers::VectorBase::has_ghosts
private

Whether this vector has ghost elements. This is true on all processors even if only one of them has any ghost elements.

Definition at line 930 of file trilinos_vector_base.h.

◆ vector

std_cxx11::shared_ptr<Epetra_FEVector> TrilinosWrappers::VectorBase::vector
private

Pointer to the actual Epetra vector object. This may represent a vector that is in fact distributed among multiple processors. The object requires an existing Epetra_Map for storing data when setting it up.

Definition at line 937 of file trilinos_vector_base.h.

◆ nonlocal_vector

std_cxx11::shared_ptr<Epetra_MultiVector> TrilinosWrappers::VectorBase::nonlocal_vector
private

A vector object in Trilinos to be used for collecting the non-local elements if the vector was constructed with an additional IndexSet describing ghost elements.

Definition at line 944 of file trilinos_vector_base.h.

◆ owned_elements

IndexSet TrilinosWrappers::VectorBase::owned_elements
private

An IndexSet storing the indices this vector owns exclusively.

Definition at line 949 of file trilinos_vector_base.h.


The documentation for this class was generated from the following files: