Reference documentation for deal.II version GIT relicensing-245-g36f19064f7 2024-03-29 07:20: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
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Private Member Functions | Friends | Related Symbols | List of all members
VectorizedArray< Number, width > Class Template Reference

#include <deal.II/base/vectorization.h>

Inheritance diagram for VectorizedArray< Number, width >:
Inheritance graph
[legend]

Public Types

using value_type = Number
 

Public Member Functions

 VectorizedArray ()=default
 
 VectorizedArray (const Number scalar)
 
template<typename U >
 VectorizedArray (const std::initializer_list< U > &list)
 
VectorizedArrayoperator= (const Number scalar) &
 
VectorizedArrayoperator= (const Number scalar) &&=delete
 
Number & operator[] (const unsigned int comp)
 
const Number & operator[] (const unsigned int comp) const
 
VectorizedArrayoperator+= (const VectorizedArray &vec)
 
VectorizedArrayoperator-= (const VectorizedArray &vec)
 
VectorizedArrayoperator*= (const VectorizedArray &vec)
 
VectorizedArrayoperator/= (const VectorizedArray &vec)
 
template<typename OtherNumber >
void load (const OtherNumber *ptr)
 
template<typename OtherNumber >
void store (OtherNumber *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
 
Number sum () const
 
constexpr VectorizedArrayIterator< VectorizedArray< Number, width > > begin ()
 
constexpr VectorizedArrayIterator< const VectorizedArray< Number, width > > begin () const
 
constexpr VectorizedArrayIterator< VectorizedArray< Number, width > > end ()
 
constexpr VectorizedArrayIterator< const VectorizedArray< Number, width > > end () const
 

Static Public Member Functions

static constexpr std::size_t size ()
 

Public Attributes

Number data
 

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 , std::size_t width2>
VectorizedArray< Number2, width2 > std::sqrt (const VectorizedArray< Number2, width2 > &)
 
template<typename Number2 , std::size_t width2>
VectorizedArray< Number2, width2 > std::abs (const VectorizedArray< Number2, width2 > &)
 
template<typename Number2 , std::size_t width2>
VectorizedArray< Number2, width2 > std::max (const VectorizedArray< Number2, width2 > &, const VectorizedArray< Number2, width2 > &)
 
template<typename Number2 , std::size_t width2>
VectorizedArray< Number2, width2 > std::min (const VectorizedArray< Number2, width2 > &, const VectorizedArray< Number2, width2 > &)
 

Related Symbols

(Note that these are not member symbols.)

template<typename Number , std::size_t width>
VectorizedArray< Number, width > sin (const ::VectorizedArray< Number, width > &x)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > cos (const ::VectorizedArray< Number, width > &x)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > tan (const ::VectorizedArray< Number, width > &x)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > exp (const ::VectorizedArray< Number, width > &x)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > log (const ::VectorizedArray< Number, width > &x)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > sqrt (const ::VectorizedArray< Number, width > &x)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > pow (const ::VectorizedArray< Number, width > &x, const Number p)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > pow (const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &p)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > abs (const ::VectorizedArray< Number, width > &x)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > max (const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &y)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > min (const ::VectorizedArray< Number, width > &x, const ::VectorizedArray< Number, width > &y)
 
Packing and unpacking of a VectorizedArray
template<typename Number , std::size_t width = internal::VectorizedArrayWidthSpecifier<Number>::max_width>
VectorizedArray< Number, width > make_vectorized_array (const Number &u)
 
template<typename VectorizedArrayType >
VectorizedArrayType make_vectorized_array (const typename VectorizedArrayType::value_type &u)
 
template<typename Number , std::size_t width>
void vectorized_load_and_transpose (const unsigned int n_entries, const Number *in, const unsigned int *offsets, VectorizedArray< Number, width > *out)
 
template<typename Number , std::size_t width>
void vectorized_transpose_and_store (const bool add_into, const unsigned int n_entries, const VectorizedArray< Number, width > *in, const unsigned int *offsets, Number *out)
 
Arithmetic operations with VectorizedArray
template<typename Number , std::size_t width>
bool operator== (const VectorizedArray< Number, width > &lhs, const VectorizedArray< Number, width > &rhs)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ (const VectorizedArray< Number, width > &u, const VectorizedArray< Number, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- (const VectorizedArray< Number, width > &u, const VectorizedArray< Number, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator* (const VectorizedArray< Number, width > &u, const VectorizedArray< Number, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator/ (const VectorizedArray< Number, width > &u, const VectorizedArray< Number, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ (const Number &u, const VectorizedArray< Number, width > &v)
 
template<std::size_t width>
VectorizedArray< float, width > operator+ (const double u, const VectorizedArray< float, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ (const VectorizedArray< Number, width > &v, const Number &u)
 
template<std::size_t width>
VectorizedArray< float, width > operator+ (const VectorizedArray< float, width > &v, const double u)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- (const Number &u, const VectorizedArray< Number, width > &v)
 
template<std::size_t width>
VectorizedArray< float, width > operator- (const double u, const VectorizedArray< float, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- (const VectorizedArray< Number, width > &v, const Number &u)
 
template<std::size_t width>
VectorizedArray< float, width > operator- (const VectorizedArray< float, width > &v, const double u)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator* (const Number &u, const VectorizedArray< Number, width > &v)
 
template<std::size_t width>
VectorizedArray< float, width > operator* (const double u, const VectorizedArray< float, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator* (const VectorizedArray< Number, width > &v, const Number &u)
 
template<std::size_t width>
VectorizedArray< float, width > operator* (const VectorizedArray< float, width > &v, const double u)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator/ (const Number &u, const VectorizedArray< Number, width > &v)
 
template<std::size_t width>
VectorizedArray< float, width > operator/ (const double u, const VectorizedArray< float, width > &v)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator/ (const VectorizedArray< Number, width > &v, const Number &u)
 
template<std::size_t width>
VectorizedArray< float, width > operator/ (const VectorizedArray< float, width > &v, const double u)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ (const VectorizedArray< Number, width > &u)
 
template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- (const VectorizedArray< Number, width > &u)
 
template<typename Number , std::size_t width>
std::ostream & operator<< (std::ostream &out, const VectorizedArray< Number, width > &p)
 

Detailed Description

template<typename Number, std::size_t width>
class VectorizedArray< Number, width >

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 with the goal to be processed in a single-instruction/multiple-data (SIMD) fashion. In the SIMD context, the elements of such a short vector are often called lanes. The number of elements packed together, i.e., the number of lanes, depends 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 (SIMD) units. Most computer systems by 2010 standards will use an array of two doubles or 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 or eight floats are used when deal.II is configured using gcc with --with-cpu=native or --with-cpu=corei7-avx. On compilations with AVX-512 support (e.g., Intel Skylake Server from 2017), eight doubles or 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>() or VectorizedArray<double> a = 0., it sets all numbers in this field to zero. This class is of standard layout type according to the C++11 standard, which means that there is an equivalent C representation and the class can e.g. be safely copied with std::memcpy. (See also https://en.cppreference.com/w/cpp/named_req/StandardLayoutType.) The standard 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 = 5.;, 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.

The user can explicitly control the width of a particular instruction set architecture (ISA) extension by specifying the number of lanes via the second template parameter of this wrapper class. For example on Intel Skylake Server, you have the following options for the data type double:

and for Intel Sandy Bridge, Haswell, Broadwell, AMD Bulldozer and Zen/Ryzen:

and for processors with AltiVec support:

for older x86 processors or in case no processor-specific compilation flags were added (i.e., without -D CMAKE_CXX_FLAGS=-march=native or similar flags):

Similar considerations also apply to the data type float.

Wrongly selecting the width, e.g., width=3 or width=8 on a processor which does not support AVX-512 leads to a static assert.

Template Parameters
Numberunderlying data type
widthvector length (optional; if not set, the maximal width of the architecture is used)

Definition at line 420 of file vectorization.h.

Member Typedef Documentation

◆ value_type

template<typename Number , std::size_t width>
using VectorizedArray< Number, width >::value_type = Number

This gives the type of the array elements.

Definition at line 427 of file vectorization.h.

Constructor & Destructor Documentation

◆ VectorizedArray() [1/3]

template<typename Number , std::size_t width>
VectorizedArray< Number, width >::VectorizedArray ( )
default

Default empty constructor, leaving the data in an uninitialized state similar to float/double.

◆ VectorizedArray() [2/3]

template<typename Number , std::size_t width>
VectorizedArray< Number, width >::VectorizedArray ( const Number  scalar)
inline

Construct an array with the given scalar broadcast to all lanes.

Definition at line 441 of file vectorization.h.

◆ VectorizedArray() [3/3]

template<typename Number , std::size_t width>
template<typename U >
VectorizedArray< Number, width >::VectorizedArray ( const std::initializer_list< U > &  list)
inline

Construct an array with the given initializer list.

Definition at line 450 of file vectorization.h.

Member Function Documentation

◆ operator=() [1/2]

template<typename Number , std::size_t width>
VectorizedArray & VectorizedArray< Number, width >::operator= ( const Number  scalar) &
inline

This function assigns a scalar to the current object.

Definition at line 459 of file vectorization.h.

◆ operator=() [2/2]

template<typename Number , std::size_t width>
VectorizedArray & VectorizedArray< Number, width >::operator= ( const Number  scalar) &&
delete

Assign a scalar to the current object. This overload is used for rvalue references; because it does not make sense to assign something to a temporary, the function is deleted.

◆ operator[]() [1/2]

template<typename Number , std::size_t width>
Number & VectorizedArray< Number, width >::operator[] ( const unsigned int  comp)
inline

Access operator (only valid with component 0 in the base class without specialization).

Definition at line 479 of file vectorization.h.

◆ operator[]() [2/2]

template<typename Number , std::size_t width>
const Number & VectorizedArray< Number, width >::operator[] ( const unsigned int  comp) const
inline

Constant access operator (only valid with component 0 in the base class without specialization).

Definition at line 492 of file vectorization.h.

◆ operator+=()

template<typename Number , std::size_t width>
VectorizedArray & VectorizedArray< Number, width >::operator+= ( const VectorizedArray< Number, width > &  vec)
inline

Addition

Definition at line 504 of file vectorization.h.

◆ operator-=()

template<typename Number , std::size_t width>
VectorizedArray & VectorizedArray< Number, width >::operator-= ( const VectorizedArray< Number, width > &  vec)
inline

Subtraction

Definition at line 515 of file vectorization.h.

◆ operator*=()

template<typename Number , std::size_t width>
VectorizedArray & VectorizedArray< Number, width >::operator*= ( const VectorizedArray< Number, width > &  vec)
inline

Multiplication

Definition at line 526 of file vectorization.h.

◆ operator/=()

template<typename Number , std::size_t width>
VectorizedArray & VectorizedArray< Number, width >::operator/= ( const VectorizedArray< Number, width > &  vec)
inline

Division

Definition at line 537 of file vectorization.h.

◆ load()

template<typename Number , std::size_t width>
template<typename OtherNumber >
void VectorizedArray< Number, width >::load ( const OtherNumber *  ptr)
inline

Load size() data items from memory into the calling class, starting at the given address. The pointer ptr needs 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 551 of file vectorization.h.

◆ store()

template<typename Number , std::size_t width>
template<typename OtherNumber >
void VectorizedArray< Number, width >::store ( OtherNumber *  ptr) const
inline

Write the content of the calling class into memory in form of size() data items to the given address. The pointer ptr needs 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 564 of file vectorization.h.

◆ streaming_store()

template<typename Number , std::size_t width>
void VectorizedArray< Number, width >::streaming_store ( Number *  ptr) const
inline

Write the content of the calling class into memory in form of size() data items 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 (meaning that the memory subsystem makes sure that consecutive stores that together span a cache line are merged, and appropriately handling the case where only part of a cache line is written), so there is 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 storing large arrays 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 617 of file vectorization.h.

◆ gather()

template<typename Number , std::size_t width>
void VectorizedArray< Number, width >::gather ( const Number *  base_ptr,
const unsigned int offsets 
)
inline

Load size() data items 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):

for (unsigned int v=0; v<VectorizedArray<Number>::size(); ++v)
this->operator[](v) = base_ptr[offsets[v]];

Definition at line 636 of file vectorization.h.

◆ scatter()

template<typename Number , std::size_t width>
void VectorizedArray< Number, width >::scatter ( const unsigned int offsets,
Number *  base_ptr 
) const
inline

Write the content of the calling class into memory in form of size() data items 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):

for (unsigned int v=0; v<VectorizedArray<Number>::size(); ++v)
base_ptr[offsets[v]] = this->operator[](v);

Definition at line 655 of file vectorization.h.

◆ sum()

template<typename Number , std::size_t width>
Number VectorizedArray< Number, width >::sum ( ) const
inline

Returns sum over entries of the data field, \(\sum_{i=1}^{\text{size}()} this->data[i]\).

Definition at line 666 of file vectorization.h.

◆ get_sqrt()

template<typename Number , std::size_t width>
VectorizedArray VectorizedArray< Number, width >::get_sqrt ( ) const
inlineprivate

Return the square root of this field. Not for use in user code. Use sqrt(x) instead.

Definition at line 685 of file vectorization.h.

◆ get_abs()

template<typename Number , std::size_t width>
VectorizedArray VectorizedArray< Number, width >::get_abs ( ) const
inlineprivate

Return the absolute value of this field. Not for use in user code. Use abs(x) instead.

Definition at line 698 of file vectorization.h.

◆ get_max()

template<typename Number , std::size_t width>
VectorizedArray VectorizedArray< Number, width >::get_max ( const VectorizedArray< Number, width > &  other) const
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 711 of file vectorization.h.

◆ get_min()

template<typename Number , std::size_t width>
VectorizedArray VectorizedArray< Number, width >::get_min ( const VectorizedArray< Number, width > &  other) const
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 724 of file vectorization.h.

◆ size()

static constexpr std::size_t VectorizedArrayBase< VectorizedArray< Number, width > , width >::size ( )
inlinestaticconstexprinherited

Return the number of elements in the array.

Definition at line 288 of file vectorization.h.

◆ begin() [1/2]

constexpr VectorizedArrayIterator< VectorizedArray< Number, width > > VectorizedArrayBase< VectorizedArray< Number, width > , width >::begin ( )
inlineconstexprinherited
Returns
An iterator pointing to the beginning of the underlying data.

Definition at line 297 of file vectorization.h.

◆ begin() [2/2]

constexpr VectorizedArrayIterator< const VectorizedArray< Number, width > > VectorizedArrayBase< VectorizedArray< Number, width > , width >::begin ( ) const
inlineconstexprinherited
Returns
An iterator pointing to the beginning of the underlying data (const version).

Definition at line 307 of file vectorization.h.

◆ end() [1/2]

constexpr VectorizedArrayIterator< VectorizedArray< Number, width > > VectorizedArrayBase< VectorizedArray< Number, width > , width >::end ( )
inlineconstexprinherited
Returns
An iterator pointing to the end of the underlying data.

Definition at line 316 of file vectorization.h.

◆ end() [2/2]

constexpr VectorizedArrayIterator< const VectorizedArray< Number, width > > VectorizedArrayBase< VectorizedArray< Number, width > , width >::end ( ) const
inlineconstexprinherited
Returns
An iterator pointing to the end of the underlying data (const version).

Definition at line 326 of file vectorization.h.

Friends And Related Symbol Documentation

◆ std::sqrt

template<typename Number , std::size_t width>
template<typename Number2 , std::size_t width2>
VectorizedArray< Number2, width2 > std::sqrt ( const VectorizedArray< Number2, width2 > &  )
friend

◆ std::abs

template<typename Number , std::size_t width>
template<typename Number2 , std::size_t width2>
VectorizedArray< Number2, width2 > std::abs ( const VectorizedArray< Number2, width2 > &  )
friend

◆ std::max

template<typename Number , std::size_t width>
template<typename Number2 , std::size_t width2>
VectorizedArray< Number2, width2 > std::max ( const VectorizedArray< Number2, width2 > &  ,
const VectorizedArray< Number2, width2 > &   
)
friend

◆ std::min

template<typename Number , std::size_t width>
template<typename Number2 , std::size_t width2>
VectorizedArray< Number2, width2 > std::min ( const VectorizedArray< Number2, width2 > &  ,
const VectorizedArray< Number2, width2 > &   
)
friend

◆ make_vectorized_array() [1/2]

template<typename Number , std::size_t width = internal::VectorizedArrayWidthSpecifier<Number>::max_width>
VectorizedArray< Number, width > make_vectorized_array ( const Number &  u)
related

Create a vectorized array that sets all entries in the array to the given scalar, i.e., broadcasts the scalar to all array elements.

Definition at line 765 of file vectorization.h.

◆ make_vectorized_array() [2/2]

template<typename VectorizedArrayType >
VectorizedArrayType make_vectorized_array ( const typename VectorizedArrayType::value_type &  u)
related

Create a vectorized array of given type and broadcast the scalar value to all array elements.

Definition at line 781 of file vectorization.h.

◆ vectorized_load_and_transpose()

template<typename Number , std::size_t width>
void vectorized_load_and_transpose ( const unsigned int  n_entries,
const Number *  in,
const unsigned int offsets,
VectorizedArray< Number, width > *  out 
)
related

This method loads VectorizedArray::size() 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:

for (unsigned int i=0; i<n_entries; ++i)
for (unsigned int v=0; v<VectorizedArray<Number>::size(); ++v)
out[i][v] = in[offsets[v]+i];

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 845 of file vectorization.h.

◆ vectorized_transpose_and_store()

template<typename Number , std::size_t width>
void vectorized_transpose_and_store ( const bool  add_into,
const unsigned int  n_entries,
const VectorizedArray< Number, width > *  in,
const unsigned int offsets,
Number *  out 
)
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 (unsigned int i=0; i<n_entries; ++i)
for (unsigned int v=0; v<VectorizedArray<Number>::size(); ++v)
out[offsets[v]+i] = in[i][v];

For add_into == true, the code implements the following action:

for (unsigned int i=0; i<n_entries; ++i)
for (unsigned int v=0; v<VectorizedArray<Number>::size(); ++v)
out[offsets[v]+i] += in[i][v];

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 920 of file vectorization.h.

◆ operator==()

template<typename Number , std::size_t width>
bool operator== ( const VectorizedArray< Number, width > &  lhs,
const VectorizedArray< Number, width > &  rhs 
)
related

Relational operator == for VectorizedArray

Definition at line 5573 of file vectorization.h.

◆ operator+() [1/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ ( const VectorizedArray< Number, width > &  u,
const VectorizedArray< Number, width > &  v 
)
related

Addition of two vectorized arrays with operator +.

Definition at line 5591 of file vectorization.h.

◆ operator-() [1/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- ( const VectorizedArray< Number, width > &  u,
const VectorizedArray< Number, width > &  v 
)
related

Subtraction of two vectorized arrays with operator -.

Definition at line 5605 of file vectorization.h.

◆ operator*() [1/5]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator* ( const VectorizedArray< Number, width > &  u,
const VectorizedArray< Number, width > &  v 
)
related

Multiplication of two vectorized arrays with operator *.

Definition at line 5619 of file vectorization.h.

◆ operator/() [1/5]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator/ ( const VectorizedArray< Number, width > &  u,
const VectorizedArray< Number, width > &  v 
)
related

Division of two vectorized arrays with operator /.

Definition at line 5633 of file vectorization.h.

◆ operator+() [2/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ ( const Number &  u,
const VectorizedArray< Number, width > &  v 
)
related

Addition of a scalar (expanded to a vectorized array with size() equal entries) and a vectorized array.

Definition at line 5648 of file vectorization.h.

◆ operator+() [3/6]

template<std::size_t width>
VectorizedArray< float, width > operator+ ( const double  u,
const VectorizedArray< float, width > &  v 
)
related

Addition of a scalar (expanded to a vectorized array with size() 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 5664 of file vectorization.h.

◆ operator+() [4/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ ( const VectorizedArray< Number, width > &  v,
const Number &  u 
)
related

Addition of a vectorized array and a scalar (expanded to a vectorized array with size() equal entries).

Definition at line 5678 of file vectorization.h.

◆ operator+() [5/6]

template<std::size_t width>
VectorizedArray< float, width > operator+ ( const VectorizedArray< float, width > &  v,
const double  u 
)
related

Addition of a vectorized array and a scalar (expanded to a vectorized array with size() 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 5693 of file vectorization.h.

◆ operator-() [2/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- ( const Number &  u,
const VectorizedArray< Number, width > &  v 
)
related

Subtraction of a vectorized array from a scalar (expanded to a vectorized array with size() equal entries).

Definition at line 5706 of file vectorization.h.

◆ operator-() [3/6]

template<std::size_t width>
VectorizedArray< float, width > operator- ( const double  u,
const VectorizedArray< float, width > &  v 
)
related

Subtraction of a vectorized array from a scalar (expanded to a vectorized array with size() 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 5722 of file vectorization.h.

◆ operator-() [4/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- ( const VectorizedArray< Number, width > &  v,
const Number &  u 
)
related

Subtraction of a scalar (expanded to a vectorized array with size() equal entries) from a vectorized array.

Definition at line 5736 of file vectorization.h.

◆ operator-() [5/6]

template<std::size_t width>
VectorizedArray< float, width > operator- ( const VectorizedArray< float, width > &  v,
const double  u 
)
related

Subtraction of a scalar (expanded to a vectorized array with size() 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 5752 of file vectorization.h.

◆ operator*() [2/5]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator* ( const Number &  u,
const VectorizedArray< Number, width > &  v 
)
related

Multiplication of a scalar (expanded to a vectorized array with size() equal entries) and a vectorized array.

Definition at line 5766 of file vectorization.h.

◆ operator*() [3/5]

template<std::size_t width>
VectorizedArray< float, width > operator* ( const double  u,
const VectorizedArray< float, width > &  v 
)
related

Multiplication of a scalar (expanded to a vectorized array with size() 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 5782 of file vectorization.h.

◆ operator*() [4/5]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator* ( const VectorizedArray< Number, width > &  v,
const Number &  u 
)
related

Multiplication of a vectorized array and a scalar (expanded to a vectorized array with size() equal entries).

Definition at line 5796 of file vectorization.h.

◆ operator*() [5/5]

template<std::size_t width>
VectorizedArray< float, width > operator* ( const VectorizedArray< float, width > &  v,
const double  u 
)
related

Multiplication of a vectorized array and a scalar (expanded to a vectorized array with size() 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 5811 of file vectorization.h.

◆ operator/() [2/5]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator/ ( const Number &  u,
const VectorizedArray< Number, width > &  v 
)
related

Quotient between a scalar (expanded to a vectorized array with size() equal entries) and a vectorized array.

Definition at line 5824 of file vectorization.h.

◆ operator/() [3/5]

template<std::size_t width>
VectorizedArray< float, width > operator/ ( const double  u,
const VectorizedArray< float, width > &  v 
)
related

Quotient between a scalar (expanded to a vectorized array with size() 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 5840 of file vectorization.h.

◆ operator/() [4/5]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator/ ( const VectorizedArray< Number, width > &  v,
const Number &  u 
)
related

Quotient between a vectorized array and a scalar (expanded to a vectorized array with size() equal entries).

Definition at line 5854 of file vectorization.h.

◆ operator/() [5/5]

template<std::size_t width>
VectorizedArray< float, width > operator/ ( const VectorizedArray< float, width > &  v,
const double  u 
)
related

Quotient between a vectorized array and a scalar (expanded to a vectorized array with size() 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 5870 of file vectorization.h.

◆ operator+() [6/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator+ ( const VectorizedArray< Number, width > &  u)
related

Unary operator + on a vectorized array.

Definition at line 5883 of file vectorization.h.

◆ operator-() [6/6]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > operator- ( const VectorizedArray< Number, width > &  u)
related

Unary operator - on a vectorized array.

Definition at line 5895 of file vectorization.h.

◆ operator<<()

template<typename Number , std::size_t width>
std::ostream & operator<< ( std::ostream &  out,
const VectorizedArray< Number, width > &  p 
)
related

Output operator for vectorized array.

Definition at line 5908 of file vectorization.h.

◆ sin()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > sin ( const ::VectorizedArray< Number, width > &  x)
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[VectorizedArray::size()-1])}.

Definition at line 6493 of file vectorization.h.

◆ cos()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > cos ( const ::VectorizedArray< Number, width > &  x)
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[size()-1])}.

Definition at line 6520 of file vectorization.h.

◆ tan()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > tan ( const ::VectorizedArray< Number, width > &  x)
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[size()-1])}.

Definition at line 6542 of file vectorization.h.

◆ exp()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > exp ( const ::VectorizedArray< Number, width > &  x)
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[size()-1])}.

Definition at line 6564 of file vectorization.h.

◆ log()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > log ( const ::VectorizedArray< Number, width > &  x)
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[size()-1])}.

Definition at line 6586 of file vectorization.h.

◆ sqrt()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > sqrt ( const ::VectorizedArray< Number, width > &  x)
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[size()-1])}.

Definition at line 6608 of file vectorization.h.

◆ pow() [1/2]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > pow ( const ::VectorizedArray< Number, width > &  x,
const Number  p 
)
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[size()-1],p)}.

Definition at line 6624 of file vectorization.h.

◆ pow() [2/2]

template<typename Number , std::size_t width>
VectorizedArray< Number, width > pow ( const ::VectorizedArray< Number, width > &  x,
const ::VectorizedArray< Number, width > &  p 
)
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[0]), pow(x[1],p[1]), ..., pow(x[size()-1],p[size()-1])}.

Definition at line 6647 of file vectorization.h.

◆ abs()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > abs ( const ::VectorizedArray< Number, width > &  x)
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[size()-1])}.

Definition at line 6670 of file vectorization.h.

◆ max()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > max ( const ::VectorizedArray< Number, width > &  x,
const ::VectorizedArray< Number, width > &  y 
)
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 6686 of file vectorization.h.

◆ min()

template<typename Number , std::size_t width>
VectorizedArray< Number, width > min ( const ::VectorizedArray< Number, width > &  x,
const ::VectorizedArray< Number, width > &  y 
)
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 6703 of file vectorization.h.

Member Data Documentation

◆ data

template<typename Number , std::size_t width>
Number VectorizedArray< Number, width >::data

Actual data field. To be consistent with the standard layout type and to enable interaction with external SIMD functionality, this member is declared public.

Definition at line 676 of file vectorization.h.


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