52 use json_module,
only : json_file
58 use,
intrinsic :: iso_c_binding
79 real(kind=
rp) :: seri_small = 1.e-14
81 real(kind=
rp) :: seri_smallr = 1.e-10
83 real(kind=
rp) :: seri_smallg = 1.e-5
85 real(kind=
rp) :: seri_smalls = 0.2
87 integer :: seri_np = 4
88 integer :: seri_np_max = 4
90 integer :: seri_elr = 0
92 real(kind=
rp),
allocatable :: eind_u(:), eind_v(:), eind_w(:)
94 real(kind=
rp),
allocatable :: sig_u(:), sig_v(:), sig_w(:)
118 type(json_file),
intent(inout) :: json
119 class(
case_t),
intent(inout),
target :: case
121 character(len=:),
allocatable :: name
122 character(len=NEKO_VARNAME_LEN) :: fields(3)
131 call json%add(
"fields", fields)
133 call this%init_base(json,
case)
134 call this%writer%init(json,
case)
144 type(
coef_t),
intent(in) :: coef
145 integer :: il, jl, aa
146 character(len=NEKO_FNAME_LEN) :: fname_speri
159 allocate(this%eind_u(coef%msh%nelv))
160 allocate(this%eind_v(coef%msh%nelv))
161 allocate(this%eind_w(coef%msh%nelv))
162 allocate(this%sig_u(coef%msh%nelv))
163 allocate(this%sig_v(coef%msh%nelv))
164 allocate(this%sig_w(coef%msh%nelv))
167 associate(lx1 => coef%Xh%lx, ly1 => coef%Xh%ly, &
169 seri_small => this%SERI_SMALL, &
170 seri_smallr => this%SERI_SMALLR, &
171 seri_smallg => this%SERI_SMALLG, &
172 seri_smalls => this%SERI_SMALLS, &
173 seri_np => this%SERI_NP, &
174 seri_np_max => this%SERI_NP_MAX, &
175 seri_elr => this%SERI_ELR &
178 if (seri_np .gt. seri_np_max)
then
179 call neko_log%message(
'SETI_NP greater than SERI_NP_MAX')
181 il = seri_np + seri_elr
185 call neko_log%message(
'SERI_NP+SERI_ELR greater than L?1')
195 if (
allocated(this%eind_u))
then
196 deallocate(this%eind_u)
199 if (
allocated(this%eind_v))
then
200 deallocate(this%eind_v)
203 if (
allocated(this%eind_w))
then
204 deallocate(this%eind_w)
207 if (
allocated(this%sig_u))
then
208 deallocate(this%sig_u)
211 if (
allocated(this%sig_v))
then
212 deallocate(this%sig_v)
215 if (
allocated(this%sig_w))
then
216 deallocate(this%sig_w)
220 call this%speri_l%free()
229 call this%writer%free()
230 call this%free_base()
237 type(time_state_t),
intent(in) :: time
239 integer :: e, i, lx, ly, lz, nelv, n
242 call this%get_indicators(this%case%fluid%c_Xh)
244 lx = this%u_hat%Xh%lx
245 ly = this%u_hat%Xh%ly
246 lz = this%u_hat%Xh%lz
247 nelv = this%u_hat%msh%nelv
252 this%u_hat%x(i, 1, 1, e) = this%eind_u(e)
253 this%v_hat%x(i, 1, 1, e) = this%eind_v(e)
254 this%w_hat%x(i, 1, 1, e) = this%eind_w(e)
261 if (neko_bcknd_device .eq. 1)
then
262 call device_memcpy(this%u_hat%x, this%u_hat%x_d, lx*ly*lz*nelv, &
263 host_to_device, sync = .true.)
264 call device_memcpy(this%v_hat%x, this%v_hat%x_d, lx*ly*lz*nelv, &
265 host_to_device, sync = .true.)
266 call device_memcpy(this%w_hat%x, this%w_hat%x_d, lx*ly*lz*nelv, &
267 host_to_device, sync = .true.)
281 type(field_t),
intent(inout) :: u_hat
282 type(field_t),
intent(inout) :: u
283 type(field_t),
intent(inout) :: wk
284 type(coef_t),
intent(inout) :: coef
285 character(len=4),
intent(in) :: space
286 integer :: i, j, k, e, nxyz, nelv, n
289 nxyz = coef%Xh%lx * coef%Xh%lx * coef%Xh%lx
294 if ((neko_bcknd_hip .eq. 1) .or. (neko_bcknd_cuda .eq. 1) .or. &
295 (neko_bcknd_opencl .eq. 1))
then
296 call device_copy(wk%x_d, u%x_d, n)
298 call copy(wk%x, u%x, n)
303 call tnsr3d(u_hat%x, coef%Xh%lx, wk%x, &
304 coef%Xh%lx, coef%Xh%vinv, &
305 coef%Xh%vinvt, coef%Xh%vinvt, nelv)
307 call tnsr3d(u_hat%x, coef%Xh%lx, wk%x, &
308 coef%Xh%lx, coef%Xh%v, &
309 coef%Xh%vt, coef%Xh%vt, nelv)
313 if ((neko_bcknd_hip .eq. 1) .or. (neko_bcknd_cuda .eq. 1) .or. &
314 (neko_bcknd_opencl .eq. 1))
then
316 call device_memcpy(u_hat%x, u_hat%x_d, n, &
317 device_to_host, sync = .true.)
326 type(coef_t),
intent(inout) :: coef
336 coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
339 coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
342 coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
351 real(kind=rp),
intent(in) :: t
354 integer lx, ly, lz, nelv
359 call this%mf_speri%write(this%speri_l, t)
375 type(coef_t),
intent(in) :: coef
380 real(kind=rp) :: eind(lnelt)
381 real(kind=rp) :: sig(lnelt)
382 real(kind=rp) :: var(lx1, ly1, lz1, lnelt)
384 real(kind=rp) :: xa(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz)
385 real(kind=rp) :: xb(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz)
389 call rzero(eind, lnelt)
390 call rzero(sig, lnelt)
393 call speri_var(this, eind, sig, var, lnelt, xa, xb, lx1, ly1, lz1)
408 subroutine speri_var(this, est, sig, var, nell, xa, xb, LX1, LY1, LZ1)
414 real(kind=rp) :: est(nell)
415 real(kind=rp) :: sig(nell)
416 real(kind=rp) :: var(lx1, ly1, lz1, nell)
417 real(kind=rp) :: xa(lx1, ly1, lz1)
418 real(kind=rp) :: xb(lx1, ly1, lz1)
421 integer :: il, jl, kl, ll, j_st, j_en, ii
423 real(kind=rp) :: coeff(lx1, ly1, lz1)
425 real(kind=rp) :: coef11
427 real(kind=rp) :: coefx(this%SERI_NP_MAX, ly1, lz1), &
428 coefy(this%SERI_NP_MAX, lx1, lz1), &
429 coefz(this%SERI_NP_MAX, lx1, ly1)
431 real(kind=rp) :: estx, esty, estz
433 real(kind=rp) :: sigx, sigy, sigz
434 real(kind=rp) :: third
435 parameter(third = 1.0/3.0)
441 do ii = 1, lx1*ly1*lz1
442 coeff(ii, 1, 1) = var(ii, 1, 1, il) * var(ii, 1, 1, il)
446 coef11 = coeff(1, 1, 1)
449 if (coef11 .ge. this%SERI_SMALL)
then
454 j_st =
max(1, lx1 - this%SERI_NP + 1 - this%SERI_ELR)
455 j_en =
max(1, lx1 - this%SERI_ELR)
459 coefx(j_en - jl + 1, kl, ll) = coeff(jl, kl, ll)
465 j_st, j_en, ly1, lz1)
470 j_st =
max(1, ly1 - this%SERI_NP + 1 - this%SERI_ELR)
471 j_en =
max(1, ly1 - this%SERI_ELR)
475 coefy(j_en - kl + 1, jl, ll) = coeff(jl, kl, ll)
481 j_st, j_en, lx1, lz1)
486 j_st =
max(1, lz1 - this%SERI_NP + 1 - this%SERI_ELR)
487 j_en =
max(1, lz1 - this%SERI_ELR)
491 coefz(j_en - ll + 1, jl, kl) = coeff(jl, kl, ll)
497 j_st, j_en, lx1, ly1)
500 est(il) = sqrt(estx + esty + estz)
501 sig(il) = third*(sigx + sigy + sigz)
531 ix_st, ix_en, nyl, nzl)
535 integer :: ix_st, ix_en, nyl, nzl
537 real(kind=rp) :: coef(this%SERI_NP_MAX, nyl, nzl)
539 real(kind=rp) :: coef11
541 real(kind=rp) :: estx, sigx
544 integer :: il, jl, kl, ll
545 integer :: nsigt, pnr, nzlt
546 real(kind=rp) :: sigt, smallr, cmin, cmax, cnm, rtmp, rtmp2, rtmp3
547 real(kind=rp) :: sumtmp(4), cffl(this%SERI_NP_MAX)
548 real(kind=rp) :: stmp, estt, clog, ctmp, cave, erlog
549 logical :: cuse(this%seri_np_max)
551 associate(seri_small => this%SERI_SMALL, &
552 seri_smallr => this%SERI_SMALLR, &
553 seri_smallg => this%SERI_SMALLG, &
554 seri_smalls => this%SERI_SMALLS, &
555 seri_np => this%SERI_NP, &
556 seri_np_max => this%SERI_NP_MAX, &
557 seri_elr => this%SERI_ELR &
564 smallr = coef11*seri_smallr
567 pnr = ix_en - ix_st + 1
577 nzlt =
max(1, nzl - seri_elr)
580 rtmp3 = 1.0/(2.0*(il - 1) + 1.0)
581 do jl = 1, nyl - seri_elr
584 cffl(1) = coef(1, jl, il)
588 cffl(kl) = coef(kl, jl, il)
589 cmin = min(cmin, cffl(kl))
590 cmax =
max(cmax, cffl(kl))
594 if ((cmin .gt. 0.0) .and. (cmax .gt. smallr))
then
607 if ((cffl(1) .lt. smallr) .and. &
608 (cffl(2) .lt. smallr))
then
609 if (cffl(3) .lt. smallr)
then
612 cnm =
real(ix_en - 2)
615 cnm =
real(ix_en - 1)
619 if ((cffl(1)/cffl(2) .lt. seri_smallg) .and. &
620 (cffl(3)/cffl(4) .lt. seri_smallg))
then
623 cnm =
real(ix_en - 1)
624 elseif ((cffl(2)/cffl(1) .lt. seri_smallg) .and. &
625 (cffl(4)/cffl(3) .lt. seri_smallg))
then
641 rtmp =
real(ix_en - kl)
642 rtmp2 = log(cffl(kl))
643 sumtmp(1) = sumtmp(1) + rtmp2
644 sumtmp(2) = sumtmp(2) + rtmp
645 sumtmp(3) = sumtmp(3) + rtmp*rtmp
646 sumtmp(4) = sumtmp(4) + rtmp2*rtmp
648 cmin = min(cffl(kl), cmin)
653 stmp = (sumtmp(1)*sumtmp(2) - sumtmp(4)*cmax)/ &
654 (sumtmp(3)*cmax - sumtmp(2)*sumtmp(2))
661 if (stmp .lt. seri_smalls)
then
665 clog = (sumtmp(1)+stmp*sumtmp(2))/cmax
668 cave = sumtmp(1)/cmax
676 erlog = clog - stmp*
real(ix_en - kl)
677 sumtmp(1) = sumtmp(1) + &
678 (erlog - log(cffl(kl)))**2
679 sumtmp(2) = sumtmp(2) + &
683 rtmp = 1.0 - sumtmp(1)/sumtmp(2)
684 if (rtmp .lt. seri_smalls)
then
688 estt = ctmp/stmp*exp(-stmp*cnm)
692 estx = estx + estt/(2.0*(jl - 1) + 1.0)*rtmp3
701 estx = estx/(2.0*(ix_en - 1) + 1.0)
705 if (nsigt .gt. 0)
then
706 sigx = 0.5*sigt/nsigt
Copy data between host and device (or device and device)
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Defines a simulation case.
integer, public pe_rank
MPI rank.
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
integer, parameter, public device_to_host
Implements the field_writer_t type.
Module for file I/O operations.
subroutine file_free(this)
File operation destructor.
Utilities for retrieving parameters from the case files.
type(log_t), public neko_log
Global log stream.
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public rzero(a, n)
Zero a real vector.
integer, parameter neko_bcknd_hip
integer, parameter neko_bcknd_device
integer, parameter neko_bcknd_opencl
integer, parameter neko_bcknd_cuda
integer, parameter, public rp
Global precision used in computations.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine compute_(this, time)
Dummy compute function.
Defines a function space.
Implements type spectral_error_t.
subroutine speri_extrap(this, estx, sigx, coef11, coef, ix_st, ix_en, nyl, nzl)
Extrapolate the Legendre spectrum from the last points.
subroutine spectral_error_write(this, t)
Write error indicators in a field file.
subroutine spectral_error_free(this)
Destructor.
subroutine spectral_error_get_indicators(this, coef)
Transform and get the spectral error indicators.
subroutine spectral_error_init_from_components(this, coef)
Actual constructor.
subroutine spectral_error_init(this, json, case)
Constructor.
subroutine spectral_error_compute(this, time)
Compute the spectral error indicator.
subroutine speri_var(this, est, sig, var, nell, xa, xb, lx1, ly1, lz1)
Calculate the indicator in a specified variable.
subroutine calculate_indicators(this, coef, eind, sig, lnelt, lx1, ly1, lz1, var)
Wrapper for old fortran 77 subroutines.
subroutine transform_to_spec_or_phys(u_hat, u, wk, coef, space)
Transform a field u to u_hat into physical or spectral space the result of the transformation is in u...
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
Module with things related to the simulation time.
integer, parameter, public neko_fname_len
integer, parameter, public neko_varname_len
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
field_list_t, To be able to group fields together
A simulation component that writes a 3d field to a file.
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Base abstract class for simulation components.
Provides tools to calculate the spectral error indicator.
A struct that contains all info about the time, expand as needed.