Neko  0.8.99
A portable framework for high-order spectral element flow simulations
spectral_error_indicator Module Reference

Implements type spectral_error_indicator_t.

Data Types

type  spectral_error_indicator_t
 Provides tools to calculate the spectral error indicator. More...
 

Functions/Subroutines

subroutine spec_err_ind_init (this, u, v, w, coef)
 Constructor. More...
 
subroutine spec_err_ind_free (this)
 Destructor. More...
 
subroutine transform_to_spec_or_phys (u_hat, u, wk, coef, space)
 Transform a field u > u_hat into physical or spectral space the result of the transformation is in u_hat. More...
 
subroutine spec_err_ind_get (this, coef)
 Transform and get the spectral error indicators. More...
 
subroutine spec_err_ind_write (this, t)
 Write error indicators in a field file. More...
 
subroutine calculate_indicators (this, coef, eind, sig, lnelt, LX1, LY1, LZ1, var)
 Wrapper for old fortran 77 subroutines. More...
 
subroutine speri_var (this, est, sig, var, nell, xa, xb, LX1, LY1, LZ1)
 Calculate the indicator in a specified variable. More...
 
subroutine speri_extrap (this, estx, sigx, coef11, coef, ix_st, ix_en, nyl, nzl)
 Extrapolate the Legendre spectrum from the last points. More...
 
subroutine list_init3 (list, uu, vv, ww)
 Helper function to initialize field list to write. More...
 
subroutine list_final3 (list)
 Helper function to finalize field list to write. More...
 

Function/Subroutine Documentation

◆ calculate_indicators()

subroutine spectral_error_indicator::calculate_indicators ( type(spectral_error_indicator_t), intent(inout)  this,
type(coef_t), intent(in)  coef,
real(kind=rp), dimension(lnelt)  eind,
real(kind=rp), dimension(lnelt)  sig,
integer  lnelt,
integer  LX1,
integer  LY1,
integer  LZ1,
real(kind=rp), dimension(lx1,ly1,lz1,lnelt)  var 
)
private
Parameters
coefcoef type
eindspectral indicator
sigcoefficient of the exponential fit
lneltnumber of elements
LX1gll points in x
LY1gll points in y
LZ1gll points in z @paran var variable to calculate indicator

Definition at line 327 of file spectral_error_indicator.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ list_final3()

subroutine spectral_error_indicator::list_final3 ( type(field_list_t), intent(inout)  list)
private
Parameters
listlist to deallocate

Definition at line 687 of file spectral_error_indicator.f90.

Here is the caller graph for this function:

◆ list_init3()

subroutine spectral_error_indicator::list_init3 ( type(field_list_t), intent(inout)  list,
type(field_t), target  uu,
type(field_t), target  vv,
type(field_t), target  ww 
)
private
Parameters
listlist to allocate
uufield to add to the list
vvfield to add to the list
wwfield to add to the list

Initialize field lists

Definition at line 672 of file spectral_error_indicator.f90.

Here is the caller graph for this function:

◆ spec_err_ind_free()

subroutine spectral_error_indicator::spec_err_ind_free ( class(spectral_error_indicator_t), intent(inout)  this)
private

finalize data related to writing

Definition at line 171 of file spectral_error_indicator.f90.

Here is the call graph for this function:

◆ spec_err_ind_get()

subroutine spectral_error_indicator::spec_err_ind_get ( class(spectral_error_indicator_t), intent(inout)  this,
type(coef_t), intent(inout)  coef 
)
private
Parameters
coeftype coef for mesh parameters and space

Definition at line 264 of file spectral_error_indicator.f90.

Here is the call graph for this function:

◆ spec_err_ind_init()

subroutine spectral_error_indicator::spec_err_ind_init ( class(spectral_error_indicator_t), intent(inout)  this,
type(field_t), intent(in), target  u,
type(field_t), intent(in), target  v,
type(field_t), intent(in), target  w,
type(coef_t), intent(in)  coef 
)
Parameters
uu velocity field
vv velocity field
ww velocity field
coeftype with all geometrical variables

