Compadre
1.5.9
|
Definition of the divergence-free polynomial basis. More...
Definition of the divergence-free polynomial basis.
Let the number of basis components in the scalar Taylor polynomial basis, for the same dimension and up to the same degree of polynomial. Let the number of basis components in the scalar Taylor polynomial basis, for one dimension lower that does not use same contain the direction v.
where is the basis element number for the divergence-free polynomial basis. The next elements are
where is the basis element number for the divergence-free polynomial basis of one dimension lower using only the variable .
where is the basis element number for the divergence-free polynomial basis. The next elements are
where is the basis element number for the divergence-free polynomial basis. The next elements are
where is the basis element number for the divergence-free polynomial basis of one dimension lower using only the variables and .
Functions | |
KOKKOS_INLINE_FUNCTION int | getSize (const int degree, const int dimension) |
Returns size of basis. More... | |
KOKKOS_INLINE_FUNCTION void | evaluate (const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const int component, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0) |
Evaluates the divergence-free polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_new_value * (calculation of this function) More... | |
KOKKOS_INLINE_FUNCTION void | evaluatePartialDerivative (const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const int component, const int partial_direction, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0) |
Evaluates the first partial derivatives of the divergence-free polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_new_value * (calculation of this function) More... | |
KOKKOS_INLINE_FUNCTION void | evaluateSecondPartialDerivative (const member_type &teamMember, double *delta, double *workspace, const int dimension, const int max_degree, const int component, const int partial_direction_1, const int partial_direction_2, const double h, const double x, const double y, const double z, const int starting_order=0, const double weight_of_original_value=0.0, const double weight_of_new_value=1.0) |
Evaluates the second partial derivatives of the divergence-free polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_new_value * (calculation of this function) More... | |
KOKKOS_INLINE_FUNCTION XYZ | evaluate (int np, const double sx, const double sy, const double sz) |
Evaluates the hard-coded divergence-free polynomial basis up to order 4 for 3D. More... | |
KOKKOS_INLINE_FUNCTION XYZ | evaluate (int np, const double sx, const double sy) |
Evaluates the hard-coded divergence-free polynomial basis up to order 4 for 2D. More... | |
KOKKOS_INLINE_FUNCTION XYZ | evaluatePartialDerivative (int np, const int partial_direction, const double h, const double sx, const double sy, const double sz) |
Evaluates the hard-coded first partial derivative of the divergence-free polynomial basis up to order 4 for 3D. More... | |
KOKKOS_INLINE_FUNCTION XYZ | evaluatePartialDerivative (int np, const int partial_direction, const double h, const double sx, const double sy) |
Evaluates the hard-coded first partial derivative of the divergence-free polynomial basis up to order 4 for 2D. More... | |
KOKKOS_INLINE_FUNCTION void Compadre::DivergenceFreePolynomialBasis::evaluate | ( | const member_type & | teamMember, |
double * | delta, | ||
double * | workspace, | ||
const int | dimension, | ||
const int | max_degree, | ||
const int | component, | ||
const double | h, | ||
const double | x, | ||
const double | y, | ||
const double | z, | ||
const int | starting_order = 0 , |
||
const double | weight_of_original_value = 0.0 , |
||
const double | weight_of_new_value = 1.0 |
||
) |
Evaluates the divergence-free polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_new_value * (calculation of this function)
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _basis_multipler*the dimension of the polynomial basis. |
workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis. |
dimension | [in] - spatial dimension to evaluate |
max_degree | [in] - highest degree of polynomial |
h | [in] - epsilon/window size |
x | [in] - x coordinate (already shifted by target) |
y | [in] - y coordinate (already shifted by target) |
z | [in] - z coordinate (already shifted by target) |
starting_order | [in] - polynomial order to start with (default=0) |
weight_of_original_value | [in] - weighting to assign to original value in delta (default=0, replace, 1-increment current value) |
weight_of_new_value | [in] - weighting to assign to new contribution (default=1, add to/replace) |
Definition at line 47 of file Compadre_DivergenceFreePolynomial.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::DivergenceFreePolynomialBasis::evaluate | ( | int | np, |
const double | sx, | ||
const double | sy, | ||
const double | sz | ||
) |
Evaluates the hard-coded divergence-free polynomial basis up to order 4 for 3D.
Polynomial degree is inferred from np
np | [in] - basis component number |
sx | [in] - x coordinate (already shifted by target and scaled) |
sy | [in] - y coordinate (already shifted by target and scaled) |
sz | [in] - z coordinate (already shifted by target and scaled) |
Definition at line 656 of file Compadre_DivergenceFreePolynomial.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::DivergenceFreePolynomialBasis::evaluate | ( | int | np, |
const double | sx, | ||
const double | sy | ||
) |
Evaluates the hard-coded divergence-free polynomial basis up to order 4 for 2D.
Polynomial degree is inferred from np
np | [in] - basis component number |
sx | [in] - x coordinate (already shifted by target and scaled) |
sy | [in] - y coordinate (already shifted by target and scaled) |
Definition at line 776 of file Compadre_DivergenceFreePolynomial.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::DivergenceFreePolynomialBasis::evaluatePartialDerivative | ( | const member_type & | teamMember, |
double * | delta, | ||
double * | workspace, | ||
const int | dimension, | ||
const int | max_degree, | ||
const int | component, | ||
const int | partial_direction, | ||
const double | h, | ||
const double | x, | ||
const double | y, | ||
const double | z, | ||
const int | starting_order = 0 , |
||
const double | weight_of_original_value = 0.0 , |
||
const double | weight_of_new_value = 1.0 |
||
) |
Evaluates the first partial derivatives of the divergence-free polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_new_value * (calculation of this function)
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _basis_multipler*the dimension of the polynomial basis. |
workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis. |
dimension | [in] - spatial dimension to evaluate |
max_degree | [in] - highest degree of polynomial |
partial_direction | [in] - direction in which to take partial derivative |
h | [in] - epsilon/window size |
x | [in] - x coordinate (already shifted by target) |
y | [in] - y coordinate (already shifted by target) |
z | [in] - z coordinate (already shifted by target) |
starting_order | [in] - polynomial order to start with (default=0) |
weight_of_original_value | [in] - weighting to assign to original value in delta (default=0, replace, 1-increment current value) |
weight_of_new_value | [in] - weighting to assign to new contribution (default=1, add to/replace) |
Definition at line 219 of file Compadre_DivergenceFreePolynomial.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::DivergenceFreePolynomialBasis::evaluatePartialDerivative | ( | int | np, |
const int | partial_direction, | ||
const double | h, | ||
const double | sx, | ||
const double | sy, | ||
const double | sz | ||
) |
Evaluates the hard-coded first partial derivative of the divergence-free polynomial basis up to order 4 for 3D.
Polynomial degree is inferred from np
np | [in] - basis component number |
partial_direction | [in] - partial derivative direction |
sx | [in] - x coordinate (already shifted by target and scaled) |
sy | [in] - y coordinate (already shifted by target and scaled) |
sz | [in] - z coordinate (already shifted by target and scaled) |
Definition at line 831 of file Compadre_DivergenceFreePolynomial.hpp.
KOKKOS_INLINE_FUNCTION XYZ Compadre::DivergenceFreePolynomialBasis::evaluatePartialDerivative | ( | int | np, |
const int | partial_direction, | ||
const double | h, | ||
const double | sx, | ||
const double | sy | ||
) |
Evaluates the hard-coded first partial derivative of the divergence-free polynomial basis up to order 4 for 2D.
Polynomial degree is inferred from np
np | [in] - basis component number |
partial_direction | [in] - partial derivative direction |
sx | [in] - x coordinate (already shifted by target and scaled) |
sy | [in] - y coordinate (already shifted by target and scaled) |
Definition at line 1163 of file Compadre_DivergenceFreePolynomial.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::DivergenceFreePolynomialBasis::evaluateSecondPartialDerivative | ( | const member_type & | teamMember, |
double * | delta, | ||
double * | workspace, | ||
const int | dimension, | ||
const int | max_degree, | ||
const int | component, | ||
const int | partial_direction_1, | ||
const int | partial_direction_2, | ||
const double | h, | ||
const double | x, | ||
const double | y, | ||
const double | z, | ||
const int | starting_order = 0 , |
||
const double | weight_of_original_value = 0.0 , |
||
const double | weight_of_new_value = 1.0 |
||
) |
Evaluates the second partial derivatives of the divergence-free polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_new_value * (calculation of this function)
delta | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _basis_multipler*the dimension of the polynomial basis. |
workspace | [in/out] - scratch space that is allocated so that each thread has its own copy. Must be at least as large as the _poly_order*the spatial dimension of the polynomial basis. |
dimension | [in] - spatial dimension to evaluate |
max_degree | [in] - highest degree of polynomial |
partial_direction_1 | [in] - direction in which to take first partial derivative |
partial_direction_2 | [in] - direction in which to take second partial derivative |
h | [in] - epsilon/window size |
x | [in] - x coordinate (already shifted by target) |
y | [in] - y coordinate (already shifted by target) |
z | [in] - z coordinate (already shifted by target) |
starting_order | [in] - polynomial order to start with (default=0) |
weight_of_original_value | [in] - weighting to assign to original value in delta (default=0, replace, 1-increment current value) |
weight_of_new_value | [in] - weighting to assign to new contribution (default=1, add to/replace) |
Definition at line 438 of file Compadre_DivergenceFreePolynomial.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::DivergenceFreePolynomialBasis::getSize | ( | const int | degree, |
const int | dimension | ||
) |
Returns size of basis.
degree | [in] - highest degree of polynomial |
dimension | [in] - spatial dimension to evaluate |
Definition at line 27 of file Compadre_DivergenceFreePolynomial.hpp.