xref: /libCEED/include/ceed/jit-source/hip/hip-ref-operator-assemble.h (revision 9330daecb0fc008043eec1b94c46ef7aecbb00cd)
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 
10*9330daecSnbeams #if CEEDSIZE
11*9330daecSnbeams typedef CeedSize IndexType;
12*9330daecSnbeams #else
13*9330daecSnbeams typedef CeedInt IndexType;
14*9330daecSnbeams #endif
15*9330daecSnbeams 
1607b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
1707b31e0eSJeremy L Thompson // Matrix assembly kernel for low-order elements (2D thread block)
1807b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
192b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BLOCK_SIZE) __global__
202b730f8bSJeremy L Thompson     void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array,
2107b31e0eSJeremy L Thompson                         CeedScalar *__restrict__ values_array) {
22ea61e9acSJeremy L Thompson   // This kernel assumes B_in and B_out have the same number of quadrature points and basis points.
2307b31e0eSJeremy L Thompson   // TODO: expand to more general cases
2407b31e0eSJeremy L Thompson   const int i = threadIdx.x;  // The output row index of each B^TDB operation
2507b31e0eSJeremy L Thompson   const int l = threadIdx.y;  // The output column index of each B^TDB operation
2607b31e0eSJeremy L Thompson                               // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
2707b31e0eSJeremy L Thompson 
28ea61e9acSJeremy L Thompson   // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element,
29ea61e9acSJeremy L Thompson   // comp_in, comp_out, node_row, node_col
30*9330daecSnbeams   const IndexType comp_out_stride = NNODES * NNODES;
31*9330daecSnbeams   const IndexType comp_in_stride  = comp_out_stride * NCOMP;
32*9330daecSnbeams   const IndexType e_stride        = comp_in_stride * NCOMP;
3307b31e0eSJeremy L Thompson   // Strides for QF array, slowest --> fastest:  emode_in, comp_in, emode_out, comp_out, elem, qpt
34*9330daecSnbeams   const IndexType qe_stride         = NQPTS;
35*9330daecSnbeams   const IndexType qcomp_out_stride  = NELEM * qe_stride;
36*9330daecSnbeams   const IndexType qemode_out_stride = qcomp_out_stride * NCOMP;
37*9330daecSnbeams   const IndexType qcomp_in_stride   = qemode_out_stride * NUMEMODEOUT;
38*9330daecSnbeams   const IndexType qemode_in_stride  = qcomp_in_stride * NCOMP;
3907b31e0eSJeremy L Thompson 
4007b31e0eSJeremy L Thompson   // Loop over each element (if necessary)
41*9330daecSnbeams   for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) {
42*9330daecSnbeams     for (IndexType comp_in = 0; comp_in < NCOMP; comp_in++) {
43*9330daecSnbeams       for (IndexType comp_out = 0; comp_out < NCOMP; comp_out++) {
4407b31e0eSJeremy L Thompson         CeedScalar result        = 0.0;
45*9330daecSnbeams         IndexType  qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e;
46*9330daecSnbeams         for (IndexType emode_in = 0; emode_in < NUMEMODEIN; emode_in++) {
47*9330daecSnbeams           IndexType b_in_index = emode_in * NQPTS * NNODES;
48*9330daecSnbeams           for (IndexType emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) {
49*9330daecSnbeams             IndexType b_out_index = emode_out * NQPTS * NNODES;
50*9330daecSnbeams             IndexType qf_index    = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in;
5107b31e0eSJeremy L Thompson             // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
52*9330daecSnbeams             for (IndexType j = 0; j < NQPTS; j++) {
5307b31e0eSJeremy L Thompson               result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l];
5407b31e0eSJeremy L Thompson             }
5507b31e0eSJeremy L Thompson           }  // end of emode_out
5607b31e0eSJeremy L Thompson         }    // end of emode_in
57*9330daecSnbeams         IndexType val_index     = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l;
5807b31e0eSJeremy L Thompson         values_array[val_index] = result;
5907b31e0eSJeremy L Thompson       }  // end of out component
6007b31e0eSJeremy L Thompson     }    // end of in component
6107b31e0eSJeremy L Thompson   }      // end of element loop
6207b31e0eSJeremy L Thompson }
6307b31e0eSJeremy L Thompson 
6407b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
6507b31e0eSJeremy L Thompson // Fallback kernel for larger orders (1D thread block)
6607b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
672b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BLOCK_SIZE) __global__
682b730f8bSJeremy L Thompson     void linearAssembleFallback(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array,
6907b31e0eSJeremy L Thompson                                 CeedScalar *__restrict__ values_array) {
70ea61e9acSJeremy L Thompson   // This kernel assumes B_in and B_out have the same number of quadrature points and basis points.
7107b31e0eSJeremy L Thompson   // TODO: expand to more general cases
7207b31e0eSJeremy L Thompson   const int l = threadIdx.x;  // The output column index of each B^TDB operation
7307b31e0eSJeremy L Thompson                               // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
7407b31e0eSJeremy L Thompson 
75ea61e9acSJeremy L Thompson   // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element,
76ea61e9acSJeremy L Thompson   // comp_in, comp_out, node_row, node_col
77*9330daecSnbeams   const IndexType comp_out_stride = NNODES * NNODES;
78*9330daecSnbeams   const IndexType comp_in_stride  = comp_out_stride * NCOMP;
79*9330daecSnbeams   const IndexType e_stride        = comp_in_stride * NCOMP;
8007b31e0eSJeremy L Thompson   // Strides for QF array, slowest --> fastest:  emode_in, comp_in, emode_out, comp_out, elem, qpt
81*9330daecSnbeams   const IndexType qe_stride         = NQPTS;
82*9330daecSnbeams   const IndexType qcomp_out_stride  = NELEM * qe_stride;
83*9330daecSnbeams   const IndexType qemode_out_stride = qcomp_out_stride * NCOMP;
84*9330daecSnbeams   const IndexType qcomp_in_stride   = qemode_out_stride * NUMEMODEOUT;
85*9330daecSnbeams   const IndexType qemode_in_stride  = qcomp_in_stride * NCOMP;
8607b31e0eSJeremy L Thompson 
8707b31e0eSJeremy L Thompson   // Loop over each element (if necessary)
88*9330daecSnbeams   for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) {
89*9330daecSnbeams     for (IndexType comp_in = 0; comp_in < NCOMP; comp_in++) {
90*9330daecSnbeams       for (IndexType comp_out = 0; comp_out < NCOMP; comp_out++) {
91*9330daecSnbeams         for (IndexType i = 0; i < NNODES; i++) {
9207b31e0eSJeremy L Thompson           CeedScalar result        = 0.0;
93*9330daecSnbeams           IndexType  qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e;
94*9330daecSnbeams           for (IndexType emode_in = 0; emode_in < NUMEMODEIN; emode_in++) {
95*9330daecSnbeams             IndexType b_in_index = emode_in * NQPTS * NNODES;
96*9330daecSnbeams             for (IndexType emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) {
97*9330daecSnbeams               IndexType b_out_index = emode_out * NQPTS * NNODES;
98*9330daecSnbeams               IndexType qf_index    = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in;
9907b31e0eSJeremy L Thompson               // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
100*9330daecSnbeams               for (IndexType j = 0; j < NQPTS; j++) {
10107b31e0eSJeremy L Thompson                 result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l];
10207b31e0eSJeremy L Thompson               }
10307b31e0eSJeremy L Thompson             }  // end of emode_out
10407b31e0eSJeremy L Thompson           }    // end of emode_in
105*9330daecSnbeams           IndexType val_index     = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l;
10607b31e0eSJeremy L Thompson           values_array[val_index] = result;
10707b31e0eSJeremy L Thompson         }  // end of loop over element node index, i
10807b31e0eSJeremy L Thompson       }    // end of out component
10907b31e0eSJeremy L Thompson     }      // end of in component
11007b31e0eSJeremy L Thompson   }        // end of element loop
11107b31e0eSJeremy L Thompson }
11207b31e0eSJeremy L Thompson 
11307b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
114