Skip to content

Commit

Permalink
test - shrink sizes in t319 for non-tensor
Browse files Browse the repository at this point in the history
  • Loading branch information
jeremylt committed Jan 7, 2025
1 parent 083ac09 commit a3fb7fa
Show file tree
Hide file tree
Showing 4 changed files with 6 additions and 12 deletions.
4 changes: 0 additions & 4 deletions backends/hip-shared/ceed-hip-shared-basis.c
Original file line number Diff line number Diff line change
Expand Up @@ -702,10 +702,6 @@ int CeedBasisCreateH1_Hip_shared(CeedElemTopology topo, CeedInt dim, CeedInt num
CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
CeedCallBackend(CeedCalloc(1, &data));

// Check max sizes
CeedCheck(dim <= 3, ceed, CEED_ERROR_BACKEND, "Backend does not implement nontensor bases with dim > 3");
CeedCheck(num_nodes * num_qpts * dim < 52 * 52 * 3, ceed, CEED_ERROR_BACKEND, "Backend does not implement nontensor bases with P * Q this large");

// Copy basis data to GPU
CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp));
CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad));
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -66,8 +66,7 @@ inline __device__ void InterpTranspose1d(SharedData_Hip &data, const CeedScalar
// Derivatives at quadrature points
//------------------------------------------------------------------------------
template <int NUM_COMP, int P, int Q>
inline __device__ void Grad1d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G,
CeedScalar *__restrict__ r_V) {
inline __device__ void Grad1d(SharedData_Hip &data, const CeedScalar *__restrict__ r_U, const CeedScalar *c_G, CeedScalar *__restrict__ r_V) {
for (CeedInt dim = 0; dim < DIM; dim++) {
for (CeedInt comp = 0; comp < NUM_COMP; comp++) {
Contract1d<NUM_COMP, P, Q>(data, &r_U[comp], &c_G[dim * P * Q], &r_V[comp + dim * NUM_COMP]);
Expand Down
9 changes: 4 additions & 5 deletions include/ceed/jit-source/hip/hip-shared-basis-nontensor.h
Original file line number Diff line number Diff line change
Expand Up @@ -9,8 +9,8 @@
/// Internal header for HIP shared memory non-tensor basis
#include <ceed/types.h>

#include "hip-shared-basis-read-write-templates.h"
#include "hip-shared-basis-nontensor-templates.h"
#include "hip-shared-basis-read-write-templates.h"

//------------------------------------------------------------------------------
// Interp kernels
Expand Down Expand Up @@ -99,8 +99,8 @@ extern "C" __launch_bounds__(BASIS_INTERP_BLOCK_SIZE) __global__
//------------------------------------------------------------------------------
// Grad kernels
//------------------------------------------------------------------------------
extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__ void Grad(const CeedInt num_elem, const CeedScalar *c_G,
const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__
void Grad(const CeedInt num_elem, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
extern __shared__ CeedScalar slice[];

SharedData_Hip data;
Expand All @@ -127,8 +127,7 @@ extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__ void Grad(const C
}

extern "C" __launch_bounds__(BASIS_GRAD_BLOCK_SIZE) __global__
void GradTranspose(const CeedInt num_elem, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U,
CeedScalar *__restrict__ d_V) {
void GradTranspose(const CeedInt num_elem, const CeedScalar *c_G, const CeedScalar *__restrict__ d_U, CeedScalar *__restrict__ d_V) {
extern __shared__ CeedScalar slice[];

SharedData_Hip data;
Expand Down
2 changes: 1 addition & 1 deletion tests/t319-basis.c
Original file line number Diff line number Diff line change
Expand Up @@ -116,7 +116,7 @@ int main(int argc, char **argv) {
for (CeedInt dim = 1; dim <= 3; dim++) {
CeedVector x_corners, x_from, x_to, u_from, u_to, du_to;
CeedBasis basis_x, basis_from, basis_to, basis_project;
CeedInt p_from = 5, p_to = 6, q = 7, x_dim = CeedIntPow(2, dim), p_from_dim = CeedIntPow(p_from, dim), p_to_dim = CeedIntPow(p_to, dim);
CeedInt p_from = 3, p_to = 4, q = 4, x_dim = CeedIntPow(2, dim), p_from_dim = CeedIntPow(p_from, dim), p_to_dim = CeedIntPow(p_to, dim);

CeedVectorCreate(ceed, x_dim * dim, &x_corners);
{
Expand Down

0 comments on commit a3fb7fa

Please sign in to comment.