ORB5  v4.9.4
pdespline Module Reference

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:psplines (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...
 

Detailed Description

Module containing B-spline utilities for solvers.

Function/Subroutine Documentation

◆ analytical_rhs()

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.

Date
02.2018

\( 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 \)

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

◆ func_int_base()

real function pdespline::func_int_base ( real, intent(in)  zs,
integer, intent(in)  zj 
)
private

Flux surface average of theta basis function.

Parameters
[in]zss position
[in]zjBase index
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ impose_bc()

subroutine, public pdespline::impose_bc ( complex, dimension(:,:), intent(inout)  v)

Impose boundary conditions on a vector before backsolve.

Author
N. Ohana
Date
12.2017
Parameters
[in,out]vVector on which boundary conditions are applied, of size (ns+p)*width_f in first dimension
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ initpde()

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.

Parameters
[in]nidbasSpline 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]iquafororder of Gaussian integration
Note
It is important to use the same quadrature method as in matrix assembly.
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ rebase_c1d()

subroutine pdespline::rebase_c1d ( complex, dimension(:), intent(inout)  v,
integer, intent(in)  nsel_which_mat 
)
private

Complex version of rebase_r1d.

Author
N. Ohana
Date
12.2017
Parameters
[in]nsel_which_matRefer to rebase_r1d
[in,out]vVector to be multiplied, of size (ns+p)*width_f

◆ rebase_c2d()

subroutine pdespline::rebase_c2d ( complex, dimension(:,:), intent(inout)  v,
integer, intent(in)  nsel_which_mat 
)
private

Call rebase_c1d for each index of second dimension of v.

Author
N. Ohana
Date
12.2017
Parameters
[in]nsel_which_matRefer to rebase_c1d
[in,out]vVector to be multiplied, of size (ns+p)*width_f in first dimension

◆ rebase_r1d()

subroutine pdespline::rebase_r1d ( real, dimension(:), intent(inout)  v,
integer, intent(in)  nsel_which_mat 
)
private

Multiply v by a boundary change of basis matrix.

Author
N. Ohana
Date
12.2017
Parameters
[in]nsel_which_matWhich matrix to use at each boundary:
  • NSEL_T for \( T = \begin{pmatrix} 1/c_1 & 0 & 0 & \dots & 0 \\\\ c_2' & 1 & 0 & \dots & 0 \\\\ c_3' & 0 & 1 & \ddots & \vdots \\\\ \vdots & \vdots & \ddots & \ddots & 0 \\\\ c_p' & 0 & \dots & 0 & 1 \end{pmatrix} \)
  • NSEL_T_TRANSP for \( T^{t} = \begin{pmatrix} 1/c_1 & c_2' & c_3' & \dots & c_p' \\\\ 0 & 1 & 0 & \dots & 0 \\\\ 0 & 0 & 1 & \ddots & \vdots \\\\ \vdots & \vdots & \ddots & \ddots & 0 \\\\ 0 & 0 & \dots & 0 & 1 \end{pmatrix} \)
  • NSEL_T_INV for \( T^{-1} = \begin{pmatrix} c_1 & 0 & 0 & \dots & 0 \\\\ c_2 & 1 & 0 & \dots & 0 \\\\ c_3 & 0 & 1 & \ddots & \vdots \\\\ \vdots & \vdots & \ddots & \ddots & 0 \\\\ c_p & 0 & \dots & 0 & 1 \end{pmatrix} \)
    where \( c_i = \Lambda_i^{(p)}(0) \) and \( c_i' = -c_i/c_1 \)
[in,out]vVector to be multiplied, of size (ns+p)*nchi

◆ rebase_r2d()

subroutine pdespline::rebase_r2d ( real, dimension(:,:), intent(inout)  v,
integer, intent(in)  nsel_which_mat 
)
private

Call rebase_r1d for each index of second dimension of v.

Author
N. Ohana
Date
12.2017
Parameters
[in]nsel_which_matRefer to rebase_r1d
[in,out]vVector to be multiplied, of size (ns+p)*nchi in first dimension

◆ splines()

real function, dimension(0:p) pdespline::splines ( real, intent(in)  z)
private

Evaluate basis functions on point z in interval [0, 1].

Parameters
[in]zPoint in [0, 1] where the functions are evaluated
Returns
The values of the contributing functions
+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Variable Documentation

◆ dchi

real, save pdespline::dchi
private

Mesh size in poloidal direction.

◆ dphi

real, save pdespline::dphi
private

Mesh size in toroidal direction.

◆ ds

real, save pdespline::ds
private

Mesh size in radial direction.

◆ mmax

integer, dimension(:,:), allocatable, save pdespline::mmax
private

maximum m value for Fourier solver

◆ mmin

integer, dimension(:,:), allocatable, save pdespline::mmin
private

minimum m value for Fourier solver

◆ nbmatf

integer, save, public pdespline::nbmatf

number of matrix per proc

◆ nchi

integer, save pdespline::nchi
private

Number of intervals in poloidal direction.

◆ nlocad

integer, dimension(:), allocatable, save pdespline::nlocad
private

local 1D addresses of contributions of one spline function

◆ nphip

integer, save pdespline::nphip
private

Number of local intervals in toroidal direction.

◆ nq

integer, save pdespline::nq
private

Number of quadrature points.

◆ ns

integer, save pdespline::ns
private

Number of intervals in radial direction.

◆ nsupp

integer, save pdespline::nsupp
private

support in number of intervals of the spline function

◆ p

integer, save pdespline::p
private

Spline order.

◆ smin

real, save pdespline::smin
private

Lower bound of radial direction.

◆ splines_g

real, dimension(:), allocatable, save pdespline::splines_g
private

basis functions at grid points

◆ splines_q

real, dimension(:,:), allocatable, save pdespline::splines_q
private

basis functions at quadrature points

◆ width_f

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.

◆ zqp

real, dimension(10), save pdespline::zqp
private

Quadrature points.

◆ zqw

real, dimension(10), save pdespline::zqw
private

Quadrature weights.

◆ zwavg

real, dimension(:,:), allocatable, save, public pdespline::zwavg

Weight function surface average.