(git:6a2e663)
pw_spline_utils Module Reference

different utils that are useful to manipulate splines on the regular grid of a pw More...

Functions/Subroutines

subroutine, public pw_spline2_interpolate_values_g (spline_g)
 calculates the FFT of the coefficients of the quadratic spline that interpolates the given values More...
 
subroutine, public pw_spline3_interpolate_values_g (spline_g)
 calculates the FFT of the coefficients of the2 cubic spline that interpolates the given values More...
 
subroutine, public pw_spline_scale_deriv (deriv_vals_r, transpose, scale)
 rescales the derivatives from gridspacing=1 to the real derivatives More...
 
subroutine, public pw_spline3_deriv_g (spline_g, idir)
 calculates the FFT of the values of the x,y,z (idir=1,2,3) derivative of the cubic spline More...
 
subroutine, public pw_spline2_deriv_g (spline_g, idir)
 calculates the FFT of the values of the x,y,z (idir=1,2,3) derivative of the quadratic spline More...
 
subroutine, public pw_nn_smear_r (pw_in, pw_out, coeffs)
 calculates the values of a nearest neighbor smearing More...
 
subroutine, public pw_nn_deriv_r (pw_in, pw_out, coeffs, idir)
 calculates a nearest neighbor central derivative. for the x dir: pw_outarray(i,j,k)=( pw_in(i+1,j,k)-pw_in(i-1,j,k) )*coeff(1)+ ( pw_in(i+1,j(+-)1,k)-pw_in(i-1,j(+-)1,k)+ pw_in(i+1,j,k(+-)1)-pw_in(i-1,j,k(+-)1) )*coeff(2)+ ( pw_in(i+1,j(+-)1,k(+-)1)-pw_in(i-1,j(+-)1,k(+-)1)+ pw_in(i+1,j(+-)1,k(-+)1)-pw_in(i-1,j(+-)1,k(-+)1) )*coeff(3) periodic boundary conditions are applied More...
 
subroutine, public add_coarse2fine (coarse_coeffs_pw, fine_values_pw, weights_1d, w_border0, w_border1, pbc, safe_computation)
 low level function that adds a coarse grid to a fine grid. If pbc is true periodic boundary conditions are applied More...
 
subroutine, public add_fine2coarse (fine_values_pw, coarse_coeffs_pw, weights_1d, w_border0, w_border1, pbc, safe_computation)
 low level function that adds a coarse grid (without boundary) to a fine grid. More...
 
subroutine, public pw_spline_precond_create (preconditioner, precond_kind, pool, pbc, transpose)
 ... More...
 
subroutine, public pw_spline_precond_set_kind (preconditioner, precond_kind, pbc, transpose)
 switches the types of precoditioner to use More...
 
subroutine, public pw_spline_precond_release (preconditioner)
 releases the preconditioner More...
 
subroutine, public pw_spline_do_precond (preconditioner, in_v, out_v)
 applies the preconditioner to the system of equations to find the coefficients of the spline More...
 
logical function, public find_coeffs (values, coeffs, linOp, preconditioner, pool, eps_r, eps_x, max_iter, sumtype)
 solves iteratively (CG) a systmes of linear equations linOp(coeffs)=values (for example those needed to find the coefficients of a spline) Returns true if the it succeeded to achieve the requested accuracy More...
 
subroutine pw_nn_compose_r_no_pbc (weights_1d, pw_in, pw_out, sharpen, normalize, transpose, smooth_boundary)
 adds to pw_out pw_in composed with the weights pw_outarray(i,j,k)=pw_outarray(i,j,k)+sum(pw_inarray(i+l,j+m,k+n)* weights_1d(abs(l)+1)*weights_1d(abs(m)+1)*weights_1d(abs(n)+1), l=-1..1,m=-1..1,n=-1..1) More...
 
subroutine, public spl3_nopbc (pw_in, pw_out)
 ... More...
 
subroutine, public spl3_nopbct (pw_in, pw_out)
 ... More...
 
subroutine, public spl3_pbc (pw_in, pw_out)
 ... More...
 
real(kind=dp) function, public eval_interp_spl3_pbc (vec, pw)
 Evaluates the PBC interpolated Spline (pw) function on the generic input vector (vec) More...
 
