Reference documentation for deal.II version GIT relicensing-468-ge3ce94fd06 2024-04-23 15:40:02+00:00
\(\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
Namespaces | Classes | Typedefs | Functions | Variables

Namespaces

namespace  distributed
 
namespace  fullydistributed
 
namespace  internal
 
namespace  shared
 

Classes

class  CellWeights
 
class  DistributedTriangulationBase
 
struct  ParallelForInteger
 
class  TriangulationBase
 

Typedefs

using SyncIterators = SynchronousIterators< Iterators >
 

Functions

template<typename InputIterator , typename OutputIterator , typename Function >
 DEAL_II_CXX20_REQUIRES ((std::invocable< Function, decltype(*std::declval< InputIterator >())> &&std::assignable_from< decltype(*std::declval< OutputIterator >()), std::invoke_result_t< Function, decltype(*std::declval< InputIterator >())> >)) void transform(const InputIterator &begin_in
 
Iterators x_begin (begin_in, out)
 
Iterators x_end (end_in, OutputIterator())
 
template<typename InputIterator1 , typename InputIterator2 , typename OutputIterator , typename Function >
 DEAL_II_CXX20_REQUIRES ((std::invocable< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >())> &&std::assignable_from< decltype(*std::declval< OutputIterator >()), std::invoke_result_t< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >())> >)) void transform(const InputIterator1 &begin_in1
 
Iterators x_begin (begin_in1, in2, out)
 
Iterators x_end (end_in1, InputIterator2(), OutputIterator())
 
template<typename InputIterator1 , typename InputIterator2 , typename InputIterator3 , typename OutputIterator , typename Function >
 DEAL_II_CXX20_REQUIRES ((std::invocable< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >()), decltype(*std::declval< InputIterator3 >())> &&std::assignable_from< decltype(*std::declval< OutputIterator >()), std::invoke_result_t< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >()), decltype(*std::declval< InputIterator3 >())> >)) void transform(const InputIterator1 &begin_in1
 
Iterators x_begin (begin_in1, in2, in3, out)
 
Iterators x_end (end_in1, InputIterator2(), InputIterator3(), OutputIterator())
 
template<typename Iterator , typename Function >
void apply_to_subranges (const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, const Function &f, const unsigned int grainsize)
 
template<typename ResultType , typename Iterator , typename Function >
 DEAL_II_CXX20_REQUIRES ((std::invocable< Function, Iterator, Iterator > &&std::convertible_to< std::invoke_result_t< Function, Iterator, Iterator >, ResultType >)) ResultType accumulate_from_subranges(const Function &f
 

Variables

const InputIterator & end_in
 
const InputIterator OutputIterator out
 
const InputIterator OutputIterator const Functionfunction
 
const InputIterator OutputIterator const Function const unsigned int grainsize
 
const InputIterator1 & end_in1
 
const InputIterator1 InputIterator2 in2
 
const InputIterator1 InputIterator2 InputIterator3 in3
 
const Iterator & begin
 
const Iterator const std_cxx20::type_identity_t< Iterator > & end
 

Detailed Description

A namespace in which we define classes and algorithms that deal with running in parallel on shared memory machines when deal.II is configured to use multiple threads (see Parallel computing with multiple processors accessing shared memory), as well as running things in parallel on distributed memory machines (see Parallel computing with multiple processors using distributed memory).

Typedef Documentation

◆ SyncIterators

Definition at line 180 of file parallel.h.

Function Documentation

◆ DEAL_II_CXX20_REQUIRES() [1/4]

template<typename InputIterator , typename OutputIterator , typename Function >
parallel::DEAL_II_CXX20_REQUIRES ( (std::invocable< Function, decltype(*std::declval< InputIterator >())> && std::assignable_from< decltype(*std::declval< OutputIterator >()), std::invoke_result_t< Function, decltype(*std::declval< InputIterator >())> >)  ) const &

An algorithm that performs the action *out++ = function(*in++) where the in iterator ranges over the given input range.

This algorithm does pretty much what std::transform does. The difference is that the function can run in parallel when deal.II is configured to use multiple threads.

If running in parallel, the iterator range is split into several chunks that are each packaged up as a task and given to the Threading Building Blocks scheduler to work on as compute resources are available. The function returns once all chunks have been worked on. The last argument denotes the minimum number of elements of the iterator range per task; the number must be large enough to amortize the startup cost of new tasks, and small enough to ensure that tasks can be reasonably load balanced.

For a discussion of the kind of problems to which this function is applicable, see the Parallel computing with multiple processors module.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::invocable<Function,
decltype(*std::declval<InputIterator>())> &&
std::assignable_from<decltype(*std::declval<OutputIterator>()),
std::invoke_result_t<Function,
decltype(*std::declval<InputIterator>())>>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ x_begin() [1/3]

Iterators parallel::x_begin ( begin_in  ,
out   
)

◆ x_end() [1/3]

Iterators parallel::x_end ( end_in  ,
OutputIterator()   
)

◆ DEAL_II_CXX20_REQUIRES() [2/4]

template<typename InputIterator1 , typename InputIterator2 , typename OutputIterator , typename Function >
parallel::DEAL_II_CXX20_REQUIRES ( (std::invocable< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >())> && std::assignable_from< decltype(*std::declval< OutputIterator >()), std::invoke_result_t< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >())> >)  ) const &

