107b31e0eSJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 207b31e0eSJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 307b31e0eSJeremy L Thompson // 407b31e0eSJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 507b31e0eSJeremy L Thompson // 607b31e0eSJeremy L Thompson // This file is part of CEED: http://github.com/ceed 707b31e0eSJeremy L Thompson 8*b2165e7aSSebastian Grimberg /// @file 9*b2165e7aSSebastian Grimberg /// Internal header for HIP operator full assembly 10*b2165e7aSSebastian Grimberg #ifndef _ceed_hip_ref_operator_assemble_h 11*b2165e7aSSebastian Grimberg #define _ceed_hip_ref_operator_assemble_h 12*b2165e7aSSebastian Grimberg 13c9c2c079SJeremy L Thompson #include <ceed.h> 1407b31e0eSJeremy L Thompson 159330daecSnbeams #if CEEDSIZE 169330daecSnbeams typedef CeedSize IndexType; 179330daecSnbeams #else 189330daecSnbeams typedef CeedInt IndexType; 199330daecSnbeams #endif 209330daecSnbeams 2107b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 2207b31e0eSJeremy L Thompson // Matrix assembly kernel for low-order elements (2D thread block) 2307b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 242b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BLOCK_SIZE) __global__ 252b730f8bSJeremy L Thompson void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array, 2607b31e0eSJeremy L Thompson CeedScalar *__restrict__ values_array) { 27ea61e9acSJeremy L Thompson // This kernel assumes B_in and B_out have the same number of quadrature points and basis points. 2807b31e0eSJeremy L Thompson // TODO: expand to more general cases 2907b31e0eSJeremy L Thompson const int i = threadIdx.x; // The output row index of each B^TDB operation 3007b31e0eSJeremy L Thompson const int l = threadIdx.y; // The output column index of each B^TDB operation 3107b31e0eSJeremy L Thompson // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 3207b31e0eSJeremy L Thompson 33ea61e9acSJeremy L Thompson // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element, 34ea61e9acSJeremy L Thompson // comp_in, comp_out, node_row, node_col 359330daecSnbeams const IndexType comp_out_stride = NNODES * NNODES; 369330daecSnbeams const IndexType comp_in_stride = comp_out_stride * NCOMP; 379330daecSnbeams const IndexType e_stride = comp_in_stride * NCOMP; 3807b31e0eSJeremy L Thompson // Strides for QF array, slowest --> fastest: emode_in, comp_in, emode_out, comp_out, elem, qpt 399330daecSnbeams const IndexType qe_stride = NQPTS; 409330daecSnbeams const IndexType qcomp_out_stride = NELEM * qe_stride; 419330daecSnbeams const IndexType qemode_out_stride = qcomp_out_stride * NCOMP; 429330daecSnbeams const IndexType qcomp_in_stride = qemode_out_stride * NUMEMODEOUT; 439330daecSnbeams const IndexType qemode_in_stride = qcomp_in_stride * NCOMP; 4407b31e0eSJeremy L Thompson 4507b31e0eSJeremy L Thompson // Loop over each element (if necessary) 469330daecSnbeams for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) { 479330daecSnbeams for (IndexType comp_in = 0; comp_in < NCOMP; comp_in++) { 489330daecSnbeams for (IndexType comp_out = 0; comp_out < NCOMP; comp_out++) { 4907b31e0eSJeremy L Thompson CeedScalar result = 0.0; 509330daecSnbeams IndexType qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e; 519330daecSnbeams for (IndexType emode_in = 0; emode_in < NUMEMODEIN; emode_in++) { 529330daecSnbeams IndexType b_in_index = emode_in * NQPTS * NNODES; 539330daecSnbeams for (IndexType emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) { 549330daecSnbeams IndexType b_out_index = emode_out * NQPTS * NNODES; 559330daecSnbeams IndexType qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in; 5607b31e0eSJeremy L Thompson // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 579330daecSnbeams for (IndexType j = 0; j < NQPTS; j++) { 5807b31e0eSJeremy L Thompson result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l]; 5907b31e0eSJeremy L Thompson } 6007b31e0eSJeremy L Thompson } // end of emode_out 6107b31e0eSJeremy L Thompson } // end of emode_in 629330daecSnbeams IndexType val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l; 6307b31e0eSJeremy L Thompson values_array[val_index] = result; 6407b31e0eSJeremy L Thompson } // end of out component 6507b31e0eSJeremy L Thompson } // end of in component 6607b31e0eSJeremy L Thompson } // end of element loop 6707b31e0eSJeremy L Thompson } 6807b31e0eSJeremy L Thompson 6907b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 7007b31e0eSJeremy L Thompson // Fallback kernel for larger orders (1D thread block) 7107b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 722b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BLOCK_SIZE) __global__ 732b730f8bSJeremy L Thompson void linearAssembleFallback(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array, 7407b31e0eSJeremy L Thompson CeedScalar *__restrict__ values_array) { 75ea61e9acSJeremy L Thompson // This kernel assumes B_in and B_out have the same number of quadrature points and basis points. 7607b31e0eSJeremy L Thompson // TODO: expand to more general cases 7707b31e0eSJeremy L Thompson const int l = threadIdx.x; // The output column index of each B^TDB operation 7807b31e0eSJeremy L Thompson // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 7907b31e0eSJeremy L Thompson 80ea61e9acSJeremy L Thompson // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element, 81ea61e9acSJeremy L Thompson // comp_in, comp_out, node_row, node_col 829330daecSnbeams const IndexType comp_out_stride = NNODES * NNODES; 839330daecSnbeams const IndexType comp_in_stride = comp_out_stride * NCOMP; 849330daecSnbeams const IndexType e_stride = comp_in_stride * NCOMP; 8507b31e0eSJeremy L Thompson // Strides for QF array, slowest --> fastest: emode_in, comp_in, emode_out, comp_out, elem, qpt 869330daecSnbeams const IndexType qe_stride = NQPTS; 879330daecSnbeams const IndexType qcomp_out_stride = NELEM * qe_stride; 889330daecSnbeams const IndexType qemode_out_stride = qcomp_out_stride * NCOMP; 899330daecSnbeams const IndexType qcomp_in_stride = qemode_out_stride * NUMEMODEOUT; 909330daecSnbeams const IndexType qemode_in_stride = qcomp_in_stride * NCOMP; 9107b31e0eSJeremy L Thompson 9207b31e0eSJeremy L Thompson // Loop over each element (if necessary) 939330daecSnbeams for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) { 949330daecSnbeams for (IndexType comp_in = 0; comp_in < NCOMP; comp_in++) { 959330daecSnbeams for (IndexType comp_out = 0; comp_out < NCOMP; comp_out++) { 969330daecSnbeams for (IndexType i = 0; i < NNODES; i++) { 9707b31e0eSJeremy L Thompson CeedScalar result = 0.0; 989330daecSnbeams IndexType qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e; 999330daecSnbeams for (IndexType emode_in = 0; emode_in < NUMEMODEIN; emode_in++) { 1009330daecSnbeams IndexType b_in_index = emode_in * NQPTS * NNODES; 1019330daecSnbeams for (IndexType emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) { 1029330daecSnbeams IndexType b_out_index = emode_out * NQPTS * NNODES; 1039330daecSnbeams IndexType qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in; 10407b31e0eSJeremy L Thompson // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 1059330daecSnbeams for (IndexType j = 0; j < NQPTS; j++) { 10607b31e0eSJeremy L Thompson result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l]; 10707b31e0eSJeremy L Thompson } 10807b31e0eSJeremy L Thompson } // end of emode_out 10907b31e0eSJeremy L Thompson } // end of emode_in 1109330daecSnbeams IndexType val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l; 11107b31e0eSJeremy L Thompson values_array[val_index] = result; 11207b31e0eSJeremy L Thompson } // end of loop over element node index, i 11307b31e0eSJeremy L Thompson } // end of out component 11407b31e0eSJeremy L Thompson } // end of in component 11507b31e0eSJeremy L Thompson } // end of element loop 11607b31e0eSJeremy L Thompson } 11707b31e0eSJeremy L Thompson 11807b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 119*b2165e7aSSebastian Grimberg 120*b2165e7aSSebastian Grimberg #endif 121