V3FIT
Data Types | Functions/Subroutines | Variables
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  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.
 

Detailed Description

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

Function/Subroutine Documentation

◆ check_int()

logical function, private integration_path::check_int ( integer, intent(in)  expected,
integer, intent(in)  received,
integer, intent(in)  testNum,
character (len=*), intent(in)  name 
)
private

Check an integer value.

Checks that the expected value matches the recieved. Otherwise report an error.

Parameters
[in]expectedThe known value.
[in]receivedThe known test value.
[in]testNumThe number of the test.
[in]nameThe name of the test.
Returns
True if the check passes and false otherwise.

Definition at line 1066 of file integration_path.f.

◆ check_log()

logical function, private integration_path::check_log ( logical, intent(in)  expected,
logical, intent(in)  received,
integer, intent(in)  testNum,
character (len=*), intent(in)  name 
)
private

Check a logical value.

Checks that the expected value matches the recieved. Otherwise report an error.

Parameters
[in]expectedThe known value.
[in]receivedThe known test value.
[in]testNumThe number of the test.
[in]nameThe name of the test.
Returns
True if the check passes and false otherwise.

Definition at line 997 of file integration_path.f.

◆ check_real()

logical function, private integration_path::check_real ( real(rprec), intent(in)  expected,
real(rprec), intent(in)  received,
integer, intent(in)  testNum,
character (len=*), intent(in)  name 
)
private

Check a real value.

Checks that the expected value matches the recieved. Otherwise report an error.

Parameters
[in]expectedThe known value.
[in]receivedThe known test value.
[in]testNumThe number of the test.
[in]nameThe name of the test.
Returns
True if the check passes and false otherwise.

Definition at line 1030 of file integration_path.f.

◆ integrate()

real (rprec) function, private integration_path::integrate ( type (path_int_class), intent(in)  this,
character(len=1), dimension(:), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2,
  integration_function 
)
private

Line integrate between to points.

This chooses the specific integration method.

Parameters
[in]thisIn instance of a path_int_class instance.
[in]contextGeneric object that contains data for the integration function.
[in]vertex1Starting point.
[in]vertex2Ending point.
[in]dxinputoptional input allowing user to define dx
integration_functionFunction pointer that defines the integrand.
Returns
The path integrated value between the vertex1 and vertex2.

Definition at line 396 of file integration_path.f.

◆ integrate_add()

real (rprec) function, private integration_path::integrate_add ( type (path_int_class), intent(in)  this,
character(len=1), dimension(:), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2,
  integration_function 
)
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.

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

Definition at line 460 of file integration_path.f.

◆ integrate_gleg()

real (rprec) function, private integration_path::integrate_gleg ( type (path_int_class), intent(in)  this,
character(len=1), dimension(:), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2,
  integration_function 
)
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.

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

Definition at line 531 of file integration_path.f.

◆ integrate_hp_gleg()

real (rprec) function, private integration_path::integrate_hp_gleg ( type (path_int_class), intent(in)  this,
character(len=1), dimension(:), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2,
  integration_function 
)
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.

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

Definition at line 607 of file integration_path.f.

◆ 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.

Definition at line 253 of file integration_path.f.

◆ path_construct_int()

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.

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

Definition at line 116 of file integration_path.f.

◆ 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.

Definition at line 167 of file integration_path.f.

◆ path_destruct_int()

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.

Parameters
[in,out]thisA path_int_class instance.

Definition at line 192 of file integration_path.f.

◆ 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.

Definition at line 222 of file integration_path.f.

◆ 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.

Definition at line 770 of file integration_path.f.

◆ path_integrate()

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.

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

Definition at line 293 of file integration_path.f.

◆ path_search()

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.

Parameters
[in]pathStarting vertex to begin search.
search_functionFunction pointer that defines the search criteria.
[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.

Definition at line 344 of file integration_path.f.

◆ 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.

Definition at line 851 of file integration_path.f.

◆ search()

logical function, private integration_path::search ( character(len=1), dimension(:), intent(in)  context,
type (vertex), intent(in)  vertex1,
type (vertex), intent(in)  vertex2,
  search_function,
real (rprec), dimension(3), intent(out)  xcart 
)
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.

Parameters
[in]contextGeneric object that contains data for the search function.
[in]vertex1Starting point.
[in]vertex2Ending point.
search_functionFunction pointer that defines the search criteria.
[out]xcartPoint where the search criteria was found.
Returns
True if the criteria was met between vertex1 and vertex2.

Definition at line 691 of file integration_path.f.

◆ test_function()

real (rprec) function, private integration_path::test_function ( character (len=1), dimension(:), intent(in)  context,
real (rprec), dimension(3), intent(in)  xcart,
real (rprec), dimension(3), intent(in)  dxcart,
real (rprec), intent(in)  length,
real (rprec), intent(in)  dx 
)
private

Call back function to test the integration.

This defined a function of f(x) = 1.

Parameters
[in]contextGeneric object that contains data for the integration function.
[in]xcartCurrent point in the integration.
[in]dxcartVector direction of the current integration path.
[in]lengthLength along the current integration.
[in]dxScaler length of the current integration path.
Returns
One

Definition at line 1100 of file integration_path.f.