real(kind=dp) function, dimension(3), public eval_d_interp_spl3_pbc (vec, pw)
 Evaluates the derivatives of the PBC interpolated Spline (pw) function on the generic input vector (vec) More...
 

Variables

integer, parameter, public no_precond = 0
 
integer, parameter, public precond_spl3_aint = 1
 
integer, parameter, public precond_spl3_1 = 2
 
integer, parameter, public precond_spl3_aint2 = 3
 
integer, parameter, public precond_spl3_2 = 4
 
integer, parameter, public precond_spl3_3 = 5
 
real(kind=dp), dimension(4), parameter, public nn10_coeffs = (/125._dp/216._dp, 25._dp/432._dp, 5._dp/864._dp, 1._dp/1728._dp/)
 
real(kind=dp), dimension(4), parameter, public spline3_coeffs = (/8._dp/(27._dp), 2._dp/(27._dp), 1._dp/(27._dp*2._dp), 1._dp/(27._dp*8._dp)/)
 
real(kind=dp), dimension(4), parameter, public spline2_coeffs = (/27._dp/(64._dp), 9._dp/(64._dp*2_dp), 3._dp/(64._dp*4._dp), 1._dp/(64._dp*8._dp)/)
 
real(kind=dp), dimension(4), parameter, public nn50_coeffs = (/15625._dp/17576._dp, 625._dp/35152._dp, 25._dp/70304._dp, 1._dp/140608._dp/)
 
real(kind=dp), dimension(4), parameter, public spl3_aint_coeff = (/46._dp/27._dp, -2._dp/(27._dp), -1._dp/(27._dp*2._dp), -1._dp/(27._dp*8._dp)/)
 
real(kind=dp), dimension(4), parameter, public spl3_precond1_coeff = (/64._dp/3._dp, -8._dp/3._dp, -1._dp/3._dp, -1._dp/24._dp/)
 
real(kind=dp), dimension(4), parameter, public spl3_1d_transf_coeffs = (/2._dp/3._dp, 23._dp/48._dp, 1._dp/6._dp, 1._dp/48._dp/)
 
real(kind=dp), dimension(3), parameter, public spline3_deriv_coeffs = (/2.0_dp/9.0_dp, 1.0_dp/18.0_dp, 1.0_dp/72.0_dp/)
 
real(kind=dp), dimension(3), parameter, public spline2_deriv_coeffs = (/9.0_dp/32.0_dp, 3.0_dp/64.0_dp, 1.0_dp/128.0_dp/)
 
real(kind=dp), dimension(3), parameter, public nn10_deriv_coeffs = (/25._dp/72._dp, 5._dp/144, 1._dp/288._dp/)
 
real(kind=dp), dimension(3), parameter, public nn50_deriv_coeffs = (/625._dp/1352._dp, 25._dp/2704._dp, 1._dp/5408._dp/)
 
real(kind=dp), dimension(3), parameter, public spl3_1d_coeffs0 = (/1._dp/6_dp, 2._dp/3._dp, 1._dp/6._dp/)
 
real(kind=dp), dimension(3), parameter, public spl3_1d_transf_border1 = (/0.517977704_dp, 0.464044595_dp, 0.17977701e-1_dp/)
 

Detailed Description

different utils that are useful to manipulate splines on the regular grid of a pw

History
05.2003 created [fawzi] 08.2004 removed spline evaluation method using more than 2 read streams (pw_compose_stripe_rs3), added linear solver based spline inversion [fawzi]
Author
Fawzi Mohamed

Function/Subroutine Documentation

◆ pw_spline2_interpolate_values_g()

subroutine, public pw_spline_utils::pw_spline2_interpolate_values_g ( type(pw_c1d_gs_type), intent(in)  spline_g)

calculates the FFT of the coefficients of the quadratic spline that interpolates the given values

Parameters
spline_gon entry the FFT of the values to interpolate as cc, will contain the FFT of the coefficients of the spline
History
06.2003 created [fawzi]
Author
Fawzi Mohamed
Note
does not work with spherical cutoff

Definition at line 128 of file pw_spline_utils.F.

Here is the caller graph for this function:

◆ pw_spline3_interpolate_values_g()

subroutine, public pw_spline_utils::pw_spline3_interpolate_values_g ( type(pw_c1d_gs_type), intent(in)  spline_g)

