Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines
sll_m_sparse_grid_3d Module Reference

Description

Implementation of a 3D sparse grid with interpolation routines.

Author
Katharina Kormann, IPP
Todo:
Implement the optimized interpolation routines for option boundary=1

<DETAILED_DESCRIPTION>

Derived types and interfaces

type  sll_t_sparse_grid_3d
 Sparse grid object for 3d with interpolation routines. More...
 

Functions/Subroutines

real(kind=f64) function interpolate_from_interpolant_value (interpolator, data, eta)
 Compute the value of the sparse grid interpolant at position eta (using standard sparse grid interpolation) More...
 
subroutine interpolate_const_disp (interpolator, dorder, displacement, data_in, data_out, hiera)
 Interpolation function for interpolation at (constantly) displaced grid points; displacement only in dimension dim. It is another implementation of the base-class function "interpolate_disp". The advantage is that we can not revisit nodes as we do in the recursive dimension-independently-programmed version. More...
 
real(kind=f64) function interpolate_from_hierarchical_surplus (interpolator, data, eta)
 Implements interpolate_from_interpolant_value for periodic sparse grid. More...
 
real(kind=f64) function interpolate_from_hierarchical_surplus_boundary (interpolator, data, eta)
 implements interpolation from hierarchical surplus (interpolate_from_interpolant_value) non-periodic More...
 
subroutine interpolate_array_disp_sgfft (interpolator, dim, displacment_in, data_in, data_out)
 Compute value at displaced grid points using trigonometric interpolation (based on SG FFT) More...
 
subroutine initialize_sg3d (interpolator, levels, order, interpolation, interpolation_type, eta_min, eta_max, boundary, modified)
 Initialization function. Set up the hierarchy of the sparse grid. More...
 
subroutine set_hierarchy_info (interpolator, counter, cdim, lvecin, kvecin, novecin)
 Helfer function for initialization. Setting all the information needed for node counter of the sparse grid along dimension cdim For a given sparse grid point fill the hierarchy information (3D specific) More...
 
subroutine set_hierarchy_info_boundary (interpolator, counter, cdim, lvecin, kvecin, novecin)
 Helfer function for initialization. Setting all the information needed for node counter of the sparse grid along dimension cdim for points at the boundary along dimension dim. More...
 
subroutine fg_to_sg (interpolator, fg_values, sg_values)
 Functions to evaluate fg on sg and sg on fg. More...
 
integer(kind=i32) function, dimension(3) fg_index (interpolator, sg_index)
 Compute the index of a sparse grid node on level "level" with index "index_on_level" on full grid with of max_level. More...
 
subroutine tohierarchical (interpolator, data_in, data_out)
 
subroutine todehi (interpolator, data_array)
 
subroutine tohira (interpolator, data_array)
 
subroutine tonodal (interpolator, data_in, data_out)
 
subroutine displace (interpolator, dim, displacement, data)
 
subroutine spfft (interpolator, data_in, data_out)
 Sparse grid FFT. More...
 
subroutine ispfft (interpolator, data_in, data_out)
 Sparse grid inverse FFT. More...
 

Function/Subroutine Documentation

◆ displace()

subroutine sll_m_sparse_grid_3d::displace ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
integer(kind=i32), intent(in)  dim,
real(kind=f64), intent(in)  displacement,
complex(kind=f64), dimension(:), intent(inout)  data 
)
private

Definition at line 942 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ fg_index()

integer(kind=i32) function, dimension(3) sll_m_sparse_grid_3d::fg_index ( class(sll_t_sparse_grid_3d), intent(in)  interpolator,
integer(kind=i32), intent(in)  sg_index 
)
private

Compute the index of a sparse grid node on level "level" with index "index_on_level" on full grid with of max_level.

Definition at line 751 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ fg_to_sg()

subroutine sll_m_sparse_grid_3d::fg_to_sg ( class(sll_t_sparse_grid_3d), intent(in)  interpolator,
real(kind=f64), dimension(:, :, :), intent(in)  fg_values,
real(kind=f64), dimension(:), intent(out)  sg_values 
)
private

Functions to evaluate fg on sg and sg on fg.

Definition at line 736 of file sll_m_sparse_grid_3d.F90.

Here is the call graph for this function:

◆ initialize_sg3d()

subroutine sll_m_sparse_grid_3d::initialize_sg3d ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
integer(kind=i32), dimension(:), intent(in)  levels,
integer(kind=i32), intent(in)  order,
integer(kind=i32), intent(in)  interpolation,
integer(kind=i32), intent(in)  interpolation_type,
real(kind=f64), dimension(:), intent(in)  eta_min,
real(kind=f64), dimension(:), intent(in)  eta_max,
integer(kind=i32), intent(in)  boundary,
integer(kind=i32), intent(in)  modified 
)
private

Initialization function. Set up the hierarchy of the sparse grid.

Parameters
[in,out]interpolatorsparse grid object
[in]eta_mineta_min defines the lower bound of the domain
[in]eta_maxeta_max defines the upper bound of the domain
[in]levelslevels defines the maximum level in the sparse grid for each direction
[in]orderorder of the sparse grid functions
[in]interpolationorder of the interpolator
[in]interpolation_typeChoose spline (interpolation_type = 0) or Lagrange (interpolation_type = 1) interpolation for the 1D interpolators if not traditional sparse grid interpolation is used.
[in]modifiedmodified defines if we have a traditional sparse grid for modified = 0 (then the l_1 norm of the levels is bounded by max(levels) ) or if the boundary is sparsified for modified = 1 (then the l_1 norm of the levels is bounded by max(levels)+1 )
[in]boundaryboundary defines the boundary conditions: define 0 for periodic boundary conditions and 1 for zero inflow boundaries

Definition at line 393 of file sll_m_sparse_grid_3d.F90.

Here is the call graph for this function:

◆ interpolate_array_disp_sgfft()

subroutine sll_m_sparse_grid_3d::interpolate_array_disp_sgfft ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
integer(kind=i32), intent(in)  dim,
real(kind=f64), intent(in)  displacment_in,
complex(kind=f64), dimension(:), intent(inout)  data_in,
real(kind=f64), dimension(:), intent(out)  data_out 
)
private

Compute value at displaced grid points using trigonometric interpolation (based on SG FFT)

Parameters
[in,out]interpolatorSparse grid object
[in]dimdimension along which the points should be displaced with displacement_in
[in]displacment_indisplacement of the data points along dimension dim
[in,out]data_inFourier transformed values on the sparse grid
[out]data_outFunction values on the sparse grid after displacement

Definition at line 373 of file sll_m_sparse_grid_3d.F90.

Here is the call graph for this function:

◆ interpolate_const_disp()

subroutine sll_m_sparse_grid_3d::interpolate_const_disp ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
integer(kind=i32), dimension(:), intent(in)  dorder,
real(kind=f64), intent(in)  displacement,
real(kind=f64), dimension(:), intent(inout)  data_in,
real(kind=f64), dimension(:), intent(out)  data_out,
logical, intent(in)  hiera 
)
private

Interpolation function for interpolation at (constantly) displaced grid points; displacement only in dimension dim. It is another implementation of the base-class function "interpolate_disp". The advantage is that we can not revisit nodes as we do in the recursive dimension-independently-programmed version.

Parameters
[in,out]interpolatorsparse grid object
[in,out]data_inhierarchical surplus
[out]data_outValue of the function or the hierarchical surplus (depending on value of hiera) for the displaced data points.
[in]dorderdorder(1) gives the dimension along which we have the displacement; dorder(2:3) give the remaining dimensions
[in]displacementConstant displacement along dimension dorder(1)
[in]hieraIf the result should be the hierarchical surplus, define hiera = .TRUE.; if the result should be the function values at the data points give hiera = .FALSE.

Definition at line 65 of file sll_m_sparse_grid_3d.F90.

◆ interpolate_from_hierarchical_surplus()

real(kind=f64) function sll_m_sparse_grid_3d::interpolate_from_hierarchical_surplus ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
real(kind=f64), dimension(:), intent(in)  data,
real(kind=f64), dimension(:), intent(in)  eta 
)
private

Implements interpolate_from_interpolant_value for periodic sparse grid.

Definition at line 124 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ interpolate_from_hierarchical_surplus_boundary()