call destructior

Assign the pointers

Initialize fields and copy data from proper one

Allocate arrays (Consider moving some to coef)

The following code has been lifted from Adams implementation

Initialize the list that holds the fields to write

Definition at line 108 of file spectral_error_indicator.f90.

Here is the call graph for this function:

◆ spec_err_ind_write()

subroutine spectral_error_indicator::spec_err_ind_write ( class(spectral_error_indicator_t), intent(inout)  this,
real(kind=rp), intent(in)  t 
)
private
Parameters
tCurrent simulation time.

Copy the element indicator into all points of the field

Write the file Remember that the list is already ponting to the fields that were just modified.

Definition at line 290 of file spectral_error_indicator.f90.

◆ speri_extrap()

subroutine spectral_error_indicator::speri_extrap ( type(spectral_error_indicator_t), intent(inout)  this,
real(kind=rp)  estx,
real(kind=rp)  sigx,
real(kind=rp)  coef11,
real(kind=rp), dimension(this%seri_np_max,nyl,nzl)  coef,
integer  ix_st,
integer  ix_en,
integer  nyl,
integer  nzl 
)
private
Parameters
estxspectral indicator
sigxcoefficient of the exponential fit
coef11legendre coefficients @paran coef legendre coefficients @paran ix_st argument list @paran ix_en argument list
nylargument list
nzlargument list
ix_stargument list
ix_enargument list
nylargument list
nzlargument list
coefLegendre coefficients; last SERI_NP columns
coef11Legendre coefficients; first value coeff(1,1,1)
estxestimated error and decay rate
sigxestimated error and decay rate

local variables

Definition at line 485 of file spectral_error_indicator.f90.

Here is the caller graph for this function:

◆ speri_var()

subroutine spectral_error_indicator::speri_var ( type(spectral_error_indicator_t), intent(inout)  this,
real(kind=rp), dimension(nell)  est,
real(kind=rp), dimension(nell)  sig,
real(kind=rp), dimension(lx1,ly1,lz1,nell)  var,
integer  nell,
real(kind=rp), dimension(lx1,ly1,lz1)  xa,
real(kind=rp), dimension(lx1,ly1,lz1)  xb,
integer  LX1,
integer  LY1,
integer  LZ1 
)
private
Parameters
estspectral indicator
sigcoefficient of the exponential fit
nellnumber of elements @paran var variable to calculate indicator @paran xa work array @paran xb work array
LX1gll points in x
LY1gll points in y
LZ1gll points in z

local variables

polynomial coefficients

Legendre coefficients; first value coeff(1,1,1)

copy of last SERI_NP columns of coefficients

estimated error

estimated decay rate

loop over elements

go to Legendre space (done in two operations) and square the coefficient

lower left corner

small value; nothing to od

extrapolate coefficients X - direction copy last SERI_NP collumns (or less if NX1 is smaller) SERI_ELR allows to exclude last row

get extrapolated values

Y - direction copy last SERI_NP collumns (or less if NY1 is smaller) SERI_ELR allows to exclude last row

get extrapolated values

Z - direction copy last SERI_NP collumns (or less if NZ1 is smaller) SERI_ELR allows to exclude last row

get extrapolated values

average

for testing

for testing; end

Definition at line 362 of file spectral_error_indicator.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ transform_to_spec_or_phys()

subroutine spectral_error_indicator::transform_to_spec_or_phys ( type(field_t), intent(inout)  u_hat,
type(field_t), intent(inout)  u,
type(field_t), intent(inout)  wk,
type(coef_t), intent(inout)  coef,
character(len=4), intent(in)  space 
)
private
Parameters
u_hattransformed field
ufield to transform
wkworking field
coeftype coef for mesh parameters
spaceString that indicates which space to transform, "spec" or "phys".

Define some constants

Copy field to working array

Definition at line 220 of file spectral_error_indicator.f90.

Here is the caller graph for this function: