15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 8b2165e7aSSebastian Grimberg /// @file 9b2165e7aSSebastian Grimberg /// Internal header for CUDA operator full assembly 10b2165e7aSSebastian Grimberg 11c9c2c079SJeremy L Thompson #include <ceed.h> 1207b31e0eSJeremy L Thompson 13ca735530SJeremy L Thompson #if USE_CEEDSIZE 14f7c1b517Snbeams typedef CeedSize IndexType; 15f7c1b517Snbeams #else 16f7c1b517Snbeams typedef CeedInt IndexType; 17f7c1b517Snbeams #endif 18f7c1b517Snbeams 1907b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 20004e4986SSebastian Grimberg // Matrix assembly kernel 2107b31e0eSJeremy L Thompson //------------------------------------------------------------------------------ 222b730f8bSJeremy L Thompson extern "C" __launch_bounds__(BLOCK_SIZE) __global__ 23004e4986SSebastian Grimberg void LinearAssemble(const CeedInt num_elem, const CeedScalar *B_in, const CeedScalar *B_out, const bool *orients_in, 24004e4986SSebastian Grimberg const CeedInt8 *curl_orients_in, const bool *orients_out, const CeedInt8 *curl_orients_out, 25004e4986SSebastian Grimberg const CeedScalar *__restrict__ qf_array, CeedScalar *__restrict__ values_array) { 26004e4986SSebastian Grimberg extern __shared__ CeedScalar s_CT[]; 27*db2becc9SJeremy L Thompson CeedScalar *s_C = &s_CT[NUM_NODES_OUT * NUM_NODES_IN]; 2807b31e0eSJeremy L Thompson 2907b31e0eSJeremy L Thompson const int l = threadIdx.x; // The output column index of each B^T D B operation 3007b31e0eSJeremy L Thompson // such that we have (Bout^T)_ij D_jk Bin_kl = C_il 3107b31e0eSJeremy L Thompson 32004e4986SSebastian Grimberg // Strides for final output ordering, determined by the reference (interface) implementation of the symbolic assembly, slowest --> fastest: e, 33ea61e9acSJeremy L Thompson // comp_in, comp_out, node_row, node_col 34004e4986SSebastian Grimberg const IndexType comp_out_stride = NUM_NODES_OUT * NUM_NODES_IN; 35004e4986SSebastian Grimberg const IndexType comp_in_stride = comp_out_stride * NUM_COMP_OUT; 36004e4986SSebastian Grimberg const IndexType e_stride = comp_in_stride * NUM_COMP_IN; 37004e4986SSebastian Grimberg 38004e4986SSebastian Grimberg // Strides for QF array, slowest --> fastest: e_in, comp_in, e_out, comp_out, e, q 39ca735530SJeremy L Thompson const IndexType q_e_stride = NUM_QPTS; 40004e4986SSebastian Grimberg const IndexType q_comp_out_stride = num_elem * q_e_stride; 41004e4986SSebastian Grimberg const IndexType q_eval_mode_out_stride = q_comp_out_stride * NUM_COMP_OUT; 42004e4986SSebastian Grimberg const IndexType q_comp_in_stride = q_eval_mode_out_stride * NUM_EVAL_MODES_OUT; 43004e4986SSebastian Grimberg const IndexType q_eval_mode_in_stride = q_comp_in_stride * NUM_COMP_IN; 4407b31e0eSJeremy L Thompson 4507b31e0eSJeremy L Thompson // Loop over each element (if necessary) 46004e4986SSebastian Grimberg for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < num_elem; e += gridDim.x * blockDim.z) { 47004e4986SSebastian Grimberg for (IndexType comp_in = 0; comp_in < NUM_COMP_IN; comp_in++) { 48004e4986SSebastian Grimberg for (IndexType comp_out = 0; comp_out < NUM_COMP_OUT; comp_out++) { 49004e4986SSebastian Grimberg for (IndexType i = threadIdx.y; i < NUM_NODES_OUT; i += BLOCK_SIZE_Y) { 5007b31e0eSJeremy L Thompson CeedScalar result = 0.0; 51ca735530SJeremy L Thompson IndexType qf_index_comp = q_comp_in_stride * comp_in + q_comp_out_stride * comp_out + q_e_stride * e; 52ca735530SJeremy L Thompson 53004e4986SSebastian Grimberg for (IndexType e_in = 0; e_in < NUM_EVAL_MODES_IN; e_in++) { 54004e4986SSebastian Grimberg IndexType b_in_index = e_in * NUM_QPTS * NUM_NODES_IN; 55ca735530SJeremy L Thompson 56004e4986SSebastian Grimberg for (IndexType e_out = 0; e_out < NUM_EVAL_MODES_OUT; e_out++) { 57004e4986SSebastian Grimberg IndexType b_out_index = e_out * NUM_QPTS * NUM_NODES_OUT; 58004e4986SSebastian Grimberg IndexType qf_index = qf_index_comp + q_eval_mode_out_stride * e_out + q_eval_mode_in_stride * e_in; 59ca735530SJeremy L Thompson 6007b31e0eSJeremy L Thompson // Perform the B^T D B operation for this 'chunk' of D (the qf_array) 61ca735530SJeremy L Thompson for (IndexType j = 0; j < NUM_QPTS; j++) { 62004e4986SSebastian Grimberg result += B_out[b_out_index + j * NUM_NODES_OUT + i] * qf_array[qf_index + j] * B_in[b_in_index + j * NUM_NODES_IN + l]; 6307b31e0eSJeremy L Thompson } 64004e4986SSebastian Grimberg } // end of out eval mode 65004e4986SSebastian Grimberg } // end of in eval mode 66004e4986SSebastian Grimberg if (orients_in) { 67004e4986SSebastian Grimberg result *= orients_in[NUM_NODES_IN * e + l] ? -1.0 : 1.0; 68004e4986SSebastian Grimberg } 69004e4986SSebastian Grimberg if (orients_out) { 70004e4986SSebastian Grimberg result *= orients_out[NUM_NODES_OUT * e + i] ? -1.0 : 1.0; 71004e4986SSebastian Grimberg } 72004e4986SSebastian Grimberg if (!curl_orients_in && !curl_orients_out) { 73004e4986SSebastian Grimberg IndexType val_index = e_stride * e + comp_in_stride * comp_in + comp_out_stride * comp_out + NUM_NODES_IN * i + l; 74ca735530SJeremy L Thompson 7507b31e0eSJeremy L Thompson values_array[val_index] = result; 76004e4986SSebastian Grimberg } else if (curl_orients_in) { 77004e4986SSebastian Grimberg s_C[NUM_NODES_IN * threadIdx.y + l] = result; 78004e4986SSebastian Grimberg __syncthreads(); 79004e4986SSebastian Grimberg s_CT[NUM_NODES_IN * i + l] = 80004e4986SSebastian Grimberg (l > 0 ? s_C[NUM_NODES_IN * threadIdx.y + l - 1] * curl_orients_in[3 * NUM_NODES_IN * e + 3 * l - 1] : 0.0) + 81004e4986SSebastian Grimberg s_C[NUM_NODES_IN * threadIdx.y + l] * curl_orients_in[3 * NUM_NODES_IN * e + 3 * l + 1] + 82004e4986SSebastian Grimberg (l < (NUM_NODES_IN - 1) ? s_C[NUM_NODES_IN * threadIdx.y + l + 1] * curl_orients_in[3 * NUM_NODES_IN * e + 3 * l + 3] : 0.0); 83004e4986SSebastian Grimberg } else { 84004e4986SSebastian Grimberg s_CT[NUM_NODES_IN * i + l] = result; 85004e4986SSebastian Grimberg } 8607b31e0eSJeremy L Thompson } // end of loop over element node index, i 87004e4986SSebastian Grimberg if (curl_orients_in || curl_orients_out) { 88004e4986SSebastian Grimberg // Compute and store the final T^T (B^T D B T) using the fully computed C T product in shared memory 89004e4986SSebastian Grimberg if (curl_orients_out) __syncthreads(); 90004e4986SSebastian Grimberg for (IndexType i = threadIdx.y; i < NUM_NODES_OUT; i += BLOCK_SIZE_Y) { 91004e4986SSebastian Grimberg IndexType val_index = e_stride * e + comp_in_stride * comp_in + comp_out_stride * comp_out + NUM_NODES_IN * i + l; 92004e4986SSebastian Grimberg 93004e4986SSebastian Grimberg if (curl_orients_out) { 94004e4986SSebastian Grimberg values_array[val_index] = 95004e4986SSebastian Grimberg (i > 0 ? s_CT[NUM_NODES_IN * (i - 1) + l] * curl_orients_out[3 * NUM_NODES_OUT * e + 3 * i - 1] : 0.0) + 96004e4986SSebastian Grimberg s_CT[NUM_NODES_IN * i + l] * curl_orients_out[3 * NUM_NODES_OUT * e + 3 * i + 1] + 97004e4986SSebastian Grimberg (i < (NUM_NODES_OUT - 1) ? s_CT[NUM_NODES_IN * (i + 1) + l] * curl_orients_out[3 * NUM_NODES_OUT * e + 3 * i + 3] : 0.0); 98004e4986SSebastian Grimberg } else { 99004e4986SSebastian Grimberg values_array[val_index] = s_CT[NUM_NODES_IN * i + l]; 100004e4986SSebastian Grimberg } 101004e4986SSebastian Grimberg } 102004e4986SSebastian Grimberg } 10307b31e0eSJeremy L Thompson } // end of out component 10407b31e0eSJeremy L Thompson } // end of in component 10507b31e0eSJeremy L Thompson } // end of element loop 10607b31e0eSJeremy L Thompson } 107