Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
DiagonalMatrix< VectorType > Class Template Reference

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

Inheritance diagram for DiagonalMatrix< VectorType >:
[legend]

Public Types

using value_type = typename VectorType::value_type
 
using size_type = typename VectorType::size_type
 

Public Member Functions

 DiagonalMatrix ()=default
 
 DiagonalMatrix (const VectorType &vec)
 
void reinit (const VectorType &vec)
 
void compress (VectorOperation::values operation)
 
VectorType & get_vector ()
 
void clear ()
 
const VectorType & get_vector () const
 
size_type m () const
 
size_type n () const
 
value_type operator() (const size_type i, const size_type j) const
 
value_typeoperator() (const size_type i, const size_type j)
 
template<typename number2 >
void add (const size_type row, const size_type n_cols, const size_type *col_indices, const number2 *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
 
void add (const size_type i, const size_type j, const value_type value)
 
void vmult (VectorType &dst, const VectorType &src) const
 
void Tvmult (VectorType &dst, const VectorType &src) const
 
void vmult_add (VectorType &dst, const VectorType &src) const
 
void Tvmult_add (VectorType &dst, const VectorType &src) const
 
value_type apply (const unsigned int index, const value_type src) const
 
void apply_to_subrange (const unsigned int begin_range, const unsigned int end_range, const value_type *src_pointer_to_current_range, value_type *dst_pointer_to_current_range) const
 
void initialize_dof_vector (VectorType &dst) const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void check_no_subscribers () const noexcept
 

Private Attributes

VectorType diagonal
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

template<typename VectorType = Vector<double>>
class DiagonalMatrix< VectorType >

This class represents a n x n diagonal matrix based on a vector of size n. The matrix-vector products are realized by VectorType::scale, so the template vector class needs to provide a scale() method.

When using this class with ConstraintsMatrix::distribute_local_to_global(), the underlying vector needs to provide write access to all entries referenced by cells in an assembly process. This means that this class also needs access to ghost entries that are owned by other processors than the calling one. In practice this requires initialization of the vector as follows

diagonal_matrix.get_vector();
diagonal_vector.reinit(locally_owned_dofs,
locally_relevant_dofs,
mpi_communicator);
VectorType & get_vector()
void reinit(const size_type size, const bool omit_zeroing_entries=false)

Definition at line 62 of file diagonal_matrix.h.

Member Typedef Documentation

◆ value_type

template<typename VectorType = Vector<double>>
using DiagonalMatrix< VectorType >::value_type = typename VectorType::value_type

Definition at line 65 of file diagonal_matrix.h.

◆ size_type

template<typename VectorType = Vector<double>>
using DiagonalMatrix< VectorType >::size_type = typename VectorType::size_type

Definition at line 66 of file diagonal_matrix.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 230 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 235 of file subscriptor.h.

Constructor & Destructor Documentation

◆ DiagonalMatrix() [1/2]

template<typename VectorType = Vector<double>>
DiagonalMatrix< VectorType >::DiagonalMatrix ( )
default

Default constructor. The object needs still to be reinitialized to be usable.

◆ DiagonalMatrix() [2/2]

template<typename VectorType = Vector<double>>
DiagonalMatrix< VectorType >::DiagonalMatrix ( const VectorType &  vec)
explicit

Constructor initializing this object as a diagonal matrix of size n x n where n is the size of the vector, and with diagonal entries equal to the elements of vec.

Member Function Documentation

◆ reinit()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::reinit ( const VectorType &  vec)

Initialize with a given vector by copying the content of the vector vec.

◆ compress()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::compress ( VectorOperation::values  operation)

Compresses the data structures and allows the resulting matrix to be used in all other operations like matrix-vector products. This is a collective operation, i.e., it needs to be run on all processors when used in parallel.

◆ get_vector() [1/2]

template<typename VectorType = Vector<double>>
VectorType & DiagonalMatrix< VectorType >::get_vector ( )

Return a reference to the underlying vector for manipulation of the entries on the matrix diagonal.

◆ clear()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::clear ( )

Clear content of this object and reset to the state of default constructor.

◆ get_vector() [2/2]

template<typename VectorType = Vector<double>>
const VectorType & DiagonalMatrix< VectorType >::get_vector ( ) const

Return a read-only reference to the underlying vector.

◆ m()

template<typename VectorType = Vector<double>>
size_type DiagonalMatrix< VectorType >::m ( ) const

Number of rows of this matrix. This number corresponds to the size of the underlying vector.

◆ n()

template<typename VectorType = Vector<double>>
size_type DiagonalMatrix< VectorType >::n ( ) const

Number of columns of this matrix. This number corresponds to the size of the underlying vector.

◆ operator()() [1/2]

template<typename VectorType = Vector<double>>
value_type DiagonalMatrix< VectorType >::operator() ( const size_type  i,
const size_type  j 
) const

Read-only access to a value. This is restricted to the case where i==j due to the matrix storage.

If the vector representing the diagonal is distributed with MPI, not all of the indices i might actually be accessible. Refer to the method get_vector().locally_owned_elements() for the entries that actually are accessible.

◆ operator()() [2/2]

template<typename VectorType = Vector<double>>
value_type & DiagonalMatrix< VectorType >::operator() ( const size_type  i,
const size_type  j 
)

Read-write access to a value. This is restricted to the case where i==j due to the matrix storage.

If the vector representing the diagonal is distributed with MPI, not all of the indices i might actually be accessible. Refer to the method get_vector().locally_owned_elements() for the entries that actually are accessible.

◆ add() [1/2]

template<typename VectorType = Vector<double>>
template<typename number2 >
void DiagonalMatrix< VectorType >::add ( const size_type  row,
const size_type  n_cols,
const size_type col_indices,
const number2 *  values,
const bool  elide_zero_values = true,
const bool  col_indices_are_sorted = false 
)

Add an array of values given by values in the given global matrix row at columns specified by col_indices. Due to the storage of this matrix, entries are only added to the diagonal of the matrix. All other entries are ignored and no exception is thrown.

This function is for a consistent interface with the other matrix classes in deal.II and can be used in AffineConstraints::distribute_local_to_global to get exactly the same diagonal as when assembling into a sparse matrix.

◆ add() [2/2]

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::add ( const size_type  i,
const size_type  j,
const value_type  value 
)

Add value to the element (i,j).

Due to the storage of this matrix, entries are only added to the diagonal of the matrix. All other entries are ignored and no exception is thrown.

◆ vmult()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::vmult ( VectorType &  dst,
const VectorType &  src 
) const

Performs a matrix-vector multiplication with the given matrix.

◆ Tvmult()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::Tvmult ( VectorType &  dst,
const VectorType &  src 
) const

Performs a transpose matrix-vector multiplication with the given matrix. Since this represents a diagonal matrix, exactly the same as vmult().

◆ vmult_add()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::vmult_add ( VectorType &  dst,
const VectorType &  src 
) const

Adds the result of a matrix-vector multiplication into the destination vector dst. Needs to create a temporary vector, which makes performance slower than for vmult().

◆ Tvmult_add()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::Tvmult_add ( VectorType &  dst,
const VectorType &  src 
) const

Adds the result of a transpose matrix-vector multiplication into the destination vector dst. Needs to create a temporary vector, which makes performance slower than for Tvmult().

◆ apply()

template<typename VectorType = Vector<double>>
value_type DiagonalMatrix< VectorType >::apply ( const unsigned int  index,
const value_type  src 
) const

Apply the preconditioner to a single vector entry. Note that index of the unknown needs to be expressed by an MPI-local index as it would be accessed in the action on a vector.

◆ apply_to_subrange()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::apply_to_subrange ( const unsigned int  begin_range,
const unsigned int  end_range,
const value_type src_pointer_to_current_range,
value_type dst_pointer_to_current_range 
) const

Apply the preconditioner only to a subrange of elements in an array src, and store the result in another array dst, for compatibility with classes support vector operations on a slice of entries, such as SolverCG or PreconditionChebyshev. Note that the range indicates MPI-local indices as they would be accessed in the action on a vector. The pointers are supposed to point to the beginning of the given range.

◆ initialize_dof_vector()

template<typename VectorType = Vector<double>>
void DiagonalMatrix< VectorType >::initialize_dof_vector ( VectorType &  dst) const

Initialize vector dst to have the same size and partition as diagonal member of this class.

This is a part of the interface required by linear_operator().

◆ memory_consumption()

template<typename VectorType = Vector<double>>
std::size_t DiagonalMatrix< VectorType >::memory_consumption ( ) const

Return the memory consumption of this object.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 136 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 156 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 53 of file subscriptor.cc.

Member Data Documentation

◆ diagonal

template<typename VectorType = Vector<double>>
VectorType DiagonalMatrix< VectorType >::diagonal
private

The stored vector.

Definition at line 256 of file diagonal_matrix.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 219 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 225 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 241 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 249 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when printing out the list of subscribers.

Definition at line 271 of file subscriptor.h.


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