39  module subroutine opr_cpu_dudxyz(du, u, dr, ds, dt, coef)
 
   40    type(coef_t), 
intent(in), 
target :: coef
 
   41    real(kind=rp), 
dimension(coef%Xh%lx, coef%Xh%ly, &
         coef%Xh%lz, coef%msh%nelv), 
intent(inout) ::  du
 
   42    real(kind=rp), 
dimension(coef%Xh%lx, coef%Xh%ly, &
         coef%Xh%lz, coef%msh%nelv), 
intent(in) ::  u, dr, ds, dt
 
   44    associate(xh => coef%Xh, msh => coef%msh, dof => coef%dof)
 
   45      select case (coef%Xh%lx)
 
   47         call cpu_dudxyz_lx14(du, u, dr, ds, dt, &
 
   48              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   50         call cpu_dudxyz_lx13(du, u, dr, ds, dt, &
 
   51              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   53         call cpu_dudxyz_lx12(du, u, dr, ds, dt, &
 
   54              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   56         call cpu_dudxyz_lx11(du, u, dr, ds, dt, &
 
   57              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   59         call cpu_dudxyz_lx10(du, u, dr, ds, dt, &
 
   60              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   62         call cpu_dudxyz_lx9(du, u, dr, ds, dt, &
 
   63              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   65         call cpu_dudxyz_lx8(du, u, dr, ds, dt, &
 
   66              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   68         call cpu_dudxyz_lx7(du, u, dr, ds, dt, &
 
   69              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   71         call cpu_dudxyz_lx6(du, u, dr, ds, dt, &
 
   72              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   74         call cpu_dudxyz_lx5(du, u, dr, ds, dt, &
 
   75              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   77         call cpu_dudxyz_lx4(du, u, dr, ds, dt, &
 
   78              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   80         call cpu_dudxyz_lx3(du, u, dr, ds, dt, &
 
   81              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   83         call cpu_dudxyz_lx2(du, u, dr, ds, dt, &
 
   84              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv)
 
   86         call cpu_dudxyz_lx(du, u, dr, ds, dt, &
 
   87              xh%dx, xh%dy, xh%dz, coef%jacinv, msh%nelv, xh%lx)
 
   92  end subroutine opr_cpu_dudxyz
 
   94  subroutine cpu_dudxyz_lx(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel, lx)
 
   95    integer, 
intent(in) :: nel, lx
 
   96    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
   97    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
   98    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
   99    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  100    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  102    integer :: e, i, j, k, l
 
  109                tmp = tmp + dx(i,k) * u(k,j,1,e)
 
  115       do i = 1, lx * lx * lx
 
  116          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  124                   tmp = tmp + dy(j,l) * u(i,l,k,e)
 
  131       do i = 1, lx * lx * lx
 
  132          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  139                tmp = tmp + dz(k,l) * u(i,1,l,e)
 
  145       do i = 1, lx * lx * lx
 
  146          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  149       do i = 1, lx * lx * lx
 
  150          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  155  end subroutine cpu_dudxyz_lx
 
  157  subroutine cpu_dudxyz_lx14(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  158    integer, 
parameter :: lx = 14
 
  159    integer, 
intent(in) :: nel
 
  160    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  161    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  162    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  163    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  164    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  165    integer :: e, i, j, k
 
  170             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  171                         + dx(i,2) * u(2,j,1,e) &
 
  172                         + dx(i,3) * u(3,j,1,e) &
 
  173                         + dx(i,4) * u(4,j,1,e) &
 
  174                         + dx(i,5) * u(5,j,1,e) &
 
  175                         + dx(i,6) * u(6,j,1,e) &
 
  176                         + dx(i,7) * u(7,j,1,e) &
 
  177                         + dx(i,8) * u(8,j,1,e) &
 
  178                         + dx(i,9) * u(9,j,1,e) &
 
  179                         + dx(i,10) * u(10,j,1,e) &
 
  180                         + dx(i,11) * u(11,j,1,e) &
 
  181                         + dx(i,12) * u(12,j,1,e) &
 
  182                         + dx(i,13) * u(13,j,1,e) &
 
  183                         + dx(i,14) * u(14,j,1,e)
 
  187       do i = 1, lx * lx * lx
 
  188          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  194                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  195                            + dy(j,2) * u(i,2,k,e) &
 
  196                            + dy(j,3) * u(i,3,k,e) &
 
  197                            + dy(j,4) * u(i,4,k,e) &
 
  198                            + dy(j,5) * u(i,5,k,e) &
 
  199                            + dy(j,6) * u(i,6,k,e) &
 
  200                            + dy(j,7) * u(i,7,k,e) &
 
  201                            + dy(j,8) * u(i,8,k,e) &
 
  202                            + dy(j,9) * u(i,9,k,e) &
 
  203                            + dy(j,10) * u(i,10,k,e) &
 
  204                            + dy(j,11) * u(i,11,k,e) &
 
  205                            + dy(j,12) * u(i,12,k,e) &
 
  206                            + dy(j,13) * u(i,13,k,e) &
 
  207                            + dy(j,14) * u(i,14,k,e)
 
  212       do i = 1, lx * lx * lx
 
  213          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  218             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  219                         + dz(k,2) * u(i,1,2,e) &
 
  220                         + dz(k,3) * u(i,1,3,e) &
 
  221                         + dz(k,4) * u(i,1,4,e) &
 
  222                         + dz(k,5) * u(i,1,5,e) &
 
  223                         + dz(k,6) * u(i,1,6,e) &
 
  224                         + dz(k,7) * u(i,1,7,e) &
 
  225                         + dz(k,8) * u(i,1,8,e) &
 
  226                         + dz(k,9) * u(i,1,9,e) &
 
  227                         + dz(k,10) * u(i,1,10,e) &
 
  228                         + dz(k,11) * u(i,1,11,e) &
 
  229                         + dz(k,12) * u(i,1,12,e) &
 
  230                         + dz(k,13) * u(i,1,13,e) &
 
  231                         + dz(k,14) * u(i,1,14,e)
 
  235       do i = 1, lx * lx * lx
 
  236          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  239       do i = 1, lx * lx * lx
 
  240          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  245  end subroutine cpu_dudxyz_lx14
 
  247  subroutine cpu_dudxyz_lx13(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  248    integer, 
parameter :: lx = 13
 
  249    integer, 
intent(in) :: nel
 
  250    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  251    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  252    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  253    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  254    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  255    integer :: e, i, j, k
 
  260             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  261                         + dx(i,2) * u(2,j,1,e) &
 
  262                         + dx(i,3) * u(3,j,1,e) &
 
  263                         + dx(i,4) * u(4,j,1,e) &
 
  264                         + dx(i,5) * u(5,j,1,e) &
 
  265                         + dx(i,6) * u(6,j,1,e) &
 
  266                         + dx(i,7) * u(7,j,1,e) &
 
  267                         + dx(i,8) * u(8,j,1,e) &
 
  268                         + dx(i,9) * u(9,j,1,e) &
 
  269                         + dx(i,10) * u(10,j,1,e) &
 
  270                         + dx(i,11) * u(11,j,1,e) &
 
  271                         + dx(i,12) * u(12,j,1,e) &
 
  272                         + dx(i,13) * u(13,j,1,e)
 
  276       do i = 1, lx * lx * lx
 
  277          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  283                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  284                            + dy(j,2) * u(i,2,k,e) &
 
  285                            + dy(j,3) * u(i,3,k,e) &
 
  286                            + dy(j,4) * u(i,4,k,e) &
 
  287                            + dy(j,5) * u(i,5,k,e) &
 
  288                            + dy(j,6) * u(i,6,k,e) &
 
  289                            + dy(j,7) * u(i,7,k,e) &
 
  290                            + dy(j,8) * u(i,8,k,e) &
 
  291                            + dy(j,9) * u(i,9,k,e) &
 
  292                            + dy(j,10) * u(i,10,k,e) &
 
  293                            + dy(j,11) * u(i,11,k,e) &
 
  294                            + dy(j,12) * u(i,12,k,e) &
 
  295                            + dy(j,13) * u(i,13,k,e)
 
  300       do i = 1, lx * lx * lx
 
  301          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  306             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  307                         + dz(k,2) * u(i,1,2,e) &
 
  308                         + dz(k,3) * u(i,1,3,e) &
 
  309                         + dz(k,4) * u(i,1,4,e) &
 
  310                         + dz(k,5) * u(i,1,5,e) &
 
  311                         + dz(k,6) * u(i,1,6,e) &
 
  312                         + dz(k,7) * u(i,1,7,e) &
 
  313                         + dz(k,8) * u(i,1,8,e) &
 
  314                         + dz(k,9) * u(i,1,9,e) &
 
  315                         + dz(k,10) * u(i,1,10,e) &
 
  316                         + dz(k,11) * u(i,1,11,e) &
 
  317                         + dz(k,12) * u(i,1,12,e) &
 
  318                         + dz(k,13) * u(i,1,13,e)
 
  322       do i = 1, lx * lx * lx
 
  323          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  326       do i = 1, lx * lx * lx
 
  327          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  332  end subroutine cpu_dudxyz_lx13
 
  334  subroutine cpu_dudxyz_lx12(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  335    integer, 
parameter :: lx = 12
 
  336    integer, 
intent(in) :: nel
 
  337    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  338    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  339    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  340    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  341    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  342    integer :: e, i, j, k
 
  347             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  348                         + dx(i,2) * u(2,j,1,e) &
 
  349                         + dx(i,3) * u(3,j,1,e) &
 
  350                         + dx(i,4) * u(4,j,1,e) &
 
  351                         + dx(i,5) * u(5,j,1,e) &
 
  352                         + dx(i,6) * u(6,j,1,e) &
 
  353                         + dx(i,7) * u(7,j,1,e) &
 
  354                         + dx(i,8) * u(8,j,1,e) &
 
  355                         + dx(i,9) * u(9,j,1,e) &
 
  356                         + dx(i,10) * u(10,j,1,e) &
 
  357                         + dx(i,11) * u(11,j,1,e) &
 
  358                         + dx(i,12) * u(12,j,1,e)
 
  362       do i = 1, lx * lx * lx
 
  363          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  369                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  370                            + dy(j,2) * u(i,2,k,e) &
 
  371                            + dy(j,3) * u(i,3,k,e) &
 
  372                            + dy(j,4) * u(i,4,k,e) &
 
  373                            + dy(j,5) * u(i,5,k,e) &
 
  374                            + dy(j,6) * u(i,6,k,e) &
 
  375                            + dy(j,7) * u(i,7,k,e) &
 
  376                            + dy(j,8) * u(i,8,k,e) &
 
  377                            + dy(j,9) * u(i,9,k,e) &
 
  378                            + dy(j,10) * u(i,10,k,e) &
 
  379                            + dy(j,11) * u(i,11,k,e) &
 
  380                            + dy(j,12) * u(i,12,k,e)
 
  385       do i = 1, lx * lx * lx
 
  386          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  391             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  392                         + dz(k,2) * u(i,1,2,e) &
 
  393                         + dz(k,3) * u(i,1,3,e) &
 
  394                         + dz(k,4) * u(i,1,4,e) &
 
  395                         + dz(k,5) * u(i,1,5,e) &
 
  396                         + dz(k,6) * u(i,1,6,e) &
 
  397                         + dz(k,7) * u(i,1,7,e) &
 
  398                         + dz(k,8) * u(i,1,8,e) &
 
  399                         + dz(k,9) * u(i,1,9,e) &
 
  400                         + dz(k,10) * u(i,1,10,e) &
 
  401                         + dz(k,11) * u(i,1,11,e) &
 
  402                         + dz(k,12) * u(i,1,12,e)
 
  406       do i = 1, lx * lx * lx
 
  407          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  410       do i = 1, lx * lx * lx
 
  411          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  416  end subroutine cpu_dudxyz_lx12
 
  418  subroutine cpu_dudxyz_lx11(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  419    integer, 
parameter :: lx = 11
 
  420    integer, 
intent(in) :: nel
 
  421    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  422    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  423    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  424    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  425    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  426    integer :: e, i, j, k
 
  431             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  432                         + dx(i,2) * u(2,j,1,e) &
 
  433                         + dx(i,3) * u(3,j,1,e) &
 
  434                         + dx(i,4) * u(4,j,1,e) &
 
  435                         + dx(i,5) * u(5,j,1,e) &
 
  436                         + dx(i,6) * u(6,j,1,e) &
 
  437                         + dx(i,7) * u(7,j,1,e) &
 
  438                         + dx(i,8) * u(8,j,1,e) &
 
  439                         + dx(i,9) * u(9,j,1,e) &
 
  440                         + dx(i,10) * u(10,j,1,e) &
 
  441                         + dx(i,11) * u(11,j,1,e)
 
  445       do i = 1, lx * lx * lx
 
  446          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  452                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  453                            + dy(j,2) * u(i,2,k,e) &
 
  454                            + dy(j,3) * u(i,3,k,e) &
 
  455                            + dy(j,4) * u(i,4,k,e) &
 
  456                            + dy(j,5) * u(i,5,k,e) &
 
  457                            + dy(j,6) * u(i,6,k,e) &
 
  458                            + dy(j,7) * u(i,7,k,e) &
 
  459                            + dy(j,8) * u(i,8,k,e) &
 
  460                            + dy(j,9) * u(i,9,k,e) &
 
  461                            + dy(j,10) * u(i,10,k,e) &
 
  462                            + dy(j,11) * u(i,11,k,e)
 
  467       do i = 1, lx * lx * lx
 
  468          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  473             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  474                         + dz(k,2) * u(i,1,2,e) &
 
  475                         + dz(k,3) * u(i,1,3,e) &
 
  476                         + dz(k,4) * u(i,1,4,e) &
 
  477                         + dz(k,5) * u(i,1,5,e) &
 
  478                         + dz(k,6) * u(i,1,6,e) &
 
  479                         + dz(k,7) * u(i,1,7,e) &
 
  480                         + dz(k,8) * u(i,1,8,e) &
 
  481                         + dz(k,9) * u(i,1,9,e) &
 
  482                         + dz(k,10) * u(i,1,10,e) &
 
  483                         + dz(k,11) * u(i,1,11,e)
 
  487       do i = 1, lx * lx * lx
 
  488          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  491       do i = 1, lx * lx * lx
 
  492          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  497  end subroutine cpu_dudxyz_lx11
 
  499  subroutine cpu_dudxyz_lx10(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  500    integer, 
parameter :: lx = 10
 
  501    integer, 
intent(in) :: nel
 
  502    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  503    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  504    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  505    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  506    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  507    integer :: e, i, j, k
 
  512             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  513                         + dx(i,2) * u(2,j,1,e) &
 
  514                         + dx(i,3) * u(3,j,1,e) &
 
  515                         + dx(i,4) * u(4,j,1,e) &
 
  516                         + dx(i,5) * u(5,j,1,e) &
 
  517                         + dx(i,6) * u(6,j,1,e) &
 
  518                         + dx(i,7) * u(7,j,1,e) &
 
  519                         + dx(i,8) * u(8,j,1,e) &
 
  520                         + dx(i,9) * u(9,j,1,e) &
 
  521                         + dx(i,10) * u(10,j,1,e)
 
  525       do i = 1, lx * lx * lx
 
  526          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  532                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  533                            + dy(j,2) * u(i,2,k,e) &
 
  534                            + dy(j,3) * u(i,3,k,e) &
 
  535                            + dy(j,4) * u(i,4,k,e) &
 
  536                            + dy(j,5) * u(i,5,k,e) &
 
  537                            + dy(j,6) * u(i,6,k,e) &
 
  538                            + dy(j,7) * u(i,7,k,e) &
 
  539                            + dy(j,8) * u(i,8,k,e) &
 
  540                            + dy(j,9) * u(i,9,k,e) &
 
  541                            + dy(j,10) * u(i,10,k,e)
 
  546       do i = 1, lx * lx * lx
 
  547          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  552             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  553                         + dz(k,2) * u(i,1,2,e) &
 
  554                         + dz(k,3) * u(i,1,3,e) &
 
  555                         + dz(k,4) * u(i,1,4,e) &
 
  556                         + dz(k,5) * u(i,1,5,e) &
 
  557                         + dz(k,6) * u(i,1,6,e) &
 
  558                         + dz(k,7) * u(i,1,7,e) &
 
  559                         + dz(k,8) * u(i,1,8,e) &
 
  560                         + dz(k,9) * u(i,1,9,e) &
 
  561                         + dz(k,10) * u(i,1,10,e)
 
  565       do i = 1, lx * lx * lx
 
  566          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  569       do i = 1, lx * lx * lx
 
  570          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  575  end subroutine cpu_dudxyz_lx10
 
  577  subroutine cpu_dudxyz_lx9(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  578    integer, 
parameter :: lx = 9
 
  579    integer, 
intent(in) :: nel
 
  580    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  581    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  582    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  583    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  584    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  585    integer :: e, i, j, k
 
  590             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  591                         + dx(i,2) * u(2,j,1,e) &
 
  592                         + dx(i,3) * u(3,j,1,e) &
 
  593                         + dx(i,4) * u(4,j,1,e) &
 
  594                         + dx(i,5) * u(5,j,1,e) &
 
  595                         + dx(i,6) * u(6,j,1,e) &
 
  596                         + dx(i,7) * u(7,j,1,e) &
 
  597                         + dx(i,8) * u(8,j,1,e) &
 
  598                         + dx(i,9) * u(9,j,1,e)
 
  602       do i = 1, lx * lx * lx
 
  603          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  609                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  610                            + dy(j,2) * u(i,2,k,e) &
 
  611                            + dy(j,3) * u(i,3,k,e) &
 
  612                            + dy(j,4) * u(i,4,k,e) &
 
  613                            + dy(j,5) * u(i,5,k,e) &
 
  614                            + dy(j,6) * u(i,6,k,e) &
 
  615                            + dy(j,7) * u(i,7,k,e) &
 
  616                            + dy(j,8) * u(i,8,k,e) &
 
  617                            + dy(j,9) * u(i,9,k,e)
 
  622       do i = 1, lx * lx * lx
 
  623          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  628             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  629                         + dz(k,2) * u(i,1,2,e) &
 
  630                         + dz(k,3) * u(i,1,3,e) &
 
  631                         + dz(k,4) * u(i,1,4,e) &
 
  632                         + dz(k,5) * u(i,1,5,e) &
 
  633                         + dz(k,6) * u(i,1,6,e) &
 
  634                         + dz(k,7) * u(i,1,7,e) &
 
  635                         + dz(k,8) * u(i,1,8,e) &
 
  636                         + dz(k,9) * u(i,1,9,e)
 
  640       do i = 1, lx * lx * lx
 
  641          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  644       do i = 1, lx * lx * lx
 
  645          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  650  end subroutine cpu_dudxyz_lx9
 
  652  subroutine cpu_dudxyz_lx8(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  653    integer, 
parameter :: lx = 8
 
  654    integer, 
intent(in) :: nel
 
  655    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  656    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  657    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  658    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  659    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  660    integer :: e, i, j, k
 
  665             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  666                         + dx(i,2) * u(2,j,1,e) &
 
  667                         + dx(i,3) * u(3,j,1,e) &
 
  668                         + dx(i,4) * u(4,j,1,e) &
 
  669                         + dx(i,5) * u(5,j,1,e) &
 
  670                         + dx(i,6) * u(6,j,1,e) &
 
  671                         + dx(i,7) * u(7,j,1,e) &
 
  672                         + dx(i,8) * u(8,j,1,e)
 
  676       do i = 1, lx * lx * lx
 
  677          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  683                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  684                            + dy(j,2) * u(i,2,k,e) &
 
  685                            + dy(j,3) * u(i,3,k,e) &
 
  686                            + dy(j,4) * u(i,4,k,e) &
 
  687                            + dy(j,5) * u(i,5,k,e) &
 
  688                            + dy(j,6) * u(i,6,k,e) &
 
  689                            + dy(j,7) * u(i,7,k,e) &
 
  690                            + dy(j,8) * u(i,8,k,e)
 
  695       do i = 1, lx * lx * lx
 
  696          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  701             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  702                         + dz(k,2) * u(i,1,2,e) &
 
  703                         + dz(k,3) * u(i,1,3,e) &
 
  704                         + dz(k,4) * u(i,1,4,e) &
 
  705                         + dz(k,5) * u(i,1,5,e) &
 
  706                         + dz(k,6) * u(i,1,6,e) &
 
  707                         + dz(k,7) * u(i,1,7,e) &
 
  708                         + dz(k,8) * u(i,1,8,e)
 
  712       do i = 1, lx * lx * lx
 
  713          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  716       do i = 1, lx * lx * lx
 
  717          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  722  end subroutine cpu_dudxyz_lx8
 
  724  subroutine cpu_dudxyz_lx7(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  725    integer, 
parameter :: lx = 7
 
  726    integer, 
intent(in) :: nel
 
  727    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  728    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  729    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  730    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  731    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  732    integer :: e, i, j, k
 
  737             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  738                         + dx(i,2) * u(2,j,1,e) &
 
  739                         + dx(i,3) * u(3,j,1,e) &
 
  740                         + dx(i,4) * u(4,j,1,e) &
 
  741                         + dx(i,5) * u(5,j,1,e) &
 
  742                         + dx(i,6) * u(6,j,1,e) &
 
  743                         + dx(i,7) * u(7,j,1,e)
 
  747       do i = 1, lx * lx * lx
 
  748          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  754                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  755                            + dy(j,2) * u(i,2,k,e) &
 
  756                            + dy(j,3) * u(i,3,k,e) &
 
  757                            + dy(j,4) * u(i,4,k,e) &
 
  758                            + dy(j,5) * u(i,5,k,e) &
 
  759                            + dy(j,6) * u(i,6,k,e) &
 
  760                            + dy(j,7) * u(i,7,k,e)
 
  765       do i = 1, lx * lx * lx
 
  766          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  771             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  772                         + dz(k,2) * u(i,1,2,e) &
 
  773                         + dz(k,3) * u(i,1,3,e) &
 
  774                         + dz(k,4) * u(i,1,4,e) &
 
  775                         + dz(k,5) * u(i,1,5,e) &
 
  776                         + dz(k,6) * u(i,1,6,e) &
 
  777                         + dz(k,7) * u(i,1,7,e)
 
  781       do i = 1, lx * lx * lx
 
  782          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  785       do i = 1, lx * lx * lx
 
  786          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  791  end subroutine cpu_dudxyz_lx7
 
  793  subroutine cpu_dudxyz_lx6(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  794    integer, 
parameter :: lx = 6
 
  795    integer, 
intent(in) :: nel
 
  796    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  797    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  798    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  799    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  800    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  801    integer :: e, i, j, k
 
  806             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  807                         + dx(i,2) * u(2,j,1,e) &
 
  808                         + dx(i,3) * u(3,j,1,e) &
 
  809                         + dx(i,4) * u(4,j,1,e) &
 
  810                         + dx(i,5) * u(5,j,1,e) &
 
  811                         + dx(i,6) * u(6,j,1,e)
 
  815       do i = 1, lx * lx * lx
 
  816          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  822                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  823                            + dy(j,2) * u(i,2,k,e) &
 
  824                            + dy(j,3) * u(i,3,k,e) &
 
  825                            + dy(j,4) * u(i,4,k,e) &
 
  826                            + dy(j,5) * u(i,5,k,e) &
 
  827                            + dy(j,6) * u(i,6,k,e)
 
  832       do i = 1, lx * lx * lx
 
  833          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  838             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  839                         + dz(k,2) * u(i,1,2,e) &
 
  840                         + dz(k,3) * u(i,1,3,e) &
 
  841                         + dz(k,4) * u(i,1,4,e) &
 
  842                         + dz(k,5) * u(i,1,5,e) &
 
  843                         + dz(k,6) * u(i,1,6,e)
 
  847       do i = 1, lx * lx * lx
 
  848          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  851       do i = 1, lx * lx * lx
 
  852          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  857  end subroutine cpu_dudxyz_lx6
 
  859  subroutine cpu_dudxyz_lx5(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  860    integer, 
parameter :: lx = 5
 
  861    integer, 
intent(in) :: nel
 
  862    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  863    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  864    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  865    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  866    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  867    integer :: e, i, j, k
 
  872             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  873                         + dx(i,2) * u(2,j,1,e) &
 
  874                         + dx(i,3) * u(3,j,1,e) &
 
  875                         + dx(i,4) * u(4,j,1,e) &
 
  876                         + dx(i,5) * u(5,j,1,e)
 
  880       do i = 1, lx * lx * lx
 
  881          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  887                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  888                            + dy(j,2) * u(i,2,k,e) &
 
  889                            + dy(j,3) * u(i,3,k,e) &
 
  890                            + dy(j,4) * u(i,4,k,e) &
 
  891                            + dy(j,5) * u(i,5,k,e)
 
  896       do i = 1, lx * lx * lx
 
  897          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  902             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  903                         + dz(k,2) * u(i,1,2,e) &
 
  904                         + dz(k,3) * u(i,1,3,e) &
 
  905                         + dz(k,4) * u(i,1,4,e) &
 
  906                         + dz(k,5) * u(i,1,5,e)
 
  910       do i = 1, lx * lx * lx
 
  911          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  914       do i = 1, lx * lx * lx
 
  915          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  920  end subroutine cpu_dudxyz_lx5
 
  922  subroutine cpu_dudxyz_lx4(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  923    integer, 
parameter :: lx = 4
 
  924    integer, 
intent(in) :: nel
 
  925    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  926    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  927    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  928    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  929    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  930    integer :: e, i, j, k
 
  935             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  936                         + dx(i,2) * u(2,j,1,e) &
 
  937                         + dx(i,3) * u(3,j,1,e) &
 
  938                         + dx(i,4) * u(4,j,1,e)
 
  942       do i = 1, lx * lx * lx
 
  943          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
  949                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
  950                            + dy(j,2) * u(i,2,k,e) &
 
  951                            + dy(j,3) * u(i,3,k,e) &
 
  952                            + dy(j,4) * u(i,4,k,e)
 
  957       do i = 1, lx * lx * lx
 
  958          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
  963             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
  964                         + dz(k,2) * u(i,1,2,e) &
 
  965                         + dz(k,3) * u(i,1,3,e) &
 
  966                         + dz(k,4) * u(i,1,4,e)
 
  970       do i = 1, lx * lx * lx
 
  971          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
  974       do i = 1, lx * lx * lx
 
  975          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
  980  end subroutine cpu_dudxyz_lx4
 
  982  subroutine cpu_dudxyz_lx3(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
  983    integer, 
parameter :: lx = 3
 
  984    integer, 
intent(in) :: nel
 
  985    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
  986    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
  987    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
  988    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
  989    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
  990    integer :: e, i, j, k
 
  995             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
  996                         + dx(i,2) * u(2,j,1,e) &
 
  997                         + dx(i,3) * u(3,j,1,e)
 
 1001       do i = 1, lx * lx * lx
 
 1002          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
 1008                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
 1009                            + dy(j,2) * u(i,2,k,e) &
 
 1010                            + dy(j,3) * u(i,3,k,e)
 
 1015       do i = 1, lx * lx * lx
 
 1016          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
 1021             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
 1022                         + dz(k,2) * u(i,1,2,e) &
 
 1023                         + dz(k,3) * u(i,1,3,e)
 
 1027       do i = 1, lx * lx * lx
 
 1028          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
 1031       do i = 1, lx * lx * lx
 
 1032          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
 1037  end subroutine cpu_dudxyz_lx3
 
 1039  subroutine cpu_dudxyz_lx2(du, u, dr, ds, dt, dx, dy, dz, jacinv, nel)
 
 1040    integer, 
parameter :: lx = 2
 
 1041    integer, 
intent(in) :: nel
 
 1042    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(inout) ::  du
 
 1043    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) ::  u, dr, ds, dt
 
 1044    real(kind=rp), 
dimension(lx, lx, lx, nel), 
intent(in) :: jacinv
 
 1045    real(kind=rp), 
dimension(lx, lx), 
intent(in) :: dx, dy, dz
 
 1046    real(kind=rp), 
dimension(lx, lx, lx) :: drst
 
 1047    integer :: e, i, j, k
 
 1052             du(i,j,1,e) = dx(i,1) * u(1,j,1,e) &
 
 1053                         + dx(i,2) * u(2,j,1,e)
 
 1057       do i = 1, lx * lx * lx
 
 1058          du(i,1,1,e) = du(i,1,1,e) * dr(i,1,1,e)
 
 1064                drst(i,j,k) = dy(j,1) * u(i,1,k,e) &
 
 1065                            + dy(j,2) * u(i,2,k,e)
 
 1070       do i = 1, lx * lx * lx
 
 1071          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * ds(i,1,1,e)
 
 1076             drst(i,1,k) = dz(k,1) * u(i,1,1,e) &
 
 1077                         + dz(k,2) * u(i,1,2,e)
 
 1081       do i = 1, lx * lx * lx
 
 1082          du(i,1,1,e) = du(i,1,1,e) + drst(i,1,1) * dt(i,1,1,e)
 
 1085       do i = 1, lx * lx * lx
 
 1086          du(i,1,1,e) = du(i,1,1,e) * jacinv(i,1,1,e)
 
 1091  end subroutine cpu_dudxyz_lx2
 
 1093end submodule cpu_dudxyz