xref: /libCEED/include/ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h (revision 69bf922daa36270062bb51e5ab4f299e9d7ad1af)
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 
807b31e0eSJeremy L Thompson #include <ceed/ceed.h>
907b31e0eSJeremy L Thompson 
1007b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
1107b31e0eSJeremy L Thompson // Diagonal assembly kernels
1207b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
1307b31e0eSJeremy L Thompson 
1407b31e0eSJeremy L Thompson typedef enum {
1507b31e0eSJeremy L Thompson   /// Perform no evaluation (either because there is no data or it is already at
1607b31e0eSJeremy L Thompson   /// quadrature points)
1707b31e0eSJeremy L Thompson   CEED_EVAL_NONE   = 0,
1807b31e0eSJeremy L Thompson   /// Interpolate from nodes to quadrature points
1907b31e0eSJeremy L Thompson   CEED_EVAL_INTERP = 1,
2007b31e0eSJeremy L Thompson   /// Evaluate gradients at quadrature points from input in a nodal basis
2107b31e0eSJeremy L Thompson   CEED_EVAL_GRAD   = 2,
2207b31e0eSJeremy L Thompson   /// Evaluate divergence at quadrature points from input in a nodal basis
2307b31e0eSJeremy L Thompson   CEED_EVAL_DIV    = 4,
2407b31e0eSJeremy L Thompson   /// Evaluate curl at quadrature points from input in a nodal basis
2507b31e0eSJeremy L Thompson   CEED_EVAL_CURL   = 8,
2607b31e0eSJeremy L Thompson   /// Using no input, evaluate quadrature weights on the reference element
2707b31e0eSJeremy L Thompson   CEED_EVAL_WEIGHT = 16,
2807b31e0eSJeremy L Thompson } CeedEvalMode;
2907b31e0eSJeremy L Thompson 
3007b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
3107b31e0eSJeremy L Thompson // Get Basis Emode Pointer
3207b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
3307b31e0eSJeremy L Thompson extern "C" __device__ void CeedOperatorGetBasisPointer_Cuda(const CeedScalar **basisptr,
3407b31e0eSJeremy L Thompson     CeedEvalMode emode, const CeedScalar *identity, const CeedScalar *interp,
3507b31e0eSJeremy L Thompson     const CeedScalar *grad) {
3607b31e0eSJeremy L Thompson   switch (emode) {
3707b31e0eSJeremy L Thompson   case CEED_EVAL_NONE:
3807b31e0eSJeremy L Thompson     *basisptr = identity;
3907b31e0eSJeremy L Thompson     break;
4007b31e0eSJeremy L Thompson   case CEED_EVAL_INTERP:
4107b31e0eSJeremy L Thompson     *basisptr = interp;
4207b31e0eSJeremy L Thompson     break;
4307b31e0eSJeremy L Thompson   case CEED_EVAL_GRAD:
4407b31e0eSJeremy L Thompson     *basisptr = grad;
4507b31e0eSJeremy L Thompson     break;
4607b31e0eSJeremy L Thompson   case CEED_EVAL_WEIGHT:
4707b31e0eSJeremy L Thompson   case CEED_EVAL_DIV:
4807b31e0eSJeremy L Thompson   case CEED_EVAL_CURL:
4907b31e0eSJeremy L Thompson     break; // Caught by QF Assembly
5007b31e0eSJeremy L Thompson   }
5107b31e0eSJeremy L Thompson }
5207b31e0eSJeremy L Thompson 
5307b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
5407b31e0eSJeremy L Thompson // Core code for diagonal assembly
5507b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
5607b31e0eSJeremy L Thompson __device__ void diagonalCore(const CeedInt nelem,
57*69bf922dSJeremy L Thompson     const bool pointBlock, const CeedScalar *identity,
5807b31e0eSJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
5907b31e0eSJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
6007b31e0eSJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
6107b31e0eSJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
6207b31e0eSJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
6307b31e0eSJeremy L Thompson   const int tid = threadIdx.x; // running with P threads, tid is evec node
64*69bf922dSJeremy L Thompson   if (tid >= NNODES) return;
6507b31e0eSJeremy L Thompson 
6607b31e0eSJeremy L Thompson   // Compute the diagonal of B^T D B
6707b31e0eSJeremy L Thompson   // Each element
6807b31e0eSJeremy L Thompson   for (CeedInt e = blockIdx.x*blockDim.z + threadIdx.z; e < nelem;
6907b31e0eSJeremy L Thompson        e += gridDim.x*blockDim.z) {
7007b31e0eSJeremy L Thompson     CeedInt dout = -1;
7107b31e0eSJeremy L Thompson     // Each basis eval mode pair
7207b31e0eSJeremy L Thompson     for (CeedInt eout = 0; eout < NUMEMODEOUT; eout++) {
7307b31e0eSJeremy L Thompson       const CeedScalar *bt = NULL;
7407b31e0eSJeremy L Thompson       if (emodeout[eout] == CEED_EVAL_GRAD)
7507b31e0eSJeremy L Thompson         dout += 1;
7607b31e0eSJeremy L Thompson       CeedOperatorGetBasisPointer_Cuda(&bt, emodeout[eout], identity, interpout,
7707b31e0eSJeremy L Thompson                                       &gradout[dout*NQPTS*NNODES]);
7807b31e0eSJeremy L Thompson       CeedInt din = -1;
7907b31e0eSJeremy L Thompson       for (CeedInt ein = 0; ein < NUMEMODEIN; ein++) {
8007b31e0eSJeremy L Thompson         const CeedScalar *b = NULL;
8107b31e0eSJeremy L Thompson         if (emodein[ein] == CEED_EVAL_GRAD)
8207b31e0eSJeremy L Thompson           din += 1;
8307b31e0eSJeremy L Thompson         CeedOperatorGetBasisPointer_Cuda(&b, emodein[ein], identity, interpin,
8407b31e0eSJeremy L Thompson                                         &gradin[din*NQPTS*NNODES]);
8507b31e0eSJeremy L Thompson         // Each component
8607b31e0eSJeremy L Thompson         for (CeedInt compOut = 0; compOut < NCOMP; compOut++) {
8707b31e0eSJeremy L Thompson           // Each qpoint/node pair
8807b31e0eSJeremy L Thompson           if (pointBlock) {
8907b31e0eSJeremy L Thompson             // Point Block Diagonal
9007b31e0eSJeremy L Thompson             for (CeedInt compIn = 0; compIn < NCOMP; compIn++) {
9107b31e0eSJeremy L Thompson               CeedScalar evalue = 0.;
9207b31e0eSJeremy L Thompson               for (CeedInt q = 0; q < NQPTS; q++) {
9307b31e0eSJeremy L Thompson                 const CeedScalar qfvalue =
9407b31e0eSJeremy L Thompson                   assembledqfarray[((((ein*NCOMP+compIn)*NUMEMODEOUT+eout)*
9507b31e0eSJeremy L Thompson                                      NCOMP+compOut)*nelem+e)*NQPTS+q];
9607b31e0eSJeremy L Thompson                 evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid];
9707b31e0eSJeremy L Thompson               }
9807b31e0eSJeremy L Thompson               elemdiagarray[((compOut*NCOMP+compIn)*nelem+e)*NNODES+tid] += evalue;
9907b31e0eSJeremy L Thompson             }
10007b31e0eSJeremy L Thompson           } else {
10107b31e0eSJeremy L Thompson             // Diagonal Only
10207b31e0eSJeremy L Thompson             CeedScalar evalue = 0.;
10307b31e0eSJeremy L Thompson             for (CeedInt q = 0; q < NQPTS; q++) {
10407b31e0eSJeremy L Thompson               const CeedScalar qfvalue =
10507b31e0eSJeremy L Thompson                 assembledqfarray[((((ein*NCOMP+compOut)*NUMEMODEOUT+eout)*
10607b31e0eSJeremy L Thompson                                    NCOMP+compOut)*nelem+e)*NQPTS+q];
10707b31e0eSJeremy L Thompson               evalue += bt[q*NNODES+tid] * qfvalue * b[q*NNODES+tid];
10807b31e0eSJeremy L Thompson             }
10907b31e0eSJeremy L Thompson             elemdiagarray[(compOut*nelem+e)*NNODES+tid] += evalue;
11007b31e0eSJeremy L Thompson           }
11107b31e0eSJeremy L Thompson         }
11207b31e0eSJeremy L Thompson       }
11307b31e0eSJeremy L Thompson     }
11407b31e0eSJeremy L Thompson   }
11507b31e0eSJeremy L Thompson }
11607b31e0eSJeremy L Thompson 
11707b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
11807b31e0eSJeremy L Thompson // Linear diagonal
11907b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
12007b31e0eSJeremy L Thompson extern "C" __global__ void linearDiagonal(const CeedInt nelem,
121*69bf922dSJeremy L Thompson     const CeedScalar *identity,
12207b31e0eSJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
12307b31e0eSJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
12407b31e0eSJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
12507b31e0eSJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
12607b31e0eSJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
127*69bf922dSJeremy L Thompson   diagonalCore(nelem, false, identity, interpin, gradin, interpout,
12807b31e0eSJeremy L Thompson                gradout, emodein, emodeout, assembledqfarray, elemdiagarray);
12907b31e0eSJeremy L Thompson }
13007b31e0eSJeremy L Thompson 
13107b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
13207b31e0eSJeremy L Thompson // Linear point block diagonal
13307b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
13407b31e0eSJeremy L Thompson extern "C" __global__ void linearPointBlockDiagonal(const CeedInt nelem,
135*69bf922dSJeremy L Thompson     const CeedScalar *identity,
13607b31e0eSJeremy L Thompson     const CeedScalar *interpin, const CeedScalar *gradin,
13707b31e0eSJeremy L Thompson     const CeedScalar *interpout, const CeedScalar *gradout,
13807b31e0eSJeremy L Thompson     const CeedEvalMode *emodein, const CeedEvalMode *emodeout,
13907b31e0eSJeremy L Thompson     const CeedScalar *__restrict__ assembledqfarray,
14007b31e0eSJeremy L Thompson     CeedScalar *__restrict__ elemdiagarray) {
141*69bf922dSJeremy L Thompson   diagonalCore(nelem, true, identity, interpin, gradin, interpout,
14207b31e0eSJeremy L Thompson                gradout, emodein, emodeout, assembledqfarray, elemdiagarray);
14307b31e0eSJeremy L Thompson }
14407b31e0eSJeremy L Thompson 
14507b31e0eSJeremy L Thompson //------------------------------------------------------------------------------
146