xref: /libCEED/include/ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
1*5aed82e4SJeremy 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 diagonal assembly
1094b7b29bSJeremy L Thompson #ifndef CEED_CUDA_REF_OPERATOR_ASSEMBLE_DIAGONAL_H
1194b7b29bSJeremy L Thompson #define CEED_CUDA_REF_OPERATOR_ASSEMBLE_DIAGONAL_H
12b2165e7aSSebastian Grimberg 
13c9c2c079SJeremy L Thompson #include <ceed.h>
1407b31e0eSJeremy L Thompson 
15ca735530SJeremy L Thompson #if USE_CEEDSIZE
16f7c1b517Snbeams typedef CeedSize IndexType;
17f7c1b517Snbeams #else
18f7c1b517Snbeams typedef CeedInt IndexType;
19f7c1b517Snbeams #endif
20f7c1b517Snbeams 
2107b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
22004e4986SSebastian Grimberg // Get basis pointer
2307b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
24004e4986SSebastian Grimberg static __device__ __inline__ void GetBasisPointer(const CeedScalar **basis_ptr, CeedEvalMode eval_modes, const CeedScalar *identity,
25004e4986SSebastian Grimberg                                                   const CeedScalar *interp, const CeedScalar *grad, const CeedScalar *div, const CeedScalar *curl) {
26004e4986SSebastian Grimberg   switch (eval_modes) {
2707b31e0eSJeremy L Thompson     case CEED_EVAL_NONE:
28ca735530SJeremy L Thompson       *basis_ptr = identity;
2907b31e0eSJeremy L Thompson       break;
3007b31e0eSJeremy L Thompson     case CEED_EVAL_INTERP:
31ca735530SJeremy L Thompson       *basis_ptr = interp;
3207b31e0eSJeremy L Thompson       break;
3307b31e0eSJeremy L Thompson     case CEED_EVAL_GRAD:
34ca735530SJeremy L Thompson       *basis_ptr = grad;
3507b31e0eSJeremy L Thompson       break;
3607b31e0eSJeremy L Thompson     case CEED_EVAL_DIV:
37004e4986SSebastian Grimberg       *basis_ptr = div;
38004e4986SSebastian Grimberg       break;
3907b31e0eSJeremy L Thompson     case CEED_EVAL_CURL:
40004e4986SSebastian Grimberg       *basis_ptr = curl;
41004e4986SSebastian Grimberg       break;
42004e4986SSebastian Grimberg     case CEED_EVAL_WEIGHT:
43004e4986SSebastian Grimberg       break;  // Caught by QF assembly
4407b31e0eSJeremy L Thompson   }
4507b31e0eSJeremy L Thompson }
4607b31e0eSJeremy L Thompson 
4707b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
4807b31e0eSJeremy L Thompson // Core code for diagonal assembly
4907b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
50cbfe683aSSebastian Grimberg extern "C" __launch_bounds__(BLOCK_SIZE) __global__
51cbfe683aSSebastian Grimberg     void LinearDiagonal(const CeedInt num_elem, const CeedScalar *identity, const CeedScalar *interp_in, const CeedScalar *grad_in,
52cbfe683aSSebastian Grimberg                         const CeedScalar *div_in, const CeedScalar *curl_in, const CeedScalar *interp_out, const CeedScalar *grad_out,
53cbfe683aSSebastian Grimberg                         const CeedScalar *div_out, const CeedScalar *curl_out, const CeedEvalMode *eval_modes_in, const CeedEvalMode *eval_modes_out,
54cbfe683aSSebastian Grimberg                         const CeedScalar *__restrict__ assembled_qf_array, CeedScalar *__restrict__ elem_diag_array) {
55004e4986SSebastian Grimberg   const int tid = threadIdx.x;  // Running with P threads
56004e4986SSebastian Grimberg 
57ca735530SJeremy L Thompson   if (tid >= NUM_NODES) return;
5807b31e0eSJeremy L Thompson 
5907b31e0eSJeremy L Thompson   // Compute the diagonal of B^T D B
6007b31e0eSJeremy L Thompson   // Each element
61ca735530SJeremy L Thompson   for (IndexType e = blockIdx.x * blockDim.z + threadIdx.z; e < num_elem; e += gridDim.x * blockDim.z) {
6207b31e0eSJeremy L Thompson     // Each basis eval mode pair
63004e4986SSebastian Grimberg     IndexType    d_out               = 0;
64004e4986SSebastian Grimberg     CeedEvalMode eval_modes_out_prev = CEED_EVAL_NONE;
65004e4986SSebastian Grimberg 
66004e4986SSebastian Grimberg     for (IndexType e_out = 0; e_out < NUM_EVAL_MODES_OUT; e_out++) {
67004e4986SSebastian Grimberg       IndexType         d_in               = 0;
68004e4986SSebastian Grimberg       CeedEvalMode      eval_modes_in_prev = CEED_EVAL_NONE;
69ca735530SJeremy L Thompson       const CeedScalar *b_t                = NULL;
70ca735530SJeremy L Thompson 
71004e4986SSebastian Grimberg       GetBasisPointer(&b_t, eval_modes_out[e_out], identity, interp_out, grad_out, div_out, curl_out);
72004e4986SSebastian Grimberg       if (e_out == 0 || eval_modes_out[e_out] != eval_modes_out_prev) d_out = 0;
73004e4986SSebastian Grimberg       else b_t = &b_t[(++d_out) * NUM_QPTS * NUM_NODES];
74004e4986SSebastian Grimberg       eval_modes_out_prev = eval_modes_out[e_out];
75ca735530SJeremy L Thompson 
76004e4986SSebastian Grimberg       for (IndexType e_in = 0; e_in < NUM_EVAL_MODES_IN; e_in++) {
7707b31e0eSJeremy L Thompson         const CeedScalar *b = NULL;
78ca735530SJeremy L Thompson 
79004e4986SSebastian Grimberg         GetBasisPointer(&b, eval_modes_in[e_in], identity, interp_in, grad_in, div_in, curl_in);
80004e4986SSebastian Grimberg         if (e_in == 0 || eval_modes_in[e_in] != eval_modes_in_prev) d_in = 0;
81004e4986SSebastian Grimberg         else b = &b[(++d_in) * NUM_QPTS * NUM_NODES];
82004e4986SSebastian Grimberg         eval_modes_in_prev = eval_modes_in[e_in];
83004e4986SSebastian Grimberg 
8407b31e0eSJeremy L Thompson         // Each component
85ca735530SJeremy L Thompson         for (IndexType comp_out = 0; comp_out < NUM_COMP; comp_out++) {
86cbfe683aSSebastian Grimberg #if USE_POINT_BLOCK
87004e4986SSebastian Grimberg           // Point block diagonal
88ca735530SJeremy L Thompson           for (IndexType comp_in = 0; comp_in < NUM_COMP; comp_in++) {
89ca735530SJeremy L Thompson             CeedScalar e_value = 0.;
90ca735530SJeremy L Thompson 
91cbfe683aSSebastian Grimberg             // Each qpoint/node pair
92ca735530SJeremy L Thompson             for (IndexType q = 0; q < NUM_QPTS; q++) {
93ca735530SJeremy L Thompson               const CeedScalar qf_value =
94cbfe683aSSebastian Grimberg                   assembled_qf_array[((((e_in * NUM_COMP + comp_in) * NUM_EVAL_MODES_OUT + e_out) * NUM_COMP + comp_out) * num_elem + e) * NUM_QPTS +
95ca735530SJeremy L Thompson                                      q];
96ca735530SJeremy L Thompson 
97ca735530SJeremy L Thompson               e_value += b_t[q * NUM_NODES + tid] * qf_value * b[q * NUM_NODES + tid];
9807b31e0eSJeremy L Thompson             }
99ca735530SJeremy L Thompson             elem_diag_array[((comp_out * NUM_COMP + comp_in) * num_elem + e) * NUM_NODES + tid] += e_value;
10007b31e0eSJeremy L Thompson           }
101cbfe683aSSebastian Grimberg #else
102004e4986SSebastian Grimberg           // Diagonal only
103ca735530SJeremy L Thompson           CeedScalar e_value = 0.;
104ca735530SJeremy L Thompson 
105cbfe683aSSebastian Grimberg           // Each qpoint/node pair
106ca735530SJeremy L Thompson           for (IndexType q = 0; q < NUM_QPTS; q++) {
107ca735530SJeremy L Thompson             const CeedScalar qf_value =
108004e4986SSebastian Grimberg                 assembled_qf_array[((((e_in * NUM_COMP + comp_out) * NUM_EVAL_MODES_OUT + e_out) * NUM_COMP + comp_out) * num_elem + e) * NUM_QPTS +
109004e4986SSebastian Grimberg                                    q];
110ca735530SJeremy L Thompson 
111ca735530SJeremy L Thompson             e_value += b_t[q * NUM_NODES + tid] * qf_value * b[q * NUM_NODES + tid];
11207b31e0eSJeremy L Thompson           }
113ca735530SJeremy L Thompson           elem_diag_array[(comp_out * num_elem + e) * NUM_NODES + tid] += e_value;
114cbfe683aSSebastian Grimberg #endif
11507b31e0eSJeremy L Thompson         }
11607b31e0eSJeremy L Thompson       }
11707b31e0eSJeremy L Thompson     }
11807b31e0eSJeremy L Thompson   }
11907b31e0eSJeremy L Thompson }
12007b31e0eSJeremy L Thompson 
12107b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
122b2165e7aSSebastian Grimberg 
12394b7b29bSJeremy L Thompson #endif  // CEED_CUDA_REF_OPERATOR_ASSEMBLE_DIAGONAL_H
124