Compadre
1.5.9
|
Definition of scalar Bernstein polynomial basis. More...
Definition of scalar Bernstein polynomial basis.
For order P, the sum of the basis is defined by:
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 Bernstein 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 Bernstein 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 Bernstein polynomial basis delta[j] = weight_of_original_value * delta[j] + weight_of_new_value * (calculation of this function) More... | |
KOKKOS_INLINE_FUNCTION void Compadre::BernsteinPolynomialBasis::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 Bernstein 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 49 of file Compadre_BernsteinPolynomial.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::BernsteinPolynomialBasis::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 Bernstein 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 155 of file Compadre_BernsteinPolynomial.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::BernsteinPolynomialBasis::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 Bernstein 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 299 of file Compadre_BernsteinPolynomial.hpp.
KOKKOS_INLINE_FUNCTION int Compadre::BernsteinPolynomialBasis::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 30 of file Compadre_BernsteinPolynomial.hpp.