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
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
SolverBase< VectorType > Class Template Reference

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

Inheritance diagram for SolverBase< VectorType >:
[legend]

Classes

struct  StateCombiner
 

Public Types

using vector_type = VectorType
 

Public Member Functions

 SolverBase (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory)
 
 SolverBase (SolverControl &solver_control)
 
boost::signals2::connection connect (const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate)> &slot)
 
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)
 

Protected Attributes

GrowingVectorMemory< VectorType > static_vector_memory
 
VectorMemory< VectorType > & memory
 
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate), StateCombineriteration_status
 

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

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<class VectorType = Vector<double>>
class SolverBase< VectorType >

A base class for iterative linear solvers. This class provides interfaces to a memory pool and the objects that determine whether a solver has converged.

Requirements common to derived solver classes

In general, iterative solvers do not rely on any special structure of matrices or the format of storage. Rather, they only require that matrices and vectors define certain operations such as matrix-vector products, or scalar products between vectors. Consequently, this class as well as the derived classes and their member functions implementing concrete linear solvers are templated on the types of matrices and vectors. However, there are some common requirements a matrix or vector type must fulfill to qualify as an acceptable type for the solvers in this hierarchy. These requirements are listed below.

The classes we show below are not any concrete class. Rather, they are intended to form a "signature" which a concrete class has to conform to. Note that the matrix and vector classes within this library of course conform to this interface; therefore, SparseMatrix and Vector are good examples for these classes as they provide the necessary signatures of member functions (although they also provide many more interfaces that solvers do not in fact need – for example, element access). In addition, you may want to take a look at step-20, step-22, or a number of classes in the LinearSolvers namespace for examples of how one can define matrix-like classes that can serve as linear operators for linear solvers.

Concretely, matrix and vector classes that can be passed to a linear solver need to provide the following interfaces:

class Matrix
{
public:
// Application of matrix to vector src. Write result into dst.
void vmult (VectorType &dst,
const VectorType &src) const;
// Application of transpose to a vector. This function is,
// however, only used by some iterative methods.
void Tvmult (VectorType &dst,
const VectorType &src) const;
};
class Vector
{
public:
// Define value type of the entries
using value_type = double;
// Resize the current object to have the same size and layout as
// the model_vector argument provided. The second argument
// indicates whether to clear the current object after resizing.
// The second argument must have a default value equal to false.
void reinit (const Vector &model_vector,
const bool leave_elements_uninitialized = false);
// Inner product between the current object and the argument.
double operator * (const Vector &v) const;
// Set all the vector entries to a constant scalar.
Vector & operator = (const double a);
// Deep copy of the vector.
// Important if Vector contains pointers to data to duplicate data.
Vector & operator = (const Vector &x);
// Addition of vectors
void add (const Vector &x);
// Scaled addition of vectors
void add (const double a,
const Vector &x);
// Scaled addition of vectors
void sadd (const double a,
const double b,
const Vector &x);
// Scaled assignment of a vector
void equ (const double a,
const Vector &x);
// Combined scaled addition of vector x into the current object and
// subsequent inner product of the current object with v.
double add_and_dot (const double a,
const Vector &x,
const Vector &v);
// Multiply the elements of the current object by a fixed value.
Vector & operator *= (const double a);
// Return the l2 norm of the vector.
double l2_norm () const;
};
std::enable_if_t< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typename ProductType< std::complex< T >, std::complex< U > >::type > operator*(const std::complex< T > &left, const std::complex< U > &right)

In addition, for some solvers there has to be a global function swap(VectorType &a, VectorType &b) that exchanges the values of the two vectors.

Finally, the solvers also expect an instantiation of GrowingVectorMemory<VectorType>. These instantiations are provided by the deal.II library for the built-in vector types, but must be explicitly added for user-provided vector classes. Otherwise, the linker will complain that it cannot find the constructors and destructors of GrowingVectorMemory that happen in the SolverBase class.

// Definition and implementation of vector class
class UserVector { ... };
// Create explicit instantiation for the vector class. If your project
// consists of multiple files, including header files, this instantiation
// must be put in a <code>.cc</code> file in order to instantiate only
// once.
#include <deal.II/lac/vector_memory.templates.h>
template class VectorMemory<UserVector>;

The preconditioners used must have the same interface as matrices, i.e. in particular they have to provide a member function vmult which denotes the application of the preconditioner.

AdditionalData

