Parameters of the B-spline orbital. The intention is that the user declares this type in the main progam and specifies: A, B, C, order, l,m and ind, no_bps. This object can then be passed to bspl_functioninit, which checks the parameters and initializes the parameters in the particular bspl_function. The default values of some of the variables below have been chosen so that an error will be triggered by the member procedure check if they have not been set to meaningful values.
More...
|
| procedure | check (this) |
| | Checks that the parameters of this function are acceptable. Returns 0 on success, positive integer on fail; the integer is the number of the error message described in the routine. The parameters that are checked are: A, B, order, l,m, ind and no_bps. All the other data components are set by bspline_grid_objinit_grid and bspline_grid_objnormalize. The only value that is not checked is 'norm'.
|
| procedure | normalize (this, ind) |
| | Returns the norm of the radial B-spline with a given index.
|
| procedure | read (this, lunit, posit, pos_after_rw, version) |
| | Reads-in the BTO data from the stream file and position given.
|
| procedure | write (this, lunit, posit, pos_after_rw) |
| | Writes the BTO data into the stream file and position given.
|
| procedure | print (this) |
| | Prints the basic grid parameters to stdout.
|
| procedure | print_grid (this) |
| | Prints all grid parameters to stdout.
|
| procedure | init_grid (this, a, b, c, order, no_bsplines) |
| | This procedure sets up the grid of radial B-splines given the values: A, B, C order, no_bps. The variables initialized by this routine are: knots, n, no_knots, bcoef, work. Note that it is here where different breakpoint sequences can be implemented. At the moment only linear sequence is implemented. This routine does not have to be used, i.e. the user can set-up the grid parameters by hand to whatever grid needed. The points A and B delimit the radial range of the B-spline grid and the point C should be a point somewhere in this range. The B-spline grid is then set-up in such a way that the point C is forced to be one of the knot points for the B-spline grid. The intention was to interpret the point C as the R-matrix radius so that the B-spline grid can extend beyond it. This can be useful for codes based on methods other than the R-matrix. If fixing a gridpoint at C is not required, input C outside of the range (A,B)
|
| procedure | bspline_range (this, ind, r1, r2) |
| | For a given index of the radial B-spline it returns the values r1, r2 corresponding to its radial range.
|
| procedure | bspline_amplitude (this, r, norm, ind, der) |
| | Calculates reduced boundary amplitude of a B-spline function (i.e. without the l,m-dependent factors delta_{m,mp}*delta_{l,lp}). The last argument der is the order of the derivative requried. For a simple evaluation of the B-spline value set der = 0.
|
| procedure | get_unique_pairs (this, pairs_to_unique_pairs, n_unique_pairs) |
| | Returns an array mapping the product of a pair of B-splines to a unique pair of B-splines. If the pair of B-splines does not overlap the corresponding value in pairs_to_unique_pairs is zero.
|
| procedure | get_last_inner_bspline (this, r) |
| | Returns the index of the last B-spline which starts before a given radial point R. If R .le. 0.0_cfp then the last inner B-spline is equal to the total number of B-splines.
|
|
| real(kind=cfp) | a |
| | The range of the radial basis of which this function is a member.
|
| real(kind=cfp) | b |
| real(kind=cfp) | c = -1 |
| integer | no_bps = -1 |
| | Number of break-points in the basis.
|
| integer | order = -1 |
| | Order of the radial B-spline basis.
|
| integer | n |
| | Number of B-splines.
|
| integer | ind_0_der = 0 |
| | Index of the first radial B-spline whose first derivative at point r=A is 0.
|
| integer | ind_last_before_c = 0 |
| | Index of the last B-spline that starts before the point C. It is set to n in case C .le. 0.0_cfp.
|
| real(kind=cfp), dimension(:), allocatable | bcoef |
| | Array of coefficients in the B-spline basis that can be used to evaluate B-spline basis functions.
|
| integer | no_knots = 0 |
| | Number of knots.
|
| real(kind=cfp), dimension(:), allocatable | knots |
| | Array of knots.
|
| real(kind=cfp), dimension(:), allocatable | work |
| | Work array used for evaluation of this spline.
|
| integer | inbv = 0 |
| | Helper variable used for evaluation of B-splines.
|
| real(kind=cfp) | tol = 0.0_cfp |
| | Relative precision (tolerance) for numerical calculation of integrals involving this B-spline grid.
|
Parameters of the B-spline orbital. The intention is that the user declares this type in the main progam and specifies: A, B, C, order, l,m and ind, no_bps. This object can then be passed to bspl_functioninit, which checks the parameters and initializes the parameters in the particular bspl_function. The default values of some of the variables below have been chosen so that an error will be triggered by the member procedure check if they have not been set to meaningful values.
◆ bspline_amplitude()
| procedure bspline_grid_gbl::bspline_grid_obj::bspline_amplitude |
( |
class(bspline_grid_obj) | this, |
|
|
real(kind=cfp), intent(in) | r, |
|
|
real(kind=cfp), intent(in) | norm, |
|
|
integer, intent(in) | ind, |
|
|
integer, intent(in) | der ) |
Calculates reduced boundary amplitude of a B-spline function (i.e. without the l,m-dependent factors delta_{m,mp}*delta_{l,lp}). The last argument der is the order of the derivative requried. For a simple evaluation of the B-spline value set der = 0.
◆ bspline_range()
| procedure bspline_grid_gbl::bspline_grid_obj::bspline_range |
( |
class(bspline_grid_obj) | this, |
|
|
integer, intent(in) | ind, |
|
|
real(kind=cfp), intent(out) | r1, |
|
|
real(kind=cfp), intent(out) | r2 ) |
For a given index of the radial B-spline it returns the values r1, r2 corresponding to its radial range.
◆ check()
Checks that the parameters of this function are acceptable. Returns 0 on success, positive integer on fail; the integer is the number of the error message described in the routine. The parameters that are checked are: A, B, order, l,m, ind and no_bps. All the other data components are set by bspline_grid_objinit_grid and bspline_grid_objnormalize. The only value that is not checked is 'norm'.
◆ get_last_inner_bspline()
| procedure bspline_grid_gbl::bspline_grid_obj::get_last_inner_bspline |
( |
class(bspline_grid_obj) | this, |
|
|
real(kind=cfp), intent(in) | r ) |
Returns the index of the last B-spline which starts before a given radial point R. If R .le. 0.0_cfp then the last inner B-spline is equal to the total number of B-splines.
◆ get_unique_pairs()
| procedure bspline_grid_gbl::bspline_grid_obj::get_unique_pairs |
( |
class(bspline_grid_obj) | this, |
|
|
integer, dimension(:,:), allocatable | pairs_to_unique_pairs, |
|
|
integer, intent(out) | n_unique_pairs ) |
Returns an array mapping the product of a pair of B-splines to a unique pair of B-splines. If the pair of B-splines does not overlap the corresponding value in pairs_to_unique_pairs is zero.
◆ init_grid()
| procedure bspline_grid_gbl::bspline_grid_obj::init_grid |
( |
class(bspline_grid_obj) | this, |
|
|
real(kind=cfp), intent(in) | a, |
|
|
real(kind=cfp), intent(in) | b, |
|
|
real(kind=cfp), intent(in) | c, |
|
|
integer, intent(in) | order, |
|
|
integer, intent(in) | no_bsplines ) |
This procedure sets up the grid of radial B-splines given the values: A, B, C order, no_bps. The variables initialized by this routine are: knots, n, no_knots, bcoef, work. Note that it is here where different breakpoint sequences can be implemented. At the moment only linear sequence is implemented. This routine does not have to be used, i.e. the user can set-up the grid parameters by hand to whatever grid needed. The points A and B delimit the radial range of the B-spline grid and the point C should be a point somewhere in this range. The B-spline grid is then set-up in such a way that the point C is forced to be one of the knot points for the B-spline grid. The intention was to interpret the point C as the R-matrix radius so that the B-spline grid can extend beyond it. This can be useful for codes based on methods other than the R-matrix. If fixing a gridpoint at C is not required, input C outside of the range (A,B)
◆ normalize()
| procedure bspline_grid_gbl::bspline_grid_obj::normalize |
( |
class(bspline_grid_obj) | this, |
|
|
integer, intent(in) | ind ) |
Returns the norm of the radial B-spline with a given index.
◆ print()
| procedure bspline_grid_gbl::bspline_grid_obj::print |
( |
class(bspline_grid_obj) | this | ) |
|
Prints the basic grid parameters to stdout.
◆ print_grid()
| procedure bspline_grid_gbl::bspline_grid_obj::print_grid |
( |
class(bspline_grid_obj) | this | ) |
|
Prints all grid parameters to stdout.
◆ read()
| procedure bspline_grid_gbl::bspline_grid_obj::read |
( |
class(bspline_grid_obj) | this, |
|
|
integer, intent(in) | lunit, |
|
|
integer, intent(in) | posit, |
|
|
integer, intent(out) | pos_after_rw, |
|
|
character(len=line_len), intent(in) | version ) |
Reads-in the BTO data from the stream file and position given.
◆ write()
| procedure bspline_grid_gbl::bspline_grid_obj::write |
( |
class(bspline_grid_obj) | this, |
|
|
integer, intent(in) | lunit, |
|
|
integer, intent(in) | posit, |
|
|
integer, intent(out) | pos_after_rw ) |
Writes the BTO data into the stream file and position given.
| real(kind=cfp) bspline_grid_gbl::bspline_grid_obj::a |
The range of the radial basis of which this function is a member.
| real(kind=cfp) bspline_grid_gbl::bspline_grid_obj::b |
◆ bcoef
| real(kind=cfp), dimension(:), allocatable bspline_grid_gbl::bspline_grid_obj::bcoef |
Array of coefficients in the B-spline basis that can be used to evaluate B-spline basis functions.
| real(kind=cfp) bspline_grid_gbl::bspline_grid_obj::c = -1 |
◆ inbv
| integer bspline_grid_gbl::bspline_grid_obj::inbv = 0 |
Helper variable used for evaluation of B-splines.
◆ ind_0_der
| integer bspline_grid_gbl::bspline_grid_obj::ind_0_der = 0 |
Index of the first radial B-spline whose first derivative at point r=A is 0.
◆ ind_last_before_c
| integer bspline_grid_gbl::bspline_grid_obj::ind_last_before_c = 0 |
Index of the last B-spline that starts before the point C. It is set to n in case C .le. 0.0_cfp.
◆ knots
| real(kind=cfp), dimension(:), allocatable bspline_grid_gbl::bspline_grid_obj::knots |
| integer bspline_grid_gbl::bspline_grid_obj::n |
◆ no_bps
| integer bspline_grid_gbl::bspline_grid_obj::no_bps = -1 |
Number of break-points in the basis.
◆ no_knots
| integer bspline_grid_gbl::bspline_grid_obj::no_knots = 0 |
◆ order
| integer bspline_grid_gbl::bspline_grid_obj::order = -1 |
Order of the radial B-spline basis.
◆ tol
| real(kind=cfp) bspline_grid_gbl::bspline_grid_obj::tol = 0.0_cfp |
Relative precision (tolerance) for numerical calculation of integrals involving this B-spline grid.
◆ work
| real(kind=cfp), dimension(:), allocatable bspline_grid_gbl::bspline_grid_obj::work |
Work array used for evaluation of this spline.
The documentation for this type was generated from the following file: