Quadruple precision Math Library reference (x86)

Table of contents

Data types

Two vector registers are combined to retain a vector of QP FP numbers. One of these registers holds the upper part of QP numbers, and the other register holds the lower part. You have to use the dedicated functions to access each element in a vector variable.

Sleef_quadx2

Sleef_quadx2 is a data type for retaining two QP FP numbers. Please use Sleef_loadq2, Sleef_storeq2, Sleef_getq2 and Sleef_setq2 functions to access the data inside variables in this data type.

Sleef_quadx4

Sleef_quadx4 is a data type for retaining four QP FP numbers. Please use Sleef_loadq4_avx2, Sleef_storeq4_avx2, Sleef_getq4_avx2 and Sleef_setq4_avx2 functions to access the data inside variables in this data type.

Sleef_quadx8

Sleef_quadx8 is a data type for retaining eight QP FP numbers. Please use Sleef_loadq8_avx512f, Sleef_storeq8_avx512f, Sleef_getq8_avx512f and Sleef_setq8_avx512f functions to access the data inside variables in this data type.

Functions for accessing elements inside vector variable

#include <sleefquad.h>

Sleef_quadx2 Sleef_loadq2(Sleef_quad * ptr);
Sleef_quadx2 Sleef_loadq2_sse2(Sleef_quad * ptr);
Sleef_quadx2 Sleef_loadq2_avx2128(Sleef_quad * ptr);
Sleef_quadx4 Sleef_loadq4_avx2(Sleef_quad * ptr);
Sleef_quadx8 Sleef_loadq8_avx512f(Sleef_quad * ptr);

Link with -lsleefquad.

These functions load QP FP values from memory and store the results in a vector variable.

#include <sleefquad.h>

void Sleef_storeq2(Sleef_quad * ptr, Sleef_quadx2 v);
void Sleef_storeq2_sse2(Sleef_quad * ptr, Sleef_quadx2 v);
void Sleef_storeq2_avx2128(Sleef_quad * ptr, Sleef_quadx2 v);
void Sleef_storeq4_avx2(Sleef_quad * ptr, Sleef_quadx4 v);
void Sleef_storeq8_avx512f(Sleef_quad * ptr, Sleef_quadx8 v);

Link with -lsleefquad.

These functions store QP FP values in the given variable to memory.

#include <sleefquad.h>

Sleef_quad Sleef_getq2(Sleef_quadx2 v, int index);
Sleef_quad Sleef_getq2_sse2(Sleef_quadx2 v, int index);
Sleef_quad Sleef_getq2_avx2128(Sleef_quadx2 v, int index);
Sleef_quad Sleef_getq4_avx2(Sleef_quadx4 v, int index);
Sleef_quad Sleef_getq8_avx512f(Sleef_quadx8 v, int index);

Link with -lsleefquad.

These functions extract the index-th element from a QP FP vector variable.

#include <sleefquad.h>

Sleef_quadx2 Sleef_setq2(Sleef_quadx2 v, int index, Sleef_quad e);
Sleef_quadx2 Sleef_setq2_sse2(Sleef_quadx2 v, int index, Sleef_quad e);
Sleef_quadx2 Sleef_setq2_avx2128(Sleef_quadx2 v, int index, Sleef_quad e);
Sleef_quadx4 Sleef_setq4_avx2(Sleef_quadx4 v, int index, Sleef_quad e);
Sleef_quadx8 Sleef_setq8_avx512f(Sleef_quadx8 v, int index, Sleef_quad e);

Link with -lsleefquad.

These functions return a QP FP vector in which the index-th element is e, and other elements are the same as v.

#include <sleefquad.h>

Sleef_quadx2 Sleef_splatq2(Sleef_quad e);
Sleef_quadx2 Sleef_splatq2_sse2(Sleef_quad e);
Sleef_quadx2 Sleef_splatq2_avx2128(Sleef_quad e);
Sleef_quadx4 Sleef_splatq4_avx2(Sleef_quad e);
Sleef_quadx8 Sleef_splatq8_avx512f(Sleef_quad e);

