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 8c9c2c079SJeremy L Thompson #include <ceed.h> 907b31e0eSJeremy L Thompson 1007b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 1107b31e0eSJeremy L Thompson // Matrix assembly kernel for low-order elements (2D thread block) 1207b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 132b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BLOCK_SIZE) __global__ 142b730f8bSJeremy L Thompson void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array, 1507b31e0eSJeremy L Thompson CeedScalar *__restrict__ values_array) { 16*ea61e9acSJeremy L Thompson // This kernel assumes B_in and B_out have the same number of quadrature points and basis points. 1707b31e0eSJeremy L Thompson // TODO: expand to more general cases 1807b31e0eSJeremy L Thompson const int i = threadIdx.x; // The output row index of each B^TDB operation 1907b31e0eSJeremy L Thompson const int l = threadIdx.y; // The output column index of each B^TDB operation 2007b31e0eSJeremy L Thompson // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 2107b31e0eSJeremy L Thompson 22*ea61e9acSJeremy L Thompson // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element, 23*ea61e9acSJeremy L Thompson // comp_in, comp_out, node_row, node_col 2407b31e0eSJeremy L Thompson const CeedInt comp_out_stride = NNODES * NNODES; 2507b31e0eSJeremy L Thompson const CeedInt comp_in_stride = comp_out_stride * NCOMP; 2607b31e0eSJeremy L Thompson const CeedInt e_stride = comp_in_stride * NCOMP; 2707b31e0eSJeremy L Thompson // Strides for QF array, slowest --> fastest: emode_in, comp_in, emode_out, comp_out, elem, qpt 2807b31e0eSJeremy L Thompson const CeedInt qe_stride = NQPTS; 2907b31e0eSJeremy L Thompson const CeedInt qcomp_out_stride = NELEM * qe_stride; 3007b31e0eSJeremy L Thompson const CeedInt qemode_out_stride = qcomp_out_stride * NCOMP; 3107b31e0eSJeremy L Thompson const CeedInt qcomp_in_stride = qemode_out_stride * NUMEMODEOUT; 3207b31e0eSJeremy L Thompson const CeedInt qemode_in_stride = qcomp_in_stride * NCOMP; 3307b31e0eSJeremy L Thompson 3407b31e0eSJeremy L Thompson // Loop over each element (if necessary) 352b730f8bSJeremy L Thompson for (CeedInt e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) { 3607b31e0eSJeremy L Thompson for (CeedInt comp_in = 0; comp_in < NCOMP; comp_in++) { 3707b31e0eSJeremy L Thompson for (CeedInt comp_out = 0; comp_out < NCOMP; comp_out++) { 3807b31e0eSJeremy L Thompson CeedScalar result = 0.0; 3907b31e0eSJeremy L Thompson CeedInt qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e; 4007b31e0eSJeremy L Thompson for (CeedInt emode_in = 0; emode_in < NUMEMODEIN; emode_in++) { 4107b31e0eSJeremy L Thompson CeedInt b_in_index = emode_in * NQPTS * NNODES; 4207b31e0eSJeremy L Thompson for (CeedInt emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) { 4307b31e0eSJeremy L Thompson CeedInt b_out_index = emode_out * NQPTS * NNODES; 4407b31e0eSJeremy L Thompson CeedInt qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in; 4507b31e0eSJeremy L Thompson // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 4607b31e0eSJeremy L Thompson for (CeedInt j = 0; j < NQPTS; j++) { 4707b31e0eSJeremy L Thompson result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l]; 4807b31e0eSJeremy L Thompson } 4907b31e0eSJeremy L Thompson } // end of emode_out 5007b31e0eSJeremy L Thompson } // end of emode_in 5107b31e0eSJeremy L Thompson CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l; 5207b31e0eSJeremy L Thompson values_array[val_index] = result; 5307b31e0eSJeremy L Thompson } // end of out component 5407b31e0eSJeremy L Thompson } // end of in component 5507b31e0eSJeremy L Thompson } // end of element loop 5607b31e0eSJeremy L Thompson } 5707b31e0eSJeremy L Thompson 5807b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 5907b31e0eSJeremy L Thompson // Fallback kernel for larger orders (1D thread block) 6007b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 612b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BLOCK_SIZE) __global__ 622b730f8bSJeremy L Thompson void linearAssembleFallback(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array, 6307b31e0eSJeremy L Thompson CeedScalar *__restrict__ values_array) { 64*ea61e9acSJeremy L Thompson // This kernel assumes B_in and B_out have the same number of quadrature points and basis points. 6507b31e0eSJeremy L Thompson // TODO: expand to more general cases 6607b31e0eSJeremy L Thompson const int l = threadIdx.x; // The output column index of each B^TDB operation 6707b31e0eSJeremy L Thompson // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 6807b31e0eSJeremy L Thompson 69*ea61e9acSJeremy L Thompson // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element, 70*ea61e9acSJeremy L Thompson // comp_in, comp_out, node_row, node_col 7107b31e0eSJeremy L Thompson const CeedInt comp_out_stride = NNODES * NNODES; 7207b31e0eSJeremy L Thompson const CeedInt comp_in_stride = comp_out_stride * NCOMP; 7307b31e0eSJeremy L Thompson const CeedInt e_stride = comp_in_stride * NCOMP; 7407b31e0eSJeremy L Thompson // Strides for QF array, slowest --> fastest: emode_in, comp_in, emode_out, comp_out, elem, qpt 7507b31e0eSJeremy L Thompson const CeedInt qe_stride = NQPTS; 7607b31e0eSJeremy L Thompson const CeedInt qcomp_out_stride = NELEM * qe_stride; 7707b31e0eSJeremy L Thompson const CeedInt qemode_out_stride = qcomp_out_stride * NCOMP; 7807b31e0eSJeremy L Thompson const CeedInt qcomp_in_stride = qemode_out_stride * NUMEMODEOUT; 7907b31e0eSJeremy L Thompson const CeedInt qemode_in_stride = qcomp_in_stride * NCOMP; 8007b31e0eSJeremy L Thompson 8107b31e0eSJeremy L Thompson // Loop over each element (if necessary) 822b730f8bSJeremy L Thompson for (CeedInt e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) { 8307b31e0eSJeremy L Thompson for (CeedInt comp_in = 0; comp_in < NCOMP; comp_in++) { 8407b31e0eSJeremy L Thompson for (CeedInt comp_out = 0; comp_out < NCOMP; comp_out++) { 8507b31e0eSJeremy L Thompson for (CeedInt i = 0; i < NNODES; i++) { 8607b31e0eSJeremy L Thompson CeedScalar result = 0.0; 8707b31e0eSJeremy L Thompson CeedInt qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e; 8807b31e0eSJeremy L Thompson for (CeedInt emode_in = 0; emode_in < NUMEMODEIN; emode_in++) { 8907b31e0eSJeremy L Thompson CeedInt b_in_index = emode_in * NQPTS * NNODES; 9007b31e0eSJeremy L Thompson for (CeedInt emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) { 9107b31e0eSJeremy L Thompson CeedInt b_out_index = emode_out * NQPTS * NNODES; 9207b31e0eSJeremy L Thompson CeedInt qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in; 9307b31e0eSJeremy L Thompson // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 9407b31e0eSJeremy L Thompson for (CeedInt j = 0; j < NQPTS; j++) { 9507b31e0eSJeremy L Thompson result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l]; 9607b31e0eSJeremy L Thompson } 9707b31e0eSJeremy L Thompson } // end of emode_out 9807b31e0eSJeremy L Thompson } // end of emode_in 9907b31e0eSJeremy L Thompson CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l; 10007b31e0eSJeremy L Thompson values_array[val_index] = result; 10107b31e0eSJeremy L Thompson } // end of loop over element node index, i 10207b31e0eSJeremy L Thompson } // end of out component 10307b31e0eSJeremy L Thompson } // end of in component 10407b31e0eSJeremy L Thompson } // end of element loop 10507b31e0eSJeremy L Thompson } 10607b31e0eSJeremy L Thompson 10707b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 108