calculates the FFT of the coefficients of the2 cubic spline that interpolates the given values

Parameters
spline_gon entry the FFT of the values to interpolate as cc, will contain the FFT of the coefficients of the spline
History
06.2003 created [fawzi]
Author
Fawzi Mohamed
Note
does not work with spherical cutoff stupid distribution for cos calculation, it should calculate only the needed cos, and avoid the mpi_allreduce

Definition at line 199 of file pw_spline_utils.F.

Here is the caller graph for this function:

◆ pw_spline_scale_deriv()

subroutine, public pw_spline_utils::pw_spline_scale_deriv ( type(pw_r3d_rs_type), dimension(3), intent(in)  deriv_vals_r,
logical, intent(in), optional  transpose,
real(kind=dp), intent(in), optional  scale 
)

rescales the derivatives from gridspacing=1 to the real derivatives

Parameters
deriv_vals_ran array of x,y,z derivatives
transposeif true applies the transpose of the map (defaults to false)
scalea scaling factor (defaults to 1.0)
History
06.2003 created [fawzi]
Author
Fawzi Mohamed

Definition at line 274 of file pw_spline_utils.F.

Here is the caller graph for this function:

◆ pw_spline3_deriv_g()

subroutine, public pw_spline_utils::pw_spline3_deriv_g ( type(pw_c1d_gs_type), intent(in)  spline_g,
integer, intent(in)  idir 
)

calculates the FFT of the values of the x,y,z (idir=1,2,3) derivative of the cubic spline

Parameters
spline_gon entry the FFT of the coefficients of the spline will contain the FFT of the derivative
idirdirection of the derivative
History
06.2003 created [fawzi]
Author
Fawzi Mohamed
Note
the distance between gridpoints is assumed to be 1

Definition at line 365 of file pw_spline_utils.F.

Here is the caller graph for this function:

◆ pw_spline2_deriv_g()

subroutine, public pw_spline_utils::pw_spline2_deriv_g ( type(pw_c1d_gs_type), intent(in)  spline_g,
integer, intent(in)  idir 
)

calculates the FFT of the values of the x,y,z (idir=1,2,3) derivative of the quadratic spline

Parameters
spline_gon entry the FFT of the coefficients of the spline will contain the FFT of the derivative
idirdirection of the derivative
History
06.2003 created [fawzi]
Author
Fawzi Mohamed
Note
the distance between gridpoints is assumed to be 1

Definition at line 500 of file pw_spline_utils.F.

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

◆ pw_nn_smear_r()

subroutine, public pw_spline_utils::pw_nn_smear_r ( type(pw_r3d_rs_type), intent(in)  pw_in,
type(pw_r3d_rs_type), intent(in)  pw_out,
real(kind=dp), dimension(4), intent(in)  coeffs 
)

calculates the values of a nearest neighbor smearing

Parameters
pw_inthe argument for the linear operator
pw_outplace where the smeared values should be added
coeffsarray with the coefficent of the smearing, ordered with the distance from the center: coeffs(1) the coeff of the central element, coeffs(2) the coeff of the 6 element with distance 1, coeff(3) the coeff of the 12 elements at distance sqrt(2), coeff(4) the coeff of the 8 elements at distance sqrt(3).
Author
Fawzi Mohamed
Note
does not normalize the smear to 1. with coeff=(/ 8._dp/27._dp, 2._dp/27._dp, 1._dp/54._dp, 1._dp/216._dp /) is equivalent to pw_spline3_evaluate_values_g, with coeff=(/ 27._dp/64._dp, 9._dp/128._dp, 3._dp/256._dp, 1._dp/512._dp /)

Definition at line 855 of file pw_spline_utils.F.

Here is the caller graph for this function:

◆ pw_nn_deriv_r()

subroutine, public pw_spline_utils::pw_nn_deriv_r ( type(pw_r3d_rs_type), intent(in)  pw_in,
type(pw_r3d_rs_type), intent(in)  pw_out,
real(kind=dp), dimension(3), intent(in)  coeffs,
integer  idir 
)

