Reference documentation for deal.II version 8.5.1
Public Member Functions | Public Attributes | Static Public Attributes | Private Member Functions | Friends | Related Functions | List of all members
VectorizedArray< Number > Class Template Reference

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

Public Member Functions

DEAL_II_ALWAYS_INLINE VectorizedArrayoperator= (const Number scalar)
 
DEAL_II_ALWAYS_INLINE Number & operator[] (const unsigned int comp)
 
DEAL_II_ALWAYS_INLINE const Number & operator[] (const unsigned int comp) const
 
DEAL_II_ALWAYS_INLINE VectorizedArrayoperator+= (const VectorizedArray< Number > &vec)
 
DEAL_II_ALWAYS_INLINE VectorizedArrayoperator-= (const VectorizedArray< Number > &vec)
 
DEAL_II_ALWAYS_INLINE VectorizedArrayoperator*= (const VectorizedArray< Number > &vec)
 
DEAL_II_ALWAYS_INLINE VectorizedArrayoperator/= (const VectorizedArray< Number > &vec)
 
DEAL_II_ALWAYS_INLINE void load (const Number *ptr)
 
DEAL_II_ALWAYS_INLINE void store (Number *ptr) const
 
DEAL_II_ALWAYS_INLINE void gather (const Number *base_ptr, const unsigned int *offsets)
 
DEAL_II_ALWAYS_INLINE 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

DEAL_II_ALWAYS_INLINE VectorizedArray get_sqrt () const
 
DEAL_II_ALWAYS_INLINE VectorizedArray get_abs () const
 
DEAL_II_ALWAYS_INLINE VectorizedArray get_max (const VectorizedArray &other) const
 
