Stellarator-Tools
integration_path Module Reference

Module is part of the LIBSTELL. This modules contains code to define and integrate along an arbitray path. More...

Data Types

interface  integration_path_class
 Base class containing the parameters of the integration method to use. More...
 
interface  integration_path_gleg_class
 Subclass to use Gauss Legendre Quadrature. More...
 
interface  integration_path_hp_glep_class
 Subclass to use hp approch. More...
 
type  vertex
 A single point in space defined by an z, y, z coordinate. A vertex is structured as a singly linked list. More...
 
type  test_context
 A type for the test cases. More...
 
type  test_search_context
 A type for the test cases. More...
 
interface  path_construct
 Construction interface using path_construct_vertex. More...
 
interface  path_destruct
 Destruct interface using path_destruct_vertex. More...
 
interface  check
 Interface for checking the results of the unit tests. More...
 

Functions/Subroutines

class(integration_path_class) function, pointer make_integrator (method, npoints, length)
 Factory method to construct a path integrator using a method. More...
 
class(integration_path_class) function, pointer integration_path_construct ()
 Construct an integration_path_class object. More...
 
class(integration_path_gleg_class) function, pointer integration_path_gleg_construct (npoints)
 Construct an integration_path_class object. More...
 
class(integration_path_hp_glep_class) function, pointer integration_path_hp_glep_construct (npoints, length)
 Construct an integration_path_class object. More...
 
type(vertex) function, pointer path_construct_vertex (position)
 Construct a single vertex. More...
 
subroutine integration_path_gleg_destruct (this)
 Deconstruct a integration_path_gleg_class object. More...
 
recursive subroutine path_destruct_vertex (this)
 Deconstruct a vertex object. More...
 
recursive subroutine path_append_vertex (this, position)
 Append a vertex to a path. More...
 
recursive real(rprec) function integration_path_integrate_paths (this, path, context)
 Integrate along the paths. More...
 
recursive real(rprec) function, dimension(3) search_paths (path, context, found)
 Search along the path. More...
 
real(rprec) function integration_path_integrate_path (this, context, vertex1, vertex2)
 Line integrate between to points. More...
 
real(rprec) function integration_path_gleg_integrate_path (this, context, vertex1, vertex2)
 Line integrate between to points. More...
 
real(rprec) function integration_path_hp_gleg_integrate_path (this, context, vertex1, vertex2)
 Line integrate between to points. More...
 
subroutine path_get_gaussqad_weights (a, b, abscissas, weights)
 @breif Calculate the weights and abscissas for Gauss-Legendre integration. More...
 
logical function path_test ()
 Path unit test function. More...
 
logical function test_search_function (context, xcart1, xcart2)
 Call back function to test the search. More...
 

Variables

real(rprec), parameter path_default_dx = 0.0025
 Default step size of the integration.
 

Detailed Description

Module is part of the LIBSTELL. This modules contains code to define and integrate along an arbitray path.

Function/Subroutine Documentation

◆ integration_path_construct()

class (integration_path_class) function, pointer integration_path::integration_path_construct

Construct an integration_path_class object.

Returns
A pointer to a constructed integration_path_class object.

◆ integration_path_gleg_construct()

class (integration_path_gleg_class) function, pointer integration_path::integration_path_gleg_construct ( integer, intent(in)  npoints)

Construct an integration_path_class object.

Parameters
[in]npointsNumber of quadrature points to use.
Returns
A pointer to a constructed integration_path_class object.

◆ integration_path_gleg_destruct()

subroutine integration_path::integration_path_gleg_destruct ( type (integration_path_gleg_class), intent(inout)  this)

Deconstruct a integration_path_gleg_class object.

Deallocates memory and uninitializes a integration_path_gleg_class object.

Parameters
[in,out]thisA integration_path_gleg_class instance.

◆ integration_path_gleg_integrate_path()

real (rprec) function integration_path::integration_path_gleg_integrate_path ( class (integration_path_gleg_class), intent(in)  this,
class (integration_path_context_class), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2 
)

Line integrate between to points.

This chooses the specific integration method.

Parameters
[in]thisIn instance of a integration_path_class instance.
[in]contextGeneric object that contains data for the integration function.
[in]vertex1Starting point.
[in]vertex2Ending point.
Returns
The path integrated value between the vertex1 and vertex2.

◆ integration_path_hp_gleg_integrate_path()

real (rprec) function integration_path::integration_path_hp_gleg_integrate_path ( class (integration_path_hp_glep_class), intent(in)  this,
class (integration_path_context_class), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2 
)

Line integrate between to points.

This chooses the specific integration method.

Parameters
[in]thisIn instance of a integration_path_class instance.
[in]contextGeneric object that contains data for the integration function.
[in]vertex1Starting point.
[in]vertex2Ending point.
Returns
The path integrated value between the vertex1 and vertex2.

◆ integration_path_hp_glep_construct()

