Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Fekete quadrature rules for a triangle.
This module contains the Fekete quadrature rule adapted to triangles. The main functions are taken from the following site but have been modified to respect Selalib's structure: http://people.sc.fsu.edu/~jburkardt/f_src/triangle_fekete_rule/
Derived types and interfaces | |
interface | function_2D |
2d real function More... | |
Functions/Subroutines | |
subroutine | fekete_degree (rule, degree) |
returns the degree of a Fekete rule for the triangle. More... | |
subroutine, public | sll_s_fekete_order_num (rule, order_num) |
returns the order of a Fekete rule for the triangle. More... | |
subroutine | fekete_rule (rule, order_num, xy, w) |
returns the points and weights of a Fekete rule. More... | |
subroutine | rearrange_fekete_rule (n, xy, w) |
Re arranges the order of the quadrature point to go from lower-left to top right. More... | |
subroutine | fekete_rule_num (rule_num) |
returns the number of Fekete rules available. More... | |
subroutine | fekete_suborder (rule, suborder_num, suborder) |
returns the suborders for a Fekete rule. More... | |
subroutine | fekete_suborder_num (rule, suborder_num) |
returns the number of suborders for a Fekete rule. More... | |
subroutine | fekete_subrule (rule, suborder_num, suborder_xyz, suborder_w) |
returns a compressed Fekete rule More... | |
integer(kind=i32) function | wrapping (ival, ilo, ihi) |
forces an I4 to lie between given limits by wrapping. More... | |
subroutine | reference_to_physical_t3 (node_xy, n, ref, phy) |
maps T3 reference points to physical points. More... | |
real(kind=f64) function, dimension(:, :), allocatable, public | sll_f_fekete_points_and_weights (node_xy2, rule) |
Gives the fekete points coordinates and associated weights for a certain rule in a given triangle. More... | |
real(kind=f64) function, public | sll_f_fekete_integral (f, pxy) |
Fekete quadrature rule over a triangle. More... | |
subroutine | triangle_area (node_xy, area) |
subroutine | write_quadrature (rule) |
Writes fekete points coordinates of a hex-mesh reference triangle. More... | |
subroutine | write_basis_values (deg, rule) |
Writes on a file values of boxsplines on fekete points. More... | |
subroutine, public | sll_s_write_all_django_files (num_cells, deg, rule, transf) |
This function is supposed to write all django input files needed for a Django/Jorek simulation. More... | |
subroutine sll_m_fekete_integration::fekete_degree | ( | integer(kind=i32), intent(in) | rule, |
integer(kind=i32), intent(out) | degree | ||
) |
returns the degree of a Fekete rule for the triangle.
This code was first written by John Burkardt and is available online under the GNU LGPL license. Reference: Mark Taylor, Beth Wingate, Rachel Vincent, An Algorithm for Computing Fekete Points in the Triangle, SIAM Journal on Numerical Analysis, Volume 38, Number 5, 2000, pages 1707-1720.
[IN] | rule integer rule of quadrature for the fekete quadrature |
[OUT] | degree integer the polynomial degree of exactness of |
Definition at line 61 of file sll_m_fekete_integration.F90.
|
private |
returns the points and weights of a Fekete rule.
This code was first written by John Burkardt and is available online under the GNU LGPL license.
[IN] | rule integer the index of the rule |
[IN] | order_num integer the order (number of points of the rule) |
[OUT] | xy real table the points of the rule |
[OUT] | w real table the weigths of the rule |
Definition at line 125 of file sll_m_fekete_integration.F90.
|
private |
returns the number of Fekete rules available.
This code was first written by John Burkardt and is available online under the GNU LGPL license.
[OUT] | rule_num integer the number of rules available |
Definition at line 237 of file sll_m_fekete_integration.F90.
|
private |
returns the suborders for a Fekete rule.
This code was first written by John Burkardt and is available online under the GNU LGPL license.
[IN] | rule integer the index of the rule |
[IN] | suborder_num integer the number of suborder of the rule |
[OUT] | suborder integer array the suborders of the rule |
Definition at line 250 of file sll_m_fekete_integration.F90.
|
private |
returns the number of suborders for a Fekete rule.
This code was first written by John Burkardt and is available online under the GNU LGPL license.
[IN] | rule integer the index of the rule |
[OUTPUT] | suborder_num the number of suborder of the rule |
Definition at line 299 of file sll_m_fekete_integration.F90.
|
private |
returns a compressed Fekete rule
This code was first written by John Burkardt and is available online under the GNU LGPL license. Discussion: The listed weights are twice what we want...since we want them to sum to 1/2, reflecting the area of a unit triangle. So we simple halve the values before exiting this routine.
[IN] | rule integer the index of the rule |
[IN] | suborder_num the number of suborders |
[OUT] | suborder_xyz(3, suborder_num) real the number of suborders of the rule |
[OUT] | suborder_w(suborder_num) real the suborder weights |
Definition at line 339 of file sll_m_fekete_integration.F90.
|
private |
Re arranges the order of the quadrature point to go from lower-left to top right.
[IN] | n integer number of quadrature points |
[INOUT] | xy real table the points of the rule |
[INOUT] | w real table the weigths of the rule |
Definition at line 196 of file sll_m_fekete_integration.F90.
|
private |
maps T3 reference points to physical points.
This code was first written by John Burkardt and is available online under the GNU LGPL license. Given the vertices of an order 3 physical triangle and a point (XSI,ETA) in the reference triangle, the routine computes the value of the corresponding image point (X,Y) in physical space.
This routine is also appropriate for an order 4 triangle, as long as the fourth node is the centroid of the triangle.
This routine may also be appropriate for an order 6 triangle, if the mapping between reference and physical space is linear. This implies, in particular, that the sides of the image triangle are straight and that the "midside" nodes in the physical triangle are literally halfway along the sides of the physical triangle.
Reference Element T3:
| 1 3 | |\ | | \ S | \ | | \ | | \ 0 1--—2 | +–0–R–1-->
[IN] | node_xy(2,3) integer the coordinates of the vertices. The vertices are assumed to be the images of (0,0), (1,0) and (0,1). |
[IN] | n integer the number of objects to transform |
[IN] | ref(2,n) points in the reference triangle |
[OUT] | phy(2,n) corresponding points in the physical triangle |
Definition at line 748 of file sll_m_fekete_integration.F90.
real(kind=f64) function, public sll_m_fekete_integration::sll_f_fekete_integral | ( | procedure(function_2d) | f, |
real(kind=f64), dimension(2, 3), intent(in) | pxy | ||
) |
Fekete quadrature rule over a triangle.
To integrate the function \( f(x) \) (real-valued and over a triangle) we use the Fekete formula
\[ \int_{\Omega} f(x)dx \approx \sum_{k=1}^{N} w_k f(x_k) \]
The only quadrature rule possible for now is 1 (10 points)
[in] | f | the function to be integrated |
[in] | pxy | array of dimesion (2,3) containg the coordinates of the edges of the triangle |
Definition at line 816 of file sll_m_fekete_integration.F90.
real(kind=f64) function, dimension(:, :), allocatable, public sll_m_fekete_integration::sll_f_fekete_points_and_weights | ( | real(kind=f64), dimension(2, 3), intent(in) | node_xy2, |
integer(kind=i32), intent(in) | rule | ||
) |
Gives the fekete points coordinates and associated weights for a certain rule in a given triangle.
This code was first written by John Burkardt and is available online under the GNU LGPL license.
[IN] | node_xy2 array of dimesion (2,3) containg the coordinates of the edges of the triangle |
[IN] | rule integer quadrature rule |
Definition at line 775 of file sll_m_fekete_integration.F90.
subroutine, public sll_m_fekete_integration::sll_s_fekete_order_num | ( | integer(kind=i32), intent(in) | rule, |
integer(kind=i32), intent(out) | order_num | ||
) |
returns the order of a Fekete rule for the triangle.
This code was first written by John Burkardt and is available online under the GNU LGPL license. Reference: Mark Taylor, Beth Wingate, Rachel Vincent, An Algorithm for Computing Fekete Points in the Triangle, SIAM Journal on Numerical Analysis, Volume 38, Number 5, 2000, pages 1707-1720.
[IN] | rule integer the index of the rule |
[OUT] | order_num integer the order (number of points) of the rule |
Definition at line 99 of file sll_m_fekete_integration.F90.
subroutine, public sll_m_fekete_integration::sll_s_write_all_django_files | ( | integer(kind=i32), intent(in) | num_cells, |
integer(kind=i32), intent(in) | deg, | ||
integer(kind=i32), intent(in) | rule, | ||
character(len=*), intent(in) | transf | ||
) |
This function is supposed to write all django input files needed for a Django/Jorek simulation.
[in] | num_cells | integer number of cells in a radius of the hexagonal mesh |
[in] | deg | integer degree of the splines that will be used for the interpolation |
Definition at line 1032 of file sll_m_fekete_integration.F90.
|
private |
Definition at line 863 of file sll_m_fekete_integration.F90.
|
private |
forces an I4 to lie between given limits by wrapping.
This code was first written by John Burkardt and is available online under the GNU LGPL license. ILO = 4, IHI = 8
I Value
-2 8 -1 4 0 5 1 6 2 7 3 8 4 4 5 5 6 6 7 7 8 8 9 4 10 5 11 6 12 7 13 8 14 4
[IN] | ival an integer value |
[IN] | ilo the desired lower bound |
[IN] | ihi the desired upper bound |
[OUT] | res a wrapped version of ival |
Definition at line 691 of file sll_m_fekete_integration.F90.
|
private |
Writes on a file values of boxsplines on fekete points.
Following CAID structure, we write a file with the values of the basis function (box splines) on a reference element (triangle) fekete points. Output for DJANGO. Output file : boxsplines_basis_values.txt
[in] | deg | integer with degree of splines |
Number of non null box splines on a cell
Number of derivatives to be computed
The displament vector correspond to the translation done to obtain the other non null basis functions
Definition at line 945 of file sll_m_fekete_integration.F90.
|
private |
Writes fekete points coordinates of a hex-mesh reference triangle.
Takes the reference triangle of a hexmesh and computes the fekete points on it. Then it writes the results in a file following CAID/Django nomenclature. Output file : boxsplines_quadrature.txt
[in] | rule | integer for the fekete quadrature rule |
Definition at line 884 of file sll_m_fekete_integration.F90.