calculates a nearest neighbor central derivative. for the x dir: pw_outarray(i,j,k)=( pw_in(i+1,j,k)-pw_in(i-1,j,k) )*coeff(1)+ ( pw_in(i+1,j(+-)1,k)-pw_in(i-1,j(+-)1,k)+ pw_in(i+1,j,k(+-)1)-pw_in(i-1,j,k(+-)1) )*coeff(2)+ ( pw_in(i+1,j(+-)1,k(+-)1)-pw_in(i-1,j(+-)1,k(+-)1)+ pw_in(i+1,j(+-)1,k(-+)1)-pw_in(i-1,j(+-)1,k(-+)1) )*coeff(3) periodic boundary conditions are applied

Parameters
pw_inthe argument for the linear operator
pw_outplace where the smeared values should be added
coeffsarray with the coefficent of the front (positive) plane of the central derivative, ordered with the distance from the center: coeffs(1) the coeff of the central element, coeffs(2) the coeff of the 4 element with distance 1, coeff(3) the coeff of the 4 elements at distance sqrt(2)
idir...
Author
Fawzi Mohamed
Note
with coeff=(/ 2.0_dp/9.0_dp, 1.0_dp/18.0_dp, 1.0_dp/72.0_dp /) is equivalent to pw_spline3_deriv_r, with coeff=(/ 9.0_dp/32.0_dp, 3.0_dp/64.0_dp, 1.0_dp/128.0_dp /) to pw_spline2_deriv_r coeff=(/ 25._dp/72._dp, 5._dp/144, 1._dp/288._dp /)

Definition at line 898 of file pw_spline_utils.F.

Here is the caller graph for this function:

◆ add_coarse2fine()

subroutine, public pw_spline_utils::add_coarse2fine ( type(pw_r3d_rs_type), intent(in)  coarse_coeffs_pw,
type(pw_r3d_rs_type), intent(in)  fine_values_pw,
real(kind=dp), dimension(4), intent(in)  weights_1d,
real(kind=dp), intent(in)  w_border0,
real(kind=dp), dimension(3), intent(in)  w_border1,
logical, intent(in)  pbc,
logical, intent(in), optional  safe_computation 
)

low level function that adds a coarse grid to a fine grid. If pbc is true periodic boundary conditions are applied

It will add to

fine_values(2*coarse_bounds(1,1):2*coarse_bounds(2,1), 2*coarse_bounds(1,2):2*coarse_bounds(2,2), 2*coarse_bounds(1,3):2*coarse_bounds(2,3))

using

coarse_coeffs(coarse_bounds(1,1):coarse_bounds(2,1), coarse_bounds(1,2):coarse_bounds(2,2), coarse_bounds(1,3):coarse_bounds(2,3))

composed with the weights obtained by the direct product of the 1d coefficients weights:

for i,j,k in -3..3 w(i,j,k)=weights_1d(abs(i)+1)*weights_1d(abs(j)+1)* weights_1d(abs(k)+1)

Parameters
coarse_coeffs_pwthe values of the coefficients
fine_values_pwwhere to add the values due to the coarse coeffs
weights_1dthe weights of the 1d smearing
w_border0the 1d weight at the border (when pbc is false)
w_border1the 1d weights for a point one off the border (w_border1(1) is the weight of the coefficent at the border) (used if pbc is false)
pbcif periodic boundary conditions should be applied
safe_computation...
Author
fawzi
Note
coarse looping is continuos, I did not check if keeping the fine looping contiguous is better. And I ask myself yet again why, why we use x-slice distribution, z-slice distribution would be much better performancewise (and would semplify this code enormously). fine2coarse has much more understandable parallel part (build up of send/rcv sizes,... but worse if you have really a lot of processors, probabily irrelevant because it is not critical) [fawzi].

Definition at line 975 of file pw_spline_utils.F.

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

◆ add_fine2coarse()

subroutine, public pw_spline_utils::add_fine2coarse ( type(pw_r3d_rs_type), intent(in)  fine_values_pw,
type(pw_r3d_rs_type), intent(in)  coarse_coeffs_pw,
real(kind=dp), dimension(4), intent(in)  weights_1d,
real(kind=dp), intent(in)  w_border0,
real(kind=dp), dimension(3), intent(in)  w_border1,
logical, intent(in)  pbc,
logical, intent(in), optional  safe_computation 
)

low level function that adds a coarse grid (without boundary) to a fine grid.

It will add to

