38 use json_module,
only : json_file
55 character(len=:),
allocatable :: delta_type
59 type(
coef_t),
pointer :: coef => null()
82 real(kind=
rp),
intent(in) :: t
83 integer,
intent(in) :: tstep
95 type(
coef_t),
intent(in) :: coef
97 type(json_file),
intent(inout) :: json
116 module subroutine les_model_factory(object, type_name,
dofmap, coef, json)
117 class(les_model_t),
allocatable,
intent(inout) :: object
118 character(len=*),
intent(in) :: type_name
119 type(
dofmap_t),
intent(in) :: dofmap
120 type(
coef_t),
intent(in) :: coef
121 type(json_file),
intent(inout) :: json
122 end subroutine les_model_factory
125 public :: les_model_factory
132 subroutine les_model_init_base(this, dofmap, coef, nut_name, delta_type)
133 class(les_model_t),
intent(inout) :: this
134 type(
dofmap_t),
intent(in) :: dofmap
135 type(
coef_t),
target,
intent(in) :: coef
136 character(len=*),
intent(in) :: nut_name
137 character(len=*),
intent(in) :: delta_type
148 this%delta_type = delta_type
150 call this%compute_delta()
151 end subroutine les_model_init_base
154 subroutine les_model_free_base(this)
155 class(les_model_t),
intent(inout) :: this
160 end subroutine les_model_free_base
168 subroutine les_model_compute_delta(this)
169 class(les_model_t),
intent(inout) :: this
170 integer :: e, i, j, k
171 integer :: im, ip, jm, jp, km, kp
172 real(kind=
rp) :: di, dj, dk, ndim_inv
173 integer :: lx_half, ly_half, lz_half
175 lx_half = this%coef%Xh%lx / 2
176 ly_half = this%coef%Xh%ly / 2
177 lz_half = this%coef%Xh%lz / 2
179 if (this%delta_type .eq.
"elementwise")
then
182 do e = 1, this%coef%msh%nelv
183 di = (this%coef%dof%x(lx_half, 1, 1, e) &
184 - this%coef%dof%x(lx_half + 1, 1, 1, e))**2 &
185 + (this%coef%dof%y(lx_half, 1, 1, e) &
186 - this%coef%dof%y(lx_half + 1, 1, 1, e))**2 &
187 + (this%coef%dof%z(lx_half, 1, 1, e) &
188 - this%coef%dof%z(lx_half + 1, 1, 1, e))**2
190 dj = (this%coef%dof%x(1, ly_half, 1, e) &
191 - this%coef%dof%x(1, ly_half + 1, 1, e))**2 &
192 + (this%coef%dof%y(1, ly_half, 1, e) &
193 - this%coef%dof%y(1, ly_half + 1, 1, e))**2 &
194 + (this%coef%dof%z(1, ly_half, 1, e) &
195 - this%coef%dof%z(1, ly_half + 1, 1, e))**2
197 dk = (this%coef%dof%x(1, 1, lz_half, e) &
198 - this%coef%dof%x(1, 1, lz_half + 1, e))**2 &
199 + (this%coef%dof%y(1, 1, lz_half, e) &
200 - this%coef%dof%y(1, 1, lz_half + 1, e))**2 &
201 + (this%coef%dof%z(1, 1, lz_half, e) &
202 - this%coef%dof%z(1, 1, lz_half + 1, e))**2
206 this%delta%x(:,:,:,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
209 else if (this%delta_type .eq.
"pointwise")
then
210 do e = 1, this%coef%msh%nelv
211 do k = 1, this%coef%Xh%lz
213 kp = min(this%coef%Xh%lz, k+1)
215 do j = 1, this%coef%Xh%ly
217 jp = min(this%coef%Xh%ly, j+1)
219 do i = 1, this%coef%Xh%lx
221 ip = min(this%coef%Xh%lx, i+1)
223 di = (this%coef%dof%x(ip, j, k, e) - &
224 this%coef%dof%x(im, j, k, e))**2 &
225 + (this%coef%dof%y(ip, j, k, e) - &
226 this%coef%dof%y(im, j, k, e))**2 &
227 + (this%coef%dof%z(ip, j, k, e) - &
228 this%coef%dof%z(im, j, k, e))**2
230 dj = (this%coef%dof%x(i, jp, k, e) - &
231 this%coef%dof%x(i, jm, k, e))**2 &
232 + (this%coef%dof%y(i, jp, k, e) - &
233 this%coef%dof%y(i, jm, k, e))**2 &
234 + (this%coef%dof%z(i, jp, k, e) - &
235 this%coef%dof%z(i, jm, k, e))**2
237 dk = (this%coef%dof%x(i, j, kp, e) - &
238 this%coef%dof%x(i, j, km, e))**2 &
239 + (this%coef%dof%y(i, j, kp, e) - &
240 this%coef%dof%y(i, j, km, e))**2 &
241 + (this%coef%dof%z(i, j, kp, e) - &
242 this%coef%dof%z(i, j, km, e))**2
244 di = sqrt(di) / (ip - im)
245 dj = sqrt(dj) / (jp - jm)
246 dk = sqrt(dk) / (kp - km)
247 this%delta%x(i,j,k,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
256 call device_memcpy(this%delta%x, this%delta%x_d, this%delta%dof%size(),&
258 call this%coef%gs_h%op(this%delta%x, this%delta%dof%size(),
gs_op_add)
259 call device_col2(this%delta%x_d, this%coef%mult_d, this%delta%dof%size())
261 call this%coef%gs_h%op(this%delta%x, this%delta%dof%size(),
gs_op_add)
262 call col2(this%delta%x, this%coef%mult, this%delta%dof%size())
265 end subroutine les_model_compute_delta
Copy data between host and device (or device and device)
subroutine, public device_col2(a_d, b_d, n)
Vector multiplication .
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
Defines a mapping of the degrees of freedom.
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines Gather-scatter operations.
integer, parameter, public gs_op_add
subroutine les_model_init_base(this, dofmap, coef, nut_name, delta_type)
LES model factory. Both constructs and initializes the object.
subroutine les_model_compute_delta(this)
Compute the LES lengthscale. For each GLL point, we take the distance between its neighbours in all 3...
subroutine les_model_free_base(this)
Destructor for the les_model_t (base) class.
subroutine, public col2(a, b, n)
Vector multiplication .
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
field_ptr_t, To easily obtain a pointer to a field
Base abstract type for LES models based on the Boussinesq approximation.