Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/vectorization.h>
Public Member Functions | |
VectorizedArray & | operator= (const Number scalar) |
Number & | operator[] (const unsigned int comp) |
const Number & | operator[] (const unsigned int comp) const |
VectorizedArray & | operator+= (const VectorizedArray< Number > &vec) |
VectorizedArray & | operator-= (const VectorizedArray< Number > &vec) |
VectorizedArray & | operator*= (const VectorizedArray< Number > &vec) |
VectorizedArray & | operator/= (const VectorizedArray< Number > &vec) |
void | load (const Number *ptr) |
void | store (Number *ptr) const |
void | streaming_store (Number *ptr) const |
void | gather (const Number *base_ptr, const unsigned int *offsets) |
void | scatter (const unsigned int *offsets, Number *base_ptr) const |
Public Attributes | |
Number | data |
Static Public Attributes | |
static const unsigned int | n_array_elements = 1 |
Private Member Functions | |
VectorizedArray | get_sqrt () const |
VectorizedArray | get_abs () const |
VectorizedArray | get_max (const VectorizedArray &other) const |
VectorizedArray | get_min (const VectorizedArray &other) const |
Friends | |
template<typename Number2 > | |
VectorizedArray< Number2 > | std::sqrt (const VectorizedArray< Number2 > &) |
Related Functions | |
(Note that these are not member functions.) | |
template<typename Number > | |
VectorizedArray< Number > | make_vectorized_array (const Number &u) |
template<typename Number > | |
void | vectorized_load_and_transpose (const unsigned int n_entries, const Number *in, const unsigned int *offsets, VectorizedArray< Number > *out) |
template<typename Number > | |
void | vectorized_transpose_and_store (const bool add_into, const unsigned int n_entries, const VectorizedArray< Number > *in, const unsigned int *offsets, Number *out) |
template<typename Number > | |
bool | operator== (const VectorizedArray< Number > &lhs, const VectorizedArray< Number > &rhs) |
template<typename Number > | |
VectorizedArray< Number > | operator+ (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator- (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator* (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator/ (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator+ (const Number &u, const VectorizedArray< Number > &v) |
VectorizedArray< float > | operator+ (const double &u, const VectorizedArray< float > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator+ (const VectorizedArray< Number > &v, const Number &u) |
VectorizedArray< float > | operator+ (const VectorizedArray< float > &v, const double &u) |
template<typename Number > | |
VectorizedArray< Number > | operator- (const Number &u, const VectorizedArray< Number > &v) |
VectorizedArray< float > | operator- (const double &u, const VectorizedArray< float > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator- (const VectorizedArray< Number > &v, const Number &u) |
VectorizedArray< float > | operator- (const VectorizedArray< float > &v, const double &u) |
template<typename Number > | |
VectorizedArray< Number > | operator* (const Number &u, const VectorizedArray< Number > &v) |
VectorizedArray< float > | operator* (const double &u, const VectorizedArray< float > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator* (const VectorizedArray< Number > &v, const Number &u) |
VectorizedArray< float > | operator* (const VectorizedArray< float > &v, const double &u) |
template<typename Number > | |
VectorizedArray< Number > | operator/ (const Number &u, const VectorizedArray< Number > &v) |
VectorizedArray< float > | operator/ (const double &u, const VectorizedArray< float > &v) |
template<typename Number > | |
VectorizedArray< Number > | operator/ (const VectorizedArray< Number > &v, const Number &u) |
VectorizedArray< float > | operator/ (const VectorizedArray< float > &v, const double &u) |
template<typename Number > | |
VectorizedArray< Number > | operator+ (const VectorizedArray< Number > &u) |
template<typename Number > | |
VectorizedArray< Number > | operator- (const VectorizedArray< Number > &u) |
template<typename Number > | |
VectorizedArray< Number > | sin (const ::VectorizedArray< Number > &x) |
template<typename Number > | |
VectorizedArray< Number > | cos (const ::VectorizedArray< Number > &x) |
template<typename Number > | |
VectorizedArray< Number > | tan (const ::VectorizedArray< Number > &x) |
template<typename Number > | |
VectorizedArray< Number > | exp (const ::VectorizedArray< Number > &x) |
template<typename Number > | |
VectorizedArray< Number > | log (const ::VectorizedArray< Number > &x) |
template<typename Number > | |
VectorizedArray< Number > | sqrt (const ::VectorizedArray< Number > &x) |
template<typename Number > | |
VectorizedArray< Number > | pow (const ::VectorizedArray< Number > &x, const Number p) |
template<typename Number > | |
VectorizedArray< Number > | abs (const ::VectorizedArray< Number > &x) |
template<typename Number > | |
VectorizedArray< Number > | max (const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y) |
template<typename Number > | |
VectorizedArray< Number > | min (const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y) |
This generic class defines a unified interface to a vectorized data type. For general template arguments, this class simply corresponds to the template argument. For example, VectorizedArray<long double> is nothing else but a wrapper around long double
with exactly one data field of type long double
and overloaded arithmetic operations. This means that VectorizedArray<ComplicatedType>
has a similar layout as ComplicatedType, provided that ComplicatedType defines basic arithmetic operations. For floats and doubles, an array of numbers are packed together, though. The number of elements packed together depend on the computer system and compiler flags that are used for compilation of deal.II. The fundamental idea of these packed data types is to use one single CPU instruction to perform arithmetic operations on the whole array using the processor's vector units. Most computer systems by 2010 standards will use an array of two doubles and four floats, respectively (this corresponds to the SSE/SSE2 data sets) when compiling deal.II on 64-bit operating systems. On Intel Sandy Bridge processors and newer or AMD Bulldozer processors and newer, four doubles and eight floats are used when deal.II is configured e.g. using gcc with –with-cpu=native or –with- cpu=corei7-avx. On compilations with AVX-512 support, eight doubles and sixteen floats are used.
This behavior of this class is made similar to the basic data types double and float. The definition of a vectorized array does not initialize the data field but rather leaves it undefined, as is the case for double and float. However, when calling something like VectorizedArray<double> a = VectorizedArray<double>(), it sets all numbers in this field to zero. In other words, this class is a plain old data (POD) type which has an equivalent C representation and can e.g. be safely copied with std::memcpy. This POD layout is also necessary for ensuring correct alignment of data with address boundaries when collected in a vector (i.e., when the first element in a vector is properly aligned, all subsequent elements will be correctly aligned, too).
Note that for proper functioning of this class, certain data alignment rules must be respected. This is because the computer expects the starting address of a VectorizedArray<double> field at specific addresses in memory (usually, the address of the vectorized array should be a multiple of the length of the array in bytes). Otherwise, a segmentation fault or a severe loss of performance might occur. When creating a single data field on the stack like VectorizedArray<double> a = VectorizedArray<double>()
, the compiler will take care of data alignment automatically. However, when allocating a long vector of VectorizedArray<double> data, one needs to respect these rules. Use the class AlignedVector or data containers based on AlignedVector (such as Table) for this purpose. It is a class very similar to std::vector otherwise but always makes sure that data is correctly aligned.
Definition at line 35 of file memory_consumption.h.
|
inline |
This function assigns a scalar to this class.
Definition at line 167 of file vectorization.h.
|
inline |
Access operator (only valid with component 0)
Definition at line 178 of file vectorization.h.
|
inline |
Constant access operator (only valid with component 0)
Definition at line 190 of file vectorization.h.
|
inline |
Addition
Definition at line 202 of file vectorization.h.
|
inline |
Subtraction
Definition at line 213 of file vectorization.h.
|
inline |
Multiplication
Definition at line 224 of file vectorization.h.
|
inline |
Division
Definition at line 235 of file vectorization.h.
|
inline |
Load n_array_elements
from memory into the calling class, starting at the given address. The memory need not be aligned by the amount of bytes in the vectorized array, as opposed to casting a double address to VectorizedArray<double>*.
Definition at line 248 of file vectorization.h.
|
inline |
Write the content of the calling class into memory in form of n_array_elements
to the given address. The memory need not be aligned by the amount of bytes in the vectorized array, as opposed to casting a double address to VectorizedArray<double>*.
Definition at line 260 of file vectorization.h.
|
inline |
Write the content of the calling class into memory in form of n_array_elements
to the given address using non-temporal stores that bypass the processor's caches, using _mm_stream_pd
store intrinsics on supported CPUs. The destination of the store ptr
must be aligned by the amount of bytes in the vectorized array.
This store operation can be faster than usual store operations in case the store is streaming because it avoids the read-for-ownership transfer typically invoked in standard stores. This approximately works as follows (see the literature on computer architecture for details): When an algorithm stores some results to a memory address, a processor typically wants to move it into some of its caches as it expects the data to be re-used again at some point. Since caches are organized in lines of sizes either 64 byte or 128 byte but writes are usually smaller, a processor must first load in the destination cache line upon a write because only part of the cache line is overwritten initially. If a series of stores write data in a chunk bigger than any of its caches could handle, the data finally has to be moved out from the caches to main memory. But since all addressed have first been read, this doubles the load on main memory, which can incur a performance penalty. Furthermore, the organization of caches in a multicore context also requires reading an address before something can be written to cache to that address, see e.g. the Wikipedia article on the MESI protocol for details. The instruction underlying this function call signals to the processor that these two prerequisites on a store are relaxed: Firstly, one expects the whole cache line to be overwritten (that the memory subsystem then handles appropriately), so no need to first read the "remainder" of the cache line. Secondly, the data behind that particular memory will not be subject to cache coherency protocol as it will be in main memory both when the same processor wants to access it again as well as any other processors in a multicore chip. Due to this particular setup, any subsequent access to the data written by this function will need to query main memory, which is slower than an access from a cache both latency-wise and throughput-wise. Thus, this command should only be used for large stores that will collectively not fit into caches, as performance will be degraded otherwise. For a typical use case, see also this blog article.
Note that streaming stores are only available in the specialized SSE/AVX classes of VectorizedArray of type double
or float
, not in the generic base class.
Definition at line 310 of file vectorization.h.
|
inline |
Load n_array_elements
from memory into the calling class, starting at the given address and with given offsets, each entry from the offset providing one element of the vectorized array.
This operation corresponds to the following code (but uses a more efficient implementation in case the hardware allows for that):
Definition at line 328 of file vectorization.h.
|
inline |
Write the content of the calling class into memory in form of n_array_elements
to the given address and the given offsets, filling the elements of the vectorized array into each offset.
This operation corresponds to the following code (but uses a more efficient implementation in case the hardware allows for that):
Definition at line 347 of file vectorization.h.
|
inlineprivate |
Return the square root of this field. Not for use in user code. Use sqrt(x) instead.
Definition at line 366 of file vectorization.h.
|
inlineprivate |
Return the absolute value of this field. Not for use in user code. Use abs(x) instead.
Definition at line 379 of file vectorization.h.
|
inlineprivate |
Return the component-wise maximum of this field and another one. Not for use in user code. Use max(x,y) instead.
Definition at line 392 of file vectorization.h.
|
inlineprivate |
Return the component-wise minimum of this field and another one. Not for use in user code. Use min(x,y) instead.
Definition at line 405 of file vectorization.h.
|
friend |
Make a few functions friends.
|
related |
Create a vectorized array that sets all entries in the array to the given scalar.
Definition at line 436 of file vectorization.h.
|
related |
This method loads VectorizedArray::n_array_elements data streams from the given array in
. The offsets to the input array are given by the array offsets
. From each stream, n_entries are read. The data is then transposed and stored it into an array of VectorizedArray type. The output array out
is expected to be an array of size n_entries
. This method operates on plain arrays, so no checks for valid data access are made. It is the user's responsibility to ensure that the given arrays are valid according to the access layout below.
This operation corresponds to a transformation of an array-of-struct (input) into a struct-of-array (output) according to the following formula:
A more optimized version of this code will be used for supported types.
This is the inverse operation to vectorized_transpose_and_store().
Definition at line 473 of file vectorization.h.
|
related |
This method stores the vectorized arrays in transposed form into the given output array out
with the given offsets offsets
. This operation corresponds to a transformation of a struct-of-array (input) into an array- of-struct (output). This method operates on plain array, so no checks for valid data access are made. It is the user's responsibility to ensure that the given arrays are valid according to the access layout below.
This method assumes that the specified offsets do not overlap. Otherwise, the behavior is undefined in the vectorized case. It is the user's responsibility to make sure that the access does not overlap and avoid undefined behavior.
The argument add_into
selects where the entries should only be written into the output arrays or the result should be added into the existing entries in the output. For add_into == false
, the following code is assumed:
For add_into == true
, the code implements the following action:
A more optimized version of this code will be used for supported types.
This is the inverse operation to vectorized_load_and_transpose().
Definition at line 526 of file vectorization.h.
|
related |
Relational operator == for VectorizedArray
Definition at line 2829 of file vectorization.h.
|
related |
Addition of two vectorized arrays with operator +.
Definition at line 2848 of file vectorization.h.
|
related |
Subtraction of two vectorized arrays with operator -.
Definition at line 2863 of file vectorization.h.
|
related |
Multiplication of two vectorized arrays with operator *.
Definition at line 2878 of file vectorization.h.
|
related |
Division of two vectorized arrays with operator /.
Definition at line 2893 of file vectorization.h.
|
related |
Addition of a scalar (expanded to a vectorized array with n_array_elements
equal entries) and a vectorized array.
Definition at line 2909 of file vectorization.h.
|
related |
Addition of a scalar (expanded to a vectorized array with n_array_elements
equal entries) and a vectorized array in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 2927 of file vectorization.h.
|
related |
Addition of a vectorized array and a scalar (expanded to a vectorized array with n_array_elements
equal entries).
Definition at line 2944 of file vectorization.h.
|
related |
Addition of a vectorized array and a scalar (expanded to a vectorized array with n_array_elements
equal entries) in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 2960 of file vectorization.h.
|
related |
Subtraction of a vectorized array from a scalar (expanded to a vectorized array with n_array_elements
equal entries).
Definition at line 2975 of file vectorization.h.
|
related |
Subtraction of a vectorized array from a scalar (expanded to a vectorized array with n_array_elements
equal entries) in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 2993 of file vectorization.h.
|
related |
Subtraction of a scalar (expanded to a vectorized array with n_array_elements
equal entries) from a vectorized array.
Definition at line 3010 of file vectorization.h.
|
related |
Subtraction of a scalar (expanded to a vectorized array with n_array_elements
equal entries) from a vectorized array in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 3028 of file vectorization.h.
|
related |
Multiplication of a scalar (expanded to a vectorized array with n_array_elements
equal entries) and a vectorized array.
Definition at line 3045 of file vectorization.h.
|
related |
Multiplication of a scalar (expanded to a vectorized array with n_array_elements
equal entries) and a vectorized array in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 3063 of file vectorization.h.
|
related |
Multiplication of a vectorized array and a scalar (expanded to a vectorized array with n_array_elements
equal entries).
Definition at line 3080 of file vectorization.h.
|
related |
Multiplication of a vectorized array and a scalar (expanded to a vectorized array with n_array_elements
equal entries) in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 3096 of file vectorization.h.
|
related |
Quotient between a scalar (expanded to a vectorized array with n_array_elements
equal entries) and a vectorized array.
Definition at line 3111 of file vectorization.h.
|
related |
Quotient between a scalar (expanded to a vectorized array with n_array_elements
equal entries) and a vectorized array in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 3129 of file vectorization.h.
|
related |
Quotient between a vectorized array and a scalar (expanded to a vectorized array with n_array_elements
equal entries).
Definition at line 3146 of file vectorization.h.
|
related |
Quotient between a vectorized array and a scalar (expanded to a vectorized array with n_array_elements
equal entries) in case the scalar is a double (needed in order to be able to write simple code with constants that are usually double numbers).
Definition at line 3164 of file vectorization.h.
|
related |
Unary operator + on a vectorized array.
Definition at line 3180 of file vectorization.h.
|
related |
Unary operator - on a vectorized array.
Definition at line 3193 of file vectorization.h.
|
related |
Compute the sine of a vectorized data field. The result is returned as vectorized array in the form {sin(x[0]), sin(x[1]), ..., sin(x[n_array_elements-1])}
.
Definition at line 3222 of file vectorization.h.
|
related |
Compute the cosine of a vectorized data field. The result is returned as vectorized array in the form {cos(x[0]), cos(x[1]), ..., cos(x[n_array_elements-1])}
.
Definition at line 3249 of file vectorization.h.
|
related |
Compute the tangent of a vectorized data field. The result is returned as vectorized array in the form {tan(x[0]), tan(x[1]), ..., tan(x[n_array_elements-1])}
.
Definition at line 3271 of file vectorization.h.
|
related |
Compute the exponential of a vectorized data field. The result is returned as vectorized array in the form {exp(x[0]), exp(x[1]), ..., exp(x[n_array_elements-1])}
.
Definition at line 3293 of file vectorization.h.
|
related |
Compute the natural logarithm of a vectorized data field. The result is returned as vectorized array in the form {log(x[0]), log(x[1]), ..., log(x[n_array_elements-1])}
.
Definition at line 3315 of file vectorization.h.
|
related |
Compute the square root of a vectorized data field. The result is returned as vectorized array in the form {sqrt(x[0]), sqrt(x[1]), ..., sqrt(x[n_array_elements-1])}
.
Definition at line 3337 of file vectorization.h.
|
related |
Raises the given number x
to the power p
for a vectorized data field. The result is returned as vectorized array in the form {pow(x[0],p), pow(x[1],p), ..., pow(x[n_array_elements-1],p)}
.
Definition at line 3354 of file vectorization.h.
|
related |
Compute the absolute value (modulus) of a vectorized data field. The result is returned as vectorized array in the form {abs(x[0]), abs(x[1]), ..., abs(x[n_array_elements-1])}
.
Definition at line 3377 of file vectorization.h.
|
related |
Compute the componentwise maximum of two vectorized data fields. The result is returned as vectorized array in the form {max(x[0],y[0]), max(x[1],y[1]), ...}
.
Definition at line 3394 of file vectorization.h.
|
related |
Compute the componentwise minimum of two vectorized data fields. The result is returned as vectorized array in the form {min(x[0],y[0]), min(x[1],y[1]), ...}
.
Definition at line 3412 of file vectorization.h.
|
static |
This gives the number of vectors collected in this class.
Definition at line 157 of file vectorization.h.
Number VectorizedArray< Number >::data |
Actual data field. Since this class represents a POD data type, it is declared public.
Definition at line 357 of file vectorization.h.