coarse_coeffs(coarse_bounds(1,1):coarse_bounds(2,1), coarse_bounds(1,2):coarse_bounds(2,2), coarse_bounds(1,3):coarse_bounds(2,3))

using

fine_values(2*coarse_bounds(1,1):2*coarse_bounds(2,1), 2*coarse_bounds(1,2):2*coarse_bounds(2,2), 2*coarse_bounds(1,3):2*coarse_bounds(2,3))

composed with the weights obtained by the direct product of the 1d coefficients weights:

for i,j,k in -3..3 w(i,j,k)=weights_1d(abs(i)+1)*weights_1d(abs(j)+1)* weights_1d(abs(k)+1)

Parameters
fine_values_pw3d array where to add the values due to the coarse coeffs
coarse_coeffs_pw3d array with boundary of size 1 with the values of the coefficients
weights_1dthe weights of the 1d smearing
w_border0the 1d weight at the border
w_border1the 1d weights for a point one off the border (w_border1(1) is the weight of the coefficent at the border)
pbc...
safe_computation...
Author
fawzi
Note
see coarse2fine for some relevant notes

Definition at line 1770 of file pw_spline_utils.F.

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

◆ pw_spline_precond_create()

subroutine, public pw_spline_utils::pw_spline_precond_create ( type(pw_spline_precond_type), intent(out)  preconditioner,
integer, intent(in)  precond_kind,
type(pw_pool_type), intent(in), pointer  pool,
logical, intent(in)  pbc,
logical, intent(in)  transpose 
)

...

Parameters
preconditionerthe preconditioner to create
precond_kindthe kind of preconditioner to use
poola pool with grids of the same type as the elements to precondition
pbcif periodic boundary conditions should be applied
transpose...
Author
fawzi

Definition at line 2473 of file pw_spline_utils.F.

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

◆ pw_spline_precond_set_kind()

subroutine, public pw_spline_utils::pw_spline_precond_set_kind ( type(pw_spline_precond_type), intent(inout)  preconditioner,
integer, intent(in)  precond_kind,
logical, intent(in), optional  pbc,
logical, intent(in), optional  transpose 
)

switches the types of precoditioner to use

Parameters
preconditionerthe preconditioner to be changed
precond_kindthe new kind of preconditioner to use
pbc...
transpose...
Author
fawzi

Definition at line 2496 of file pw_spline_utils.F.

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

◆ pw_spline_precond_release()

subroutine, public pw_spline_utils::pw_spline_precond_release ( type(pw_spline_precond_type), intent(inout)  preconditioner)

releases the preconditioner

Parameters
preconditionerthe preconditioner to release
Author
fawzi

Definition at line 2580 of file pw_spline_utils.F.

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

◆ pw_spline_do_precond()

subroutine, public pw_spline_utils::pw_spline_do_precond ( type(pw_spline_precond_type), intent(in)  preconditioner,
type(pw_r3d_rs_type), intent(in)  in_v,
type(pw_r3d_rs_type), intent(inout)  out_v 
)

applies the preconditioner to the system of equations to find the coefficients of the spline

Parameters
preconditionerthe preconditioner to apply
in_vthe grid on which the preconditioner should be applied
out_vplace to store the preconditioner applied on v_out
Author
fawzi

Definition at line 2594 of file pw_spline_utils.F.

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

◆ find_coeffs()

logical function, public pw_spline_utils::find_coeffs ( type(pw_r3d_rs_type), intent(in)  values,
type(pw_r3d_rs_type), intent(inout)  coeffs,
  linOp,
type(pw_spline_precond_type), intent(in)  preconditioner,
type(pw_pool_type), pointer  pool,
real(kind=dp), intent(in)  eps_r,
real(kind=dp), intent(in)  eps_x,
integer, intent(in)  max_iter,
integer, intent(in), optional  sumtype 
)

solves iteratively (CG) a systmes of linear equations linOp(coeffs)=values (for example those needed to find the coefficients of a spline) Returns true if the it succeeded to achieve the requested accuracy

Parameters
valuesthe right hand side of the system
coeffswill contain the solution of the system (and on entry it contains the starting point)
linOpthe linear operator to be inverted
preconditionerthe preconditioner to apply
poola pool of grids (for the temporary objects)
eps_rthe requested precision on the residual
eps_xthe requested precision on the solution
max_itermaximum number of iteration allowed
sumtype...
Returns
...
Author
fawzi