Several solvers need additional data, like the damping parameter omega of the SolverRichardson class or the maximum number of temporary vectors of SolverGMRES. To have a standardized way of constructing solvers, each solver class has a struct AdditionalData as a member, and constructors of all solver classes take such an argument. Some solvers need no additional data, or may not at the current time. For these solvers the struct AdditionalData is empty and calling the constructor may be done without giving the additional structure as an argument as a default AdditionalData is set by default.

With this, creating a solver looks like one of the following blocks:

// GMRES with restart every 50 iterations
SolverGMRES solver_gmres (solver_control, vector_memory,
// Richardson with omega=0.8
SolverRichardson solver_richardson (solver_control, vector_memory,
// CG with default AdditionalData
SolverCG solver_cg (solver_control, vector_memory);

Using a unified constructor parameter list for all solvers supports the SolverSelector class; the unified interface enables us to use this class unchanged even if the number of types of parameters to a certain solver changes and it is still possible in a simple way to give these additional data to the SolverSelector object for each solver which it may use.

Observing the progress of linear solver iterations

The SolverBase class, being the base class for all of the iterative solvers such as SolverCG, SolverGMRES, etc, provides the facilities by which actual solver implementations determine whether the iteration is converged, not yet converged, or has failed. Typically, this is done using an object of type SolverControl that is passed to the solver classes's constructors and from them down to the constructor of this base class. Every one of the tutorial programs that solves a linear problem (starting with step-3) uses this method and it is described in detail there. However, the underlying mechanism is more general and allows for many other uses to observe how the linear solver iterations progress.

The basic approach is that the iterative solvers invoke a signal at the end of each iteration to determine whether the solution is converged. A signal is a class that has, conceptually, a list of pointers to functions and every time the signal is invoked, each of these functions are called. In the language of signals, the functions called are called slots and one can attach any number of slots to a signal. (The implementation of signals and slots we use here is the one from the BOOST.signals2 library.) A number of details may clarify what is happening underneath:

Given all this, it should now be obvious how the SolverControl object fits into this scheme: when a SolverControl object is passed to the constructor of the current class, we simply connect the SolverControl::check() function of that object as a slot to the signal we maintain here. In other words, since a SolverBase object is always constructed using a SolverControl object, there is always at least one slot associated with the signal, namely the one that determines convergence.

On the other hand, using the connect() member function, it is possible to connect any number of other slots to the signal to observe whatever it is you want to observe. The connect() function also returns an object that describes the connection from the signal to the slot, and the corresponding BOOST functions then allow you to disconnect the slot if you want.

An example may illuminate these issues. In the step-3 tutorial program, let us add a member function as follows to the main class:

Step3::write_intermediate_solution (
const unsigned int iteration,
const double , //check_value
const Vector<double> &current_iterate) const
{
DataOut<2> data_out;
data_out.attach_dof_handler (dof_handler);
data_out.add_data_vector (current_iterate, "solution");
data_out.build_patches ();
std::ofstream output ("solution-"
+ Utilities::int_to_string(iteration,4)
+ ".vtu");
data_out.write_vtu (output);
}
void write_vtu(std::ostream &out) const
void attach_dof_handler(const DoFHandler< dim, spacedim > &)
void add_data_vector(const VectorType &data, const std::vector< std::string > &names, const DataVectorType type=type_automatic, const std::vector< DataComponentInterpretation::DataComponentInterpretation > &data_component_interpretation={})
virtual void build_patches(const unsigned int n_subdivisions=0)
Definition data_out.cc:1063
@ success
Stop iteration, goal reached.
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:471

The function satisfies the signature necessary to be a slot for the signal discussed above, with the exception that it is a member function and consequently requires a this pointer. What the function does is to take the vector given as last argument and write it into a file in VTU format with a file name derived from the number of the iteration.

This function can then be hooked into the CG solver by modifying the Step3::solve() function as follows:

void Step3::solve ()
{
SolverControl solver_control (1000, 1e-12);
SolverCG<> solver (solver_control);
solver.connect ([this](const unsigned int iteration,
const double check_value,
const Vector<double> *current_iterate){
this->write_intermediate_solution(
iteration, check_value, current_iterate);
});
solver.solve (system_matrix, solution, system_rhs,
}

The use of a lambda function here ensures that we convert the member function with its three arguments plus the this pointer, to a function that only takes three arguments, by fixing the implicit this argument of the function to the this pointer in the current function.

It is well understood that the CG method is a smoothing iteration (in the same way as the more commonly used Jacobi or SSOR iterations are smoothers). The code above therefore allows to observe how the solution becomes smoother and smoother in every iteration. This is best observed by initializing the solution vector with randomly distributed numbers in \([-1,1]\), using code such as

for (unsigned int i=0; i<solution.size(); ++i)
solution(i) = 2.*rand()/RAND_MAX-1;

Using this, the slot will then generate files that when visualized look like this over the course of iterations zero to five:

Definition at line 340 of file solver.h.

Member Typedef Documentation

◆ vector_type

template<class VectorType = Vector<double>>
using SolverBase< VectorType >::vector_type = VectorType

An alias for the underlying vector type

Definition at line 346 of file solver.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

◆ SolverBase() [1/2]

template<class VectorType = Vector<double>>
SolverBase< VectorType >::SolverBase ( SolverControl solver_control,
VectorMemory< VectorType > &  vector_memory 
)

Constructor. Takes a control object which evaluates the conditions for convergence, and an object that allows solvers to allocate memory for temporary objects.

Of both objects, a reference is stored, so it is the user's responsibility to guarantee that the lifetime of the two arguments is at least as long as that of the solver object.

◆ SolverBase() [2/2]

template<class VectorType = Vector<double>>
SolverBase< VectorType >::SolverBase ( SolverControl solver_control)

Constructor. Takes a control object which evaluates the conditions for convergence. In contrast to the other constructor, this constructor designates an internal object of type GrowingVectorMemory to allocate memory.

A reference to the control object is stored, so it is the user's responsibility to guarantee that the lifetime of the argument is at least as long as that of the solver object.

Member Function Documentation

◆ connect()

template<class VectorType = Vector<double>>
boost::signals2::connection SolverBase< VectorType >::connect ( const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType &current_iterate)> &  slot)

Connect a function object that will be called periodically within iterative solvers. This function is used to attach monitors to iterative solvers, either to determine when convergence has happened, or simply to observe the progress of an iteration. See the documentation of this class for more information.

Parameters
slotA function object specified here will, with each call, receive the number of the current iteration, the value that is used to check for convergence (typically the residual of the current iterate with respect to the linear system to be solved) and the currently best available guess for the current iterate. Note that some solvers do not update the approximate solution in every iteration but only after convergence or failure has been determined (GMRES is an example); in such cases, the vector passed as the last argument to the signal is simply the best approximate at the time the signal is called, but not the vector that will be returned if the signal's return value indicates that the iteration should be terminated. The function object must return a SolverControl::State value that indicates whether the iteration should continue, has failed, or has succeeded. The results of all connected functions will then be combined to determine what should happen with the iteration.
Returns
A connection object that represents the connection from the signal to the function object. It can be used to disconnect the function object again from the signal. See the documentation of the BOOST Signals2 library for more information on connection management.

◆ 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

◆ static_vector_memory

template<class VectorType = Vector<double>>
GrowingVectorMemory<VectorType> SolverBase< VectorType >::static_vector_memory
mutableprotected

A static vector memory object to be used whenever no such object has been given to the constructor.

Definition at line 414 of file solver.h.

◆ memory

template<class VectorType = Vector<double>>
VectorMemory<VectorType>& SolverBase< VectorType >::memory
protected

A reference to an object that provides memory for auxiliary vectors.

Definition at line 419 of file solver.h.

◆ iteration_status

template<class VectorType = Vector<double>>
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType & current_iterate), StateCombiner> SolverBase< VectorType >::iteration_status
protected

A signal that iterative solvers can execute at the end of every iteration (or in an otherwise periodic fashion) to find out whether we should continue iterating or not. The signal may call one or more slots that each will make this determination by themselves, and the result over all slots (function calls) will be determined by the StateCombiner object.

The arguments passed to the signal are (i) the number of the current iteration; (ii) the value that is used to determine convergence (oftentimes the residual, but in other cases other quantities may be used as long as they converge to zero as the iterate approaches the solution of the linear system); and (iii) a vector that corresponds to the current best guess for the solution at the point where the signal is called. Note that some solvers do not update the approximate solution in every iteration but only after convergence or failure has been determined (GMRES is an example); in such cases, the vector passed as the last argument to the signal is simply the best approximate at the time the signal is called, but not the vector that will be returned if the signal's return value indicates that the iteration should be terminated.

Definition at line 470 of file solver.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: