xref: /libCEED/include/ceed/jit-source/hip/hip-ref-operator-assemble.h (revision edf0491998c1d524f2f70fdd683669b8904cb3b6)
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 #if CEEDSIZE
11 typedef CeedSize IndexType;
12 #else
13 typedef CeedInt IndexType;
14 #endif
15 
16 //------------------------------------------------------------------------------
17 // Matrix assembly kernel for low-order elements (2D thread block)
18 //------------------------------------------------------------------------------
19 extern "C" __launch_bounds__(BLOCK_SIZE) __global__
20     void linearAssemble(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array,
21                         CeedScalar *__restrict__ values_array) {
22   // This kernel assumes B_in and B_out have the same number of quadrature points and basis points.
23   // TODO: expand to more general cases
24   const int i = threadIdx.x;  // The output row index of each B^TDB operation
25   const int l = threadIdx.y;  // The output column index of each B^TDB operation
26                               // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
27 
28   // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element,
29   // comp_in, comp_out, node_row, node_col
30   const IndexType comp_out_stride = NNODES * NNODES;
31   const IndexType comp_in_stride  = comp_out_stride * NCOMP;
32   const IndexType e_stride        = comp_in_stride * NCOMP;
33   // Strides for QF array, slowest --> fastest:  emode_in, comp_in, emode_out, comp_out, elem, qpt
34   const IndexType qe_stride         = NQPTS;
35   const IndexType qcomp_out_stride  = NELEM * qe_stride;
36   const IndexType qemode_out_stride = qcomp_out_stride * NCOMP;
37   const IndexType qcomp_in_stride   = qemode_out_stride * NUMEMODEOUT;
38   const IndexType qemode_in_stride  = qcomp_in_stride * NCOMP;
39 
40   // Loop over each element (if necessary)
41   for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) {
42     for (IndexType comp_in = 0; comp_in < NCOMP; comp_in++) {
43       for (IndexType comp_out = 0; comp_out < NCOMP; comp_out++) {
44         CeedScalar result        = 0.0;
45         IndexType  qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e;
46         for (IndexType emode_in = 0; emode_in < NUMEMODEIN; emode_in++) {
47           IndexType b_in_index = emode_in * NQPTS * NNODES;
48           for (IndexType emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) {
49             IndexType b_out_index = emode_out * NQPTS * NNODES;
50             IndexType qf_index    = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in;
51             // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
52             for (IndexType j = 0; j < NQPTS; j++) {
53               result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l];
54             }
55           }  // end of emode_out
56         }    // end of emode_in
57         IndexType val_index     = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l;
58         values_array[val_index] = result;
59       }  // end of out component
60     }    // end of in component
61   }      // end of element loop
62 }
63 
64 //------------------------------------------------------------------------------
65 // Fallback kernel for larger orders (1D thread block)
66 //------------------------------------------------------------------------------
67 extern "C" __launch_bounds__(BLOCK_SIZE) __global__
68     void linearAssembleFallback(const CeedScalar *B_in, const CeedScalar *B_out, const CeedScalar *__restrict__ qf_array,
69                                 CeedScalar *__restrict__ values_array) {
70   // This kernel assumes B_in and B_out have the same number of quadrature points and basis points.
71   // TODO: expand to more general cases
72   const int l = threadIdx.x;  // The output column index of each B^TDB operation
73                               // such that we have (Bout^T)_ij D_jk Bin_kl = C_il
74 
75   // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: element,
76   // comp_in, comp_out, node_row, node_col
77   const IndexType comp_out_stride = NNODES * NNODES;
78   const IndexType comp_in_stride  = comp_out_stride * NCOMP;
79   const IndexType e_stride        = comp_in_stride * NCOMP;
80   // Strides for QF array, slowest --> fastest:  emode_in, comp_in, emode_out, comp_out, elem, qpt
81   const IndexType qe_stride         = NQPTS;
82   const IndexType qcomp_out_stride  = NELEM * qe_stride;
83   const IndexType qemode_out_stride = qcomp_out_stride * NCOMP;
84   const IndexType qcomp_in_stride   = qemode_out_stride * NUMEMODEOUT;
85   const IndexType qemode_in_stride  = qcomp_in_stride * NCOMP;
86 
87   // Loop over each element (if necessary)
88   for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < NELEM; e += gridDim.x * blockDim.z) {
89     for (IndexType comp_in = 0; comp_in < NCOMP; comp_in++) {
90       for (IndexType comp_out = 0; comp_out < NCOMP; comp_out++) {
91         for (IndexType i = 0; i < NNODES; i++) {
92           CeedScalar result        = 0.0;
93           IndexType  qf_index_comp = qcomp_in_stride * comp_in + qcomp_out_stride * comp_out + qe_stride * e;
94           for (IndexType emode_in = 0; emode_in < NUMEMODEIN; emode_in++) {
95             IndexType b_in_index = emode_in * NQPTS * NNODES;
96             for (IndexType emode_out = 0; emode_out < NUMEMODEOUT; emode_out++) {
97               IndexType b_out_index = emode_out * NQPTS * NNODES;
98               IndexType qf_index    = qf_index_comp + qemode_out_stride * emode_out + qemode_in_stride * emode_in;
99               // Perform the B^T D B operation for this 'chunk' of D (the qf_array)
100               for (IndexType j = 0; j < NQPTS; j++) {
101                 result += B_out[b_out_index + j * NNODES + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NNODES + l];
102               }
103             }  // end of emode_out
104           }    // end of emode_in
105           IndexType val_index     = comp_in_stride * comp_in + comp_out_stride * comp_out + e_stride * e + NNODES * i + l;
106           values_array[val_index] = result;
107         }  // end of loop over element node index, i
108       }    // end of out component
109     }      // end of in component
110   }        // end of element loop
111 }
112 
113 //------------------------------------------------------------------------------
114