Definition at line 2647 of file pw_spline_utils.F.

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

◆ pw_nn_compose_r_no_pbc()

subroutine pw_spline_utils::pw_nn_compose_r_no_pbc ( real(kind=dp), dimension(-1:1)  weights_1d,
type(pw_r3d_rs_type), intent(in)  pw_in,
type(pw_r3d_rs_type), intent(in)  pw_out,
logical, intent(in), optional  sharpen,
logical, intent(in), optional  normalize,
logical, intent(in), optional  transpose,
logical, intent(in), optional  smooth_boundary 
)

adds to pw_out pw_in composed with the weights pw_outarray(i,j,k)=pw_outarray(i,j,k)+sum(pw_inarray(i+l,j+m,k+n)* weights_1d(abs(l)+1)*weights_1d(abs(m)+1)*weights_1d(abs(n)+1), l=-1..1,m=-1..1,n=-1..1)

Parameters
weights_1d...
pw_in...
pw_out...
sharpen...
normalize...
transpose...
smooth_boundary...
Author
fawzi

Definition at line 2756 of file pw_spline_utils.F.

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

◆ spl3_nopbc()

subroutine, public pw_spline_utils::spl3_nopbc ( type(pw_r3d_rs_type), intent(in)  pw_in,
type(pw_r3d_rs_type), intent(inout)  pw_out 
)

...

Parameters
pw_in...
pw_out...

Definition at line 3096 of file pw_spline_utils.F.

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

◆ spl3_nopbct()

subroutine, public pw_spline_utils::spl3_nopbct ( type(pw_r3d_rs_type), intent(in)  pw_in,
type(pw_r3d_rs_type), intent(inout)  pw_out 
)

...

Parameters
pw_in...
pw_out...

Definition at line 3112 of file pw_spline_utils.F.

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

◆ spl3_pbc()

subroutine, public pw_spline_utils::spl3_pbc ( type(pw_r3d_rs_type), intent(in)  pw_in,
type(pw_r3d_rs_type), intent(inout)  pw_out 
)

...

Parameters
pw_in...
pw_out...

Definition at line 3128 of file pw_spline_utils.F.

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

◆ eval_interp_spl3_pbc()

real(kind=dp) function, public pw_spline_utils::eval_interp_spl3_pbc ( real(kind=dp), dimension(3), intent(in)  vec,
type(pw_r3d_rs_type), intent(in)  pw 
)

Evaluates the PBC interpolated Spline (pw) function on the generic input vector (vec)

Parameters
vec...
pw...
Returns
...
History
12.2007 Adapted for use with distributed grids [rdeclerck]
Author
Teodoro Laino 12/2005 [tlaino]
Note
Requires the Spline coefficients to be computed with PBC

Definition at line 3150 of file pw_spline_utils.F.

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

◆ eval_d_interp_spl3_pbc()

real(kind=dp) function, dimension(3), public pw_spline_utils::eval_d_interp_spl3_pbc ( real(kind=dp), dimension(3), intent(in)  vec,
type(pw_r3d_rs_type), intent(in)  pw 
)

Evaluates the derivatives of the PBC interpolated Spline (pw) function on the generic input vector (vec)

Parameters
vec...
pw...
Returns
...
History
12.2007 Adapted for use with distributed grids [rdeclerck]
Author
Teodoro Laino 12/2005 [tlaino]
Note
Requires the Spline coefficients to be computed with PBC

Definition at line 3299 of file pw_spline_utils.F.

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

Variable Documentation

◆ no_precond

integer, parameter, public pw_spline_utils::no_precond = 0

Definition at line 44 of file pw_spline_utils.F.

◆ precond_spl3_aint

integer, parameter, public pw_spline_utils::precond_spl3_aint = 1

Definition at line 44 of file pw_spline_utils.F.

◆ precond_spl3_1

integer, parameter, public pw_spline_utils::precond_spl3_1 = 2

Definition at line 44 of file pw_spline_utils.F.

◆ precond_spl3_aint2

integer, parameter, public pw_spline_utils::precond_spl3_aint2 = 3

Definition at line 44 of file pw_spline_utils.F.

◆ precond_spl3_2

integer, parameter, public pw_spline_utils::precond_spl3_2 = 4