Link with -lsleefquad.

These functions return a QP FP vector in which all elements are set to e.

Conversion functions

Convert QP number to double-precision number

#include <sleefquad.h>

__m128d Sleef_cast_to_doubleq2(Sleef_quadx2 a);
__m128d Sleef_cast_to_doubleq2_sse2(Sleef_quadx2 a);
__m128d Sleef_cast_to_doubleq2_avx2128(Sleef_quadx2 a);
__m256d Sleef_cast_to_doubleq4_avx2(Sleef_quadx4 a);
__m512d Sleef_cast_to_doubleq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These functions convert a QP FP value to a double-precision value.

Convert double-precision number to QP number

#include <sleefquad.h>

Sleef_quadx2 Sleef_cast_from_doubleq2(__m128d a);
Sleef_quadx2 Sleef_cast_from_doubleq2_sse2(__m128d a);
Sleef_quadx2 Sleef_cast_from_doubleq2_avx2128(__m128d a);
Sleef_quadx4 Sleef_cast_from_doubleq4_avx2(__m256d a);
Sleef_quadx8 Sleef_cast_from_doubleq8_avx512f(__m512d a);

Link with -lsleefquad.

These functions convert a double-precision value to a QP FP value.

Convert QP number to 64-bit signed integer

#include <sleefquad.h>

__m128i Sleef_cast_to_int64q2(Sleef_quadx2 a);
__m128i Sleef_cast_to_int64q2_sse2(Sleef_quadx2 a);
__m128i Sleef_cast_to_int64q2_avx2128(Sleef_quadx2 a);
__m256i Sleef_cast_to_int64q4_avx2(Sleef_quadx4 a);
__m512i Sleef_cast_to_int64q8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These functions convert a QP FP value to a 64-bit signed integer.

Convert 64-bit signed integer to QP number

#include <sleefquad.h>

Sleef_quadx2 Sleef_cast_from_int64q2(__m128i a);
Sleef_quadx2 Sleef_cast_from_int64q2_sse2(__m128i a);
Sleef_quadx2 Sleef_cast_from_int64q2_avx2128(__m128i a);
Sleef_quadx4 Sleef_cast_from_int64q4_avx2(__m256i a);
Sleef_quadx8 Sleef_cast_from_int64q8_avx512f(__m512i a);

Link with -lsleefquad.

These functions convert a 64-bit signed integer to a QP FP value.

Convert QP number to 64-bit unsigned integer

#include <sleefquad.h>

__m128i Sleef_cast_to_uint64q2(Sleef_quadx2 a);
__m128i Sleef_cast_to_uint64q2_sse2(Sleef_quadx2 a);
__m128i Sleef_cast_to_uint64q2_avx2128(Sleef_quadx2 a);
__m256i Sleef_cast_to_uint64q4_avx2(Sleef_quadx4 a);
__m512i Sleef_cast_to_uint64q8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These functions convert a QP FP value to a 64-bit signed integer.

Convert 64-bit unsigned integer to QP number

#include <sleefquad.h>

Sleef_quadx2 Sleef_cast_from_uint64q2(__m128i a);
Sleef_quadx2 Sleef_cast_from_uint64q2_sse2(__m128i a);
Sleef_quadx2 Sleef_cast_from_uint64q2_avx2128(__m128i a);
Sleef_quadx4 Sleef_cast_from_uint64q4_avx2(__m256i a);
Sleef_quadx8 Sleef_cast_from_uint64q8_avx512f(__m512i a);

Link with -lsleefquad.

These functions convert a 64-bit unsigned integer to a QP FP value.

Comparison and selection functions

QP comparison functions

#include <sleefquad.h>

