Skip to content

Commit

Permalink
basis - add CeedBasisApplyAddAtPoints + default impl
Browse files Browse the repository at this point in the history
  • Loading branch information
jeremylt committed Aug 12, 2024
1 parent 18d93ab commit 3b2de1f
Show file tree
Hide file tree
Showing 6 changed files with 334 additions and 15 deletions.
1 change: 1 addition & 0 deletions include/ceed-impl.h
Original file line number Diff line number Diff line change
Expand Up @@ -188,6 +188,7 @@ struct CeedBasis_private {
int (*Apply)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, CeedVector, CeedVector);
int (*ApplyAdd)(CeedBasis, CeedInt, CeedTransposeMode, CeedEvalMode, CeedVector, CeedVector);
int (*ApplyAtPoints)(CeedBasis, CeedInt, const CeedInt *, CeedTransposeMode, CeedEvalMode, CeedVector, CeedVector, CeedVector);
int (*ApplyAddAtPoints)(CeedBasis, CeedInt, const CeedInt *, CeedTransposeMode, CeedEvalMode, CeedVector, CeedVector, CeedVector);
int (*Destroy)(CeedBasis);
int ref_count;
bool is_tensor_basis; /* flag for tensor basis */
Expand Down
2 changes: 2 additions & 0 deletions include/ceed/ceed.h
Original file line number Diff line number Diff line change
Expand Up @@ -308,6 +308,8 @@ CEED_EXTERN int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeM
CEED_EXTERN int CeedBasisApplyAdd(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, CeedVector v);
CEED_EXTERN int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
CeedVector x_ref, CeedVector u, CeedVector v);
CEED_EXTERN int CeedBasisApplyAddAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v);
CEED_EXTERN int CeedBasisGetCeed(CeedBasis basis, Ceed *ceed);
CEED_EXTERN Ceed CeedBasisReturnCeed(CeedBasis basis);
CEED_EXTERN int CeedBasisGetDimension(CeedBasis basis, CeedInt *dim);
Expand Down
124 changes: 109 additions & 15 deletions interface/ceed-basis.c
Original file line number Diff line number Diff line change
Expand Up @@ -1579,8 +1579,8 @@ int CeedBasisApply(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mode,
@param[in] basis `CeedBasis` to evaluate
@param[in] num_elem The number of elements to apply the basis evaluation to;
the backend will specify the ordering in @ref CeedElemRestrictionCreate()
@param[in] t_mode @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes.
@ref CEED_NOTRANSPOSE is not valid for CeedBasisApplyAdd.
@param[in] t_mode @ref CEED_TRANSPOSE to apply the transpose, mapping from quadrature points to nodes;
@ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAdd()`
@param[in] eval_mode @ref CEED_EVAL_NONE to use values directly,
@ref CEED_EVAL_INTERP to use interpolated values,
@ref CEED_EVAL_GRAD to use gradients,
Expand All @@ -1603,7 +1603,7 @@ int CeedBasisApplyAdd(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mod
}

/**
@brief Apply basis evaluation from nodes to arbitrary points
@brief Check input vector dimensions for CeedBasisApply[Add]AtPoints
@param[in] basis `CeedBasis` to evaluate
@param[in] num_elem The number of elements to apply the basis evaluation to;
Expand All @@ -1620,11 +1620,10 @@ int CeedBasisApplyAdd(CeedBasis basis, CeedInt num_elem, CeedTransposeMode t_mod
@return An error code: 0 - success, otherwise - failure
@ref User
@ref Developer
**/
int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
CeedVector x_ref, CeedVector u, CeedVector v) {
bool is_tensor_basis;
static int CeedBasisApplyAtPointsCheckDims(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
CeedInt dim, num_comp, num_q_comp, num_nodes, P_1d = 1, Q_1d = 1, total_num_points = 0;
CeedSize x_length = 0, u_length = 0, v_length;
Ceed ceed;
Expand Down Expand Up @@ -1677,16 +1676,48 @@ int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num
// LCOV_EXCL_STOP
}
CeedCheck(has_good_dims, ceed, CEED_ERROR_DIMENSION, "Input/output vectors too short for basis and evaluation mode");
return CEED_ERROR_SUCCESS;
}

// Backend method
if (basis->ApplyAtPoints) {
CeedCall(basis->ApplyAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
return CEED_ERROR_SUCCESS;
}
/**
@brief Default implimentation to apply basis evaluation from nodes to arbitrary points
@param[in] basis `CeedBasis` to evaluate
@param[in] apply_add Sum result into target vector or overwrite
@param[in] num_elem The number of elements to apply the basis evaluation to;
the backend will specify the ordering in @ref CeedElemRestrictionCreate()
@param[in] num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
@param[in] t_mode @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
@ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
@param[in] eval_mode @ref CEED_EVAL_INTERP to use interpolated values,
@ref CEED_EVAL_GRAD to use gradients,
@ref CEED_EVAL_WEIGHT to use quadrature weights
@param[in] x_ref `CeedVector` holding reference coordinates of each point
@param[in] u Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
@param[out] v Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
@return An error code: 0 - success, otherwise - failure
@ref Developer
**/
static int CeedBasisApplyAtPoints_Core(CeedBasis basis, bool apply_add, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode,
CeedEvalMode eval_mode, CeedVector x_ref, CeedVector u, CeedVector v) {
CeedInt dim, num_comp, P_1d = 1, Q_1d = 1, total_num_points = num_points[0];
Ceed ceed;

CeedCall(CeedBasisGetCeed(basis, &ceed));
CeedCall(CeedBasisGetDimension(basis, &dim));
CeedCall(CeedBasisGetNumNodes1D(basis, &P_1d));
CeedCall(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
CeedCall(CeedBasisGetNumComponents(basis, &num_comp));

// Default implementation
CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
CeedCheck(is_tensor_basis, ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases");
{
bool is_tensor_basis;

CeedCall(CeedBasisIsTensor(basis, &is_tensor_basis));
CeedCheck(is_tensor_basis, ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for tensor product bases");
}
CeedCheck(num_elem == 1, ceed, CEED_ERROR_UNSUPPORTED, "Evaluation at arbitrary points only supported for a single element at a time");
if (eval_mode == CEED_EVAL_WEIGHT) {
CeedCall(CeedVectorSetValue(v, 1.0));
Expand Down Expand Up @@ -1858,13 +1889,76 @@ int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num
CeedCall(CeedVectorRestoreArrayRead(u, &u_array));

// -- Interpolate transpose from Chebyshev coefficients
CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
if (apply_add) CeedCall(CeedBasisApplyAdd(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
else CeedCall(CeedBasisApply(basis->basis_chebyshev, 1, CEED_TRANSPOSE, CEED_EVAL_INTERP, basis->vec_chebyshev, v));
break;
}
}
return CEED_ERROR_SUCCESS;
}

/**
@brief Apply basis evaluation from nodes to arbitrary points
@param[in] basis `CeedBasis` to evaluate
@param[in] num_elem The number of elements to apply the basis evaluation to;
the backend will specify the ordering in @ref CeedElemRestrictionCreate()
@param[in] num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
@param[in] t_mode @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
@ref CEED_TRANSPOSE to apply the transpose, mapping from points to nodes
@param[in] eval_mode @ref CEED_EVAL_INTERP to use interpolated values,
@ref CEED_EVAL_GRAD to use gradients,
@ref CEED_EVAL_WEIGHT to use quadrature weights
@param[in] x_ref `CeedVector` holding reference coordinates of each point
@param[in] u Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
@param[out] v Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
@return An error code: 0 - success, otherwise - failure
@ref User
**/
int CeedBasisApplyAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
CeedVector x_ref, CeedVector u, CeedVector v) {
CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
if (basis->ApplyAtPoints) {
CeedCall(basis->ApplyAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
} else {
CeedCall(CeedBasisApplyAtPoints_Core(basis, false, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
}
return CEED_ERROR_SUCCESS;
}

/**
@brief Apply basis evaluation from nodes to arbitrary points and sum into target vector
@param[in] basis `CeedBasis` to evaluate
@param[in] num_elem The number of elements to apply the basis evaluation to;
the backend will specify the ordering in @ref CeedElemRestrictionCreate()
@param[in] num_points Array of the number of points to apply the basis evaluation to in each element, size `num_elem`
@param[in] t_mode @ref CEED_NOTRANSPOSE to evaluate from nodes to points;
@ref CEED_NOTRANSPOSE is not valid for `CeedBasisApplyAddAtPoints()`
@param[in] eval_mode @ref CEED_EVAL_INTERP to use interpolated values,
@ref CEED_EVAL_GRAD to use gradients,
@ref CEED_EVAL_WEIGHT to use quadrature weights
@param[in] x_ref `CeedVector` holding reference coordinates of each point
@param[in] u Input `CeedVector`, of length `num_nodes * num_comp` for @ref CEED_NOTRANSPOSE
@param[out] v Output `CeedVector`, of length `num_points * num_q_comp` for @ref CEED_NOTRANSPOSE with @ref CEED_EVAL_INTERP
@return An error code: 0 - success, otherwise - failure
@ref User
**/
int CeedBasisApplyAddAtPoints(CeedBasis basis, CeedInt num_elem, const CeedInt *num_points, CeedTransposeMode t_mode, CeedEvalMode eval_mode,
CeedVector x_ref, CeedVector u, CeedVector v) {
CeedCall(CeedBasisApplyAtPointsCheckDims(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
if (basis->ApplyAddAtPoints) {
CeedCall(basis->ApplyAddAtPoints(basis, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
} else {
CeedCall(CeedBasisApplyAtPoints_Core(basis, true, num_elem, num_points, t_mode, eval_mode, x_ref, u, v));
}
return CEED_ERROR_SUCCESS;
}

/**
@brief Get the `Ceed` associated with a `CeedBasis`
Expand Down
1 change: 1 addition & 0 deletions interface/ceed.c
Original file line number Diff line number Diff line change
Expand Up @@ -954,6 +954,7 @@ int CeedInit(const char *resource, Ceed *ceed) {
CEED_FTABLE_ENTRY(CeedBasis, Apply),
CEED_FTABLE_ENTRY(CeedBasis, ApplyAdd),
CEED_FTABLE_ENTRY(CeedBasis, ApplyAtPoints),
CEED_FTABLE_ENTRY(CeedBasis, ApplyAddAtPoints),
CEED_FTABLE_ENTRY(CeedBasis, Destroy),
CEED_FTABLE_ENTRY(CeedTensorContract, Apply),
CEED_FTABLE_ENTRY(CeedTensorContract, Destroy),
Expand Down
98 changes: 98 additions & 0 deletions tests/t364-basis.c
Original file line number Diff line number Diff line change
@@ -0,0 +1,98 @@
/// @file
/// Test polynomial interpolation transpose ApplyAdd from arbitrary points in 1D
/// \test Test polynomial interpolation transpose ApplyAdd from arbitrary points in 1D
#include <ceed.h>
#include <math.h>
#include <stdio.h>

#define ALEN(a) (sizeof(a) / sizeof((a)[0]))

static CeedScalar Eval(CeedScalar x, CeedInt n, const CeedScalar *c) {
CeedScalar y = c[n - 1];
for (CeedInt i = n - 2; i >= 0; i--) y = y * x + c[i];
return y;
}

int main(int argc, char **argv) {
Ceed ceed;
CeedVector x, x_nodes, x_points, x_point, u, v, u_point, v_point;
CeedBasis basis_x, basis_u;
const CeedInt p = 5, q = 5, num_points = 4;
const CeedScalar c[4] = {1, 2, 3, 4}; // 1 + 2x + 3x^2 + ...

CeedInit(argv[1], &ceed);

CeedVectorCreate(ceed, 2, &x);
CeedVectorCreate(ceed, p, &x_nodes);
CeedVectorCreate(ceed, num_points, &x_points);
CeedVectorCreate(ceed, 1, &x_point);
CeedVectorCreate(ceed, p, &u);
CeedVectorCreate(ceed, num_points, &v);
CeedVectorCreate(ceed, p, &u_point);
CeedVectorCreate(ceed, 1, &v_point);
CeedVectorSetValue(v_point, 1.0);

// Get nodal coordinates
CeedBasisCreateTensorH1Lagrange(ceed, 1, 1, 2, p, CEED_GAUSS_LOBATTO, &basis_x);
{
CeedScalar x_array[2];

for (CeedInt i = 0; i < 2; i++) x_array[i] = CeedIntPow(-1, i + 1);
CeedVectorSetArray(x, CEED_MEM_HOST, CEED_COPY_VALUES, x_array);
}
CeedBasisApply(basis_x, 1, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, x, x_nodes);

// Set values of u at nodes
{
const CeedScalar *x_array;
CeedScalar u_array[p];

CeedVectorGetArrayRead(x_nodes, CEED_MEM_HOST, &x_array);
for (CeedInt i = 0; i < p; i++) u_array[i] = Eval(x_array[i], ALEN(c), c);
CeedVectorRestoreArrayRead(x_nodes, &x_array);
CeedVectorSetArray(u, CEED_MEM_HOST, CEED_COPY_VALUES, (CeedScalar *)&u_array);
}

// Interpolate to arbitrary points
CeedBasisCreateTensorH1Lagrange(ceed, 1, 1, p, q, CEED_GAUSS, &basis_u);
{
CeedScalar x_array[4] = {-0.33, -0.65, 0.16, 0.99};

CeedVectorSetArray(x_points, CEED_MEM_HOST, CEED_COPY_VALUES, x_array);
}
CeedBasisApplyAtPoints(basis_u, 1, &num_points, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, x_points, u, v);

for (CeedInt i = 0; i < num_points; i++) {
const CeedInt num_point[1] = {1};
CeedScalar fx = 0.0;
const CeedScalar *x_array, *u_array, *v_array, *u_point_array;

CeedVectorGetArrayRead(x_points, CEED_MEM_HOST, &x_array);
CeedVectorGetArrayRead(u, CEED_MEM_HOST, &u_array);
CeedVectorGetArrayRead(v, CEED_MEM_HOST, &v_array);
CeedVectorSetValue(x_point, x_array[i]);
CeedBasisApplyAtPoints(basis_u, 1, num_point, CEED_TRANSPOSE, CEED_EVAL_INTERP, x_point, v_point, u_point);
// Double it
CeedBasisApplyAddAtPoints(basis_u, 1, num_point, CEED_TRANSPOSE, CEED_EVAL_INTERP, x_point, v_point, u_point);
CeedVectorGetArrayRead(u_point, CEED_MEM_HOST, &u_point_array);
for (CeedInt j = 0; j < p; j++) fx += u_array[j] * u_point_array[j];
if (fabs(v_array[i] * 2.0 - fx) > 100. * CEED_EPSILON) printf("%f != %f = f(%f)\n", v_array[i] * 2.0, fx, x_array[i]);
CeedVectorRestoreArrayRead(u_point, &u_point_array);
CeedVectorRestoreArrayRead(x_points, &x_array);
CeedVectorRestoreArrayRead(u, &u_array);
CeedVectorRestoreArrayRead(v, &v_array);
}

CeedVectorDestroy(&x);
CeedVectorDestroy(&x_nodes);
CeedVectorDestroy(&x_points);
CeedVectorDestroy(&x_point);
CeedVectorDestroy(&u);
CeedVectorDestroy(&v);
CeedVectorDestroy(&u_point);
CeedVectorDestroy(&v_point);
CeedBasisDestroy(&basis_x);
CeedBasisDestroy(&basis_u);
CeedDestroy(&ceed);
return 0;
}
Loading

0 comments on commit 3b2de1f

Please sign in to comment.