Neko  0.9.99
A portable framework for high-order spectral element flow simulations
rhs_maker_cpu.f90
Go to the documentation of this file.
4  use field_series, only : field_series_t
5  use field, only : field_t
6  use num_types, only : rp, c_rp
8  implicit none
9  private
10 
11  type, public, extends(rhs_maker_sumab_t) :: rhs_maker_sumab_cpu_t
12  contains
13  procedure, nopass :: compute_fluid => rhs_maker_sumab_cpu
14  end type rhs_maker_sumab_cpu_t
15 
16  type, public, extends(rhs_maker_ext_t) :: rhs_maker_ext_cpu_t
17  contains
18  procedure, nopass :: compute_fluid => rhs_maker_ext_cpu
19  procedure, nopass :: compute_scalar => scalar_rhs_maker_ext_cpu
20  end type rhs_maker_ext_cpu_t
21 
22  type, public, extends(rhs_maker_bdf_t) :: rhs_maker_bdf_cpu_t
23  contains
24  procedure, nopass :: compute_fluid => rhs_maker_bdf_cpu
25  procedure, nopass :: compute_scalar => scalar_rhs_maker_bdf_cpu
26  end type rhs_maker_bdf_cpu_t
27 
28  type, public, extends(rhs_maker_oifs_t) :: rhs_maker_oifs_cpu_t
29  contains
30  procedure, nopass :: compute_fluid => rhs_maker_oifs_cpu
31  procedure, nopass :: compute_scalar => scalar_rhs_maker_oifs_cpu
32  end type rhs_maker_oifs_cpu_t
33 
34 contains
35 
36  subroutine rhs_maker_sumab_cpu(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, &
37  ab, nab)
38  type(field_t), intent(inout) :: u,v, w
39  type(field_t), intent(inout) :: uu, vv, ww
40  type(field_series_t), intent(inout) :: uulag, vvlag, wwlag
41  real(kind=rp), dimension(3), intent(in) :: ab
42  integer, intent(in) :: nab
43  integer :: i, n
44 
45  n = uu%dof%size()
46 
47  do concurrent(i = 1:n)
48  u%x(i,1,1,1) = ab(1) * uu%x(i,1,1,1) + ab(2) * uulag%lf(1)%x(i,1,1,1)
49  v%x(i,1,1,1) = ab(1) * vv%x(i,1,1,1) + ab(2) * vvlag%lf(1)%x(i,1,1,1)
50  w%x(i,1,1,1) = ab(1) * ww%x(i,1,1,1) + ab(2) * wwlag%lf(1)%x(i,1,1,1)
51  end do
52 
53  if (nab .eq. 3) then
54  do concurrent(i = 1:n)
55  u%x(i,1,1,1) = u%x(i,1,1,1) + ab(3) * uulag%lf(2)%x(i,1,1,1)
56  v%x(i,1,1,1) = v%x(i,1,1,1) + ab(3) * vvlag%lf(2)%x(i,1,1,1)
57  w%x(i,1,1,1) = w%x(i,1,1,1) + ab(3) * wwlag%lf(2)%x(i,1,1,1)
58  end do
59  end if
60 
61  end subroutine rhs_maker_sumab_cpu
62 
63  subroutine rhs_maker_ext_cpu(fx_lag, fy_lag, fz_lag, &
64  fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, &
65  rho, ext_coeffs, n)
66  type(field_t), intent(inout) :: fx_lag, fy_lag, fz_lag
67  type(field_t), intent(inout) :: fx_laglag, fy_laglag, fz_laglag
68  real(kind=rp), intent(inout) :: rho, ext_coeffs(4)
69  integer, intent(in) :: n
70  real(kind=rp), intent(inout) :: fx(n), fy(n), fz(n)
71  integer :: i
72  type(field_t), pointer :: temp1, temp2, temp3
73  integer :: temp_indices(3)
74 
75  call neko_scratch_registry%request_field(temp1, temp_indices(1))
76  call neko_scratch_registry%request_field(temp2, temp_indices(2))
77  call neko_scratch_registry%request_field(temp3, temp_indices(3))
78 
79  do concurrent(i = 1:n)
80  temp1%x(i,1,1,1) = ext_coeffs(2) * fx_lag%x(i,1,1,1) + &
81  ext_coeffs(3) * fx_laglag%x(i,1,1,1)
82  temp2%x(i,1,1,1) = ext_coeffs(2) * fy_lag%x(i,1,1,1) + &
83  ext_coeffs(3) * fy_laglag%x(i,1,1,1)
84  temp3%x(i,1,1,1) = ext_coeffs(2) * fz_lag%x(i,1,1,1) + &
85  ext_coeffs(3) * fz_laglag%x(i,1,1,1)
86  end do
87 
88  do concurrent(i = 1:n)
89  fx_laglag%x(i,1,1,1) = fx_lag%x(i,1,1,1)
90  fy_laglag%x(i,1,1,1) = fy_lag%x(i,1,1,1)
91  fz_laglag%x(i,1,1,1) = fz_lag%x(i,1,1,1)
92  fx_lag%x(i,1,1,1) = fx(i)
93  fy_lag%x(i,1,1,1) = fy(i)
94  fz_lag%x(i,1,1,1) = fz(i)
95  end do
96 
97  do concurrent(i = 1:n)
98  fx(i) = (ext_coeffs(1) * fx(i) + temp1%x(i,1,1,1)) * rho
99  fy(i) = (ext_coeffs(1) * fy(i) + temp2%x(i,1,1,1)) * rho
100  fz(i) = (ext_coeffs(1) * fz(i) + temp3%x(i,1,1,1)) * rho
101  end do
102 
103  call neko_scratch_registry%relinquish_field(temp_indices)
104 
105  end subroutine rhs_maker_ext_cpu
106 
107  subroutine scalar_rhs_maker_ext_cpu(fs_lag, fs_laglag, fs, rho, &
108  ext_coeffs, n)
109  type(field_t), intent(inout) :: fs_lag
110  type(field_t), intent(inout) :: fs_laglag
111  real(kind=rp), intent(inout) :: rho, ext_coeffs(4)
112  integer, intent(in) :: n
113  real(kind=rp), intent(inout) :: fs(n)
114  integer :: i
115  type(field_t), pointer :: temp1
116  integer :: temp_index
117 
118  call neko_scratch_registry%request_field(temp1, temp_index)
119 
120  do concurrent(i = 1:n)
121  temp1%x(i,1,1,1) = ext_coeffs(2) * fs_lag%x(i,1,1,1) + &
122  ext_coeffs(3) * fs_laglag%x(i,1,1,1)
123  end do
124 
125  do concurrent(i = 1:n)
126  fs_laglag%x(i,1,1,1) = fs_lag%x(i,1,1,1)
127  fs_lag%x(i,1,1,1) = fs(i)
128  end do
129 
130  do concurrent(i = 1:n)
131  fs(i) = (ext_coeffs(1) * fs(i) + temp1%x(i,1,1,1)) * rho
132  end do
133 
134  call neko_scratch_registry%relinquish_field(temp_index)
135  end subroutine scalar_rhs_maker_ext_cpu
136 
137  subroutine rhs_maker_bdf_cpu(ulag, vlag, wlag, bfx, bfy, bfz, &
138  u, v, w, B, rho, dt, bd, nbd, n)
139  integer, intent(in) :: n, nbd
140  type(field_t), intent(in) :: u, v, w
141  type(field_series_t), intent(in) :: ulag, vlag, wlag
142  real(kind=rp), intent(inout) :: bfx(n), bfy(n), bfz(n)
143  real(kind=rp), intent(in) :: b(n)
144  real(kind=rp), intent(in) :: dt, rho, bd(4)
145  type(field_t), pointer :: tb1, tb2, tb3
146  type(field_t), pointer :: ta1, ta2, ta3
147  integer :: temp_indices(6)
148  integer :: i, ilag
149 
150  call neko_scratch_registry%request_field(ta1, temp_indices(1))
151  call neko_scratch_registry%request_field(ta2, temp_indices(2))
152  call neko_scratch_registry%request_field(ta3, temp_indices(3))
153  call neko_scratch_registry%request_field(tb1, temp_indices(4))
154  call neko_scratch_registry%request_field(tb2, temp_indices(5))
155  call neko_scratch_registry%request_field(tb3, temp_indices(6))
156 
157  do concurrent(i = 1:n)
158  tb1%x(i,1,1,1) = u%x(i,1,1,1) * b(i) * bd(2)
159  tb2%x(i,1,1,1) = v%x(i,1,1,1) * b(i) * bd(2)
160  tb3%x(i,1,1,1) = w%x(i,1,1,1) * b(i) * bd(2)
161  end do
162 
163  do ilag = 2, nbd
164  do concurrent(i = 1:n)
165  ta1%x(i,1,1,1) = ulag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
166  ta2%x(i,1,1,1) = vlag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
167  ta3%x(i,1,1,1) = wlag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
168  end do
169 
170  do concurrent(i = 1:n)
171  tb1%x(i,1,1,1) = tb1%x(i,1,1,1) + ta1%x(i,1,1,1)
172  tb2%x(i,1,1,1) = tb2%x(i,1,1,1) + ta2%x(i,1,1,1)
173  tb3%x(i,1,1,1) = tb3%x(i,1,1,1) + ta3%x(i,1,1,1)
174  end do
175  end do
176 
177  do concurrent(i = 1:n)
178  bfx(i) = bfx(i) + tb1%x(i,1,1,1) * (rho / dt)
179  bfy(i) = bfy(i) + tb2%x(i,1,1,1) * (rho / dt)
180  bfz(i) = bfz(i) + tb3%x(i,1,1,1) * (rho / dt)
181  end do
182 
183  call neko_scratch_registry%relinquish_field(temp_indices)
184 
185  end subroutine rhs_maker_bdf_cpu
186 
187  subroutine scalar_rhs_maker_bdf_cpu(s_lag, fs, s, B, rho, dt, bd, nbd, n)
188  integer, intent(in) :: n, nbd
189  type(field_t), intent(in) :: s
190  type(field_series_t), intent(in) :: s_lag
191  real(kind=rp), intent(inout) :: fs(n)
192  real(kind=rp), intent(in) :: b(n)
193  real(kind=rp), intent(in) :: dt, rho, bd(4)
194  integer :: i, ilag
195  type(field_t), pointer :: temp1, temp2
196  integer :: temp_indices(2)
197 
198  call neko_scratch_registry%request_field(temp1, temp_indices(1))
199  call neko_scratch_registry%request_field(temp2, temp_indices(2))
200 
201  do concurrent(i = 1:n)
202  temp2%x(i,1,1,1) = s%x(i,1,1,1) * b(i) * bd(2)
203  end do
204 
205  do ilag = 2, nbd
206  do concurrent(i = 1:n)
207  temp1%x(i,1,1,1) = s_lag%lf(ilag-1)%x(i,1,1,1) * b(i) * bd(ilag+1)
208  end do
209 
210  do concurrent(i = 1:n)
211  temp2%x(i,1,1,1) = temp2%x(i,1,1,1) + temp1%x(i,1,1,1)
212  end do
213  end do
214 
215  do concurrent(i = 1:n)
216  fs(i) = fs(i) + temp2%x(i,1,1,1) * (rho / dt)
217  end do
218 
219  call neko_scratch_registry%relinquish_field(temp_indices)
220  end subroutine scalar_rhs_maker_bdf_cpu
221 
222  subroutine rhs_maker_oifs_cpu(phi_x, phi_y, phi_z, bf_x, bf_y, bf_z, &
223  rho, dt, n)
224  real(kind=rp), intent(in) :: rho, dt
225  integer, intent(in) :: n
226  real(kind=rp), intent(inout) :: bf_x(n), bf_y(n), bf_z(n)
227  real(kind=rp), intent(inout) :: phi_x(n), phi_y(n), phi_z(n)
228  integer :: i
229 
230  do concurrent(i = 1:n)
231  bf_x(i) = bf_x(i) + phi_x(i) * (rho / dt)
232  bf_y(i) = bf_y(i) + phi_y(i) * (rho / dt)
233  bf_z(i) = bf_z(i) + phi_z(i) * (rho / dt)
234  end do
235 
236  end subroutine rhs_maker_oifs_cpu
237 
238  subroutine scalar_rhs_maker_oifs_cpu(phi_s, bf_s, rho, dt, n)
239  real(kind=rp), intent(in) :: rho, dt
240  integer, intent(in) :: n
241  real(kind=rp), intent(inout) :: bf_s(n)
242  real(kind=rp), intent(inout) :: phi_s(n)
243  integer :: i
244 
245  do concurrent(i = 1:n)
246  bf_s(i) = bf_s(i) + phi_s(i) * (rho / dt)
247  end do
248 
249  end subroutine scalar_rhs_maker_oifs_cpu
250 
251 end module rhs_maker_cpu
252 
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 scalar_rhs_maker_ext_cpu(fs_lag, fs_laglag, fs, rho, ext_coeffs, n)
subroutine rhs_maker_sumab_cpu(u, v, w, uu, vv, ww, uulag, vvlag, wwlag, ab, nab)
subroutine scalar_rhs_maker_oifs_cpu(phi_s, bf_s, rho, dt, n)
subroutine rhs_maker_bdf_cpu(ulag, vlag, wlag, bfx, bfy, bfz, u, v, w, B, rho, dt, bd, nbd, n)
subroutine scalar_rhs_maker_bdf_cpu(s_lag, fs, s, B, rho, dt, bd, nbd, n)
subroutine rhs_maker_oifs_cpu(phi_x, phi_y, phi_z, bf_x, bf_y, bf_z, rho, dt, n)
subroutine rhs_maker_ext_cpu(fx_lag, fy_lag, fz_lag, fx_laglag, fy_laglag, fz_laglag, fx, fy, fz, 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 add contributions of kth order OIFS scheme.
Definition: rhs_maker.f90:66
Abstract type to compute extrapolated velocity field for the pressure equation.
Definition: rhs_maker.f90:46