Neko  0.9.99
A portable framework for high-order spectral element flow simulations
spectral_error.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
35  use num_types, only: rp
36  use field, only: field_t
37  use coefs, only: coef_t
38  use field_list, only: field_list_t
39  use math, only: rzero, copy
40  use file, only: file_t, file_free
41  use tensor, only: tnsr3d
42  use device_math, only: device_copy
43  use gather_scatter
44  use neko_config
45  use logger, only: neko_log
47  use comm, only: pe_rank
48  use utils, only: neko_fname_len, neko_error
49  use field_writer, only: field_writer_t
51  use json_module, only: json_file
53  use case, only: case_t
55 
56  use, intrinsic :: iso_c_binding
57  implicit none
58  private
59 
65  type, public, extends(simulation_component_t) :: spectral_error_t
67  type(field_t), pointer :: u => null()
68  type(field_t), pointer :: v => null()
69  type(field_t), pointer :: w => null()
71  type(field_t), pointer :: u_hat => null()
72  type(field_t), pointer :: v_hat => null()
73  type(field_t), pointer :: w_hat => null()
75  type(field_t) :: wk
77  real(kind=rp) :: seri_small = 1.e-14
79  real(kind=rp) :: seri_smallr = 1.e-10
81  real(kind=rp) :: seri_smallg = 1.e-5
83  real(kind=rp) :: seri_smalls = 0.2
85  integer :: seri_np = 4
86  integer :: seri_np_max = 4
88  integer :: seri_elr = 0
90  real(kind=rp), allocatable :: eind_u(:), eind_v(:), eind_w(:)
92  real(kind=rp), allocatable :: sig_u(:), sig_v(:), sig_w(:)
94  type(field_list_t) :: speri_l
96  type(file_t) :: mf_speri
98  type(field_writer_t) :: writer
99  contains
101  procedure, pass(this) :: init => spectral_error_init
103  procedure, pass(this) :: free => spectral_error_free
105  procedure, pass(this) :: compute_ => spectral_error_compute
107  procedure, pass(this) :: get_indicators => spectral_error_get_indicators
108 
109  end type spectral_error_t
110 
111 contains
112 
114  subroutine spectral_error_init(this, json, case)
115  class(spectral_error_t), intent(inout) :: this
116  type(json_file), intent(inout) :: json
117  class(case_t), intent(inout), target :: case
118 
119  character(len=20) :: fields(3)
120 
122  ! picks it up. Will also add those fields to the registry.
123  fields(1) = "u_hat"
124  fields(2) = "v_hat"
125  fields(3) = "w_hat"
126  call json%add("fields", fields)
127 
128  call this%init_base(json, case)
129  call this%writer%init(json, case)
130 
131  call spectral_error_init_from_attributes(this, case%fluid%c_Xh)
132 
133  end subroutine spectral_error_init
134 
137  subroutine spectral_error_init_from_attributes(this, coef)
138  class(spectral_error_t), intent(inout) :: this
139  type(coef_t), intent(in) :: coef
140  integer :: il, jl, aa
141  character(len=NEKO_FNAME_LEN) :: fname_speri
142 
143  this%u => neko_field_registry%get_field("u")
144  this%v => neko_field_registry%get_field("v")
145  this%w => neko_field_registry%get_field("w")
146  this%u_hat => neko_field_registry%get_field("u_hat")
147  this%v_hat => neko_field_registry%get_field("v_hat")
148  this%w_hat => neko_field_registry%get_field("w_hat")
149 
151  this%wk = this%u
152 
154  allocate(this%eind_u(coef%msh%nelv))
155  allocate(this%eind_v(coef%msh%nelv))
156  allocate(this%eind_w(coef%msh%nelv))
157  allocate(this%sig_u(coef%msh%nelv))
158  allocate(this%sig_v(coef%msh%nelv))
159  allocate(this%sig_w(coef%msh%nelv))
160 
162  associate(lx1 => coef%Xh%lx, ly1 => coef%Xh%ly, &
163  lz1 => coef%Xh%lz, &
164  seri_small => this%SERI_SMALL, &
165  seri_smallr => this%SERI_SMALLR, &
166  seri_smallg => this%SERI_SMALLG, &
167  seri_smalls => this%SERI_SMALLS, &
168  seri_np => this%SERI_NP, &
169  seri_np_max => this%SERI_NP_MAX, &
170  seri_elr => this%SERI_ELR &
171  )
172  ! correctness check
173  if (seri_np.gt.seri_np_max) then
174  call neko_log%message('SETI_NP greater than SERI_NP_MAX')
175  endif
176  il = seri_np+seri_elr
177  jl = min(lx1,ly1)
178  jl = min(jl,lz1)
179  if (il.gt.jl) then
180  call neko_log%message('SERI_NP+SERI_ELR greater than L?1')
181  endif
182  end associate
183 
185 
187  subroutine spectral_error_free(this)
188  class(spectral_error_t), intent(inout) :: this
189 
190  if(allocated(this%eind_u)) then
191  deallocate(this%eind_u)
192  end if
193 
194  if(allocated(this%eind_v)) then
195  deallocate(this%eind_v)
196  end if
197 
198  if(allocated(this%eind_w)) then
199  deallocate(this%eind_w)
200  end if
201 
202  if(allocated(this%sig_u)) then
203  deallocate(this%sig_u)
204  end if
205 
206  if(allocated(this%sig_w)) then
207  deallocate(this%sig_w)
208  end if
209 
210  if(allocated(this%sig_w)) then
211  deallocate(this%sig_w)
212  end if
213 
214  call this%wk%free()
215 
216  nullify(this%u)
217  nullify(this%v)
218  nullify(this%w)
219  nullify(this%u_hat)
220  nullify(this%v_hat)
221  nullify(this%w_hat)
222 
223  call this%writer%free()
224  call this%free_base()
225 
226  end subroutine spectral_error_free
227 
229  subroutine spectral_error_compute(this, t, tstep)
230  class(spectral_error_t), intent(inout) :: this
231  real(kind=rp), intent(in) :: t
232  integer, intent(in) :: tstep
233 
234  integer :: e, i, lx, ly, lz, nelv, n
235 
236 
237  call this%get_indicators(this%case%fluid%c_Xh)
238 
239  lx = this%u_hat%Xh%lx
240  ly = this%u_hat%Xh%ly
241  lz = this%u_hat%Xh%lz
242  nelv = this%u_hat%msh%nelv
243 
245  do e = 1,nelv
246  do i = 1,lx*ly*lx
247  this%u_hat%x(i,1,1,e) = this%eind_u(e)
248  this%v_hat%x(i,1,1,e) = this%eind_v(e)
249  this%w_hat%x(i,1,1,e) = this%eind_w(e)
250  end do
251  end do
252 
253  ! We need this copy to GPU since the field writer does an internal copy
254  ! GPU -> CPU before writing the field files. This overwrites the *_hat
255  ! host arrays that contain the values.
256  if (neko_bcknd_device .eq. 1) then
257  call device_memcpy(this%u_hat%x, this%u_hat%x_d, lx*ly*lz*nelv, &
258  host_to_device, sync = .true.)
259  call device_memcpy(this%v_hat%x, this%v_hat%x_d, lx*ly*lz*nelv, &
260  host_to_device, sync = .true.)
261  call device_memcpy(this%w_hat%x, this%w_hat%x_d, lx*ly*lz*nelv, &
262  host_to_device, sync = .true.)
263  end if
264 
265  end subroutine spectral_error_compute
266 
274  subroutine transform_to_spec_or_phys(u_hat, u, wk, coef, space)
275  type(field_t), intent(inout) :: u_hat
276  type(field_t), intent(inout) :: u
277  type(field_t), intent(inout) :: wk
278  type(coef_t), intent(inout) :: coef
279  character(len=4), intent(in) :: space
280  integer :: i, j, k, e, nxyz, nelv, n
281 
283  nxyz = coef%Xh%lx*coef%Xh%lx*coef%Xh%lx
284  nelv = coef%msh%nelv
285  n = nxyz*nelv
286 
288  if ((neko_bcknd_hip .eq. 1) .or. (neko_bcknd_cuda .eq. 1) .or. &
289  (neko_bcknd_opencl .eq. 1)) then
290  call device_copy(wk%x_d, u%x_d, n)
291  else
292  call copy(wk%x,u%x,n)
293  end if
294 
295  select case(space)
296  case('spec')
297  call tnsr3d(u_hat%x, coef%Xh%lx, wk%x, &
298  coef%Xh%lx,coef%Xh%vinv, &
299  coef%Xh%vinvt, coef%Xh%vinvt, nelv)
300  case('phys')
301  call tnsr3d(u_hat%x, coef%Xh%lx, wk%x, &
302  coef%Xh%lx,coef%Xh%v, &
303  coef%Xh%vt, coef%Xh%vt, nelv)
304  end select
305 
306  ! Synchronize
307  if ((neko_bcknd_hip .eq. 1) .or. (neko_bcknd_cuda .eq. 1) .or. &
308  (neko_bcknd_opencl .eq. 1)) then
309 
310  call device_memcpy(u_hat%x,u_hat%x_d, n, &
311  device_to_host, sync=.true.)
312  end if
313 
314  end subroutine transform_to_spec_or_phys
315 
318  subroutine spectral_error_get_indicators(this,coef)
319  class(spectral_error_t), intent(inout) :: this
320  type(coef_t), intent(inout) :: coef
321  integer :: i
322 
323  ! Generate the uvwhat field (legendre coeff)
324  call transform_to_spec_or_phys(this%u_hat, this%u, this%wk, coef, 'spec')
325  call transform_to_spec_or_phys(this%v_hat, this%v, this%wk, coef, 'spec')
326  call transform_to_spec_or_phys(this%w_hat, this%w, this%wk, coef, 'spec')
327 
328  ! Get the spectral error indicator
329  call calculate_indicators(this, coef, this%eind_u, this%sig_u, &
330  coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
331  this%u_hat%x)
332  call calculate_indicators(this, coef, this%eind_v, this%sig_v, &
333  coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
334  this%v_hat%x)
335  call calculate_indicators(this, coef, this%eind_w, this%sig_w, &
336  coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
337  this%w_hat%x)
338 
339  end subroutine spectral_error_get_indicators
340 
343  subroutine spectral_error_write(this, t)
344  class(spectral_error_t), intent(inout) :: this
345  real(kind=rp), intent(in) :: t
346 
347  integer i, e
348  integer lx, ly, lz, nelv
349 
353  call this%mf_speri%write(this%speri_l,t)
354 
355  end subroutine spectral_error_write
356 
366  subroutine calculate_indicators(this, coef, eind, sig, lnelt, LX1, LY1, LZ1, var)
367  type(spectral_error_t), intent(inout) :: this
368  type(coef_t), intent(in) :: coef
369  integer :: lnelt
370  integer :: LX1
371  integer :: LY1
372  integer :: LZ1
373  real(kind=rp) :: eind(lnelt)
374  real(kind=rp) :: sig(lnelt)
375  real(kind=rp) :: var(lx1,ly1,lz1,lnelt)
376 
377  real(kind=rp) :: xa(coef%Xh%lx,coef%Xh%ly,coef%Xh%lz)
378  real(kind=rp) :: xb(coef%Xh%lx,coef%Xh%ly,coef%Xh%lz)
379  integer :: i, e
380 
381  ! zero arrays
382  call rzero(eind,lnelt)
383  call rzero(sig,lnelt)
384 
385  ! Get the indicator
386  call speri_var(this, eind,sig,var,lnelt,xa,xb, lx1, ly1, lz1)
387 
388  end subroutine calculate_indicators
389 
390 
401  subroutine speri_var(this, est,sig,var,nell,xa,xb,LX1,LY1,LZ1)
402  type(spectral_error_t), intent(inout) :: this
403  integer :: nell
404  integer :: LX1
405  integer :: LY1
406  integer :: LZ1
407  real(kind=rp) :: est(nell)
408  real(kind=rp) :: sig(nell)
409  real(kind=rp) :: var(lx1,ly1,lz1,nell)
410  real(kind=rp) :: xa(lx1,ly1,lz1)
411  real(kind=rp) :: xb(lx1,ly1,lz1)
412 
414  integer :: il, jl, kl, ll, j_st, j_en, ii
416  real(kind=rp) :: coeff(lx1,ly1,lz1)
418  real(kind=rp) :: coef11
420  real(kind=rp) :: coefx(this%SERI_NP_MAX,ly1,lz1), &
421  coefy(this%SERI_NP_MAX,lx1,lz1), &
422  coefz(this%SERI_NP_MAX,lx1,ly1)
424  real(kind=rp) :: estx, esty, estz
426  real(kind=rp) :: sigx, sigy, sigz
427  real(kind=rp) :: third
428  parameter(third = 1.0/3.0)
429 
431  do il = 1,nell
434  do ii = 1, lx1*ly1*lz1
435  coeff(ii,1,1) = var(ii,1,1,il) * var(ii,1,1,il)
436  end do
437 
439  coef11 = coeff(1,1,1)
440 
442  if (coef11.ge.this%SERI_SMALL) then
447  j_st = max(1,lx1-this%SERI_NP+1-this%SERI_ELR)
448  j_en = max(1,lx1-this%SERI_ELR)
449  do ll = 1,lz1
450  do kl = 1,ly1
451  do jl = j_st,j_en
452  coefx(j_en-jl+1,kl,ll) = coeff(jl,kl,ll)
453  enddo
454  enddo
455  enddo
457  call speri_extrap(this,estx,sigx,coef11,coefx, &
458  j_st,j_en,ly1,lz1)
459 
463  j_st = max(1,ly1-this%SERI_NP+1-this%SERI_ELR)
464  j_en = max(1,ly1-this%SERI_ELR)
465  do ll = 1,lz1
466  do kl = j_st,j_en
467  do jl = 1,lx1
468  coefy(j_en-kl+1,jl,ll) = coeff(jl,kl,ll)
469  enddo
470  enddo
471  enddo
473  call speri_extrap(this, esty,sigy,coef11,coefy, &
474  j_st,j_en,lx1,lz1)
475 
479  j_st = max(1,lz1-this%SERI_NP+1-this%SERI_ELR)
480  j_en = max(1,lz1-this%SERI_ELR)
481  do ll = j_st,j_en
482  do kl = 1,ly1
483  do jl = 1,lx1
484  coefz(j_en-ll+1,jl,kl) = coeff(jl,kl,ll)
485  enddo
486  enddo
487  enddo
489  call speri_extrap(this, estz,sigz,coef11,coefz, &
490  j_st,j_en,lx1,ly1)
491 
493  est(il) = sqrt(estx + esty + estz)
494  sig(il) = third*(sigx + sigy + sigz)
495 
496  else
498  estx = 0.0
499  esty = 0.0
500  estz = 0.0
501  sigx = -1.0
502  sigy = -1.0
503  sigz = -1.0
505 
506  est(il) = 0.0
507  sig(il) = -1.0
508  endif
509 
510  end do
511 
512  end subroutine speri_var
513 
523  subroutine speri_extrap(this,estx,sigx,coef11,coef, &
524  ix_st,ix_en,nyl,nzl)
525  implicit none
526  type(spectral_error_t), intent(inout) :: this
528  integer :: ix_st,ix_en,nyl,nzl
530  real(kind=rp) :: coef(this%SERI_NP_MAX,nyl,nzl)
532  real(kind=rp) :: coef11
534  real(kind=rp) :: estx, sigx
535 
537  integer :: il, jl, kl, ll ! loop index
538  integer :: nsigt, pnr, nzlt
539  real(kind=rp) :: sigt, smallr, cmin, cmax, cnm, rtmp, rtmp2, rtmp3
540  real(kind=rp) :: sumtmp(4), cffl(this%SERI_NP_MAX)
541  real(kind=rp) :: stmp, estt, clog, ctmp, cave, erlog
542  logical :: cuse(this%seri_np_max)
543 
544  associate(seri_small => this%SERI_SMALL, &
545  seri_smallr => this%SERI_SMALLR, &
546  seri_smallg => this%SERI_SMALLG, &
547  seri_smalls => this%SERI_SMALLS, &
548  seri_np => this%SERI_NP, &
549  seri_np_max => this%SERI_NP_MAX, &
550  seri_elr => this%SERI_ELR &
551  )
552  ! initial values
553  estx = 0.0
554  sigx = -1.0
555 
556  ! relative cutoff
557  smallr = coef11*seri_smallr
558 
559  ! number of points
560  pnr = ix_en - ix_st +1
561 
562  ! to few points to interpolate
563  !if ((ix_en - ix_st).le.1) return
564 
565  ! for averaging, initial values
566  sigt = 0.0
567  nsigt = 0
568 
569  ! loop over all face points
570  nzlt = max(1,nzl - seri_elr) ! for 2D runs
571  do il=1,nzlt
572  ! weight
573  rtmp3 = 1.0/(2.0*(il-1)+1.0)
574  do jl=1,nyl - seri_elr
575 
576  ! find min and max coef along single row
577  cffl(1) = coef(1,jl,il)
578  cmin = cffl(1)
579  cmax = cmin
580  do kl =2,pnr
581  cffl(kl) = coef(kl,jl,il)
582  cmin = min(cmin,cffl(kl))
583  cmax = max(cmax,cffl(kl))
584  enddo
585 
586  ! are coefficients sufficiently big
587  if((cmin.gt.0.0).and.(cmax.gt.smallr)) then
588  ! mark array position we use in iterpolation
589  do kl =1,pnr
590  cuse(kl) = .true.
591  enddo
592  ! max n for polynomial order
593  cnm = real(ix_en)
594 
595  ! check if all the points should be taken into account
596  ! in original code by Catherine Mavriplis this part is written
597  ! for 4 points, so I place if statement first
598  if (pnr.eq.4) then
599  ! should we neglect last values
600  if ((cffl(1).lt.smallr).and. &
601  (cffl(2).lt.smallr)) then
602  if (cffl(3).lt.smallr) then
603  cuse(1) = .false.
604  cuse(2) = .false.
605  cnm = real(ix_en-2)
606  else
607  cuse(1) = .false.
608  cnm = real(ix_en-1)
609  endif
610  else
611  ! should we take stronger gradient
612  if ((cffl(1)/cffl(2).lt.seri_smallg).and. &
613  (cffl(3)/cffl(4).lt.seri_smallg)) then
614  cuse(1) = .false.
615  cuse(3) = .false.
616  cnm = real(ix_en-1)
617  elseif ((cffl(2)/cffl(1).lt.seri_smallg).and. &
618  (cffl(4)/cffl(3).lt.seri_smallg)) then
619  cuse(2) = .false.
620  cuse(4) = .false.
621  endif
622  endif
623  endif
624 
625  ! get sigma for given face point
626  do kl =1,4
627  sumtmp(kl) = 0.0
628  enddo
629  ! find new min and count number of points
630  cmin = cmax
631  cmax = 0.0
632  do kl =1,pnr
633  if(cuse(kl)) then
634  rtmp = real(ix_en-kl)
635  rtmp2 = log(cffl(kl))
636  sumtmp(1) = sumtmp(1) +rtmp2
637  sumtmp(2) = sumtmp(2) +rtmp
638  sumtmp(3) = sumtmp(3) +rtmp*rtmp
639  sumtmp(4) = sumtmp(4) +rtmp2*rtmp
640  ! find new min and count used points
641  cmin = min(cffl(kl),cmin)
642  cmax = cmax + 1.0
643  endif
644  enddo
645  ! decay rate along single row
646  stmp = (sumtmp(1)*sumtmp(2) - sumtmp(4)*cmax)/ &
647  (sumtmp(3)*cmax - sumtmp(2)*sumtmp(2))
648  ! for averaging
649  sigt = sigt + stmp
650  nsigt = nsigt + 1
651 
652  ! get error estimator depending on calculated decay rate
653  estt = 0.0
654  if (stmp.lt.seri_smalls) then
655  estt = cmin
656  else
657  ! get averaged constant in front of c*exp(-sig*n)
658  clog = (sumtmp(1)+stmp*sumtmp(2))/cmax
659  ctmp = exp(clog)
660  ! average exponent
661  cave = sumtmp(1)/cmax
662  ! check quality of approximation comparing is to the constant cave
663  do kl =1,2
664  sumtmp(kl) = 0.0
665  enddo
666  do kl =1,pnr
667  if(cuse(kl)) then
668  erlog = clog - stmp*real(ix_en-kl)
669  sumtmp(1) = sumtmp(1)+ &
670  (erlog-log(cffl(kl)))**2
671  sumtmp(2) = sumtmp(2)+ &
672  (erlog-cave)**2
673  endif
674  enddo
675  rtmp = 1.0 - sumtmp(1)/sumtmp(2)
676  if (rtmp.lt.seri_smalls) then
677  estt = cmin
678  else
679  ! last coefficient is not included in error estimator
680  estt = ctmp/stmp*exp(-stmp*cnm)
681  endif
682  endif
683  ! add contribution to error estimator; variable weight
684  estx = estx + estt/(2.0*(jl-1)+1.0)*rtmp3
685  endif ! if((cmin.gt.0.0).and.(cmax.gt.smallr))
686  enddo
687  enddo
688  ! constant weight
689  ! Multiplication by 4 in 2D / 8 in 3D
690  ! Normalization of the error by the volume of the reference element
691  ! which is equal to 4 in 2D / 8 in 3D
692  ! ==> Both operations cancel each other
693  estx = estx/(2.0*(ix_en-1)+1.0)
694 
695  ! final everaging
696  ! sigt = 2*sigma so we divide by 2
697  if (nsigt.gt.0) then
698  sigx = 0.5*sigt/nsigt
699  endif
700 
701  end associate
702 
703  end subroutine speri_extrap
704 
705 end module spectral_error
double real
Definition: device_config.h:12
Copy data between host and device (or device and device)
Definition: device.F90:51
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:54
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:45
Defines a simulation case.
Definition: case.f90:34
Coefficients.
Definition: coef.f90:34
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:28
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
Definition: device_math.F90:76
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public host_to_device
Definition: device.F90:47
integer, parameter, public device_to_host
Definition: device.F90:47
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Implements the field_writer_t type.
Defines a field.
Definition: field.f90:34
Module for file I/O operations.
Definition: file.f90:34
subroutine file_free(this)
File operation destructor.
Definition: file.f90:134
Gather-scatter.
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Logging routines.
Definition: log.f90:34
type(log_t), public neko_log
Global log stream.
Definition: log.f90:65
Definition: math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:239
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:195
Build configurations.
Definition: neko_config.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
Defines a function space.
Definition: space.f90:34
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(this, json, case)
Constructor.
subroutine spectral_error_compute(this, t, tstep)
Compute the spectral error indicator.
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 calculate_indicators(this, coef, eind, sig, lnelt, LX1, LY1, LZ1, var)
Wrapper for old fortran 77 subroutines.
subroutine spectral_error_init_from_attributes(this, coef)
Actual constructor.
subroutine speri_var(this, est, sig, var, nell, xa, xb, LX1, LY1, LZ1)
Calculate the indicator in a specified variable.
Tensor operations.
Definition: tensor.f90:61
subroutine, public tnsr3d(v, nv, u, nu, A, Bt, Ct, nelv)
Tensor product performed on nelv elements.
Definition: tensor.f90:223
Utilities.
Definition: utils.f90:35
integer, parameter, public neko_fname_len
Definition: utils.f90:40
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
field_list_t, To be able to group fields together
Definition: field_list.f90:13
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...
Definition: file.f90:54
Base abstract class for simulation components.
Provides tools to calculate the spectral error indicator.
#define max(a, b)
Definition: tensor.cu:40