Neko  0.8.1
A portable framework for high-order spectral element flow simulations
scalar_residual_cpu.f90
Go to the documentation of this file.
1 
5  use operators
6  implicit none
7  private
8 
10  type, public, extends(scalar_residual_t) :: scalar_residual_cpu_t
11  contains
13  procedure, nopass :: compute => scalar_residual_cpu_compute
14  end type scalar_residual_cpu_t
15 
16 contains
17 
31  subroutine scalar_residual_cpu_compute(Ax, s, s_res, f_Xh, c_Xh, msh, Xh, &
32  lambda, rhocp, bd, dt, n)
33  class(ax_t), intent(in) :: Ax
34  type(mesh_t), intent(inout) :: msh
35  type(space_t), intent(inout) :: Xh
36  type(field_t), intent(inout) :: s
37  type(field_t), intent(inout) :: s_res
38  type(field_t), intent(inout) :: f_Xh
39  type(coef_t), intent(inout) :: c_Xh
40  real(kind=rp), intent(in) :: lambda
41  real(kind=rp), intent(in) :: rhocp
42  real(kind=rp), intent(in) :: bd
43  real(kind=rp), intent(in) :: dt
44  integer, intent(in) :: n
45  integer :: i
46 
47  do i = 1, n
48  c_xh%h1(i,1,1,1) = lambda
49  ! todo :should not be just rho here.
50  ! Tim M. 2023-12-19: What is this todo?
51  c_xh%h2(i,1,1,1) = rhocp * (bd / dt)
52  end do
53  c_xh%ifh2 = .true.
54 
55  call ax%compute(s_res%x, s%x, c_xh, msh, xh)
56 
57  do i = 1, n
58  s_res%x(i,1,1,1) = (-s_res%x(i,1,1,1)) + f_xh%x(i,1,1,1)
59  end do
60 
61  end subroutine scalar_residual_cpu_compute
62 
63 end module scalar_residual_cpu
Gather-scatter.
Operators.
Definition: operators.f90:34
Residuals in the scalar equation (CPU version).
subroutine scalar_residual_cpu_compute(Ax, s, s_res, f_Xh, c_Xh, msh, Xh, lambda, rhocp, bd, dt, n)
Compute the residual.
Defines the residual for the scalar transport equation.
Abstract type to compute scalar residual.
Wrapper type for the routine to compute the scalar residual on the CPU.