An algorithm that performs the action *out++ = function(*in1++, *in2++) where the in1 iterator ranges over the given input range, using the parallel for operator of tbb.

This algorithm does pretty much what std::transform does. The difference is that the function can run in parallel when deal.II is configured to use multiple threads.

If running in parallel, the iterator range is split into several chunks that are each packaged up as a task and given to the Threading Building Blocks scheduler to work on as compute resources are available. The function returns once all chunks have been worked on. The last argument denotes the minimum number of elements of the iterator range per task; the number must be large enough to amortize the startup cost of new tasks, and small enough to ensure that tasks can be reasonably load balanced.

For a discussion of the kind of problems to which this function is applicable, see the Parallel computing with multiple processors module.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::invocable<Function,
decltype(*std::declval<InputIterator1>()),
decltype(*std::declval<InputIterator2>())> &&
std::assignable_from<decltype(*std::declval<OutputIterator>()),
std::invoke_result_t<Function,
decltype(*std::declval<InputIterator1>()),
decltype(*std::declval<InputIterator2>())>>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ x_begin() [2/3]

Iterators parallel::x_begin ( begin_in1  ,
in2  ,
out   
)

◆ x_end() [2/3]

Iterators parallel::x_end ( end_in1  ,
InputIterator2()  ,
OutputIterator()   
)

◆ DEAL_II_CXX20_REQUIRES() [3/4]

template<typename InputIterator1 , typename InputIterator2 , typename InputIterator3 , typename OutputIterator , typename Function >
parallel::DEAL_II_CXX20_REQUIRES ( (std::invocable< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >()), decltype(*std::declval< InputIterator3 >())> && std::assignable_from< decltype(*std::declval< OutputIterator >()), std::invoke_result_t< Function, decltype(*std::declval< InputIterator1 >()), decltype(*std::declval< InputIterator2 >()), decltype(*std::declval< InputIterator3 >())> >)  ) const &

An algorithm that performs the action *out++ = function(*in1++, *in2++, *in3++) where the in1 iterator ranges over the given input range.

This algorithm does pretty much what std::transform does. The difference is that the function can run in parallel when deal.II is configured to use multiple threads.

If running in parallel, the iterator range is split into several chunks that are each packaged up as a task and given to the Threading Building Blocks scheduler to work on as compute resources are available. The function returns once all chunks have been worked on. The last argument denotes the minimum number of elements of the iterator range per task; the number must be large enough to amortize the startup cost of new tasks, and small enough to ensure that tasks can be reasonably load balanced.

For a discussion of the kind of problems to which this function is applicable, see the Parallel computing with multiple processors module.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::invocable<Function,
decltype(*std::declval<InputIterator1>()),
decltype(*std::declval<InputIterator2>()),
decltype(*std::declval<InputIterator3>())> &&
std::assignable_from<decltype(*std::declval<OutputIterator>()),
std::invoke_result_t<Function,
decltype(*std::declval<InputIterator1>()),
decltype(*std::declval<InputIterator2>()),
decltype(*std::declval<InputIterator3>())>>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

◆ x_begin() [3/3]

Iterators parallel::x_begin ( begin_in1  ,
in2  ,
in3  ,
out   
)

◆ x_end() [3/3]

Iterators parallel::x_end ( end_in1  ,
InputIterator2()  ,
InputIterator3()  ,
OutputIterator()   
)

◆ apply_to_subranges()

template<typename Iterator , typename Function >
void parallel::apply_to_subranges ( const Iterator &  begin,
const std_cxx20::type_identity_t< Iterator > &  end,
const Function f,
const unsigned int  grainsize 
)

This function applies the given function argument f to all elements in the range [begin,end) and may do so in parallel. An example of its use is given in step-69.

However, in many cases it is not efficient to call a function on each element, so this function calls the given function object on sub-ranges. In other words: if the given range [begin,end) is smaller than grainsize or if multithreading is not enabled, then we call f(begin,end); otherwise, we may execute, possibly in parallel, a sequence of calls f(b,e) where [b,e) are subintervals of [begin,end) and the collection of calls we do to f(.,.) will happen on disjoint subintervals that collectively cover the original interval [begin,end).

Oftentimes, the called function will of course have to get additional information, such as the object to work on for a given value of the iterator argument. This can be achieved by binding certain arguments. For example, here is an implementation of a matrix-vector multiplication \(y=Ax\) for a full matrix \(A\) and vectors \(x,y\):

