![]() |
ORB5
v4.9.4
|
Module containing B-spline utilities for solvers. More...
Data Types | |
interface | rebase |
Functions/Subroutines | |
subroutine, public | initpde (nidbas, ns_, ds_, smin_, nchi_, dchi_, nphip_, dphi_, mmin_, mmax_, iquafor) |
Initialize the module. More... | |
real function, dimension(0:p) | splines (z) |
Evaluate basis functions on point z in interval [0, 1]. More... | |
subroutine, public | analytical_rhs (rhs) |
Project an analytical function on B-splines. More... | |
subroutine | rebase_r1d (v, nsel_which_mat) |
Multiply v by a boundary change of basis matrix. More... | |
subroutine | rebase_r2d (v, nsel_which_mat) |
Call rebase_r1d for each index of second dimension of v. More... | |
subroutine | rebase_c1d (v, nsel_which_mat) |
Complex version of rebase_r1d. More... | |
subroutine | rebase_c2d (v, nsel_which_mat) |
Call rebase_c1d for each index of second dimension of v. More... | |
subroutine, public | impose_bc (v) |
Impose boundary conditions on a vector before backsolve. More... | |
real function | func_int_base (zs, zj) |
Flux surface average of theta basis function. More... | |
Variables | |
integer, save, public | width_f |
Total width of diagonal filter as defined using mmin_f and mmax_f in the def_matrix_filter subroutine. Note that it may differ from width_m defined using the filter array. More... | |
integer, save, public | nbmatf |
number of matrix per proc More... | |
integer, save | p |
Spline order. More... | |
integer, save | nsupp |
support in number of intervals of the spline function More... | |
integer, dimension(:), allocatable, save | nlocad |
local 1D addresses of contributions of one spline function More... | |
integer, save | ns |
Number of intervals in radial direction. More... | |
real, save | ds |
Mesh size in radial direction. More... | |
real, save | smin |
Lower bound of radial direction. More... | |
integer, save | nchi |
Number of intervals in poloidal direction. More... | |
real, save | dchi |
Mesh size in poloidal direction. More... | |
integer, save | nphip |
Number of local intervals in toroidal direction. More... | |
real, save | dphi |
Mesh size in toroidal direction. More... | |
integer, dimension(:,:), allocatable, save | mmin |
minimum m value for Fourier solver More... | |
integer, dimension(:,:), allocatable, save | mmax |
maximum m value for Fourier solver More... | |
real, dimension(:), allocatable, save | splines_g |
basis functions at grid points More... | |
real, dimension(:,:), allocatable, save | splines_q |
basis functions at quadrature points More... | |
real, dimension(10), save | zqp |
Quadrature points. More... | |
real, dimension(10), save | zqw |
Quadrature weights. More... | |
integer, save | nq |
Number of quadrature points. More... | |
real, dimension(:,:), allocatable, save, public | zwavg |
Weight function surface average. More... | |
Module containing B-spline utilities for solvers.
subroutine, public pdespline::analytical_rhs | ( | real, dimension(0:ns+p-1, 0:nchi+p-1, 0:nphip+p-1), intent(out) | rhs | ) |
Project an analytical function on B-splines.
\( b_{ijk} = \int J_{\theta^\star s\varphi}\Lambda_i(s)\Lambda_j(\theta^\star)\Lambda_k(\varphi)f(s,\theta^\star,\varphi) ds d\theta^\star d\varphi \)
|
private |
Flux surface average of theta basis function.
[in] | zs | s position |
[in] | zj | Base index |
subroutine, public pdespline::impose_bc | ( | complex, dimension(:,:), intent(inout) | v | ) |
Impose boundary conditions on a vector before backsolve.
[in,out] | v | Vector on which boundary conditions are applied, of size (ns+p)*width_f in first dimension |
subroutine, public pdespline::initpde | ( | integer, intent(in) | nidbas, |
integer, intent(in) | ns_, | ||
real, intent(in) | ds_, | ||
real, intent(in) | smin_, | ||
integer, intent(in) | nchi_, | ||
real, intent(in) | dchi_, | ||
integer, intent(in) | nphip_, | ||
real, intent(in) | dphi_, | ||
integer, dimension(:,:), intent(in) | mmin_, | ||
integer, dimension(:,:), intent(in) | mmax_, | ||
integer, intent(in) | iquafor | ||
) |
Initialize the module.
[in] | nidbas | Spline order |
[in] | ns_ | Number of intervals in radial direction |
[in] | ds_ | Mesh size in radial direction |
[in] | smin_ | Lower bound of radial direction |
[in] | nchi_ | Number of intervals in poloidal direction |
[in] | dchi_ | Mesh size in poloidal direction |
[in] | nphip_ | Number of local intervals in toroidal direction |
[in] | dphi_ | Mesh size in toroidal direction |
[in] | mmin_ | minimum m value for Fourier solver |
[in] | mmax_ | maximum m value for Fourier solver |
[in] | iquafor | order of Gaussian integration |
|
private |
Complex version of rebase_r1d.
[in] | nsel_which_mat | Refer to rebase_r1d |
[in,out] | v | Vector to be multiplied, of size (ns+p)*width_f |
|
private |
Call rebase_c1d for each index of second dimension of v.
[in] | nsel_which_mat | Refer to rebase_c1d |
[in,out] | v | Vector to be multiplied, of size (ns+p)*width_f in first dimension |
|
private |
Multiply v by a boundary change of basis matrix.
[in] | nsel_which_mat | Which matrix to use at each boundary:
|
[in,out] | v | Vector to be multiplied, of size (ns+p)*nchi |
|
private |
Call rebase_r1d for each index of second dimension of v.
[in] | nsel_which_mat | Refer to rebase_r1d |
[in,out] | v | Vector to be multiplied, of size (ns+p)*nchi in first dimension |
|
private |
Evaluate basis functions on point z in interval [0, 1].
[in] | z | Point in [0, 1] where the functions are evaluated |
|
private |
Mesh size in poloidal direction.
|
private |
Mesh size in toroidal direction.
|
private |
Mesh size in radial direction.
|
private |
maximum m value for Fourier solver
|
private |
minimum m value for Fourier solver
integer, save, public pdespline::nbmatf |
number of matrix per proc
|
private |
Number of intervals in poloidal direction.
|
private |
local 1D addresses of contributions of one spline function
|
private |
Number of local intervals in toroidal direction.
|
private |
Number of quadrature points.
|
private |
Number of intervals in radial direction.
|
private |
support in number of intervals of the spline function
|
private |
Spline order.
|
private |
Lower bound of radial direction.
|
private |
basis functions at grid points
|
private |
basis functions at quadrature points
integer, save, public pdespline::width_f |
Total width of diagonal filter as defined using mmin_f
and mmax_f
in the def_matrix_filter
subroutine. Note that it may differ from width_m
defined using the filter array.
|
private |
Quadrature points.
|
private |
Quadrature weights.
real, dimension(:,:), allocatable, save, public pdespline::zwavg |
Weight function surface average.