ALL 0.9.3
A Loadbalacing Library
|
The function API for ALL. More...
Data Types | |
interface | all_set_communicator |
type | all_t |
The object oriented API is this object. It contains all relevant functions. More... | |
Functions/Subroutines | |
subroutine, public | all_balance (this) |
Run loadbalancer and calculate new vertices. | |
integer function, public | all_error () |
character(len=all_error_length) function, public | all_error_description () |
subroutine, public | all_finalize (this) |
Delete the loadbalance object. | |
subroutine, public | all_get_dimension (this, dim) |
Get current dimension from loadbalancer. | |
subroutine, public | all_get_gamma (this, gamma) |
Retrieve currently set gamma value of balancer. | |
subroutine, public | all_get_length_of_work (this, length) |
Retrieve length of work array. | |
subroutine, public | all_get_neighbors (this, neighbors) |
Retrieve list of neighboring ranks (must be correct size already) | |
subroutine, public | all_get_number_of_neighbors (this, count) |
Retrieve number of neighbors (i.e. length of neighbors list) | |
subroutine, public | all_get_number_of_vertices (this, n) |
Retrieve number of vertices for current vertices. | |
subroutine, public | all_get_prev_vertices (this, vertices) |
Retrieve vertices from before loadbalancing. | |
subroutine, public | all_get_vertices (this, vertices) |
Retrieve current vertices. | |
subroutine, public | all_get_vertices_alloc (this, vertices) |
Same function as get_vertices, but takes an allocatable array, and resizes automatically. | |
subroutine, public | all_get_work (this, work) |
Retrieve first element of work array. | |
subroutine, public | all_get_work_array (this, work) |
Retrieve work array, which must already be the correct size. | |
subroutine, public | all_init (this, method, dim, gamma) |
Initialises the loadbalancer. | |
subroutine, public | all_print_vtk_outlines (this, step) |
Print VTK outlines (must be enabled in build step) | |
subroutine, public | all_print_vtk_vertices (this, step) |
Print VTK domain vertices (must be enabled in build step) | |
subroutine, public | all_reset_error () |
subroutine, public | all_set_gamma (this, gamma) |
Set gamma value for load balancer. | |
subroutine, public | all_set_method_data_histgram (this, nbins) |
Set number of bins for histogram method. | |
subroutine, public | all_set_min_domain_size (this, domain_size) |
Set a minimum domain size. | |
subroutine, public | all_set_proc_grid_params (this, loc, ranks) |
Set the grid parameters for this process. | |
subroutine, public | all_set_proc_tag (this, tag) |
Set process identifying tag for output. | |
subroutine, public | all_set_sys_size (this, syssize) |
Set size of system, which is required for the histogram method. | |
subroutine, public | all_set_vertices (this, vertices) |
Set new vertices. | |
subroutine, public | all_set_work (this, work) |
Set the work of this process. | |
subroutine, public | all_set_work_multi (this, work) |
Set multi dimensional work of this process. | |
subroutine, public | all_setup (this) |
Set up the loadbalancer. | |
Variables | |
integer(c_int), parameter, public | all_error_filesystem = 6 |
integer(c_int), parameter, public | all_error_generic = 1 |
integer(c_int), parameter, public | all_error_internal = 6 |
integer(c_int), parameter, public | all_error_invalidargument = 4 |
integer(c_int), parameter, public | all_error_invalidcommtype = 3 |
integer, parameter, public | all_error_length = 1024 |
integer(c_int), parameter, public | all_error_outofbounds = 5 |
integer(c_int), parameter, public | all_error_pointdimensionmissmatch = 2 |
integer(c_int), parameter, public | all_forcebased = 2 |
integer(c_int), parameter, public | all_histogram = 4 |
integer(c_int), parameter, public | all_staggered = 0 |
integer(c_int), parameter, public | all_tensor = 1 |
integer(c_int), parameter, public | all_tensor_max = 5 |
integer(c_int), parameter, public | all_unimplemented = 128 |
integer(c_int), parameter, public | all_voronoi = 3 |
The function API for ALL.
The *_c
functions should never be called directly, but their non suffixed counterparts should. Safety checks and conveniences provided by Fortran are available there. For a simple example on how to use this see examples/
directory.
subroutine, public all_module::all_balance | ( | class(all_t), intent(in) | this | ) |
Run loadbalancer and calculate new vertices.
Definition at line 412 of file ALL_module.F90.
integer function, public all_module::all_error |
Definition at line 532 of file ALL_module.F90.
character(len=all_error_length) function, public all_module::all_error_description |
Definition at line 539 of file ALL_module.F90.
subroutine, public all_module::all_finalize | ( | class(all_t), intent(in) | this | ) |
Delete the loadbalance object.
Definition at line 312 of file ALL_module.F90.
subroutine, public all_module::all_get_dimension | ( | class(all_t), intent(in) | this, |
integer(c_int), intent(out) | dim ) |
Get current dimension from loadbalancer.
Definition at line 470 of file ALL_module.F90.
subroutine, public all_module::all_get_gamma | ( | class(all_t), intent(in) | this, |
real(c_double), intent(out) | gamma ) |
Retrieve currently set gamma value of balancer.
Definition at line 417 of file ALL_module.F90.
subroutine, public all_module::all_get_length_of_work | ( | class(all_t), intent(in) | this, |
integer(c_int), intent(out) | length ) |
Retrieve length of work array.
Definition at line 476 of file ALL_module.F90.
subroutine, public all_module::all_get_neighbors | ( | class(all_t), intent(in) | this, |
integer(c_int), dimension(:), intent(out) | neighbors ) |
Retrieve list of neighboring ranks (must be correct size already)
Definition at line 506 of file ALL_module.F90.
subroutine, public all_module::all_get_number_of_neighbors | ( | class(all_t), intent(in) | this, |
integer(c_int), intent(out) | count ) |
Retrieve number of neighbors (i.e. length of neighbors list)
Definition at line 500 of file ALL_module.F90.
subroutine, public all_module::all_get_number_of_vertices | ( | class(all_t), intent(in) | this, |
integer(c_int), intent(out) | n ) |
Retrieve number of vertices for current vertices.
[out] | n | set to number of new vertices |
Definition at line 423 of file ALL_module.F90.
subroutine, public all_module::all_get_prev_vertices | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(:, :), intent(out) | vertices ) |
Retrieve vertices from before loadbalancing.
[out] | vertices | set to prev vertices, must be exact size (dim,n), unchecked! |
Definition at line 457 of file ALL_module.F90.
subroutine, public all_module::all_get_vertices | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(:, :), intent(out) | vertices ) |
Retrieve current vertices.
[out] | vertices | set to new vertices, must be exact size (dim,n) |
Definition at line 429 of file ALL_module.F90.
subroutine, public all_module::all_get_vertices_alloc | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(:, :), intent(inout), allocatable | vertices ) |
Same function as get_vertices, but takes an allocatable array, and resizes automatically.
[in,out] | vertices | set to new vertices, may be reallocated to fit |
Definition at line 443 of file ALL_module.F90.
subroutine, public all_module::all_get_work | ( | class(all_t), intent(in) | this, |
real(c_double), intent(out) | work ) |
Retrieve first element of work array.
Definition at line 482 of file ALL_module.F90.
subroutine, public all_module::all_get_work_array | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(:), intent(out) | work ) |
Retrieve work array, which must already be the correct size.
Definition at line 488 of file ALL_module.F90.
subroutine, public all_module::all_init | ( | class(all_t), intent(out) | this, |
integer, intent(in) | method, | ||
integer, intent(in) | dim, | ||
real(c_double), intent(in) | gamma ) |
Initialises the loadbalancer.
[out] | this | teh ALL object is returned |
[in] | method | Must be one of the ALL_* method values |
[in] | dim | dimensionality of system |
[in] | gamma | gamma value for load balancer (ignored for TENSOR and STAGGERED) |
Definition at line 303 of file ALL_module.F90.
subroutine, public all_module::all_print_vtk_outlines | ( | class(all_t), intent(in) | this, |
integer(c_int), intent(in) | step ) |
Print VTK outlines (must be enabled in build step)
[in] | step | current step, used for filename |
Definition at line 520 of file ALL_module.F90.
subroutine, public all_module::all_print_vtk_vertices | ( | class(all_t), intent(in) | this, |
integer(c_int), intent(in) | step ) |
Print VTK domain vertices (must be enabled in build step)
[in] | step | current step, used for filename |
Definition at line 526 of file ALL_module.F90.
subroutine, public all_module::all_reset_error |
Definition at line 536 of file ALL_module.F90.
subroutine, public all_module::all_set_gamma | ( | class(all_t), intent(in) | this, |
real(c_double), intent(in) | gamma ) |
Set gamma value for load balancer.
Definition at line 317 of file ALL_module.F90.
subroutine, public all_module::all_set_method_data_histgram | ( | class(all_t), intent(in) | this, |
integer(c_int), dimension(3), intent(in) | nbins ) |
Set number of bins for histogram method.
[in] | nbins | Number of bins per dimension |
Definition at line 401 of file ALL_module.F90.
subroutine, public all_module::all_set_min_domain_size | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(this%dim), intent(in) | domain_size ) |
Set a minimum domain size.
[in] | domain_size | minimum domain size |
Definition at line 336 of file ALL_module.F90.
subroutine, public all_module::all_set_proc_grid_params | ( | class(all_t), intent(in) | this, |
integer, dimension(this%dim), intent(in) | loc, | ||
integer, dimension(this%dim), intent(in) | ranks ) |
Set the grid parameters for this process.
[in] | loc | index of domain in dim directions (0-indexed!) |
[in] | ranks | total number of domains in dim directions |
Definition at line 323 of file ALL_module.F90.
subroutine, public all_module::all_set_proc_tag | ( | class(all_t), intent(in) | this, |
integer, intent(in) | tag ) |
Set process identifying tag for output.
[in] | tag | tag of local process, only output in VTK outlines |
Definition at line 330 of file ALL_module.F90.
subroutine, public all_module::all_set_sys_size | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(:), intent(in) | syssize ) |
Set size of system, which is required for the histogram method.
[in] | syssize | system size |
Definition at line 389 of file ALL_module.F90.
subroutine, public all_module::all_set_vertices | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(:, :), intent(in) | vertices ) |
Set new vertices.
[in] | vertices | vertices of domain, for n domains must have shape: (dim,n) |
Definition at line 362 of file ALL_module.F90.
subroutine, public all_module::all_set_work | ( | class(all_t), intent(in) | this, |
real(c_double), intent(in) | work ) |
Set the work of this process.
[in] | work | work of this domain |
Definition at line 342 of file ALL_module.F90.
subroutine, public all_module::all_set_work_multi | ( | class(all_t), intent(in) | this, |
real(c_double), dimension(:), intent(in) | work ) |
Set multi dimensional work of this process.
[in] | work | multi dimensional work of this domain |
Definition at line 348 of file ALL_module.F90.
subroutine, public all_module::all_setup | ( | class(all_t), intent(in) | this | ) |
Set up the loadbalancer.
Definition at line 407 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_error_filesystem = 6 |
Definition at line 59 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_error_generic = 1 |
Definition at line 53 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_error_internal = 6 |
Definition at line 58 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_error_invalidargument = 4 |
Definition at line 56 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_error_invalidcommtype = 3 |
Definition at line 55 of file ALL_module.F90.
integer, parameter, public all_module::all_error_length = 1024 |
Definition at line 61 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_error_outofbounds = 5 |
Definition at line 57 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_error_pointdimensionmissmatch = 2 |
Definition at line 54 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_forcebased = 2 |
Definition at line 45 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_histogram = 4 |
Definition at line 49 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_staggered = 0 |
Definition at line 43 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_tensor = 1 |
Definition at line 44 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_tensor_max = 5 |
Definition at line 50 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_unimplemented = 128 |
Definition at line 51 of file ALL_module.F90.
integer(c_int), parameter, public all_module::all_voronoi = 3 |
Definition at line 47 of file ALL_module.F90.