V3FIT
|
Module is part of the LIBSTELL. This modules contains code to define and integrate along an arbitray path. More...
Data Types | |
interface | check |
Interface for checking the results of the unit tests. More... | |
interface | path_construct |
Construction interface using either path_construct_int or path_construct_vertex. More... | |
interface | path_destruct |
Destruct interface using either path_destruct_int or path_destruct_vertex. More... | |
type | path_int_class |
Class containing the parameters of the integration method to use. 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... | |
Functions/Subroutines | |
type(path_int_class) function, pointer | path_construct_int (method, npoints, length) |
Construct a single path_int_class. More... | |
type(vertex) function, pointer | path_construct_vertex (position) |
Construct a single vertex. More... | |
subroutine | path_destruct_int (this) |
Deconstruct a path_int_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 | path_integrate (this, path, integration_function, context) |
Integrate along the path. More... | |
recursive real(rprec) function, dimension(3) | path_search (path, search_function, context, found) |
Search along the path. More... | |
real(rprec) function, private | integrate (this, context, vertex1, vertex2, integration_function) |
Line integrate between to points. More... | |
real(rprec) function, private | integrate_add (this, context, vertex1, vertex2, integration_function) |
Line integrate between to points. More... | |
real(rprec) function, private | integrate_gleg (this, context, vertex1, vertex2, integration_function) |
Integrate between to points using Guass-Legendre quadrature. More... | |
real(rprec) function, private | integrate_hp_gleg (this, context, vertex1, vertex2, integration_function) |
Integrate between to points using hp Guass-Legendre quadrature method. More... | |
logical function, private | search (context, vertex1, vertex2, search_function, xcart) |
Search line 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, private | check_log (expected, received, testNum, name) |
Check a logical value. More... | |
logical function, private | check_real (expected, received, testNum, name) |
Check a real value. More... | |
logical function, private | check_int (expected, received, testNum, name) |
Check an integer value. More... | |
real(rprec) function, private | test_function (context, xcart, dxcart, length, dx) |
Call back function to test the integration. More... | |
Variables | |
integer, parameter | path_none_type = -1 |
No qaudrature type. | |
integer, parameter | path_add_type = 0 |
Integrate the length in addition. | |
integer, parameter | path_gleg_type = 1 |
Integrate the path using Gauss Legendre Quadrature. | |
integer, parameter | path_hp_glep_type = 2 |
Integrate the path using a hp approch. | |
real(rprec), parameter | path_default_dx = 0.0025 |
Default step size of the integration. | |
Module is part of the LIBSTELL. This modules contains code to define and integrate along an arbitray path.
|
private |
Check an integer value.
Checks that the expected value matches the recieved. Otherwise report an error.
[in] | expected | The known value. |
[in] | received | The known test value. |
[in] | testNum | The number of the test. |
[in] | name | The name of the test. |
Definition at line 1066 of file integration_path.f.
|
private |
Check a logical value.
Checks that the expected value matches the recieved. Otherwise report an error.
[in] | expected | The known value. |
[in] | received | The known test value. |
[in] | testNum | The number of the test. |
[in] | name | The name of the test. |
Definition at line 997 of file integration_path.f.
|
private |
Check a real value.
Checks that the expected value matches the recieved. Otherwise report an error.
[in] | expected | The known value. |
[in] | received | The known test value. |
[in] | testNum | The number of the test. |
[in] | name | The name of the test. |
Definition at line 1030 of file integration_path.f.
|
private |
Line integrate between to points.
This chooses the specific integration method.
[in] | this | In instance of a path_int_class instance. |
[in] | context | Generic object that contains data for the integration function. |
[in] | vertex1 | Starting point. |
[in] | vertex2 | Ending point. |
[in] | dxinput | optional input allowing user to define dx |
integration_function | Function pointer that defines the integrand. |
Definition at line 396 of file integration_path.f.
|
private |
Line integrate between to points.
This divids the straight line path defined by two vertices and line integrates the function. The integrand is proveded by means of call back function.
[in] | this | In instance of a path_int_class instance. |
[in] | context | Generic object that contains data for the integration function. |
[in] | vertex1 | Starting point. |
[in] | vertex2 | Ending point. |
integration_function | Function pointer that defines the integrand. |
Definition at line 460 of file integration_path.f.
|
private |
Integrate between to points using Guass-Legendre quadrature.
This divids the straight line path defined by two vertices and line integrates the function. The integrand is proveded by means of call back function.
[in] | this | In instance of a path_int_class instance. |
[in] | context | Generic object that contains data for the integration function. |
[in] | vertex1 | Starting point. |
[in] | vertex2 | Ending point. |
integration_function | Function pointer that defines the integrand. |
Definition at line 531 of file integration_path.f.
|
private |
Integrate between to points using hp Guass-Legendre quadrature method.
This divids the straight line path defined by two vertices into a set of intergval of length len. Then Gauss-Legendre quadrature is used to integrate the function across each interval. The integrand is provided by means of call back function.
[in] | this | In instance of a path_int_class instance. |
[in] | context | Generic object that contains data for the integration function. |
[in] | vertex1 | Starting point. |
[in] | vertex2 | Ending point. |
integration_function | Function pointer that defines the integrand. |
Definition at line 607 of file integration_path.f.
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.
[in,out] | this | Vertex to append path to. |
[in] | position | Cartesian position of the vertex object. |
Definition at line 253 of file integration_path.f.
type (path_int_class) function, pointer integration_path::path_construct_int | ( | character (len=*), intent(in) | method, |
integer, intent(in) | npoints, | ||
real (rprec), intent(in) | length | ||
) |
Construct a single path_int_class.
Allocates memory and initializes a path_int_class object.
[in] | method | Integartion method to use. |
[in] | npoints | Number of quadrature points to use. |
[in] | length | Length of the interval. |
Definition at line 116 of file integration_path.f.
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.
[in] | position | Cartesian position of the vertex object. |
Definition at line 167 of file integration_path.f.
subroutine integration_path::path_destruct_int | ( | type (path_int_class), pointer | this | ) |
Deconstruct a path_int_class object.
Deallocates memory and uninitializes a path_int_class object.
[in,out] | this | A path_int_class instance. |
Definition at line 192 of file integration_path.f.
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.
[in,out] | this | A vertex instance. |
Definition at line 222 of file integration_path.f.
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.
[in] | a | Start point for integration. |
[in] | b | End point for integration. |
[out] | abscissas | Array of abscissas. |
[out] | weights | Array of weights. |
Definition at line 770 of file integration_path.f.
recursive real (rprec) function integration_path::path_integrate | ( | type (path_int_class), intent(in) | this, |
type (vertex), intent(in) | path, | ||
integration_function, | |||
character (len=1), dimension(:), intent(in) | context | ||
) |
Integrate along the path.
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.
[in] | this | In instance of a path_int_class instance. |
[in] | path | Starting vertex to of the path to integrate. |
integration_function | Function pointer that defines the integrand. | |
[in] | context | Generic object that contains data for the integration function. |
Definition at line 293 of file integration_path.f.
recursive real (rprec) function, dimension(3) integration_path::path_search | ( | type (vertex), intent(in) | path, |
search_function, | |||
character (len=1), dimension(:), 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.
[in] | path | Starting vertex to begin search. |
search_function | Function pointer that defines the search criteria. | |
[in] | context | Generic object that contains data for the integration function. |
[out] | found | Signals if the condition was met. |
Definition at line 344 of file integration_path.f.
logical function integration_path::path_test |
Path unit test function.
This runs the associated unit tests and returns the result.
Definition at line 851 of file integration_path.f.
|
private |
Search line between to points.
This divides the straight line path defined by two vertices and searched for a condition. The search criteria is proveded by means of call back function.
[in] | context | Generic object that contains data for the search function. |
[in] | vertex1 | Starting point. |
[in] | vertex2 | Ending point. |
search_function | Function pointer that defines the search criteria. | |
[out] | xcart | Point where the search criteria was found. |
Definition at line 691 of file integration_path.f.
|
private |
Call back function to test the integration.
This defined a function of f(x) = 1.
[in] | context | Generic object that contains data for the integration function. |
[in] | xcart | Current point in the integration. |
[in] | dxcart | Vector direction of the current integration path. |
[in] | length | Length along the current integration. |
[in] | dx | Scaler length of the current integration path. |
Definition at line 1100 of file integration_path.f.