1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include <ceed.h> 9 10 //------------------------------------------------------------------------------ 11 // Matrix assembly kernel for low-order elements (2D thread block) 12 //------------------------------------------------------------------------------ 13 extern "C" __launch_bounds__(BLOCK_SIZE) __global__ 14 void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array, 15 CeedScalar *__restrict__ values_array) { 16 // This kernel assumes B_in and B_out have the same number of quadrature points and basis points. 17 // TODO: expand to more general cases 18 const int i = threadIdx.x; // The output row index of each B^TDB operation 19 const int l = threadIdx.y; // The output column index of each B^TDB operation 20 // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 21 22 // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element, 23 // comp_in, comp_out, node_row, node_col 24 const CeedInt comp_out_stride = NNODES * NNODES; 25 const CeedInt comp_in_stride = comp_out_stride * NCOMP; 26 const CeedInt e_stride = comp_in_stride * NCOMP; 27 // Strides for QF array, slowest --> fastest: emode_in, comp_in, emode_out, comp_out, elem, qpt 28 const CeedInt qe_stride = NQPTS; 29 const CeedInt qcomp_out_stride = NELEM * qe_stride; 30 const CeedInt qemode_out_stride = qcomp_out_stride * NCOMP; 31 const CeedInt qcomp_in_stride = qemode_out_stride * NUMEMODEOUT; 32 const CeedInt qemode_in_stride = qcomp_in_stride * NCOMP; 33 34 // Loop over each element (if necessary) 35 for (CeedInt e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) { 36 for (CeedInt comp_in = 0; comp_in < NCOMP; comp_in++) { 37 for (CeedInt comp_out = 0; comp_out < NCOMP; comp_out++) { 38 CeedScalar result = 0.0; 39 CeedInt qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e; 40 for (CeedInt emode_in = 0; emode_in < NUMEMODEIN; emode_in++) { 41 CeedInt b_in_index = emode_in * NQPTS * NNODES; 42 for (CeedInt emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) { 43 CeedInt b_out_index = emode_out * NQPTS * NNODES; 44 CeedInt qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in; 45 // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 46 for (CeedInt j = 0; j < NQPTS; j++) { 47 result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l]; 48 } 49 } // end of emode_out 50 } // end of emode_in 51 CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l; 52 values_array[val_index] = result; 53 } // end of out component 54 } // end of in component 55 } // end of element loop 56 } 57 58 //------------------------------------------------------------------------------ 59 // Fallback kernel for larger orders (1D thread block) 60 //------------------------------------------------------------------------------ 61 extern "C" __launch_bounds__(BLOCK_SIZE) __global__ 62 void linearAssembleFallback(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array, 63 CeedScalar *__restrict__ values_array) { 64 // This kernel assumes B_in and B_out have the same number of quadrature points and basis points. 65 // TODO: expand to more general cases 66 const int l = threadIdx.x; // The output column index of each B^TDB operation 67 // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 68 69 // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element, 70 // comp_in, comp_out, node_row, node_col 71 const CeedInt comp_out_stride = NNODES * NNODES; 72 const CeedInt comp_in_stride = comp_out_stride * NCOMP; 73 const CeedInt e_stride = comp_in_stride * NCOMP; 74 // Strides for QF array, slowest --> fastest: emode_in, comp_in, emode_out, comp_out, elem, qpt 75 const CeedInt qe_stride = NQPTS; 76 const CeedInt qcomp_out_stride = NELEM * qe_stride; 77 const CeedInt qemode_out_stride = qcomp_out_stride * NCOMP; 78 const CeedInt qcomp_in_stride = qemode_out_stride * NUMEMODEOUT; 79 const CeedInt qemode_in_stride = qcomp_in_stride * NCOMP; 80 81 // Loop over each element (if necessary) 82 for (CeedInt e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) { 83 for (CeedInt comp_in = 0; comp_in < NCOMP; comp_in++) { 84 for (CeedInt comp_out = 0; comp_out < NCOMP; comp_out++) { 85 for (CeedInt i = 0; i < NNODES; i++) { 86 CeedScalar result = 0.0; 87 CeedInt qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e; 88 for (CeedInt emode_in = 0; emode_in < NUMEMODEIN; emode_in++) { 89 CeedInt b_in_index = emode_in * NQPTS * NNODES; 90 for (CeedInt emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) { 91 CeedInt b_out_index = emode_out * NQPTS * NNODES; 92 CeedInt qf_index = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in; 93 // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 94 for (CeedInt j = 0; j < NQPTS; j++) { 95 result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l]; 96 } 97 } // end of emode_out 98 } // end of emode_in 99 CeedInt val_index = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l; 100 values_array[val_index] = result; 101 } // end of loop over element node index, i 102 } // end of out component 103 } // end of in component 104 } // end of element loop 105 } 106 107 //------------------------------------------------------------------------------ 108