Neko  0.8.1
A portable framework for high-order spectral element flow simulations
rhs_maker_sx.f90
Go to the documentation of this file.
2  use rhs_maker
3  use field_series, only : field_series_t
4  use field, only : field_t
5  use num_types, only : rp, c_rp
7  implicit none
8  private
9 
10  type, public, extends(rhs_maker_sumab_t) :: rhs_maker_sumab_sx_t
11  contains
12  procedure, nopass :: compute_fluid => rhs_maker_sumab_sx
13  end type rhs_maker_sumab_sx_t
14 
15  type, public, extends(rhs_maker_ext_t) :: rhs_maker_ext_sx_t
16  contains
17  procedure, nopass :: compute_fluid => rhs_maker_ext_sx
18  procedure, nopass :: compute_scalar => scalar_rhs_maker_ext_sx
19  end type rhs_maker_ext_sx_t
20 
21  type, public, extends(rhs_maker_bdf_t) :: rhs_maker_bdf_sx_t
22  contains
23  procedure, nopass :: compute_fluid => rhs_maker_bdf_sx
24  procedure, nopass :: compute_scalar => scalar_rhs_maker_bdf_sx
25  end type rhs_maker_bdf_sx_t
26 
27 contains
28 
29  subroutine rhs_maker_sumab_sx(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, ab, nab)
30  type(field_t), intent(inout) :: u,v, w
31  type(field_t), intent(inout) :: uu, vv, ww
32  type(field_series_t), intent(inout) :: uulag, vvlag, wwlag
33  real(kind=rp), dimension(3), intent(in) :: ab
34  integer, intent(in) :: nab
35  integer :: i, n
36 
37  n = uu%dof%size()
38 
39  do i = 1, n
40  u%x(i,1,1,1) = ab(1) * uu%x(i,1,1,1) + ab(2) * uulag%lf(1)%x(i,1,1,1)
41  v%x(i,1,1,1) = ab(1) * vv%x(i,1,1,1) + ab(2) * vvlag%lf(1)%x(i,1,1,1)
42  w%x(i,1,1,1) = ab(1) * ww%x(i,1,1,1) + ab(2) * wwlag%lf(1)%x(i,1,1,1)
43  end do
44 
45  if (nab .eq. 3) then
46  do i = 1, n
47  u%x(i,1,1,1) = u%x(i,1,1,1) + ab(3) * uulag%lf(2)%x(i,1,1,1)
48  v%x(i,1,1,1) = v%x(i,1,1,1) + ab(3) * vvlag%lf(2)%x(i,1,1,1)
49  w%x(i,1,1,1) = w%x(i,1,1,1) + ab(3) * wwlag%lf(2)%x(i,1,1,1)
50  end do
51  end if
52 
53  end subroutine rhs_maker_sumab_sx
54 
55  subroutine rhs_maker_ext_sx(fx_lag, fy_lag, fz_lag, &
56  fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, &
57  rho, ext_coeffs, n)
58  type(field_t), intent(inout) :: fx_lag, fy_lag, fz_lag
59  type(field_t), intent(inout) :: fx_laglag, fy_laglag, fz_laglag
60  real(kind=rp), intent(inout) :: rho, ext_coeffs(4)
61  integer, intent(in) :: n
62  real(kind=rp), intent(inout) :: fx(n), fy(n), fz(n)
63  integer :: i
64  type(field_t), pointer :: temp1, temp2, temp3
65  integer :: temp_indices(3)
66 
67  call neko_scratch_registry%request_field(temp1, temp_indices(1))
68  call neko_scratch_registry%request_field(temp2, temp_indices(2))
69  call neko_scratch_registry%request_field(temp3, temp_indices(3))
70 
71  do i = 1, n
72  temp1%x(i,1,1,1) = ext_coeffs(2) * fx_lag%x(i,1,1,1) + &
73  ext_coeffs(3) * fx_laglag%x(i,1,1,1)
74  temp2%x(i,1,1,1) = ext_coeffs(2) * fy_lag%x(i,1,1,1) + &
75  ext_coeffs(3) * fy_laglag%x(i,1,1,1)
76  temp3%x(i,1,1,1) = ext_coeffs(2) * fz_lag%x(i,1,1,1) + &
77  ext_coeffs(3) * fz_laglag%x(i,1,1,1)
78  end do
79 
80  do i = 1, n
81  fx_laglag%x(i,1,1,1) = fx_lag%x(i,1,1,1)
82  fy_laglag%x(i,1,1,1) = fy_lag%x(i,1,1,1)
83  fz_laglag%x(i,1,1,1) = fz_lag%x(i,1,1,1)
84  fx_lag%x(i,1,1,1) = fx(i)
85  fy_lag%x(i,1,1,1) = fy(i)
86  fz_lag%x(i,1,1,1) = fz(i)
87  end do
88 
89  do i = 1, n
90  fx(i) = (ext_coeffs(1) * fx(i) + temp1%x(i,1,1,1)) * rho
91  fy(i) = (ext_coeffs(1) * fy(i) + temp2%x(i,1,1,1)) * rho
92  fz(i) = (ext_coeffs(1) * fz(i) + temp3%x(i,1,1,1)) * rho
93  end do
94 
95  call neko_scratch_registry%relinquish_field(temp_indices)
96 
97  end subroutine rhs_maker_ext_sx
98 
99  subroutine scalar_rhs_maker_ext_sx(fs_lag, fs_laglag, fs, rho, ext_coeffs, n)
100  type(field_t), intent(inout) :: fs_lag
101  type(field_t), intent(inout) :: fs_laglag
102  real(kind=rp), intent(inout) :: rho, ext_coeffs(4)
103  integer, intent(in) :: n
104  real(kind=rp), intent(inout) :: fs(n)
105  integer :: i
106  type(field_t), pointer :: temp1
107  integer :: temp_index
108 
109  call neko_scratch_registry%request_field(temp1, temp_index)
110 
111  do i = 1, n
112  temp1%x(i,1,1,1) = ext_coeffs(2) * fs_lag%x(i,1,1,1) + &
113  ext_coeffs(3) * fs_laglag%x(i,1,1,1)
114  end do
115 
116  do i = 1, n
117  fs_laglag%x(i,1,1,1) = fs_lag%x(i,1,1,1)
118  fs_lag%x(i,1,1,1) = fs(i)
119  end do
120 
121  do i = 1, n
122  fs(i) = (ext_coeffs(1) * fs(i) + temp1%x(i,1,1,1)) * rho
123  end do
124 
125  call neko_scratch_registry%relinquish_field(temp_index)
126  end subroutine scalar_rhs_maker_ext_sx
127 
128  subroutine rhs_maker_bdf_sx(ulag, vlag, wlag, bfx, bfy, bfz, &
129  u, v, w, B, rho, dt, bd, nbd, n)
130  integer, intent(in) :: n, nbd
131  type(field_t), intent(in) :: u, v, w
132  type(field_series_t), intent(in) :: ulag, vlag, wlag
133  real(kind=rp), intent(inout) :: bfx(n), bfy(n), bfz(n)
134  real(kind=rp), intent(in) :: b(n)
135  real(kind=rp), intent(in) :: dt, rho, bd(4)
136  type(field_t), pointer :: tb1, tb2, tb3
137  type(field_t), pointer :: ta1, ta2, ta3
138  integer :: temp_indices(6)
139  integer :: i, ilag
140 
141  call neko_scratch_registry%request_field(ta1, temp_indices(1))
142  call neko_scratch_registry%request_field(ta2, temp_indices(2))
143  call neko_scratch_registry%request_field(ta3, temp_indices(3))
144  call neko_scratch_registry%request_field(tb1, temp_indices(4))
145  call neko_scratch_registry%request_field(tb2, temp_indices(5))
146  call neko_scratch_registry%request_field(tb3, temp_indices(6))
147 
148  do i = 1, n
149  tb1%x(i,1,1,1) = u%x(i,1,1,1) * b(i) * bd(2)
150  tb2%x(i,1,1,1) = v%x(i,1,1,1) * b(i) * bd(2)
151  tb3%x(i,1,1,1) = w%x(i,1,1,1) * b(i) * bd(2)
152  end do
153 
154  do ilag = 2, nbd
155  do i = 1, n
156  ta1%x(i,1,1,1) = ulag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
157  ta2%x(i,1,1,1) = vlag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
158  ta3%x(i,1,1,1) = wlag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
159  end do
160 
161  do i = 1, n
162  tb1%x(i,1,1,1) = tb1%x(i,1,1,1) + ta1%x(i,1,1,1)
163  tb2%x(i,1,1,1) = tb2%x(i,1,1,1) + ta2%x(i,1,1,1)
164  tb3%x(i,1,1,1) = tb3%x(i,1,1,1) + ta3%x(i,1,1,1)
165  end do
166  end do
167 
168  do i = 1, n
169  bfx(i) = bfx(i) + tb1%x(i,1,1,1) * (rho / dt)
170  bfy(i) = bfy(i) + tb2%x(i,1,1,1) * (rho / dt)
171  bfz(i) = bfz(i) + tb3%x(i,1,1,1) * (rho / dt)
172  end do
173 
174  call neko_scratch_registry%relinquish_field(temp_indices)
175 
176  end subroutine rhs_maker_bdf_sx
177 
178  subroutine scalar_rhs_maker_bdf_sx(s_lag, fs, s, B, rho, dt, bd, nbd, n)
179  integer, intent(in) :: n, nbd
180  type(field_t), intent(in) :: s
181  type(field_series_t), intent(in) :: s_lag
182  real(kind=rp), intent(inout) :: fs(n)
183  real(kind=rp), intent(in) :: b(n)
184  real(kind=rp), intent(in) :: dt, rho, bd(4)
185  integer :: i, ilag
186  type(field_t), pointer :: temp1, temp2
187  integer :: temp_indices(2)
188 
189  call neko_scratch_registry%request_field(temp1, temp_indices(1))
190  call neko_scratch_registry%request_field(temp2, temp_indices(2))
191 
192  do i = 1, n
193  temp2%x(i,1,1,1) = s%x(i,1,1,1) * b(i) * bd(2)
194  end do
195 
196  do ilag = 2, nbd
197  do i = 1, n
198  temp1%x(i,1,1,1) = s_lag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
199  end do
200 
201  do i = 1, n
202  temp2%x(i,1,1,1) = temp2%x(i,1,1,1) + temp1%x(i,1,1,1)
203  end do
204  end do
205 
206  do i = 1, n
207  fs(i) = fs(i) + temp2%x(i,1,1,1) * (rho / dt)
208  end do
209 
210  call neko_scratch_registry%relinquish_field(temp_indices)
211  end subroutine scalar_rhs_maker_bdf_sx
212 
213 end module rhs_maker_sx
214 
Stores a series fields.
Defines a field.
Definition: field.f90:34
integer, parameter, public c_rp
Definition: num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
subroutine rhs_maker_bdf_sx(ulag, vlag, wlag, bfx, bfy, bfz, u, v, w, B, rho, dt, bd, nbd, n)
subroutine scalar_rhs_maker_bdf_sx(s_lag, fs, s, B, rho, dt, bd, nbd, n)
subroutine rhs_maker_ext_sx(fx_lag, fy_lag, fz_lag, fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, rho, ext_coeffs, n)
subroutine rhs_maker_sumab_sx(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, ab, nab)
subroutine scalar_rhs_maker_ext_sx(fs_lag, fs_laglag, fs, rho, ext_coeffs, n)
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Definition: rhs_maker.f90:38
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Abstract type to add contributions to F from lagged BD terms.
Definition: rhs_maker.f90:59
Abstract type to sum up contributions to kth order extrapolation scheme.
Definition: rhs_maker.f90:52
Abstract type to compute extrapolated velocity field for the pressure equation.
Definition: rhs_maker.f90:46