Reference documentation for deal.II version 9.2.0
|
The classes in this module are wrappers around functionality provided by the Trilinos library. They provide a modern object-oriented interface that is compatible with the interfaces of the other linear algebra classes in deal.II. All classes and functions in this group reside in a namespace TrilinosWrappers
.
More...
Namespaces | |
TrilinosWrappers::internal | |
TrilinosWrappers::internal::BlockLinearOperatorImplementation | |
TrilinosWrappers | |
TrilinosWrappers::MPI | |
internal | |
internal::LinearOperatorImplementation | |
Access to underlying Trilinos data | |
Epetra_Operator & | TrilinosWrappers::PreconditionBase::trilinos_operator () const |
Partitioners | |
class | TrilinosWrappers::PreconditionBase::SolverBase |
Teuchos::RCP< Epetra_Operator > | TrilinosWrappers::PreconditionBase::preconditioner |
Epetra_MpiComm | TrilinosWrappers::PreconditionBase::communicator |
std::shared_ptr< Epetra_Map > | TrilinosWrappers::PreconditionBase::vector_distributor |
IndexSet | TrilinosWrappers::PreconditionBase::locally_owned_domain_indices () const |
IndexSet | TrilinosWrappers::PreconditionBase::locally_owned_range_indices () const |
static ::ExceptionBase & | TrilinosWrappers::PreconditionBase::ExcNonMatchingMaps (std::string arg1) |
2: Data-Access | |
reference | TrilinosWrappers::MPI::Vector::operator() (const size_type index) |
TrilinosScalar | TrilinosWrappers::MPI::Vector::operator() (const size_type index) const |
reference | TrilinosWrappers::MPI::Vector::operator[] (const size_type index) |
TrilinosScalar | TrilinosWrappers::MPI::Vector::operator[] (const size_type index) const |
void | TrilinosWrappers::MPI::Vector::extract_subvector_to (const std::vector< size_type > &indices, std::vector< TrilinosScalar > &values) const |
template<typename ForwardIterator , typename OutputIterator > | |
void | TrilinosWrappers::MPI::Vector::extract_subvector_to (ForwardIterator indices_begin, const ForwardIterator indices_end, OutputIterator values_begin) const |
iterator | TrilinosWrappers::MPI::Vector::begin () |
const_iterator | TrilinosWrappers::MPI::Vector::begin () const |
iterator | TrilinosWrappers::MPI::Vector::end () |
const_iterator | TrilinosWrappers::MPI::Vector::end () const |
The classes in this module are wrappers around functionality provided by the Trilinos library. They provide a modern object-oriented interface that is compatible with the interfaces of the other linear algebra classes in deal.II. All classes and functions in this group reside in a namespace TrilinosWrappers
.
These classes are only available if a Trilinos installation was detected during configuration of deal.II. Refer to the README file for more details about this.
Typedef the base class for simpler access to its own alias.
Definition at line 78 of file trilinos_block_sparse_matrix.h.
Typedef the type of the underlying matrix.
Definition at line 83 of file trilinos_block_sparse_matrix.h.
Import the alias from the base class.
Definition at line 88 of file trilinos_block_sparse_matrix.h.
Definition at line 89 of file trilinos_block_sparse_matrix.h.
Definition at line 90 of file trilinos_block_sparse_matrix.h.
Definition at line 91 of file trilinos_block_sparse_matrix.h.
Definition at line 92 of file trilinos_block_sparse_matrix.h.
Definition at line 93 of file trilinos_block_sparse_matrix.h.
Definition at line 94 of file trilinos_block_sparse_matrix.h.
Definition at line 95 of file trilinos_block_sparse_matrix.h.
using TrilinosWrappers::internal::BlockLinearOperatorImplementation::TrilinosBlockPayload< PayloadBlockType >::BlockType = PayloadBlockType |
Type of payload held by each subblock
Definition at line 590 of file trilinos_block_sparse_matrix.h.
Typedef the base class for simpler access to its own alias.
Definition at line 81 of file trilinos_parallel_block_vector.h.
Typedef the type of the underlying vector.
Definition at line 86 of file trilinos_parallel_block_vector.h.
Import the alias from the base class.
Definition at line 91 of file trilinos_parallel_block_vector.h.
Definition at line 92 of file trilinos_parallel_block_vector.h.
Definition at line 93 of file trilinos_parallel_block_vector.h.
Definition at line 94 of file trilinos_parallel_block_vector.h.
Definition at line 95 of file trilinos_parallel_block_vector.h.
Definition at line 96 of file trilinos_parallel_block_vector.h.
Definition at line 97 of file trilinos_parallel_block_vector.h.
Definition at line 98 of file trilinos_parallel_block_vector.h.
Declare the type for container size.
Definition at line 88 of file trilinos_precondition.h.
Declare some of the standard types used in all containers. These types parallel those in the C
standard libraries vector<...>
class.
Definition at line 408 of file trilinos_vector.h.
Definition at line 409 of file trilinos_vector.h.
Definition at line 410 of file trilinos_vector.h.
Definition at line 411 of file trilinos_vector.h.
using TrilinosWrappers::MPI::Vector::const_iterator = const value_type * |
Definition at line 412 of file trilinos_vector.h.
using TrilinosWrappers::MPI::Vector::reference = internal::VectorReference |
Definition at line 413 of file trilinos_vector.h.
using TrilinosWrappers::MPI::Vector::const_reference = const internal::VectorReference |
Definition at line 414 of file trilinos_vector.h.
|
default |
Initialize the matrix empty, that is with no memory allocated. This is useful if you want such objects as member variables in other classes. You can make the structure usable by calling the reinit() function.
Initialize the matrix with the given number of block rows and columns. The blocks themselves are still empty, and you have to call collect_sizes() after you assign them sizes.
Definition at line 511 of file block_sparsity_pattern.cc.
BlockSparsityPattern::BlockSparsityPattern | ( | const std::vector< size_type > & | row_block_sizes, |
const std::vector< size_type > & | col_block_sizes | ||
) |
Initialize the pattern with two BlockIndices for the block structures of matrix rows and columns. This function is equivalent to calling the previous constructor with the length of the two index vector and then entering the index values.
Definition at line 518 of file block_sparsity_pattern.cc.
BlockSparsityPattern::BlockSparsityPattern | ( | const std::vector< IndexSet > & | parallel_partitioning, |
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Initialize the pattern with an array of index sets that specifies both rows and columns of the matrix (so the final matrix will be a square matrix), where the size() of the IndexSets specifies the size of the blocks and the values in each IndexSet denotes the rows that are going to be saved in each block.
Definition at line 532 of file block_sparsity_pattern.cc.
BlockSparsityPattern::BlockSparsityPattern | ( | const std::vector< IndexSet > & | row_parallel_partitioning, |
const std::vector< IndexSet > & | column_parallel_partitioning, | ||
const std::vector< IndexSet > & | writeable_rows, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Initialize the pattern with two arrays of index sets that specify rows and columns of the matrix, where the size() of the IndexSets specifies the size of the blocks and the values in each IndexSet denotes the rows that are going to be saved in each block. The additional index set writable_rows is used to set all rows that we allow to write locally. This constructor is used to create matrices that allow several threads to write simultaneously into the matrix (to different rows, of course), see the method TrilinosWrappers::SparsityPattern::reinit method with three index set arguments for more details.
Definition at line 548 of file block_sparsity_pattern.cc.
void BlockSparsityPattern::reinit | ( | const std::vector< size_type > & | row_block_sizes, |
const std::vector< size_type > & | col_block_sizes | ||
) |
Resize the matrix to a tensor product of matrices with dimensions defined by the arguments.
The matrix will have as many block rows and columns as there are entries in the two arguments. The block at position (i,j) will have the dimensions row_block_sizes[i]
times col_block_sizes[j]
.
Definition at line 569 of file block_sparsity_pattern.cc.
void BlockSparsityPattern::reinit | ( | const std::vector< IndexSet > & | parallel_partitioning, |
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Resize the matrix to a square tensor product of matrices. See the constructor that takes a vector of IndexSets for details.
Definition at line 583 of file block_sparsity_pattern.cc.
void BlockSparsityPattern::reinit | ( | const std::vector< IndexSet > & | row_parallel_partitioning, |
const std::vector< IndexSet > & | column_parallel_partitioning, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Resize the matrix to a rectangular block matrices. This method allows rows and columns to be different, both in the outer block structure and within the blocks.
Definition at line 600 of file block_sparsity_pattern.cc.
void BlockSparsityPattern::reinit | ( | const std::vector< IndexSet > & | row_parallel_partitioning, |
const std::vector< IndexSet > & | column_parallel_partitioning, | ||
const std::vector< IndexSet > & | writeable_rows, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Resize the matrix to a rectangular block matrices that furthermore explicitly specify the writable rows in each of the blocks. This method is used to create matrices that allow several threads to write simultaneously into the matrix (to different rows, of course), see the method TrilinosWrappers::SparsityPattern::reinit method with three index set arguments for more details.
Definition at line 618 of file block_sparsity_pattern.cc.
|
default |
Constructor; initializes the matrix to be empty, without any structure, i.e. the matrix is not usable at all. This constructor is therefore only useful for matrices which are members of a class. All other matrices should be created at a point in the data flow where all necessary information is available.
You have to initialize the matrix before usage with reinit(BlockSparsityPattern). The number of blocks per row and column are then determined by that function.
|
override |
Destructor.
Definition at line 27 of file trilinos_block_sparse_matrix.cc.
|
default |
Pseudo copy operator only copying empty objects. The sizes of the block matrices need to be the same.
|
inline |
This operator assigns a scalar to a matrix. Since this does usually not make much sense (should we set all matrix entries to this value? Only the nonzero entries of the sparsity pattern?), this operation is only allowed if the actual value to be assigned is zero. This operator only exists to allow for the obvious notation matrix=0
, which sets all elements of the matrix to zero, but keep the sparsity pattern previously used.
Definition at line 412 of file trilinos_block_sparse_matrix.h.
void BlockSparseMatrix< number >::reinit | ( | const size_type | n_block_rows, |
const size_type | n_block_columns | ||
) |
Resize the matrix, by setting the number of block rows and columns. This deletes all blocks and replaces them with uninitialized ones, i.e. ones for which also the sizes are not yet set. You have to do that by calling the reinit
functions of the blocks themselves. Do not forget to call collect_sizes() after that on this object.
The reason that you have to set sizes of the blocks yourself is that the sizes may be varying, the maximum number of elements per row may be varying, etc. It is simpler not to reproduce the interface of the SparsityPattern
class here but rather let the user call whatever function they desire.
Definition at line 36 of file petsc_parallel_block_sparse_matrix.cc.
void BlockSparseMatrix< BlockSparsityPatternType >::reinit | ( | const std::vector< IndexSet > & | input_maps, |
const BlockSparsityPatternType & | block_sparsity_pattern, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const bool | exchange_data = false |
||
) |
Resize the matrix, by using an array of index sets to determine the parallel distribution of the individual matrices. This function assumes that a quadratic block matrix is generated.
Definition at line 73 of file trilinos_block_sparse_matrix.cc.
void BlockSparseMatrix< BlockSparsityPatternType >::reinit | ( | const BlockSparsityPatternType & | block_sparsity_pattern | ) |
Resize the matrix and initialize it by the given sparsity pattern. Since no distribution map is given, the result is a block matrix for which all elements are stored locally.
Definition at line 127 of file trilinos_block_sparse_matrix.cc.
void BlockSparseMatrix< number >::reinit | ( | const std::vector< IndexSet > & | parallel_partitioning, |
const ::BlockSparseMatrix< double > & | dealii_block_sparse_matrix, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const double | drop_tolerance = 1e-13 |
||
) |
This function initializes the Trilinos matrix using the deal.II sparse matrix and the entries stored therein. It uses a threshold to copy only elements whose modulus is larger than the threshold (so zeros in the deal.II matrix can be filtered away).
Definition at line 164 of file trilinos_block_sparse_matrix.cc.
void BlockSparseMatrix< number >::reinit | ( | const ::BlockSparseMatrix< double > & | deal_ii_sparse_matrix, |
const double | drop_tolerance = 1e-13 |
||
) |
This function initializes the Trilinos matrix using the deal.II sparse matrix and the entries stored therein. It uses a threshold to copy only elements whose modulus is larger than the threshold (so zeros in the deal.II matrix can be filtered away). Since no Epetra_Map is given, all the elements will be locally stored.
Definition at line 201 of file trilinos_block_sparse_matrix.cc.
|
inline |
Return the state of the matrix, i.e., whether compress() needs to be called after an operation requiring data exchange. Does only return non-true values when used in debug
mode, since it is quite expensive to keep track of all operations that lead to the need for compress().
Definition at line 426 of file trilinos_block_sparse_matrix.h.
void BlockSparseMatrix< number >::collect_sizes | ( | ) |
This function collects the sizes of the sub-objects and stores them in internal arrays, in order to be able to relay global indices into the matrix to indices into the subobjects. You must call this function each time after you have changed the size of the sub-objects. Note that this is a collective operation, i.e., it needs to be called on all MPI processes. This command internally calls the method compress()
, so you don't need to call that function in case you use collect_sizes()
.
Definition at line 227 of file trilinos_block_sparse_matrix.cc.
BlockSparseMatrix::size_type BlockSparseMatrix< number >::n_nonzero_elements | ( | ) | const |
Return the total number of nonzero elements of this matrix (summed over all MPI processes).
Definition at line 236 of file trilinos_block_sparse_matrix.cc.
MPI_Comm BlockSparseMatrix< number >::get_mpi_communicator | ( | ) | const |
Return the MPI communicator object in use with this matrix.
Definition at line 309 of file trilinos_block_sparse_matrix.cc.
|
inline |
Return the partitioning of the domain space for the individual blocks of this matrix, i.e., the partitioning of the block vectors this matrix has to be multiplied with.
Definition at line 533 of file trilinos_block_sparse_matrix.h.
|
inline |
Return the partitioning of the range space for the individual blocks of this matrix, i.e., the partitioning of the block vectors that result from matrix-vector products.
Definition at line 549 of file trilinos_block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication: let \(dst = M*src\) with \(M\) being this matrix. The vector types can be block vectors or non-block vectors (only if the matrix has only one row or column, respectively), and need to define TrilinosWrappers::SparseMatrix::vmult.
Definition at line 444 of file trilinos_block_sparse_matrix.h.
|
inline |
Matrix-vector multiplication: let \(dst = M^T*src\) with \(M\) being this matrix. This function does the same as vmult() but takes the transposed matrix.
Definition at line 457 of file trilinos_block_sparse_matrix.h.
TrilinosScalar BlockSparseMatrix< number >::residual | ( | MPI::BlockVector & | dst, |
const MPI::BlockVector & | x, | ||
const MPI::BlockVector & | b | ||
) | const |
Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst
. The l2 norm of the residual vector is returned.
Source x and destination dst must not be the same vector.
Note that both vectors have to be distributed vectors generated using the same Map as was used for the matrix.
This function only applicable if the matrix only has one block row.
Definition at line 249 of file trilinos_block_sparse_matrix.cc.
TrilinosScalar BlockSparseMatrix< number >::residual | ( | MPI::BlockVector & | dst, |
const MPI::Vector & | x, | ||
const MPI::BlockVector & | b | ||
) | const |
Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst
. The l2 norm of the residual vector is returned.
This function is only applicable if the matrix only has one block row.
Definition at line 267 of file trilinos_block_sparse_matrix.cc.
TrilinosScalar BlockSparseMatrix< number >::residual | ( | MPI::Vector & | dst, |
const MPI::BlockVector & | x, | ||
const MPI::Vector & | b | ||
) | const |
Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst
. The l2 norm of the residual vector is returned.
This function is only applicable if the matrix only has one block column.
Definition at line 281 of file trilinos_block_sparse_matrix.cc.
TrilinosScalar BlockSparseMatrix< number >::residual | ( | MPI::Vector & | dst, |
const MPI::Vector & | x, | ||
const MPI::Vector & | b | ||
) | const |
Compute the residual of an equation Mx=b, where the residual is defined to be r=b-Mx. Write the residual into dst
. The l2 norm of the residual vector is returned.
This function is only applicable if the matrix only has one block.
Definition at line 295 of file trilinos_block_sparse_matrix.cc.
|
inlineprivate |
Internal version of (T)vmult with two block vectors
Definition at line 470 of file trilinos_block_sparse_matrix.h.
|
inlineprivate |
Internal version of (T)vmult where the source vector is a block vector but the destination vector is a non-block vector
Definition at line 486 of file trilinos_block_sparse_matrix.h.
|
inlineprivate |
Internal version of (T)vmult where the source vector is a non-block vector but the destination vector is a block vector
Definition at line 502 of file trilinos_block_sparse_matrix.h.
|
inlineprivate |
Internal version of (T)vmult where both source vector and the destination vector are non-block vectors (only defined if the matrix consists of only one block)
Definition at line 518 of file trilinos_block_sparse_matrix.h.
|
inline |
Default constructor
This simply checks that the payload for each block has been chosen correctly (i.e. is of type TrilinosPayload). Apart from this, this class does not do anything in particular and needs no special configuration, we have only one generic constructor that can be called under any conditions.
Definition at line 602 of file trilinos_block_sparse_matrix.h.
|
inline |
A function that encapsulates generic matrix
objects, based on an operator_exemplar
, that act on a compatible Vector type into a LinearOperator.
This function is the equivalent of the linear_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 77 of file trilinos_linear_operator.h.
|
inline |
A function that encapsulates generic matrix
objects that act on a compatible Vector type into a LinearOperator.
This function is the equivalent of the linear_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 108 of file trilinos_linear_operator.h.
|
inline |
A function that encapsulates a block_matrix
into a BlockLinearOperator.
This function is the equivalent of the block_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 145 of file trilinos_linear_operator.h.
|
inline |
A variant of above function that builds up a block diagonal linear operator from an array ops
of diagonal elements (off-diagonal blocks are assumed to be 0).
This function is the equivalent of the block_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 182 of file trilinos_linear_operator.h.
|
inline |
This function extracts the diagonal blocks of block_matrix
(either a block matrix type or a BlockLinearOperator) and creates a BlockLinearOperator with the diagonal. Off-diagonal elements are initialized as null_operator (with correct reinit_range_vector and reinit_domain_vector methods).
This function is the equivalent of the block_diagonal_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 224 of file trilinos_linear_operator.h.
|
inline |
A variant of above function that builds up a block diagonal linear operator from an array ops
of diagonal elements (off-diagonal blocks are assumed to be 0).
This function is the equivalent of the block_diagonal_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 260 of file trilinos_linear_operator.h.
|
default |
Default constructor. Generate an empty vector without any blocks.
|
inlineexplicit |
Constructor. Generate a block vector with as many blocks as there are entries in partitioning
. Each IndexSet together with the MPI communicator contains the layout of the distribution of data among the MPI processes.
Definition at line 318 of file trilinos_parallel_block_vector.h.
|
inline |
Creates a BlockVector with ghost elements. See the respective reinit() method for more details. ghost_values
may contain any elements in parallel_partitioning
, they will be ignored.
Definition at line 327 of file trilinos_parallel_block_vector.h.
|
inline |
Copy-Constructor. Set all the properties of the parallel vector to those of the given argument and copy the elements.
Definition at line 348 of file trilinos_parallel_block_vector.h.
|
inlinenoexcept |
Move constructor. Creates a new vector by stealing the internal data of the vector v
.
Definition at line 360 of file trilinos_parallel_block_vector.h.
|
inlineexplicit |
Creates a block vector consisting of num_blocks
components, but there is no content in the individual components and the user has to fill appropriate data using a reinit of the blocks.
Definition at line 341 of file trilinos_parallel_block_vector.h.
|
overridedefault |
Destructor. Clears memory
BlockVector & BlockVector< Number >::operator= | ( | const value_type | s | ) |
Copy operator: fill all components of the vector that are locally stored with the given scalar value.
Definition at line 31 of file trilinos_block_vector.cc.
BlockVector & BlockVector< Number >::operator= | ( | const BlockVector & | v | ) |
Copy operator for arguments of the same type.
Definition at line 40 of file trilinos_block_vector.cc.
|
noexcept |
Move the given vector. This operator replaces the present vector with v
by efficiently swapping the internal data structures.
Definition at line 61 of file trilinos_block_vector.cc.
BlockVector & BlockVector< Number >::operator= | ( | const ::BlockVector< Number > & | v | ) |
Another copy function. This one takes a deal.II block vector and copies it into a TrilinosWrappers block vector. Note that the number of blocks has to be the same in the vector as in the input vector. Use the reinit() command for resizing the BlockVector or for changing the internal structure of the block components.
Since Trilinos only works on doubles, this function is limited to accept only one possible number type in the deal.II vector.
Definition at line 371 of file trilinos_parallel_block_vector.h.
void BlockVector< Number >::reinit | ( | const std::vector< IndexSet > & | parallel_partitioning, |
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const bool | omit_zeroing_entries = false |
||
) |
Reinitialize the BlockVector to contain as many blocks as there are index sets given in the input argument, according to the parallel distribution of the individual components described in the maps.
If omit_zeroing_entries==false
, the vector is filled with zeros.
Definition at line 70 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const std::vector< IndexSet > & | partitioning, |
const std::vector< IndexSet > & | ghost_values, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const bool | vector_writable = false |
||
) |
Reinit functionality. This function destroys the old vector content and generates a new one based on the input partitioning. In addition to just specifying one index set as in all the other methods above, this method allows to supply an additional set of ghost entries. There are two different versions of a vector that can be created. If the flag vector_writable
is set to false
, the vector only allows read access to the joint set of parallel_partitioning
and ghost_entries
. The effect of the reinit method is then equivalent to calling the other reinit method with an index set containing both the locally owned entries and the ghost entries.
If the flag vector_writable
is set to true, this creates an alternative storage scheme for ghost elements that allows multiple threads to write into the vector (for the other reinit methods, only one thread is allowed to write into the ghost entries at a time).
Definition at line 95 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const BlockVector & | V, |
const bool | omit_zeroing_entries = false |
||
) |
Change the dimension to that of the vector V
. The same applies as for the other reinit() function.
The elements of V
are not copied, i.e. this function is the same as calling reinit (V.size(), omit_zeroing_entries)
.
Note that you must call this (or the other reinit() functions) function, rather than calling the reinit() functions of an individual block, to allow the block vector to update its caches of vector sizes. If you call reinit() on one of the blocks, then subsequent actions on this object may yield unpredictable results since they may be routed to the wrong block.
Definition at line 123 of file trilinos_block_vector.cc.
void BlockVector< Number >::reinit | ( | const size_type | num_blocks | ) |
Change the number of blocks to num_blocks
. The individual blocks will get initialized with zero size, so it is assumed that the user resizes the individual blocks by herself in an appropriate way, and calls collect_sizes
afterwards.
Definition at line 138 of file trilinos_block_vector.cc.
void BlockVector< Number >::import_nonlocal_data_for_fe | ( | const TrilinosWrappers::BlockSparseMatrix & | m, |
const BlockVector & | v | ||
) |
This reinit function is meant to be used for parallel calculations where some non-local data has to be used. The typical situation where one needs this function is the call of the FEValues<dim>::get_function_values function (or of some derivatives) in parallel. Since it is usually faster to retrieve the data in advance, this function can be called before the assembly forks out to the different processors. What this function does is the following: It takes the information in the columns of the given matrix and looks which data couples between the different processors. That data is then queried from the input vector. Note that you should not write to the resulting vector any more, since the some data can be stored several times on different processors, leading to unpredictable results. In particular, such a vector cannot be used for matrix- vector products as for example done during the solution of linear systems.
Definition at line 154 of file trilinos_block_vector.cc.
|
inline |
Return if this Vector contains ghost elements.
Definition at line 392 of file trilinos_parallel_block_vector.h.
|
inline |
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.
Limitation: right now this function only works if both vectors have the same number of blocks. If needed, the numbers of blocks should be exchanged, too.
This function is analogous to 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 405 of file trilinos_parallel_block_vector.h.
void BlockVector< Number >::print | ( | std::ostream & | out, |
const unsigned int | precision = 3 , |
||
const bool | scientific = true , |
||
const bool | across = true |
||
) | const |
Print to a stream.
Definition at line 178 of file trilinos_block_vector.cc.
|
inline |
Global function 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 423 of file trilinos_parallel_block_vector.h.
|
inlinestatic |
Definition at line 452 of file trilinos_parallel_block_vector.h.
|
inlinestatic |
Definition at line 463 of file trilinos_parallel_block_vector.h.
TrilinosWrappers::PreconditionBase::PreconditionBase | ( | ) |
Constructor. Does not do anything. The initialize
function of the derived classes will have to create the preconditioner from a given sparse matrix.
Definition at line 34 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionBase::PreconditionBase | ( | const PreconditionBase & | base | ) |
Copy constructor.
Definition at line 42 of file trilinos_precondition.cc.
|
overridedefault |
Destructor.
void TrilinosWrappers::PreconditionBase::clear | ( | ) |
Destroys the preconditioner, leaving an object like just after having called the constructor.
Definition at line 56 of file trilinos_precondition.cc.
MPI_Comm TrilinosWrappers::PreconditionBase::get_mpi_communicator | ( | ) | const |
Return the MPI communicator object in use with this matrix.
Definition at line 67 of file trilinos_precondition.cc.
void TrilinosWrappers::PreconditionBase::transpose | ( | ) |
Sets an internal flag so that all operations performed by the matrix, i.e., multiplications, are done in transposed order. However, this does not reshape the matrix to transposed form directly, so care should be taken when using this flag.
|
virtual |
Apply the preconditioner.
Reimplemented in TrilinosWrappers::PreconditionIdentity.
|
virtual |
Apply the transpose preconditioner.
Reimplemented in TrilinosWrappers::PreconditionIdentity.
|
virtual |
Apply the preconditioner on deal.II data structures instead of the ones provided in the Trilinos wrapper class.
Reimplemented in TrilinosWrappers::PreconditionIdentity.
|
virtual |
Apply the transpose preconditioner on deal.II data structures instead of the ones provided in the Trilinos wrapper class.
Reimplemented in TrilinosWrappers::PreconditionIdentity.
|
virtual |
Apply the preconditioner on deal.II parallel data structures instead of the ones provided in the Trilinos wrapper class.
|
virtual |
Apply the transpose preconditioner on deal.II parallel data structures instead of the ones provided in the Trilinos wrapper class.
Epetra_Operator & TrilinosWrappers::PreconditionBase::trilinos_operator | ( | ) | const |
Calling this function from an uninitialized object will cause an exception.
Definition at line 78 of file trilinos_precondition.cc.
IndexSet TrilinosWrappers::PreconditionBase::locally_owned_domain_indices | ( | ) | const |
Return the partitioning of the domain space of this matrix, i.e., the partitioning of the vectors this matrix has to be multiplied with.
Definition at line 87 of file trilinos_precondition.cc.
IndexSet TrilinosWrappers::PreconditionBase::locally_owned_range_indices | ( | ) | const |
Return the partitioning of the range space of this matrix, i.e., the partitioning of the vectors that are result from matrix-vector products.
Definition at line 94 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionJacobi::AdditionalData::AdditionalData | ( | const double | omega = 1 , |
const double | min_diagonal = 0 , |
||
const unsigned int | n_sweeps = 1 |
||
) |
Constructor. By default, set the damping parameter to one, and do not modify the diagonal.
Definition at line 101 of file trilinos_precondition.cc.
void PreconditionJacobi< MatrixType >::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Take the sparse matrix the preconditioner object should be built of, and additional flags (damping parameter, etc.) if there are any.
Definition at line 113 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionSSOR::AdditionalData::AdditionalData | ( | const double | omega = 1 , |
const double | min_diagonal = 0 , |
||
const unsigned int | overlap = 0 , |
||
const unsigned int | n_sweeps = 1 |
||
) |
Constructor. By default, set the damping parameter to one, we do not modify the diagonal, and there is no overlap (i.e. in parallel, we run a BlockJacobi preconditioner, where each block is inverted approximately by an SSOR).
Definition at line 153 of file trilinos_precondition.cc.
void PreconditionSSOR< MatrixType >::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Take the sparse matrix the preconditioner object should be built of, and additional flags (damping parameter, overlap in parallel computations, etc.) if there are any.
Definition at line 166 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionSOR::AdditionalData::AdditionalData | ( | const double | omega = 1 , |
const double | min_diagonal = 0 , |
||
const unsigned int | overlap = 0 , |
||
const unsigned int | n_sweeps = 1 |
||
) |
Constructor. By default, set the damping parameter to one, we do not modify the diagonal, and there is no overlap (i.e. in parallel, we run a BlockJacobi preconditioner, where each block is inverted approximately by an SOR.
Definition at line 206 of file trilinos_precondition.cc.
void PreconditionSOR< MatrixType >::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Take the sparse matrix the preconditioner object should be built of, and additional flags (damping parameter, overlap in parallel computations etc.) if there are any.
Definition at line 219 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionBlockJacobi::AdditionalData::AdditionalData | ( | const unsigned int | block_size = 1 , |
const std::string & | block_creation_type = "linear" , |
||
const double | omega = 1 , |
||
const double | min_diagonal = 0 , |
||
const unsigned int | n_sweeps = 1 |
||
) |
Constructor. By default, use a block size of 1, use linear subdivision of the rows, set the damping parameter to one, and do not modify the diagonal.
Definition at line 259 of file trilinos_precondition.cc.
void PreconditionBlockJacobi< MatrixType, inverse_type >::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Take the sparse matrix the preconditioner object should be built of, and additional flags (damping parameter, etc.) if there are any.
Definition at line 275 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionBlockSSOR::AdditionalData::AdditionalData | ( | const unsigned int | block_size = 1 , |
const std::string & | block_creation_type = "linear" , |
||
const double | omega = 1 , |
||
const double | min_diagonal = 0 , |
||
const unsigned int | overlap = 0 , |
||
const unsigned int | n_sweeps = 1 |
||
) |
Constructor. By default, use a block size of 1, use linear subdivision of the rows, set the damping parameter to one, we do not modify the diagonal, and there is no overlap (i.e. in parallel, we run a BlockJacobi preconditioner, where each block is inverted approximately by a block SOR).
Definition at line 326 of file trilinos_precondition.cc.
void PreconditionBlockSSOR< MatrixType, inverse_type >::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Take the sparse matrix the preconditioner object should be built of, and additional flags (damping parameter, overlap in parallel computations, etc.) if there are any.
Definition at line 344 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionBlockSOR::AdditionalData::AdditionalData | ( | const unsigned int | block_size = 1 , |
const std::string & | block_creation_type = "linear" , |
||
const double | omega = 1 , |
||
const double | min_diagonal = 0 , |
||
const unsigned int | overlap = 0 , |
||
const unsigned int | n_sweeps = 1 |
||
) |
Constructor. By default, use a block size of 1, use linear subdivision of the rows, set the damping parameter to one, we do not modify the diagonal, and there is no overlap (i.e. in parallel, we run a BlockJacobi preconditioner, where each block is inverted approximately by a block SOR).
Definition at line 395 of file trilinos_precondition.cc.
void PreconditionBlockSOR< MatrixType, inverse_type >::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Take the sparse matrix the preconditioner object should be built of, and additional flags (damping parameter, overlap in parallel computations etc.) if there are any.
Definition at line 413 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionIC::AdditionalData::AdditionalData | ( | const unsigned int | ic_fill = 0 , |
const double | ic_atol = 0. , |
||
const double | ic_rtol = 1. , |
||
const unsigned int | overlap = 0 |
||
) |
Constructor. By default, set the drop tolerance to 0, the level of extra fill-ins is set to be zero (just use the matrix structure, do not generate any additional fill-in), the tolerance level are 0 and 1, respectively, and the overlap in case of a parallel execution is zero. This overlap in a block-application of the IC in the parallel case makes the preconditioner a so-called additive Schwarz preconditioner.
Definition at line 464 of file trilinos_precondition.cc.
void TrilinosWrappers::PreconditionIC::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initialize function. Takes the matrix the preconditioner should be computed of, and additional flags if there are any.
Definition at line 477 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionILU::AdditionalData::AdditionalData | ( | const unsigned int | ilu_fill = 0 , |
const double | ilu_atol = 0. , |
||
const double | ilu_rtol = 1. , |
||
const unsigned int | overlap = 0 |
||
) |
Constructor with default values for all parameters.
Definition at line 514 of file trilinos_precondition.cc.
void TrilinosWrappers::PreconditionILU::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initialize function. Takes the matrix which is used to form the preconditioner, and additional flags if there are any.
Definition at line 527 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionILUT::AdditionalData::AdditionalData | ( | const double | ilut_drop = 0. , |
const unsigned int | ilut_fill = 0 , |
||
const double | ilut_atol = 0. , |
||
const double | ilut_rtol = 1. , |
||
const unsigned int | overlap = 0 |
||
) |
Constructor. By default, no element will be dropped, the level of extra fill-ins is set to be zero (just use the matrix structure, do not generate any additional fill-in except the one that results from non-dropping large elements), the tolerance level are 0 and 1, respectively, and the overlap in case of a parallel execution is zero. This overlap in a block-application of the ILU in the parallel case makes the preconditioner a so-called additive Schwarz preconditioner.
Definition at line 565 of file trilinos_precondition.cc.
void TrilinosWrappers::PreconditionILUT::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initialize function. Takes the matrix which is used to form the preconditioner, and additional flags if there are any.
Definition at line 580 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionBlockwiseDirect::AdditionalData::AdditionalData | ( | const unsigned int | overlap = 0 | ) |
Constructor.
Definition at line 619 of file trilinos_precondition.cc.
void TrilinosWrappers::PreconditionBlockwiseDirect::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initialize function. Takes the matrix which is used to form the preconditioner, and additional flags if there are any.
Definition at line 627 of file trilinos_precondition.cc.
PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::AdditionalData::AdditionalData | ( | const unsigned int | degree = 1 , |
const double | max_eigenvalue = 10. , |
||
const double | eigenvalue_ratio = 30. , |
||
const double | min_eigenvalue = 1. , |
||
const double | min_diagonal = 1e-12 , |
||
const bool | nonzero_starting = false |
||
) |
Constructor.
Definition at line 661 of file trilinos_precondition.cc.
void PreconditionChebyshev< MatrixType, VectorType, PreconditionerType >::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Initialize function. Takes the matrix which is used to form the preconditioner, and additional flags if there are any.
Definition at line 679 of file trilinos_precondition.cc.
TrilinosWrappers::PreconditionAMG::AdditionalData::AdditionalData | ( | const bool | elliptic = true , |
const bool | higher_order_elements = false , |
||
const unsigned int | n_cycles = 1 , |
||
const bool | w_cyle = false , |
||
const double | aggregation_threshold = 1e-4 , |
||
const std::vector< std::vector< bool >> & | constant_modes = std::vector<std::vector<bool>>(0) , |
||
const unsigned int | smoother_sweeps = 2 , |
||
const unsigned int | smoother_overlap = 0 , |
||
const bool | output_details = false , |
||
const char * | smoother_type = "Chebyshev" , |
||
const char * | coarse_type = "Amesos-KLU" |
||
) |
Constructor. By default, we pretend to work on elliptic problems with linear finite elements on a scalar equation.
Making use of the DoFTools::extract_constant_modes() function, the constant_modes
vector can be initialized for a given field in the following manner:
Definition at line 39 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_parameters | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const Epetra_RowMatrix & | matrix | ||
) | const |
Fill in a parameter_list
that can be used to initialize the AMG preconditioner.
The matrix
is used in conjunction with the constant_modes
to configure the null space settings for the preconditioner. The distributed_constant_modes
are initialized by this function, and must remain in scope until PreconditionAMG::initialize() has been called.
See the documentation for the Trilinos ML package for details on what options are available for modification.
Definition at line 67 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_parameters | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const SparseMatrix & | matrix | ||
) | const |
Fill in a parameter list that can be used to initialize the AMG preconditioner.
Definition at line 188 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_operator_null_space | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const Epetra_RowMatrix & | matrix | ||
) | const |
Configure the null space setting in the parameter_list
for the input matrix
based on the state of the constant_modes
variable.
Definition at line 127 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::AdditionalData::set_operator_null_space | ( | Teuchos::ParameterList & | parameter_list, |
std::unique_ptr< Epetra_MultiVector > & | distributed_constant_modes, | ||
const SparseMatrix & | matrix | ||
) | const |
Configure the null space setting in the parameter_list
for the input matrix
based on the state of the constant_modes
variable.
Definition at line 201 of file trilinos_precondition_ml.cc.
|
override |
Destructor.
Definition at line 213 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. The function uses the matrix format specified in TrilinosWrappers::SparseMatrix.
Definition at line 222 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::initialize | ( | const Epetra_RowMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. As opposed to the other initialize function above, this function uses an abstract interface to an object of type Epetra_RowMatrix which allows a user to pass quite general objects to the ML preconditioner.
This initialization routine is useful in cases where the operator to be preconditioned is not a TrilinosWrappers::SparseMatrix object but still allows getting a copy of the entries in each of the locally owned matrix rows (method ExtractMyRowCopy) and implements a matrix-vector product (methods Multiply or Apply). An example are operators which provide faster matrix-vector multiplications than possible with matrix entries (matrix-free methods). These implementations can be beneficially combined with Chebyshev smoothers that only perform matrix-vector products. The interface class Epetra_RowMatrix is very flexible to enable this kind of implementation.
Definition at line 231 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::initialize | ( | const SparseMatrix & | matrix, |
const Teuchos::ParameterList & | ml_parameters | ||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. The function uses the matrix format specified in TrilinosWrappers::SparseMatrix.
This function is similar to the one above, but allows the user to set all the options of the Trilinos ML preconditioner. In order to find out about all the options for ML, we refer to the ML documentation page. In particular, users need to follow the ML instructions in case a vector-valued problem ought to be solved.
Definition at line 257 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::initialize | ( | const Epetra_RowMatrix & | matrix, |
const Teuchos::ParameterList & | ml_parameters | ||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. As opposed to the other initialize function above, this function uses an abstract interface to an object of type Epetra_RowMatrix which allows a user to pass quite general objects to the ML preconditioner.
Definition at line 266 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::initialize | ( | const ::SparseMatrix< number > & | deal_ii_sparse_matrix, |
const AdditionalData & | additional_data = AdditionalData() , |
||
const double | drop_tolerance = 1e-13 , |
||
const ::SparsityPattern * | use_this_sparsity = nullptr |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. This function takes a deal.II matrix and copies the content into a Trilinos matrix, so the function can be considered rather inefficient.
Definition at line 277 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::reinit | ( | ) |
This function can be used for a faster recalculation of the preconditioner construction when the matrix entries underlying the preconditioner have changed, but the matrix sparsity pattern has remained the same. What this function does is taking the already generated coarsening structure, computing the AMG prolongation and restriction according to a smoothed aggregation strategy and then building the whole multilevel hierarchy. This function can be considerably faster than the initialize function, since the coarsening pattern is usually the most difficult thing to do when setting up the AMG ML preconditioner.
Definition at line 311 of file trilinos_precondition_ml.cc.
void TrilinosWrappers::PreconditionAMG::clear | ( | ) |
Destroys the preconditioner, leaving an object like just after having called the constructor.
Definition at line 321 of file trilinos_precondition_ml.cc.
PreconditionAMG::size_type TrilinosWrappers::PreconditionAMG::memory_consumption | ( | ) | const |
Prints an estimate of the memory consumption of this class.
Definition at line 330 of file trilinos_precondition_ml.cc.
TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::AdditionalData | ( | const bool | elliptic = true , |
const unsigned int | n_cycles = 1 , |
||
const bool | w_cyle = false , |
||
const double | aggregation_threshold = 1e-4 , |
||
const std::vector< std::vector< bool >> & | constant_modes = std::vector<std::vector<bool>>(0) , |
||
const unsigned int | smoother_sweeps = 2 , |
||
const unsigned int | smoother_overlap = 0 , |
||
const bool | output_details = false , |
||
const char * | smoother_type = "Chebyshev" , |
||
const char * | coarse_type = "Amesos-KLU" |
||
) |
Constructor. By default, we pretend to work on elliptic problems with linear finite elements on a scalar equation.
Definition at line 31 of file trilinos_precondition_muelu.cc.
TrilinosWrappers::PreconditionAMGMueLu::PreconditionAMGMueLu | ( | ) |
Constructor.
Definition at line 56 of file trilinos_precondition_muelu.cc.
|
overridevirtualdefault |
Destructor.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. The function uses the matrix format specified in TrilinosWrappers::SparseMatrix.
Definition at line 73 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const Epetra_CrsMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. As opposed to the other initialize function above, this function uses an object of type Epetra_CrsMatrixCrs.
Definition at line 82 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const SparseMatrix & | matrix, |
Teuchos::ParameterList & | muelu_parameters | ||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. The function uses the matrix format specified in TrilinosWrappers::SparseMatrix.
This function is similar to the one above, but allows the user to set most of the options of the Trilinos ML preconditioner. In order to find out about all the options for ML, we refer to the ML documentation page. Not all ML options have a corresponding MueLu option.
Definition at line 186 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const Epetra_CrsMatrix & | matrix, |
Teuchos::ParameterList & | muelu_parameters | ||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. As opposed to the other initialize function above, this function uses an object of type Epetra_CrsMatrix.
Definition at line 195 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::initialize | ( | const ::SparseMatrix< number > & | deal_ii_sparse_matrix, |
const AdditionalData & | additional_data = AdditionalData() , |
||
const double | drop_tolerance = 1e-13 , |
||
const ::SparsityPattern * | use_this_sparsity = nullptr |
||
) |
Let Trilinos compute a multilevel hierarchy for the solution of a linear system with the given matrix. This function takes a deal.ii matrix and copies the content into a Trilinos matrix, so the function can be considered rather inefficient.
Definition at line 208 of file trilinos_precondition_muelu.cc.
void TrilinosWrappers::PreconditionAMGMueLu::clear | ( | ) |
Destroys the preconditioner, leaving an object like just after having called the constructor.
Definition at line 242 of file trilinos_precondition_muelu.cc.
PreconditionAMGMueLu::size_type TrilinosWrappers::PreconditionAMGMueLu::memory_consumption | ( | ) | const |
Prints an estimate of the memory consumption of this class.
Definition at line 251 of file trilinos_precondition_muelu.cc.
void PreconditionIdentity::initialize | ( | const SparseMatrix & | matrix, |
const AdditionalData & | additional_data = AdditionalData() |
||
) |
The matrix argument is ignored and here just for compatibility with more complex preconditioners.
Definition at line 722 of file trilinos_precondition.cc.
|
overridevirtual |
Apply the preconditioner, i.e., dst = src.
Reimplemented from TrilinosWrappers::PreconditionBase.
Definition at line 766 of file trilinos_precondition.cc.
|
overridevirtual |
Apply the transport conditioner, i.e., dst = src.
Reimplemented from TrilinosWrappers::PreconditionBase.
Definition at line 772 of file trilinos_precondition.cc.
|
overridevirtual |
Apply the preconditioner on deal.II data structures instead of the ones provided in the Trilinos wrapper class, i.e., dst = src.
Reimplemented from TrilinosWrappers::PreconditionBase.
Definition at line 778 of file trilinos_precondition.cc.
|
overridevirtual |
Apply the transpose preconditioner on deal.II data structures instead of the ones provided in the Trilinos wrapper class, i.e. dst = src.
Reimplemented from TrilinosWrappers::PreconditionBase.
Definition at line 785 of file trilinos_precondition.cc.
|
override |
Apply the preconditioner on deal.II parallel data structures instead of the ones provided in the Trilinos wrapper class, i.e., dst = src.
|
override |
Apply the transpose preconditioner on deal.II parallel data structures instead of the ones provided in the Trilinos wrapper class, i.e., dst = src.
Definition at line 199 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 and distribution among processes in case of an MPI run.
Definition at line 70 of file trilinos_vector.cc.
Copy constructor using the given vector.
Definition at line 90 of file trilinos_vector.cc.
|
explicit |
This constructor takes an IndexSet that defines how to distribute the individual components among the MPI processors. Since it also includes information about the size of the vector, this is all we need to generate a parallel vector.
Depending on whether the parallel_partitioning
argument uniquely subdivides elements among processors or not, the resulting vector may or may not have ghost elements. See the general documentation of this class for more information.
In case the provided IndexSet forms an overlapping partitioning, it is not clear which elements are owned by which process and locally_owned_elements() will return an IndexSet of size zero.
Definition at line 81 of file trilinos_vector.cc.
Vector< Number >::Vector | ( | const IndexSet & | local, |
const IndexSet & | ghost, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Creates a ghosted parallel vector.
Depending on whether the ghost
argument uniquely subdivides elements among processors or not, the resulting vector may or may not have ghost elements. See the general documentation of this class for more information.
Definition at line 128 of file trilinos_vector.cc.
Vector< Number >::Vector | ( | const IndexSet & | parallel_partitioning, |
const Vector & | v, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Copy constructor from the TrilinosWrappers vector class. Since a vector of this class does not necessarily need to be distributed among processes, the user needs to supply us with an IndexSet and an MPI communicator that set the partitioning details.
Depending on whether the parallel_partitioning
argument uniquely subdivides elements among processors or not, the resulting vector may or may not have ghost elements. See the general documentation of this class for more information.
Definition at line 109 of file trilinos_vector.cc.
TrilinosWrappers::MPI::Vector::Vector | ( | const IndexSet & | parallel_partitioning, |
const ::Vector< Number > & | v, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD |
||
) |
Copy-constructor from deal.II vectors. Sets the dimension to that of the given vector, and copies all the elements.
Depending on whether the parallel_partitioning
argument uniquely subdivides elements among processors or not, the resulting vector may or may not have ghost elements. See the general documentation of this class for more information.
Move constructor. Creates a new vector by stealing the internal data of the vector v
.
Definition at line 100 of file trilinos_vector.cc.
|
overridedefault |
Destructor.
void Vector< Number >::clear | ( | ) |
Release all memory and return to a state just like after having called the default constructor.
Definition at line 139 of file trilinos_vector.cc.
void Vector< Number >::reinit | ( | const Vector & | v, |
const bool | omit_zeroing_entries = false , |
||
const bool | allow_different_maps = false |
||
) |
Reinit functionality. This function sets the calling vector to the dimension and the parallel distribution of the input vector, but does not copy the elements in v
. If 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.
This function has a third argument, allow_different_maps
, that allows for an exchange of data between two equal-sized vectors (but being distributed differently among the processors). A trivial application of this function is to generate a replication of a whole vector on each machine, when the calling vector is built with a map consisting of all indices on each process, and v
is a distributed vector. In this case, the variable omit_zeroing_entries
needs to be set to false
, since it does not make sense to exchange data between differently parallelized vectors without touching the elements.
Definition at line 198 of file trilinos_vector.cc.
void Vector< Number >::reinit | ( | const IndexSet & | parallel_partitioning, |
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const bool | omit_zeroing_entries = false |
||
) |
Reinit functionality. This function destroys the old vector content and generates a new one based on the input partitioning. The flag omit_zeroing_entries
determines whether the vector should be filled with zero (false). If the flag 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.
Depending on whether the parallel_partitioning
argument uniquely subdivides elements among processors or not, the resulting vector may or may not have ghost elements. See the general documentation of this class for more information.
In case parallel_partitioning
is overlapping, it is not clear which process should own which elements. Hence, locally_owned_elements() returns an empty IndexSet in this case.
Definition at line 157 of file trilinos_vector.cc.
void Vector< Number >::reinit | ( | const IndexSet & | locally_owned_entries, |
const IndexSet & | ghost_entries, | ||
const MPI_Comm & | communicator = MPI_COMM_WORLD , |
||
const bool | vector_writable = false |
||
) |
Reinit functionality. This function destroys the old vector content and generates a new one based on the input partitioning. In addition to just specifying one index set as in all the other methods above, this method allows to supply an additional set of ghost entries. There are two different versions of a vector that can be created. If the flag vector_writable
is set to false
, the vector only allows read access to the joint set of parallel_partitioning
and ghost_entries
. The effect of the reinit method is then equivalent to calling the other reinit method with an index set containing both the locally owned entries and the ghost entries.
If the flag vector_writable
is set to true, this creates an alternative storage scheme for ghost elements that allows multiple threads to write into the vector (for the other reinit methods, only one thread is allowed to write into the ghost entries at a time).
Depending on whether the ghost_entries
argument uniquely subdivides elements among processors or not, the resulting vector may or may not have ghost elements. See the general documentation of this class for more information.
Definition at line 360 of file trilinos_vector.cc.
void Vector< Number >::reinit | ( | const BlockVector & | v, |
const bool | import_data = false |
||
) |
Create vector by merging components from a block vector.
Definition at line 283 of file trilinos_vector.cc.
void Vector< Number >::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 583 of file trilinos_vector.cc.
Vector& TrilinosWrappers::MPI::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. the constructor explicit work.
Since the semantics of assigning a scalar to a vector are not immediately clear, this operator can only be used if you want to set the entire vector to zero. This allows the intuitive notation v=0
.
Copy the given vector. Resize the present vector if necessary. In this case, also the Epetra_Map that designs the parallel partitioning is taken from the input vector.
Definition at line 418 of file trilinos_vector.cc.
Move the given vector. This operator replaces the present vector with v
by efficiently swapping the internal data structures.
Definition at line 499 of file trilinos_vector.cc.
Vector& TrilinosWrappers::MPI::Vector::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.
void Vector< Number >::import_nonlocal_data_for_fe | ( | const ::TrilinosWrappers::SparseMatrix & | matrix, |
const Vector & | vector | ||
) |
This reinit function is meant to be used for parallel calculations where some non-local data has to be used. The typical situation where one needs this function is the call of the FEValues<dim>::get_function_values function (or of some derivatives) in parallel. Since it is usually faster to retrieve the data in advance, this function can be called before the assembly forks out to the different processors. What this function does is the following: It takes the information in the columns of the given matrix and looks which data couples between the different processors. That data is then queried from the input vector. Note that you should not write to the resulting vector any more, since the some data can be stored several times on different processors, leading to unpredictable results. In particular, such a vector cannot be used for matrix- vector products as for example done during the solution of linear systems.
Definition at line 528 of file trilinos_vector.cc.
void Vector< Number >::import | ( | const LinearAlgebra::ReadWriteVector< double > & | rwv, |
const VectorOperation::values | operation | ||
) |
Imports all the elements present in the vector's IndexSet from the input vector rwv
. VectorOperation::values operation
is used to decide if the elements in rwv
should be added to the current vector or replace the current elements.
Definition at line 552 of file trilinos_vector.cc.
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 716 of file trilinos_vector.cc.
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 733 of file trilinos_vector.cc.
size_type TrilinosWrappers::MPI::Vector::size | ( | ) | const |
Return the global dimension of the vector.
size_type TrilinosWrappers::MPI::Vector::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.
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)\)
Return whether index
is in the local range or not, see also local_range().
IndexSet TrilinosWrappers::MPI::Vector::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
bool TrilinosWrappers::MPI::Vector::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.
void TrilinosWrappers::MPI::Vector::update_ghost_values | ( | ) | const |
This function only exists for compatibility with the LinearAlgebra::distributed::Vector
class and does nothing: this class implements ghost value updates in a different way that is a better fit with the underlying Trilinos vector object.
TrilinosScalar TrilinosWrappers::MPI::Vector::operator* | ( | const Vector & | vec | ) | const |
Return the scalar (inner) product of two vectors. The vectors must have the same size.
real_type TrilinosWrappers::MPI::Vector::norm_sqr | ( | ) | const |
Return the square of the \(l_2\)-norm.
TrilinosScalar TrilinosWrappers::MPI::Vector::mean_value | ( | ) | const |
Mean value of the elements of this vector.
TrilinosScalar TrilinosWrappers::MPI::Vector::min | ( | ) | const |
Compute the minimal value of the elements of this vector.
TrilinosScalar TrilinosWrappers::MPI::Vector::max | ( | ) | const |
Compute the maximal value of the elements of this vector.
real_type TrilinosWrappers::MPI::Vector::l1_norm | ( | ) | const |
\(l_1\)-norm of the vector. The sum of the absolute values.
real_type TrilinosWrappers::MPI::Vector::l2_norm | ( | ) | const |
\(l_2\)-norm of the vector. The square root of the sum of the squares of the elements.
real_type TrilinosWrappers::MPI::Vector::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.
real_type TrilinosWrappers::MPI::Vector::linfty_norm | ( | ) | const |
Maximum absolute value of the elements.
TrilinosScalar TrilinosWrappers::MPI::Vector::add_and_dot | ( | const TrilinosScalar | a, |
const Vector & | V, | ||
const Vector & | 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
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.
For complex-valued vectors, the scalar product in the second step is implemented as \(\left<v,w\right>=\sum_i v_i \bar{w_i}\).
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 743 of file trilinos_vector.cc.
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 776 of file trilinos_vector.cc.
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.
TrilinosScalar Vector< Number >::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.
Definition at line 655 of file trilinos_vector.cc.
Provide access to a given element, both read and write.
Exactly the same as operator().
TrilinosScalar TrilinosWrappers::MPI::Vector::operator[] | ( | const size_type | index | ) | const |
Provide read-only access to an element.
Exactly the same as operator().
void TrilinosWrappers::MPI::Vector::extract_subvector_to | ( | const std::vector< size_type > & | indices, |
std::vector< TrilinosScalar > & | values | ||
) | const |
Instead of getting individual elements of a vector via operator(), this function allows getting 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.
If the current vector is called v
, then this function is the equivalent to the code
indices
and values
arrays must be identical. void TrilinosWrappers::MPI::Vector::extract_subvector_to | ( | ForwardIterator | indices_begin, |
const ForwardIterator | indices_end, | ||
OutputIterator | values_begin | ||
) | const |
Instead of getting individual elements of a vector via operator(), this function allows getting a whole set of elements at once. In contrast to the previous function, this function obtains the indices of the elements by dereferencing all elements of the iterator range provided by the first two arguments, and puts the vector values into memory locations obtained by dereferencing a range of iterators starting at the location pointed to by the third argument.
If the current vector is called v
, then this function is the equivalent to the code
values_begin
as there are iterators between indices_begin
and indices_end
. iterator TrilinosWrappers::MPI::Vector::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().
const_iterator TrilinosWrappers::MPI::Vector::begin | ( | ) | const |
Return constant iterator to the start of the locally owned elements of the vector.
iterator TrilinosWrappers::MPI::Vector::end | ( | ) |
Return an iterator pointing to the element past the end of the array of locally owned entries.
const_iterator TrilinosWrappers::MPI::Vector::end | ( | ) | const |
Return a constant iterator pointing to the element past the end of the array of the locally owned entries.
void TrilinosWrappers::MPI::Vector::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.
void TrilinosWrappers::MPI::Vector::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.
void TrilinosWrappers::MPI::Vector::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.
void TrilinosWrappers::MPI::Vector::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
.
void TrilinosWrappers::MPI::Vector::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.
void TrilinosWrappers::MPI::Vector::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.
Vector& TrilinosWrappers::MPI::Vector::operator*= | ( | const TrilinosScalar | factor | ) |
Multiply the entire vector by a fixed factor.
Vector& TrilinosWrappers::MPI::Vector::operator/= | ( | const TrilinosScalar | factor | ) |
Divide the entire vector by a fixed factor.
Add the given vector to the present one.
Subtract the given vector from the present one.
void TrilinosWrappers::MPI::Vector::add | ( | const TrilinosScalar | s | ) |
Addition of s
to all components. Note that s
is a scalar and not a vector.
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 681 of file trilinos_vector.cc.
void TrilinosWrappers::MPI::Vector::add | ( | const TrilinosScalar | a, |
const Vector & | V | ||
) |
Simple addition of a multiple of a vector, i.e. *this += a*V
.
void TrilinosWrappers::MPI::Vector::add | ( | const TrilinosScalar | a, |
const Vector & | V, | ||
const TrilinosScalar | b, | ||
const Vector & | W | ||
) |
Multiple addition of scaled vectors, i.e. *this += a*V + b*W
.
void TrilinosWrappers::MPI::Vector::sadd | ( | const TrilinosScalar | s, |
const Vector & | V | ||
) |
Scaling and simple vector addition, i.e. *this = s*(*this) + V
.
void TrilinosWrappers::MPI::Vector::sadd | ( | const TrilinosScalar | s, |
const TrilinosScalar | a, | ||
const Vector & | V | ||
) |
Scaling and simple addition, i.e. *this = s*(*this) + a*V
.
void TrilinosWrappers::MPI::Vector::scale | ( | const Vector & | 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.
void TrilinosWrappers::MPI::Vector::equ | ( | const TrilinosScalar | a, |
const Vector & | V | ||
) |
Assignment *this = a*V
.
const Epetra_MultiVector& TrilinosWrappers::MPI::Vector::trilinos_vector | ( | ) | const |
Return a const reference to the underlying Trilinos Epetra_MultiVector class.
Epetra_FEVector& TrilinosWrappers::MPI::Vector::trilinos_vector | ( | ) |
Return a (modifiable) reference to the underlying Trilinos Epetra_FEVector class.
const Epetra_BlockMap& TrilinosWrappers::MPI::Vector::trilinos_partitioner | ( | ) | const |
Return a const reference to the underlying Trilinos Epetra_BlockMap that sets the parallel partitioning of the vector.
void Vector< Number >::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 815 of file trilinos_vector.cc.
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 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 863 of file trilinos_vector.cc.
std::size_t Vector< Number >::memory_consumption | ( | ) | const |
Estimate for the memory consumption in bytes.
Definition at line 876 of file trilinos_vector.cc.
const MPI_Comm& TrilinosWrappers::MPI::Vector::get_mpi_communicator | ( | ) | const |
Return a reference to the MPI communicator object in use with this object.
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 1347 of file trilinos_vector.h.
|
inlinestatic |
Definition at line 2218 of file trilinos_vector.h.
|
inlinestatic |
Definition at line 2229 of file trilinos_vector.h.
|
related |
A function that encapsulates generic matrix
objects, based on an operator_exemplar
, that act on a compatible Vector type into a LinearOperator.
This function is the equivalent of the linear_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 77 of file trilinos_linear_operator.h.
|
related |
A function that encapsulates generic matrix
objects that act on a compatible Vector type into a LinearOperator.
This function is the equivalent of the linear_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 108 of file trilinos_linear_operator.h.
|
related |
A function that encapsulates a block_matrix
into a BlockLinearOperator.
This function is the equivalent of the block_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 145 of file trilinos_linear_operator.h.
|
related |
A variant of above function that builds up a block diagonal linear operator from an array ops
of diagonal elements (off-diagonal blocks are assumed to be 0).
This function is the equivalent of the block_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 182 of file trilinos_linear_operator.h.
|
related |
This function extracts the diagonal blocks of block_matrix
(either a block matrix type or a BlockLinearOperator) and creates a BlockLinearOperator with the diagonal. Off-diagonal elements are initialized as null_operator (with correct reinit_range_vector and reinit_domain_vector methods).
This function is the equivalent of the block_diagonal_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 224 of file trilinos_linear_operator.h.
|
related |
A variant of above function that builds up a block diagonal linear operator from an array ops
of diagonal elements (off-diagonal blocks are assumed to be 0).
This function is the equivalent of the block_diagonal_operator, but ensures full compatibility with Trilinos operations by preselecting the appropriate template parameters.
Definition at line 260 of file trilinos_linear_operator.h.
|
related |
Global function 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 423 of file trilinos_parallel_block_vector.h.
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 1347 of file trilinos_vector.h.
|
static |
Exception.
|
static |
Exception
|
static |
Exception
|
static |
Exception
|
protected |
This is a pointer to the preconditioner object that is used when applying the preconditioner.
Definition at line 239 of file trilinos_precondition.h.
|
protected |
Internal communication pattern in case the matrix needs to be copied from deal.II format.
Definition at line 246 of file trilinos_precondition.h.
|
protected |
Internal Trilinos map in case the matrix needs to be copied from deal.II format.
Definition at line 255 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionJacobi::AdditionalData::omega |
This specifies the relaxation parameter in the Jacobi preconditioner.
Definition at line 302 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionJacobi::AdditionalData::min_diagonal |
This specifies the minimum value the diagonal elements should have. This might be necessary when the Jacobi preconditioner is used on matrices with zero diagonal elements. In that case, a straight- forward application would not be possible since we would divide by zero.
Definition at line 311 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionJacobi::AdditionalData::n_sweeps |
Sets how many times the given operation should be applied during the vmult() operation.
Definition at line 317 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionSSOR::AdditionalData::omega |
This specifies the (over-) relaxation parameter in the SSOR preconditioner.
Definition at line 390 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionSSOR::AdditionalData::min_diagonal |
This specifies the minimum value the diagonal elements should have. This might be necessary when the SSOR preconditioner is used on matrices with zero diagonal elements. In that case, a straight- forward application would not be possible since we divide by the diagonal element.
Definition at line 399 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionSSOR::AdditionalData::overlap |
This determines how large the overlap of the local matrix portions on each processor in a parallel application should be.
Definition at line 405 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionSSOR::AdditionalData::n_sweeps |
Sets how many times the given operation should be applied during the vmult() operation.
Definition at line 411 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionSOR::AdditionalData::omega |
This specifies the (over-) relaxation parameter in the SOR preconditioner.
Definition at line 485 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionSOR::AdditionalData::min_diagonal |
This specifies the minimum value the diagonal elements should have. This might be necessary when the SOR preconditioner is used on matrices with zero diagonal elements. In that case, a straight- forward application would not be possible since we divide by the diagonal element.
Definition at line 494 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionSOR::AdditionalData::overlap |
This determines how large the overlap of the local matrix portions on each processor in a parallel application should be.
Definition at line 500 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionSOR::AdditionalData::n_sweeps |
Sets how many times the given operation should be applied during the vmult() operation.
Definition at line 506 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockJacobi::AdditionalData::block_size |
This specifies the size of blocks.
Definition at line 574 of file trilinos_precondition.h.
std::string TrilinosWrappers::PreconditionBlockJacobi::AdditionalData::block_creation_type |
Strategy for creation of blocks passed on to Ifpack block relaxation (variable 'partitioner: type') with this string as the given value. Available types in Ifpack include "linear" (i.e., divide the local range of the matrix in slices of the block size), "greedy" "metis". For a full list, see the documentation of Ifpack.
Definition at line 583 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionBlockJacobi::AdditionalData::omega |
This specifies the (over-) relaxation parameter in the Jacobi preconditioner.
Definition at line 589 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionBlockJacobi::AdditionalData::min_diagonal |
This specifies the minimum value the diagonal elements should have. This might be necessary when the block Jacobi preconditioner is used on matrices with zero diagonal elements. In that case, a straight- forward application would not be possible since we divide by the diagonal element.
Definition at line 598 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockJacobi::AdditionalData::n_sweeps |
Sets how many times the given operation should be applied during the vmult() operation.
Definition at line 604 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockSSOR::AdditionalData::block_size |
This specifies the size of blocks.
Definition at line 679 of file trilinos_precondition.h.
std::string TrilinosWrappers::PreconditionBlockSSOR::AdditionalData::block_creation_type |
Strategy for creation of blocks passed on to Ifpack block relaxation (variable 'partitioner: type') with this string as the given value. Available types in Ifpack include "linear" (i.e., divide the local range of the matrix in slices of the block size), "greedy" "metis". For a full list, see the documentation of Ifpack.
Definition at line 688 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionBlockSSOR::AdditionalData::omega |
This specifies the (over-) relaxation parameter in the SOR preconditioner.
Definition at line 694 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionBlockSSOR::AdditionalData::min_diagonal |
This specifies the minimum value the diagonal elements should have. This might be necessary when the SSOR preconditioner is used on matrices with zero diagonal elements. In that case, a straight- forward application would not be possible since we divide by the diagonal element.
Definition at line 703 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockSSOR::AdditionalData::overlap |
This determines how large the overlap of the local matrix portions on each processor in a parallel application should be.
Definition at line 709 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockSSOR::AdditionalData::n_sweeps |
Sets how many times the given operation should be applied during the vmult() operation.
Definition at line 715 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockSOR::AdditionalData::block_size |
This specifies the size of blocks.
Definition at line 791 of file trilinos_precondition.h.
std::string TrilinosWrappers::PreconditionBlockSOR::AdditionalData::block_creation_type |
Strategy for creation of blocks passed on to Ifpack block relaxation (variable 'partitioner: type') with this string as the given value. Available types in Ifpack include "linear" (i.e., divide the local range of the matrix in slices of the block size), "greedy" "metis". For a full list, see the documentation of Ifpack.
Definition at line 800 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionBlockSOR::AdditionalData::omega |
This specifies the (over-) relaxation parameter in the SOR preconditioner.
Definition at line 806 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionBlockSOR::AdditionalData::min_diagonal |
This specifies the minimum value the diagonal elements should have. This might be necessary when the SOR preconditioner is used on matrices with zero diagonal elements. In that case, a straight- forward application would not be possible since we divide by the diagonal element.
Definition at line 815 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockSOR::AdditionalData::overlap |
This determines how large the overlap of the local matrix portions on each processor in a parallel application should be.
Definition at line 821 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockSOR::AdditionalData::n_sweeps |
Sets how many times the given operation should be applied during the vmult() operation.
Definition at line 827 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionIC::AdditionalData::ic_fill |
This specifies the amount of additional fill-in elements besides the sparse matrix structure. When ic_fill
is large, this means that many fill-ins will be added, so that the IC preconditioner comes closer to a direct sparse Cholesky decomposition. Note, however, that this will drastically increase the memory requirement, especially when the preconditioner is used in 3D.
Definition at line 920 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionIC::AdditionalData::ic_atol |
This specifies the amount of an absolute perturbation that will be added to the diagonal of the matrix, which sometimes can help to get better preconditioners.
Definition at line 927 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionIC::AdditionalData::ic_rtol |
This specifies the factor by which the diagonal of the matrix will be scaled, which sometimes can help to get better preconditioners.
Definition at line 933 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionIC::AdditionalData::overlap |
This determines how large the overlap of the local matrix portions on each processor in a parallel application should be.
Definition at line 939 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionILU::AdditionalData::ilu_fill |
Additional fill-in, see class documentation above.
Definition at line 1035 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionILU::AdditionalData::ilu_atol |
The amount of perturbation to add to diagonal entries. See the class documentation above for details.
Definition at line 1041 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionILU::AdditionalData::ilu_rtol |
Scaling actor for diagonal entries. See the class documentation above for details.
Definition at line 1047 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionILU::AdditionalData::overlap |
Overlap between processors. See the class documentation for details.
Definition at line 1052 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionILUT::AdditionalData::ilut_drop |
This specifies the relative size of elements which should be dropped when forming an incomplete LU decomposition with threshold.
Definition at line 1143 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionILUT::AdditionalData::ilut_fill |
This specifies the amount of additional fill-in elements besides the sparse matrix structure. When ilu_fill
is large, this means that many fill-ins will be added, so that the ILU preconditioner comes closer to a (direct) sparse LU decomposition. Note, however, that this will drastically increase the memory requirement, especially when the preconditioner is used in 3D.
Definition at line 1153 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionILUT::AdditionalData::ilut_atol |
This specifies the amount of an absolute perturbation that will be added to the diagonal of the matrix, which sometimes can help to get better preconditioners.
Definition at line 1160 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionILUT::AdditionalData::ilut_rtol |
This specifies the factor by which the diagonal of the matrix will be scaled, which sometimes can help to get better preconditioners.
Definition at line 1166 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionILUT::AdditionalData::overlap |
This determines how large the overlap of the local matrix portions on each processor in a parallel application should be.
Definition at line 1172 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionBlockwiseDirect::AdditionalData::overlap |
This determines how large the overlap of the local matrix portions on each processor in a parallel application should be.
Definition at line 1223 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionChebyshev::AdditionalData::degree |
This determines the degree of the Chebyshev polynomial. The degree of the polynomial gives the number of matrix-vector products to be performed for one application of the vmult() operation.
Definition at line 1270 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionChebyshev::AdditionalData::max_eigenvalue |
This sets the maximum eigenvalue of the matrix, which needs to be set properly for appropriate performance of the Chebyshev preconditioner.
Definition at line 1276 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionChebyshev::AdditionalData::eigenvalue_ratio |
This sets the ratio between the maximum and the minimum eigenvalue.
Definition at line 1281 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionChebyshev::AdditionalData::min_eigenvalue |
This sets the minimum eigenvalue, which is an optional parameter only used internally for checking whether we use an identity matrix.
Definition at line 1287 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionChebyshev::AdditionalData::min_diagonal |
This sets a threshold below which the diagonal element will not be inverted in the Chebyshev algorithm.
Definition at line 1293 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionChebyshev::AdditionalData::nonzero_starting |
When this flag is set to true
, it enables the method vmult(dst, src)
to use non-zero data in the vector dst
, appending to it the Chebyshev corrections. This can be useful in some situations (e.g. when used for high-frequency error smoothing), but not the way the solver classes expect a preconditioner to work (where one ignores the content in dst
for the preconditioner application). The user should really know what they are doing when touching this flag.
Definition at line 1305 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMG::AdditionalData::elliptic |
Determines whether the AMG preconditioner should be optimized for elliptic problems (ML option smoothed aggregation SA, using a Chebyshev smoother) or for non-elliptic problems (ML option non- symmetric smoothed aggregation NSSA, smoother is SSOR with underrelaxation).
Definition at line 1486 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMG::AdditionalData::higher_order_elements |
Determines whether the matrix that the preconditioner is built upon is generated from linear or higher-order elements.
Definition at line 1492 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMG::AdditionalData::n_cycles |
Defines how many multigrid cycles should be performed by the preconditioner.
Definition at line 1498 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMG::AdditionalData::w_cycle |
Defines whether a w-cycle should be used instead of the standard setting of a v-cycle.
Definition at line 1504 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionAMG::AdditionalData::aggregation_threshold |
This threshold tells the AMG setup how the coarsening should be performed. In the AMG used by ML, all points that strongly couple with the tentative coarse-level point form one aggregate. The term strong coupling is controlled by the variable aggregation_threshold
, meaning that all elements that are not smaller than aggregation_threshold
times the diagonal element do couple strongly.
Definition at line 1515 of file trilinos_precondition.h.
std::vector<std::vector<bool> > TrilinosWrappers::PreconditionAMG::AdditionalData::constant_modes |
Specifies the constant modes (near null space) of the matrix. This parameter tells AMG whether we work on a scalar equation (where the near null space only consists of ones, and default value is OK) or on a vector-valued equation. For vector-valued equation problem with n_component
, the provided constant_modes
should fulfill the following requirements:
n_component
ic
][id
] == "id
th DoF is corresponding to component ic
Definition at line 1533 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMG::AdditionalData::smoother_sweeps |
Determines how many sweeps of the smoother should be performed. When the flag elliptic
is set to true
, i.e., for elliptic or almost elliptic problems, the polynomial degree of the Chebyshev smoother is set to smoother_sweeps
. The term sweeps refers to the number of matrix-vector products performed in the Chebyshev case. In the non-elliptic case, smoother_sweeps
sets the number of SSOR relaxation sweeps for post-smoothing to be performed.
Definition at line 1545 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMG::AdditionalData::smoother_overlap |
Determines the overlap in the SSOR/Chebyshev error smoother when run in parallel.
Definition at line 1551 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMG::AdditionalData::output_details |
If this flag is set to true
, then internal information from the ML preconditioner is printed to screen. This can be useful when debugging the preconditioner.
Definition at line 1558 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMG::AdditionalData::smoother_type |
Determines which smoother to use for the AMG cycle. Possibilities for smoother_type are the following:
Definition at line 1594 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMG::AdditionalData::coarse_type |
Determines which solver to use on the coarsest level. The same settings as for the smoother type are possible.
Definition at line 1600 of file trilinos_precondition.h.
|
private |
A copy of the deal.II matrix into Trilinos format.
Definition at line 1714 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::elliptic |
Determines whether the AMG preconditioner should be optimized for elliptic problems (MueLu option smoothed aggregation SA, using a Chebyshev smoother) or for non-elliptic problems (MueLu option non- symmetric smoothed aggregation NSSA, smoother is SSOR with underrelaxation).
Definition at line 1773 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::n_cycles |
Defines how many multigrid cycles should be performed by the preconditioner.
Definition at line 1779 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::w_cycle |
Defines whether a w-cycle should be used instead of the standard setting of a v-cycle.
Definition at line 1785 of file trilinos_precondition.h.
double TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::aggregation_threshold |
This threshold tells the AMG setup how the coarsening should be performed. In the AMG used by MueLu, all points that strongly couple with the tentative coarse-level point form one aggregate. The term strong coupling is controlled by the variable aggregation_threshold
, meaning that all elements that are not smaller than aggregation_threshold
times the diagonal element do couple strongly.
Definition at line 1796 of file trilinos_precondition.h.
std::vector<std::vector<bool> > TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::constant_modes |
Specifies the constant modes (near null space) of the matrix. This parameter tells AMG whether we work on a scalar equation (where the near null space only consists of ones) or on a vector-valued equation.
Definition at line 1804 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::smoother_sweeps |
Determines how many sweeps of the smoother should be performed. When the flag elliptic
is set to true
, i.e., for elliptic or almost elliptic problems, the polynomial degree of the Chebyshev smoother is set to smoother_sweeps
. The term sweeps refers to the number of matrix-vector products performed in the Chebyshev case. In the non-elliptic case, smoother_sweeps
sets the number of SSOR relaxation sweeps for post-smoothing to be performed.
Definition at line 1816 of file trilinos_precondition.h.
unsigned int TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::smoother_overlap |
Determines the overlap in the SSOR/Chebyshev error smoother when run in parallel.
Definition at line 1822 of file trilinos_precondition.h.
bool TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::output_details |
If this flag is set to true
, then internal information from the ML preconditioner is printed to screen. This can be useful when debugging the preconditioner.
Definition at line 1829 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::smoother_type |
Determines which smoother to use for the AMG cycle. Possibilities for smoother_type are the following:
Definition at line 1865 of file trilinos_precondition.h.
const char* TrilinosWrappers::PreconditionAMGMueLu::AdditionalData::coarse_type |
Determines which solver to use on the coarsest level. The same settings as for the smoother type are possible.
Definition at line 1871 of file trilinos_precondition.h.
|
private |
A copy of the deal.II matrix into Trilinos format.
Definition at line 1959 of file trilinos_precondition.h.
|
private |
Trilinos doesn't allow to mix additions to matrix entries and overwriting them (to make synchronization 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 1296 of file trilinos_vector.h.
|
private |
A boolean variable to hold information on whether the vector is compressed or not.
Definition at line 1302 of file trilinos_vector.h.
|
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 1308 of file trilinos_vector.h.
|
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 1315 of file trilinos_vector.h.
|
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 1322 of file trilinos_vector.h.
|
private |
An IndexSet storing the indices this vector owns exclusively.
Definition at line 1327 of file trilinos_vector.h.
|
friend |
Definition at line 232 of file trilinos_precondition.h.
|
friend |
Definition at line 1330 of file trilinos_vector.h.