__m128i Sleef_icmpltq2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpleq2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpgtq2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpgeq2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpeqq2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpneq2(Sleef_quadx2 a, Sleef_quadx2 b);

__m128i Sleef_icmpltq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpleq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpgtq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpgeq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpeqq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpneq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);

__m128i Sleef_icmpltq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpleq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpgtq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpgeq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpeqq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpneq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);

__m128i Sleef_icmpltq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);
__m128i Sleef_icmpleq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);
__m128i Sleef_icmpgtq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);
__m128i Sleef_icmpgeq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);
__m128i Sleef_icmpeqq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);
__m128i Sleef_icmpneq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);

__m256i Sleef_icmpltq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
__m256i Sleef_icmpleq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
__m256i Sleef_icmpgtq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
__m256i Sleef_icmpgeq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
__m256i Sleef_icmpeqq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
__m256i Sleef_icmpneq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);

Link with -lsleefquad.

These are the vectorized functions of comparison functions.

QP comparison functions of the second kind

#include <sleefquad.h>

__m128i Sleef_icmpq2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_icmpq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);
__m256i Sleef_icmpq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);

Link with -lsleefquad.

These are the vectorized functions of Sleef_icmpq1_purec.

Check orderedness

#include <sleefquad.h>

__m128i Sleef_iunordq2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_iunordq2_sse2(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_iunordq2_avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
__m128i Sleef_iunordq4_avx2(Sleef_quadx4 a, Sleef_quadx4 b);
__m256i Sleef_iunordq8_avx512f(Sleef_quadx8 a, Sleef_quadx8 b);

Link with -lsleefquad.

These are the vectorized functions of Sleef_iunordq1_purec.

Select elements

#include <sleefquad.h>

Sleef_quadx2 Sleef_iselectq2(__m128i c, Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_iselectq2_sse2(__m128i c, Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_iselectq2_avx2128(__m128i c, Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx4 Sleef_iselectq4_avx2(__m128i c, Sleef_quadx4 a, Sleef_quadx4 b);
Sleef_quadx8 Sleef_iselectq8_avx512f(__m256i c, Sleef_quadx8 a, Sleef_quadx8 b);

Link with -lsleefquad.

These are the vectorized functions that operate in the same way as the ternary operator.

Math functions

QP functions for basic arithmetic operations

#include <sleefquad.h>

Sleef_quadx2 Sleef_addq2_u05(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_subq2_u05(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_mulq2_u05(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_divq2_u05(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_negq2(Sleef_quadx2 a);

Sleef_quadx2 Sleef_addq2_u05sse2(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_subq2_u05sse2(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_mulq2_u05sse2(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_divq2_u05sse2(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_negq2_sse2(Sleef_quadx2 a);

Sleef_quadx2 Sleef_addq2_u05avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_subq2_u05avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_mulq2_u05avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_divq2_u05avx2128(Sleef_quadx2 a, Sleef_quadx2 b);
Sleef_quadx2 Sleef_negq2_avx2128(Sleef_quadx2 a);

Sleef_quadx4 Sleef_addq4_u05avx2(Sleef_quadx4 a, Sleef_quadx4 b);
Sleef_quadx4 Sleef_subq4_u05avx2(Sleef_quadx4 a, Sleef_quadx4 b);
Sleef_quadx4 Sleef_mulq4_u05avx2(Sleef_quadx4 a, Sleef_quadx4 b);
Sleef_quadx4 Sleef_divq4_u05avx2(Sleef_quadx4 a, Sleef_quadx4 b);
Sleef_quadx4 Sleef_negq4_avx2(Sleef_quadx4 a);

Sleef_quadx8 Sleef_addq8_u05avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
Sleef_quadx8 Sleef_subq8_u05avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
Sleef_quadx8 Sleef_mulq8_u05avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
Sleef_quadx8 Sleef_divq8_u05avx512f(Sleef_quadx8 a, Sleef_quadx8 b);
Sleef_quadx8 Sleef_negq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of the basic arithmetic operations.

Square root functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_sqrtq2_u05(Sleef_quadx2 a);
Sleef_quadx2 Sleef_sqrtq2_u05sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_sqrtq2_u05avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_sqrtq4_u05avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_sqrtq8_u05avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_sqrtq1_u05purec.

Sine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_sinq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_sinq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_sinq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_sinq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_sinq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_sinq1_u10purec.

Cosine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_cosq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_cosq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_cosq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_cosq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_cosq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_cosq1_u10purec.

Tangent functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_tanq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_tanq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_tanq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_tanq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_tanq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_tanq1_u10purec.

Arc sine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_asinq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_asinq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_asinq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_asinq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_asinq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_asinq1_u10purec.

Arc cosine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_acosq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_acosq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_acosq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_acosq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_acosq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_acosq1_u10purec.

Arc tangent functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_atanq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_atanq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_atanq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_atanq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_atanq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_atanq1_u10purec.

Arc tangent functions of two variables

#include <sleefquad.h>

Sleef_quadx2 Sleef_atanq2_u10(Sleef_quadx2 y, Sleef_quadx2 x);
Sleef_quadx2 Sleef_atanq2_u10sse2(Sleef_quadx2 y, Sleef_quadx2 x);
Sleef_quadx2 Sleef_atanq2_u10avx2128(Sleef_quadx2 y, Sleef_quadx2 x);
Sleef_quadx4 Sleef_atanq4_u10avx2(Sleef_quadx4 y, Sleef_quadx4 x);
Sleef_quadx8 Sleef_atanq8_u10avx512f(Sleef_quadx8 y, Sleef_quadx8 x);

Link with -lsleefquad.

These are the vectorized functions of Sleef_atan2q1_u10purec.

Base-e exponential functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_expq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_expq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_expq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_expq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_expq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_expq1_u10purec.

Base-2 exponential functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_exp2q2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_exp2q2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_exp2q2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_exp2q4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_exp2q8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_exp2q1_u10purec.

Base-10 exponentail

#include <sleefquad.h>

Sleef_quadx2 Sleef_exp10q2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_exp10q2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_exp10q2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_exp10q4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_exp10q8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_exp10q1_u10purec.

Base-e exponential functions minus 1

#include <sleefquad.h>

Sleef_quadx2 Sleef_expm1q2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_expm1q2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_expm1q2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_expm1q4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_expm1q8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_expm1q1_u10purec.

Natural logarithmic functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_logq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_logq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_logq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_logq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_logq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_logq1_u10purec.

Base-2 logarithmic functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_log2q2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_log2q2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_log2q2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_log2q4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_log2q8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_log2q1_u10purec.

Base-10 logarithmic functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_log10q2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_log10q2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_log10q2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_log10q4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_log10q8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_log10q1_u10purec.

Logarithm of one plus argument

#include <sleefquad.h>

Sleef_quadx2 Sleef_log1pq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_log1pq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_log1pq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_log1pq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_log1pq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_log1pq1_u10purec.

Power functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_powq2_u10(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_powq2_u10sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_powq2_u10avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_powq4_u10avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_powq8_u10avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_powq1_u10purec.

Hyperbolic sine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_sinhq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_sinhq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_sinhq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_sinhq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_sinhq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_sinhq1_u10purec.

Hyperbolic cosine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_coshq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_coshq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_coshq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_coshq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_coshq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_coshq1_u10purec.

Hyperbolic tangent functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_tanhq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_tanhq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_tanhq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_tanhq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_tanhq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_tanhq1_u10purec.

Inverse hyperbolic sine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_asinhq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_asinhq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_asinhq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_asinhq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_asinhq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_asinhq1_u10purec.

Inverse hyperbolic cosine functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_acoshq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_acoshq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_acoshq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_acoshq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_acoshq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_acoshq1_u10purec.

Inverse hyperbolic tangent functions

#include <sleefquad.h>

Sleef_quadx2 Sleef_atanhq2_u10(Sleef_quadx2 a);
Sleef_quadx2 Sleef_atanhq2_u10sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_atanhq2_u10avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_atanhq4_u10avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_atanhq8_u10avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_atanhq1_u10purec.

Round to integer towards zero

#include <sleefquad.h>

Sleef_quadx2 Sleef_truncq2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_truncq2_sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_truncq2_avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_truncq4_avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_truncq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_truncq1_purec.

Round to integer towards minus infinity

#include <sleefquad.h>

Sleef_quadx2 Sleef_floorq2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_floorq2_sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_floorq2_avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_floorq4_avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_floorq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_floorq1_purec.

Round to integer towards plus infinity

#include <sleefquad.h>

Sleef_quadx2 Sleef_ceilq2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_ceilq2_sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_ceilq2_avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_ceilq4_avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_ceilq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_ceilq1_purec.

Round to integer away from zero

#include <sleefquad.h>

Sleef_quadx2 Sleef_roundq2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_roundq2_sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_roundq2_avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_roundq4_avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_roundq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_roundq1_purec.

Round to integer, ties round to even

#include <sleefquad.h>

Sleef_quadx2 Sleef_rintq2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_rintq2_sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_rintq2_avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_rintq4_avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_rintq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_rintq1_purec.

Absolute value

#include <sleefquad.h>

Sleef_quadx2 Sleef_fabsq2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_fabsq2_sse2(Sleef_quadx2 a);
Sleef_quadx2 Sleef_fabsq2_avx2128(Sleef_quadx2 a);
Sleef_quadx4 Sleef_fabsq4_avx2(Sleef_quadx4 a);
Sleef_quadx8 Sleef_fabsq8_avx512f(Sleef_quadx8 a);

Link with -lsleefquad.

These are the vectorized functions of Sleef_fabsq1_purec.

Copy sign of a number

#include <sleefquad.h>

Sleef_quadx2 Sleef_copysignq2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_copysignq2_sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_copysignq2_avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_copysignq4_avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_copysignq8_avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_copysignq1_purec.

Maximum of two numbers

#include <sleefquad.h>

Sleef_quadx2 Sleef_fmaxq2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fmaxq2_sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fmaxq2_avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_fmaxq4_avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_fmaxq8_avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_fmaxq1_purec.

Minimum of two numbers

#include <sleefquad.h>

Sleef_quadx2 Sleef_fminq2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fminq2_sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fminq2_avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_fminq4_avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_fminq8_avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_fminq1_purec.

Positive difference

#include <sleefquad.h>

Sleef_quadx2 Sleef_fdimq2_u05(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fdimq2_u05sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fdimq2_u05avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_fdimq4_u05avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_fdimq8_u05avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_fdimq1_u05purec.

Floating point remainder

#include <sleefquad.h>

Sleef_quadx2 Sleef_fmodq2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fmodq2_sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_fmodq2_avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_fmodq4_avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_fmodq8_avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_fmodq1_purec.

Floating point remainder

#include <sleefquad.h>

Sleef_quadx2 Sleef_remainderq2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_remainderq2_sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_remainderq2_avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_remainderq4_avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_remainderq8_avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_remainderq1_purec.

Split a number to fractional and integral components

#include <sleefquad.h>

Sleef_quadx2 Sleef_frexpq2(Sleef_quadx2 x, __m128i * ptr);
Sleef_quadx2 Sleef_frexpq2_sse2(Sleef_quadx2 x, __m128i * ptr);
Sleef_quadx2 Sleef_frexpq2_avx2128(Sleef_quadx2 x, __m128i * ptr);
Sleef_quadx4 Sleef_frexpq4_avx2(Sleef_quadx4 x, __m128i * ptr);
Sleef_quadx8 Sleef_frexpq8_avx512f(Sleef_quadx8 x, __m256i * ptr);

Link with -lsleefquad.

These are the vectorized functions of Sleef_frexpq1_purec.

Break a number into integral and fractional parts

#include <sleefquad.h>

Sleef_quadx2 Sleef_modfq2(Sleef_quadx2 x, Sleef_quadx2 * ptr);
Sleef_quadx2 Sleef_modfq2_sse2(Sleef_quadx2 x, Sleef_quadx2 * ptr);
Sleef_quadx2 Sleef_modfq2_avx2128(Sleef_quadx2 x, Sleef_quadx2 * ptr);
Sleef_quadx4 Sleef_modfq4_avx2(Sleef_quadx4 x, Sleef_quadx4 * ptr);
Sleef_quadx8 Sleef_modfq8_avx512f(Sleef_quadx8 x, Sleef_quadx8 * ptr);

Link with -lsleefquad.

These are the vectorized functions of Sleef_modfq1_purec.

2D Euclidian distance

#include <sleefquad.h>

Sleef_quadx2 Sleef_hypotq2_u05(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_hypotq2_u05sse2(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx2 Sleef_hypotq2_u05avx2128(Sleef_quadx2 x, Sleef_quadx2 y);
Sleef_quadx4 Sleef_hypotq4_u05avx2(Sleef_quadx4 x, Sleef_quadx4 y);
Sleef_quadx8 Sleef_hypotq8_u05avx512f(Sleef_quadx8 x, Sleef_quadx8 y);

Link with -lsleefquad.

These are the vectorized functions of Sleef_hypotq1_u05purec.

Fused multiply and accumulate

#include <sleefquad.h>

Sleef_quadx2 Sleef_fmaq2_u05(Sleef_quadx2 x, Sleef_quadx2 y, Sleef_quadx2 z);
Sleef_quadx2 Sleef_fmaq2_u05sse2(Sleef_quadx2 x, Sleef_quadx2 y, Sleef_quadx2 z);
Sleef_quadx2 Sleef_fmaq2_u05avx2128(Sleef_quadx2 x, Sleef_quadx2 y, Sleef_quadx2 z);
Sleef_quadx4 Sleef_fmaq4_u05avx2(Sleef_quadx4 x, Sleef_quadx4 y, Sleef_quadx4 z);
Sleef_quadx8 Sleef_fmaq8_u05avx512f(Sleef_quadx8 x, Sleef_quadx8 y, Sleef_quadx8 z);

Link with -lsleefquad.

These are the vectorized functions of Sleef_fmaq1_u05purec.

Multiply by integral power of 2

#include <sleefquad.h>

Sleef_quadx2 Sleef_ldexpq2(Sleef_quadx2 x, __m128i e);
Sleef_quadx2 Sleef_ldexpq2_sse2(Sleef_quadx2 x, __m128i e);
Sleef_quadx2 Sleef_ldexpq2_avx2128(Sleef_quadx2 x, __m128i e);
Sleef_quadx4 Sleef_ldexpq4_avx2(Sleef_quadx4 x, __m128i e);
Sleef_quadx8 Sleef_ldexpq8_avx512f(Sleef_quadx8 x, __m256i e);

Link with -lsleefquad.

These are the vectorized functions of Sleef_ldexpq1_purec.

Integer exponent of an FP number

#include <sleefquad.h>

__m128i Sleef_ilogbq2(Sleef_quadx2 x);
__m128i Sleef_ilogbq2_sse2(Sleef_quadx2 x);
__m128i Sleef_ilogbq2_avx2128(Sleef_quadx2 x);
__m128i Sleef_ilogbq4_avx2(Sleef_quadx4 x);
__m256i Sleef_ilogbq8_avx512f(Sleef_quadx8 x);

Link with -lsleefquad.

These are the vectorized functions of Sleef_ilogbq1_purec.