real(kind=f64) function sll_m_sparse_grid_3d::interpolate_from_hierarchical_surplus_boundary ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
real(kind=f64), dimension(:), intent(in)  data,
real(kind=f64), dimension(:), intent(in)  eta 
)
private

implements interpolation from hierarchical surplus (interpolate_from_interpolant_value) non-periodic

Definition at line 181 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ interpolate_from_interpolant_value()

real(kind=f64) function sll_m_sparse_grid_3d::interpolate_from_interpolant_value ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
real(kind=f64), dimension(:), intent(in)  data,
real(kind=f64), dimension(:), intent(in)  eta 
)
private

Compute the value of the sparse grid interpolant at position eta (using standard sparse grid interpolation)

Parameters
[in,out]interpolatorsparse grid object
Returns
interpolated value
Parameters
[in]dataValue of hierarchical surplus
[in]etaCoordinates of the point where to interpolate

Definition at line 47 of file sll_m_sparse_grid_3d.F90.

Here is the call graph for this function:

◆ ispfft()

subroutine sll_m_sparse_grid_3d::ispfft ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
complex(kind=f64), dimension(:), intent(inout)  data_in,
real(kind=f64), dimension(:), intent(out)  data_out 
)
private

Sparse grid inverse FFT.

Definition at line 995 of file sll_m_sparse_grid_3d.F90.

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

◆ set_hierarchy_info()

subroutine sll_m_sparse_grid_3d::set_hierarchy_info ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
integer(kind=i32), intent(in)  counter,
integer(kind=i32), intent(in)  cdim,
integer(kind=i32), dimension(:), intent(in)  lvecin,
integer(kind=i32), dimension(:), intent(in)  kvecin,
integer(kind=i32), dimension(:), intent(in)  novecin 
)
private

Helfer function for initialization. Setting all the information needed for node counter of the sparse grid along dimension cdim For a given sparse grid point fill the hierarchy information (3D specific)

Parameters
[in]cdimcurrent dimension
[in]countercounter for node
[in]lvecinvector of current levels
[in]kvecinvector of current index within level
[in]novecinvector with number of points on the current level along each dimension

Definition at line 542 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ set_hierarchy_info_boundary()

subroutine sll_m_sparse_grid_3d::set_hierarchy_info_boundary ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
integer(kind=i32), intent(in)  counter,
integer(kind=i32), intent(in)  cdim,
integer(kind=i32), dimension(:), intent(in)  lvecin,
integer(kind=i32), dimension(:), intent(in)  kvecin,
integer(kind=i32), dimension(:), intent(in)  novecin 
)
private

Helfer function for initialization. Setting all the information needed for node counter of the sparse grid along dimension cdim for points at the boundary along dimension dim.

Parameters
[in]cdimcurrent dimension
[in]countercounter for node
[in]lvecinvector of current levels
[in]kvecinvector of current index within level
[in]novecinvector with number of points on the current level along each dimension

Definition at line 621 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ spfft()

subroutine sll_m_sparse_grid_3d::spfft ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
real(kind=f64), dimension(:), intent(in)  data_in,
complex(kind=f64), dimension(:), intent(out)  data_out 
)
private

Sparse grid FFT.

Definition at line 984 of file sll_m_sparse_grid_3d.F90.

Here is the call graph for this function:

◆ todehi()

subroutine sll_m_sparse_grid_3d::todehi ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
complex(kind=f64), dimension(:), intent(inout)  data_array 
)
private

Definition at line 815 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ tohierarchical()

subroutine sll_m_sparse_grid_3d::tohierarchical ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
real(kind=f64), dimension(:), intent(in)  data_in,
complex(kind=f64), dimension(:), intent(out)  data_out 
)
private

Definition at line 771 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ tohira()

subroutine sll_m_sparse_grid_3d::tohira ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
complex(kind=f64), dimension(:), intent(inout)  data_array 
)
private

Definition at line 857 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:

◆ tonodal()

subroutine sll_m_sparse_grid_3d::tonodal ( class(sll_t_sparse_grid_3d), intent(inout)  interpolator,
complex(kind=f64), dimension(:), intent(inout)  data_in,
real(kind=f64), dimension(:), intent(out)  data_out 
)
private

Definition at line 899 of file sll_m_sparse_grid_3d.F90.

Here is the caller graph for this function:
    Report Typos and Errors