class (integration_path_hp_glep_class) function, pointer integration_path::integration_path_hp_glep_construct ( integer, intent(in)  npoints,
real (rprec), intent(in)  length 
)

Construct an integration_path_class object.

Parameters
[in]npointsNumber of quadrature points to use.
[in]lengthLength of the interval.
Returns
A pointer to a constructed integration_path_class object.

◆ integration_path_integrate_path()

real (rprec) function integration_path::integration_path_integrate_path ( class (integration_path_class), intent(in)  this,
class (integration_path_context_class), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2 
)

Line integrate between to points.

This chooses the specific integration method.

Parameters
[in]thisIn instance of a integration_path_class instance.
[in]contextGeneric object that contains data for the integration function.
[in]vertex1Starting point.
[in]vertex2Ending point.

◆ integration_path_integrate_paths()

recursive real (rprec) function integration_path::integration_path_integrate_paths ( class (integration_path_class), intent(in)  this,
type (vertex), intent(in)  path,
class (integration_path_context_class), intent(in)  context 
)

Integrate along the paths.

Recursively runs through the next vertex to find the last vertex. Once the last vertex is found, integrate alone that path. The integrand is proveded by means of call back function.

Parameters
[in]thisIn instance of a integration_path_class instance.
[in]pathStarting vertex to of the path to integrate.
[in]contextGeneric object that contains data for the integration function.
Returns
The total integrated path to the end.

◆ make_integrator()

class (integration_path_class) function, pointer integration_path::make_integrator ( character (len=*), intent(in)  method,
integer, intent(in)  npoints,
real (rprec), intent(in)  length 
)

Factory method to construct a path integrator using a method.

Parameters
[in]methodIntegartion method to use.
[in]npointsNumber of quadrature points to use.
[in]lengthLength of the interval.
Returns
A pointer to a constructed integration_path_class object.

◆ path_append_vertex()

recursive subroutine integration_path::path_append_vertex ( type (vertex), pointer  this,
real (rprec), dimension(3), intent(in)  position 
)

Append a vertex to a path.

Recursively runs through the next vertex to find the last vertex. Once the last vertex is found, a new vertex is allocated and appended to the path. This allow works as a constructer and allocates the first vertex if needed.

Parameters
[in,out]thisVertex to append path to.
[in]positionCartesian position of the vertex object.

◆ path_construct_vertex()

type (vertex) function, pointer integration_path::path_construct_vertex ( real (rprec), dimension(3), intent(in)  position)

Construct a single vertex.

Allocates memory and initializes a vertex object.

Parameters
[in]positionCartesian position of the vertex object.
Returns
A pointer to a constructed vertex object.

◆ path_destruct_vertex()

recursive subroutine integration_path::path_destruct_vertex ( type (vertex), pointer  this)

Deconstruct a vertex object.

Deallocates memory and uninitializes a vertex object. This recursively deconstructed the next vertex until the last in the linked list is found.

Parameters
[in,out]thisA vertex instance.

◆ path_get_gaussqad_weights()

subroutine integration_path::path_get_gaussqad_weights ( real (rprec), intent(in)  a,
real (rprec), intent(in)  b,
real(rprec), dimension(:), intent(out)  abscissas,
real(rprec), dimension(:), intent(out)  weights 
)

@breif Calculate the weights and abscissas for Gauss-Legendre integration.

This suprogram calculated the weights and abscissas for Gauss-Legendre integration. The subroutine is adapted from NIMROD which was adapted from Numerical Recipies, 2ed., Cambridge Press.

Parameters
[in]aStart point for integration.
[in]bEnd point for integration.
[out]abscissasArray of abscissas.
[out]weightsArray of weights.

◆ path_test()

logical function integration_path::path_test

Path unit test function.

This runs the associated unit tests and returns the result.

Returns
True if the tests pass and false otherwise.

◆ search_paths()

recursive real (rprec) function, dimension(3) integration_path::search_paths ( type (vertex), intent(in)  path,
class (search_path_context_class), intent(in)  context,
logical, intent(out)  found 
)

Search along the path.

Recursively runs through the next vertex to find the last vertex. Once the last vertex is found, a new vertex is allocated and appended to the path. The integrand is proveded by means of call back function.

Parameters
[in]pathStarting vertex to begin search.
[in]contextGeneric object that contains data for the integration function.
[out]foundSignals if the condition was met.
Returns
The vertex position along the path where the search condition was found.

◆ test_search_function()

logical function integration_path::test_search_function ( class (test_search_context), intent(in)  context,
real (rprec), dimension(3), intent(in)  xcart1,
real (rprec), dimension(3), intent(in)  xcart2 
)

Call back function to test the search.

Returns true if the positions bound the loaction we are looking for.

Parameters
[in]contextGeneric object that contains data for the integration function.
[in]xcart1Lower test point.
[in]xcart2Higher test point.