268 subroutine les_model_compute_delta(this)
269 class(les_model_t),
intent(inout) :: this
270 integer :: e, i, j, k
271 integer :: im, ip, jm, jp, km, kp
272 real(kind=rp) :: di, dj, dk, ndim_inv, volume_element
273 integer :: lx_half, ly_half, lz_half
274 character(len=:),
allocatable :: type_string
276 lx_half = this%coef%Xh%lx / 2
277 ly_half = this%coef%Xh%ly / 2
278 lz_half = this%coef%Xh%lz / 2
280 if (this%delta_type .eq.
"elementwise_max")
then
283 do e = 1, this%coef%msh%nelv
284 di = (this%coef%dof%x(lx_half, 1, 1, e) &
285 - this%coef%dof%x(lx_half + 1, 1, 1, e))**2 &
286 + (this%coef%dof%y(lx_half, 1, 1, e) &
287 - this%coef%dof%y(lx_half + 1, 1, 1, e))**2 &
288 + (this%coef%dof%z(lx_half, 1, 1, e) &
289 - this%coef%dof%z(lx_half + 1, 1, 1, e))**2
291 dj = (this%coef%dof%x(1, ly_half, 1, e) &
292 - this%coef%dof%x(1, ly_half + 1, 1, e))**2 &
293 + (this%coef%dof%y(1, ly_half, 1, e) &
294 - this%coef%dof%y(1, ly_half + 1, 1, e))**2 &
295 + (this%coef%dof%z(1, ly_half, 1, e) &
296 - this%coef%dof%z(1, ly_half + 1, 1, e))**2
298 dk = (this%coef%dof%x(1, 1, lz_half, e) &
299 - this%coef%dof%x(1, 1, lz_half + 1, e))**2 &
300 + (this%coef%dof%y(1, 1, lz_half, e) &
301 - this%coef%dof%y(1, 1, lz_half + 1, e))**2 &
302 + (this%coef%dof%z(1, 1, lz_half, e) &
303 - this%coef%dof%z(1, 1, lz_half + 1, e))**2
307 this%delta%x(:,:,:,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
309 else if (this%delta_type .eq.
"elementwise_average")
then
312 do e = 1, this%coef%msh%nelv
313 volume_element = 0.0_rp
314 do k = 1, this%coef%Xh%lx * this%coef%Xh%ly * this%coef%Xh%lz
315 volume_element = volume_element + this%coef%B(k, 1, 1, e)
317 this%delta%x(:,:,:,e) = (volume_element / this%coef%Xh%lx &
318 / this%coef%Xh%ly / this%coef%Xh%lz)**(1.0_rp / 3.0_rp)
320 else if (this%delta_type .eq.
"pointwise")
then
321 do e = 1, this%coef%msh%nelv
322 do k = 1, this%coef%Xh%lz
324 kp = min(this%coef%Xh%lz, k+1)
326 do j = 1, this%coef%Xh%ly
328 jp = min(this%coef%Xh%ly, j+1)
330 do i = 1, this%coef%Xh%lx
332 ip = min(this%coef%Xh%lx, i+1)
334 di = (this%coef%dof%x(ip, j, k, e) - &
335 this%coef%dof%x(im, j, k, e))**2 &
336 + (this%coef%dof%y(ip, j, k, e) - &
337 this%coef%dof%y(im, j, k, e))**2 &
338 + (this%coef%dof%z(ip, j, k, e) - &
339 this%coef%dof%z(im, j, k, e))**2
341 dj = (this%coef%dof%x(i, jp, k, e) - &
342 this%coef%dof%x(i, jm, k, e))**2 &
343 + (this%coef%dof%y(i, jp, k, e) - &
344 this%coef%dof%y(i, jm, k, e))**2 &
345 + (this%coef%dof%z(i, jp, k, e) - &
346 this%coef%dof%z(i, jm, k, e))**2
348 dk = (this%coef%dof%x(i, j, kp, e) - &
349 this%coef%dof%x(i, j, km, e))**2 &
350 + (this%coef%dof%y(i, j, kp, e) - &
351 this%coef%dof%y(i, j, km, e))**2 &
352 + (this%coef%dof%z(i, j, kp, e) - &
353 this%coef%dof%z(i, j, km, e))**2
355 di = sqrt(di) / (ip - im)
356 dj = sqrt(dj) / (jp - jm)
357 dk = sqrt(dk) / (kp - km)
358 this%delta%x(i,j,k,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
365 call neko_type_error(
"delta_type for LES model", &
366 this%delta_type, delta_known_types)
370 if (neko_bcknd_device .eq. 1)
then
371 call device_memcpy(this%delta%x, this%delta%x_d, this%delta%dof%size(),&
372 host_to_device, sync = .false.)
373 call this%coef%gs_h%op(this%delta%x, this%delta%dof%size(), gs_op_add)
374 call device_col2(this%delta%x_d, this%coef%mult_d, this%delta%dof%size())
376 call this%coef%gs_h%op(this%delta%x, this%delta%dof%size(), gs_op_add)
377 call col2(this%delta%x, this%coef%mult, this%delta%dof%size())