void matrix_vector_product (const FullMatrix &A,
const Vector &x,
Vector &y)
{
(0, A.n_rows(),
[&](const unsigned int begin_row,
const unsigned int end_row)
{
mat_vec_on_subranges(begin_row, end_row, A, x, y);
},
50);
}
void mat_vec_on_subranges (const unsigned int begin_row,
const unsigned int end_row,
const FullMatrix &A,
const Vector &x,
Vector &y)
{
for (unsigned int row=begin_row; row!=end_row; ++row)
for (unsigned int col=0; col<x.size(); ++col)
y(row) += A(row,col) * x(col);
}
virtual size_type size() const override
void apply_to_subranges(const Iterator &begin, const std_cxx20::type_identity_t< Iterator > &end, const Function &f, const unsigned int grainsize)
Definition parallel.h:452

Note how we use the lambda function to convert mat_vec_on_subranges from a function that takes 5 arguments to one taking 2 by binding the remaining arguments. The resulting function object requires only two arguments, begin_row and end_row, with all other arguments fixed.

The code, if in single-thread mode, will call mat_vec_on_subranges on the entire range [0,n_rows) exactly once. In multi-threaded mode, however, it may be called multiple times on subranges of this interval, possibly allowing more than one CPU core to take care of part of the work.

The grainsize argument (50 in the example above) makes sure that subranges do not become too small, to avoid spending more time on scheduling subranges to CPU resources than on doing actual work.

For a discussion of the kind of problems to which this function is applicable, see also the Parallel computing with multiple processors module.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::invocable<Function, Iterator, Iterator>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 452 of file parallel.h.

◆ DEAL_II_CXX20_REQUIRES() [4/4]

template<typename ResultType , typename Iterator , typename Function >
parallel::DEAL_II_CXX20_REQUIRES ( (std::invocable< Function, Iterator, Iterator > && std::convertible_to< std::invoke_result_t< Function, Iterator, Iterator >, ResultType >)  ) const &

This function works a lot like the apply_to_subranges() function, but it allows to accumulate numerical results computed on each subrange into one number. The type of this number is given by the ResultType template argument that needs to be explicitly specified, and results are added up (i.e., the reduction of results from subranges happens by adding up these results).

An example of use of this function is to compute the value of the expression \(x^T A x\) for a square matrix \(A\) and a vector \(x\). The sum over rows can be parallelized and the whole code might look like this:

void matrix_norm (const FullMatrix &A,
const Vector &x)
{
return
(parallel::accumulate_from_subranges<double>
([&](const unsigned int begin_row,
const unsigned int end_row)
{
mat_norm_sqr_on_subranges(begin_row, end_row, A, x);
},
0, A.n_rows(),
50);
}
double
mat_norm_sqr_on_subranges (const unsigned int begin_row,
const unsigned int end_row,
const FullMatrix &A,
const Vector &x)
{
double norm_sqr = 0;
for (unsigned int row=begin_row; row!=end_row; ++row)
for (unsigned int col=0; col<x.size(); ++col)
norm_sqr += x(row) * A(row,col) * x(col);
return norm_sqr;
}
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)

Here, mat_norm_sqr_on_subranges is called on the entire range [0,A.n_rows()) if this range is less than the minimum grainsize (above chosen as 50) or if deal.II is configured to not use multithreading. Otherwise, it may be called on subsets of the given range, with results from the individual subranges accumulated internally.

Warning
If ResultType is a floating point type, then accumulation is not an associative operation. In other words, if the given function object is called three times on three subranges, returning values \(a,b,c\), then the returned result of this function is \((a+b)+c\). However, depending on how the three sub-tasks are distributed on available CPU resources, the result may also be \((a+c)+b\) or any other permutation; because floating point addition is not associative (as opposed, of course, to addition of real numbers), the result of invoking this function several times may differ on the order of round-off.

For a discussion of the kind of problems to which this function is applicable, see also the Parallel computing with multiple processors module.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true:
(std::invocable<Function, Iterator, Iterator> &&
std::convertible_to<std::invoke_result_t<Function, Iterator, Iterator>,
ResultType>)
If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Variable Documentation

◆ end_in

const InputIterator& parallel::end_in

Definition at line 166 of file parallel.h.

◆ out

const InputIterator1 InputIterator2 InputIterator3 OutputIterator parallel::out

Definition at line 167 of file parallel.h.

◆ function

const InputIterator1 InputIterator2 InputIterator3 OutputIterator const Function & parallel::function

Definition at line 168 of file parallel.h.

◆ grainsize

const Iterator const std_cxx20::type_identity_t< Iterator > const unsigned int parallel::grainsize
Initial value:
{
using Iterators = std::tuple<InputIterator, OutputIterator>

Definition at line 169 of file parallel.h.

◆ end_in1

const InputIterator1 & parallel::end_in1

Definition at line 241 of file parallel.h.

◆ in2

const InputIterator1 InputIterator2 parallel::in2

Definition at line 242 of file parallel.h.

◆ in3

const InputIterator1 InputIterator2 InputIterator3 parallel::in3

Definition at line 325 of file parallel.h.

◆ begin

const Iterator& parallel::begin

Definition at line 609 of file parallel.h.

◆ end

const Iterator const std_cxx20::type_identity_t<Iterator>& parallel::end

Definition at line 610 of file parallel.h.