Definition at line 44 of file pw_spline_utils.F.

◆ precond_spl3_3

integer, parameter, public pw_spline_utils::precond_spl3_3 = 5

Definition at line 44 of file pw_spline_utils.F.

◆ nn10_coeffs

real(kind=dp), dimension(4), parameter, public pw_spline_utils::nn10_coeffs = (/125._dp/216._dp, 25._dp/432._dp, 5._dp/864._dp, 1._dp/1728._dp/)

Definition at line 51 of file pw_spline_utils.F.

◆ spline3_coeffs

real(kind=dp), dimension(4), parameter, public pw_spline_utils::spline3_coeffs = (/8._dp/(27._dp), 2._dp/(27._dp), 1._dp/(27._dp*2._dp), 1._dp/(27._dp*8._dp)/)

Definition at line 51 of file pw_spline_utils.F.

◆ spline2_coeffs

real(kind=dp), dimension(4), parameter, public pw_spline_utils::spline2_coeffs = (/27._dp/(64._dp), 9._dp/(64._dp*2_dp), 3._dp/(64._dp*4._dp), 1._dp/(64._dp*8._dp)/)

Definition at line 51 of file pw_spline_utils.F.

◆ nn50_coeffs

real(kind=dp), dimension(4), parameter, public pw_spline_utils::nn50_coeffs = (/15625._dp/17576._dp, 625._dp/35152._dp, 25._dp/70304._dp, 1._dp/140608._dp/)

Definition at line 51 of file pw_spline_utils.F.

◆ spl3_aint_coeff

real(kind=dp), dimension(4), parameter, public pw_spline_utils::spl3_aint_coeff = (/46._dp/27._dp, -2._dp/(27._dp), -1._dp/(27._dp*2._dp), -1._dp/(27._dp*8._dp)/)

Definition at line 51 of file pw_spline_utils.F.

◆ spl3_precond1_coeff

real(kind=dp), dimension(4), parameter, public pw_spline_utils::spl3_precond1_coeff = (/64._dp/3._dp, -8._dp/3._dp, -1._dp/3._dp, -1._dp/24._dp/)

Definition at line 51 of file pw_spline_utils.F.

◆ spl3_1d_transf_coeffs

real(kind=dp), dimension(4), parameter, public pw_spline_utils::spl3_1d_transf_coeffs = (/2._dp/3._dp, 23._dp/48._dp, 1._dp/6._dp, 1._dp/48._dp/)

Definition at line 51 of file pw_spline_utils.F.

◆ spline3_deriv_coeffs

real(kind=dp), dimension(3), parameter, public pw_spline_utils::spline3_deriv_coeffs = (/2.0_dp/9.0_dp, 1.0_dp/18.0_dp, 1.0_dp/72.0_dp/)

Definition at line 70 of file pw_spline_utils.F.

◆ spline2_deriv_coeffs

real(kind=dp), dimension(3), parameter, public pw_spline_utils::spline2_deriv_coeffs = (/9.0_dp/32.0_dp, 3.0_dp/64.0_dp, 1.0_dp/128.0_dp/)

Definition at line 70 of file pw_spline_utils.F.

◆ nn10_deriv_coeffs

real(kind=dp), dimension(3), parameter, public pw_spline_utils::nn10_deriv_coeffs = (/25._dp/72._dp, 5._dp/144, 1._dp/288._dp/)

Definition at line 70 of file pw_spline_utils.F.

◆ nn50_deriv_coeffs

real(kind=dp), dimension(3), parameter, public pw_spline_utils::nn50_deriv_coeffs = (/625._dp/1352._dp, 25._dp/2704._dp, 1._dp/5408._dp/)

Definition at line 70 of file pw_spline_utils.F.

◆ spl3_1d_coeffs0

real(kind=dp), dimension(3), parameter, public pw_spline_utils::spl3_1d_coeffs0 = (/1._dp/6_dp, 2._dp/3._dp, 1._dp/6._dp/)

Definition at line 70 of file pw_spline_utils.F.

◆ spl3_1d_transf_border1

real(kind=dp), dimension(3), parameter, public pw_spline_utils::spl3_1d_transf_border1 = (/0.517977704_dp, 0.464044595_dp, 0.17977701e-1_dp/)

Definition at line 70 of file pw_spline_utils.F.