DEAL_II_ALWAYS_INLINE 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 >
DEAL_II_ALWAYS_INLINE 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 >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator- (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator* (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator/ (const VectorizedArray< Number > &u, const VectorizedArray< Number > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ (const Number &u, const VectorizedArray< Number > &v)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator+ (const double &u, const VectorizedArray< float > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ (const VectorizedArray< Number > &v, const Number &u)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator+ (const VectorizedArray< float > &v, const double &u)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator- (const Number &u, const VectorizedArray< Number > &v)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator- (const double &u, const VectorizedArray< float > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator- (const VectorizedArray< Number > &v, const Number &u)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator- (const VectorizedArray< float > &v, const double &u)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator* (const Number &u, const VectorizedArray< Number > &v)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator* (const double &u, const VectorizedArray< float > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator* (const VectorizedArray< Number > &v, const Number &u)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator* (const VectorizedArray< float > &v, const double &u)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator/ (const Number &u, const VectorizedArray< Number > &v)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator/ (const double &u, const VectorizedArray< float > &v)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator/ (const VectorizedArray< Number > &v, const Number &u)
 
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator/ (const VectorizedArray< float > &v, const double &u)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ (const VectorizedArray< Number > &u)
 
template<typename Number >
DEAL_II_ALWAYS_INLINE 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)
 

Detailed Description

template<typename Number>
class VectorizedArray< Number >

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.

Author
Katharina Kormann, Martin Kronbichler, 2010, 2011

Definition at line 36 of file memory_consumption.h.

Member Function Documentation

◆ operator=()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray& VectorizedArray< Number >::operator= ( const Number  scalar)
inline

This function assigns a scalar to this class.

Definition at line 190 of file vectorization.h.

◆ operator[]() [1/2]

template<typename Number>
DEAL_II_ALWAYS_INLINE Number& VectorizedArray< Number >::operator[] ( const unsigned int  comp)
inline

Access operator (only valid with component 0)

Definition at line 201 of file vectorization.h.

◆ operator[]() [2/2]

template<typename Number>
DEAL_II_ALWAYS_INLINE const Number& VectorizedArray< Number >::operator[] ( const unsigned int  comp) const
inline

Constant access operator (only valid with component 0)

Definition at line 213 of file vectorization.h.

◆ operator+=()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray& VectorizedArray< Number >::operator+= ( const VectorizedArray< Number > &  vec)
inline

Addition

Definition at line 225 of file vectorization.h.

◆ operator-=()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray& VectorizedArray< Number >::operator-= ( const VectorizedArray< Number > &  vec)
inline

Subtraction

Definition at line 236 of file vectorization.h.

◆ operator*=()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray& VectorizedArray< Number >::operator*= ( const VectorizedArray< Number > &  vec)
inline

Multiplication

Definition at line 247 of file vectorization.h.

◆ operator/=()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray& VectorizedArray< Number >::operator/= ( const VectorizedArray< Number > &  vec)
inline

Division

Definition at line 258 of file vectorization.h.

◆ load()

template<typename Number>
DEAL_II_ALWAYS_INLINE void VectorizedArray< Number >::load ( const Number *  ptr)
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 271 of file vectorization.h.

◆ store()

template<typename Number>
DEAL_II_ALWAYS_INLINE void VectorizedArray< Number >::store ( Number *  ptr) const
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 283 of file vectorization.h.

◆ gather()

template<typename Number>
DEAL_II_ALWAYS_INLINE void VectorizedArray< Number >::gather ( const Number *  base_ptr,
const unsigned int *  offsets 
)
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):

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

Definition at line 301 of file vectorization.h.

◆ scatter()

template<typename Number>
DEAL_II_ALWAYS_INLINE void VectorizedArray< Number >::scatter ( const unsigned int *  offsets,
Number *  base_ptr 
) const
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):

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

Definition at line 320 of file vectorization.h.

◆ get_sqrt()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray VectorizedArray< Number >::get_sqrt ( ) const
inlineprivate

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

Definition at line 339 of file vectorization.h.

◆ get_abs()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray VectorizedArray< Number >::get_abs ( ) const
inlineprivate

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

Definition at line 352 of file vectorization.h.

◆ get_max()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray VectorizedArray< Number >::get_max ( const VectorizedArray< Number > &  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 365 of file vectorization.h.

◆ get_min()

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray VectorizedArray< Number >::get_min ( const VectorizedArray< Number > &  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 378 of file vectorization.h.

Friends And Related Function Documentation

◆ std::sqrt

template<typename Number>
template<typename Number2 >
VectorizedArray<Number2> std::sqrt ( const VectorizedArray< Number2 > &  )
friend

Make a few functions friends.

◆ make_vectorized_array()

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > make_vectorized_array ( const Number &  u)
related

Create a vectorized array that sets all entries in the array to the given scalar.

Definition at line 409 of file vectorization.h.

◆ vectorized_load_and_transpose()

template<typename Number >
void vectorized_load_and_transpose ( const unsigned int  n_entries,
const Number *  in,
const unsigned int *  offsets,
VectorizedArray< Number > *  out 
)
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:

for (unsigned int i=0; i<n_entries; ++i)
for (unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++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 446 of file vectorization.h.

◆ vectorized_transpose_and_store()

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 
)
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>::n_array_elements; ++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>::n_array_elements; ++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 499 of file vectorization.h.

◆ operator+() [1/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ ( const VectorizedArray< Number > &  u,
const VectorizedArray< Number > &  v 
)
related

Addition of two vectorized arrays with operator +.

Definition at line 2764 of file vectorization.h.

◆ operator-() [1/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator- ( const VectorizedArray< Number > &  u,
const VectorizedArray< Number > &  v 
)
related

Subtraction of two vectorized arrays with operator -.

Definition at line 2779 of file vectorization.h.

◆ operator*() [1/5]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator* ( const VectorizedArray< Number > &  u,
const VectorizedArray< Number > &  v 
)
related

Multiplication of two vectorized arrays with operator *.

Definition at line 2794 of file vectorization.h.

◆ operator/() [1/5]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator/ ( const VectorizedArray< Number > &  u,
const VectorizedArray< Number > &  v 
)
related

Division of two vectorized arrays with operator /.

Definition at line 2809 of file vectorization.h.

◆ operator+() [2/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ ( const Number &  u,
const VectorizedArray< Number > &  v 
)
related

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

Definition at line 2825 of file vectorization.h.

◆ operator+() [3/6]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator+ ( const double &  u,
const VectorizedArray< float > &  v 
)
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 2843 of file vectorization.h.

◆ operator+() [4/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ ( const VectorizedArray< Number > &  v,
const Number &  u 
)
related

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

Definition at line 2860 of file vectorization.h.

◆ operator+() [5/6]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator+ ( const VectorizedArray< float > &  v,
const double &  u 
)
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 2876 of file vectorization.h.

◆ operator-() [2/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator- ( const Number &  u,
const VectorizedArray< Number > &  v 
)
related

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

Definition at line 2891 of file vectorization.h.

◆ operator-() [3/6]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator- ( const double &  u,
const VectorizedArray< float > &  v 
)
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 2909 of file vectorization.h.

◆ operator-() [4/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator- ( const VectorizedArray< Number > &  v,
const Number &  u 
)
related

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

Definition at line 2926 of file vectorization.h.

◆ operator-() [5/6]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator- ( const VectorizedArray< float > &  v,
const double &  u 
)
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 2944 of file vectorization.h.

◆ operator*() [2/5]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator* ( const Number &  u,
const VectorizedArray< Number > &  v 
)
related

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

Definition at line 2961 of file vectorization.h.

◆ operator*() [3/5]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator* ( const double &  u,
const VectorizedArray< float > &  v 
)
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 2979 of file vectorization.h.

◆ operator*() [4/5]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator* ( const VectorizedArray< Number > &  v,
const Number &  u 
)
related

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

Definition at line 2996 of file vectorization.h.

◆ operator*() [5/5]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator* ( const VectorizedArray< float > &  v,
const double &  u 
)
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 3012 of file vectorization.h.

◆ operator/() [2/5]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator/ ( const Number &  u,
const VectorizedArray< Number > &  v 
)
related

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

Definition at line 3027 of file vectorization.h.

◆ operator/() [3/5]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator/ ( const double &  u,
const VectorizedArray< float > &  v 
)
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 3045 of file vectorization.h.

◆ operator/() [4/5]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator/ ( const VectorizedArray< Number > &  v,
const Number &  u 
)
related

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

Definition at line 3062 of file vectorization.h.

◆ operator/() [5/5]

template<typename Number>
DEAL_II_ALWAYS_INLINE VectorizedArray< float > operator/ ( const VectorizedArray< float > &  v,
const double &  u 
)
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 3080 of file vectorization.h.

◆ operator+() [6/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator+ ( const VectorizedArray< Number > &  u)
related

Unary operator + on a vectorized array.

Definition at line 3096 of file vectorization.h.

◆ operator-() [6/6]

template<typename Number >
DEAL_II_ALWAYS_INLINE VectorizedArray< Number > operator- ( const VectorizedArray< Number > &  u)
related

Unary operator - on a vectorized array.

Definition at line 3109 of file vectorization.h.

◆ sin()

template<typename Number >
VectorizedArray< Number > sin ( const ::VectorizedArray< Number > &  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[n_array_elements-1])}.

Definition at line 3138 of file vectorization.h.

◆ cos()

template<typename Number >
VectorizedArray< Number > cos ( const ::VectorizedArray< Number > &  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[n_array_elements-1])}.

Definition at line 3165 of file vectorization.h.

◆ tan()

template<typename Number >
VectorizedArray< Number > tan ( const ::VectorizedArray< Number > &  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[n_array_elements-1])}.

Definition at line 3187 of file vectorization.h.

◆ exp()

template<typename Number >
VectorizedArray< Number > exp ( const ::VectorizedArray< Number > &  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[n_array_elements-1])}.

Definition at line 3209 of file vectorization.h.

◆ log()

template<typename Number >
VectorizedArray< Number > log ( const ::VectorizedArray< Number > &  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[n_array_elements-1])}.

Definition at line 3231 of file vectorization.h.

◆ sqrt()

template<typename Number >
VectorizedArray< Number > sqrt ( const ::VectorizedArray< Number > &  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[n_array_elements-1])}.

Definition at line 3253 of file vectorization.h.

◆ pow()

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

Definition at line 3270 of file vectorization.h.

◆ abs()

template<typename Number >
VectorizedArray< Number > abs ( const ::VectorizedArray< Number > &  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[n_array_elements-1])}.

Definition at line 3293 of file vectorization.h.

◆ max()

template<typename Number >
VectorizedArray< Number > max ( const ::VectorizedArray< Number > &  x,
const ::VectorizedArray< Number > &  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 3310 of file vectorization.h.

◆ min()

template<typename Number >
VectorizedArray< Number > min ( const ::VectorizedArray< Number > &  x,
const ::VectorizedArray< Number > &  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 3328 of file vectorization.h.

Member Data Documentation

◆ n_array_elements

template<typename Number>
const unsigned int VectorizedArray< Number >::n_array_elements = 1
static

This gives the number of vectors collected in this class.

Definition at line 180 of file vectorization.h.

◆ data

template<typename Number>
Number VectorizedArray< Number >::data

Actual data field. Since this class represents a POD data type, it is declared public.

Definition at line 330 of file vectorization.h.


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