| Zoltan2
    | 
Multi Jagged coordinate partitioning algorithm. More...
#include <Zoltan2_AlgMultiJagged.hpp>
| Public Member Functions | |
| AlgMJ () | |
| Multi Jagged coordinate partitioning algorithm default constructor.  More... | |
| void | multi_jagged_part (const RCP< const Environment > &env, RCP< const Comm< int > > &problemComm, double imbalance_tolerance, int num_teams, size_t num_global_parts, Kokkos::View< mj_part_t *, Kokkos::HostSpace > &part_no_array, int recursion_depth, int coord_dim, mj_lno_t num_local_coords, mj_gno_t num_global_coords, Kokkos::View< const mj_gno_t *, device_t > &initial_mj_gnos, Kokkos::View< mj_scalar_t **, Kokkos::LayoutLeft, device_t > &mj_coordinates, int num_weights_per_coord, Kokkos::View< bool *, Kokkos::HostSpace > &mj_uniform_weights, Kokkos::View< mj_scalar_t **, device_t > &mj_weights, Kokkos::View< bool *, Kokkos::HostSpace > &mj_uniform_parts, Kokkos::View< mj_part_t *, device_t > &result_assigned_part_ids, Kokkos::View< mj_gno_t *, device_t > &result_mj_gnos) | 
| Multi Jagged coordinate partitioning algorithm.  More... | |
| void | set_partitioning_parameters (bool distribute_points_on_cut_lines_, int max_concurrent_part_calculation_, int check_migrate_avoid_migration_option_, double minimum_migration_imbalance_, int migration_type_=0) | 
| Multi Jagged coordinate partitioning algorithm.  More... | |
| void | set_to_keep_part_boxes () | 
| Function call, if the part boxes are intended to be kept.  More... | |
| RCP< mj_partBox_t > | get_global_box () const | 
| DOCWORK: Documentation.  More... | |
| RCP< mj_partBoxVector_t > | get_kept_boxes () const | 
| DOCWORK: Documentation.  More... | |
| RCP< mj_partBoxVector_t > | compute_global_box_boundaries (RCP< mj_partBoxVector_t > &localPartBoxes) const | 
| DOCWORK: Documentation.  More... | |
| void | sequential_task_partitioning (const RCP< const Environment > &env, mj_lno_t num_total_coords, mj_lno_t num_selected_coords, size_t num_target_part, int coord_dim, Kokkos::View< mj_scalar_t **, Kokkos::LayoutLeft, device_t > &mj_coordinates_, Kokkos::View< mj_lno_t *, device_t > &initial_selected_coords_output_permutation, mj_lno_t *output_xadj, int recursion_depth_, const Kokkos::View< mj_part_t *, Kokkos::HostSpace > &part_no_array, bool partition_along_longest_dim, int num_ranks_per_node, bool divide_to_prime_first_, mj_part_t num_first_level_parts_=1, const Kokkos::View< mj_part_t *, Kokkos::HostSpace > &first_level_distribution_=Kokkos::View< mj_part_t *, Kokkos::HostSpace >()) | 
| Special function for partitioning for task mapping. Runs sequential, and performs deterministic partitioning for the partitioning the points along a cutline.  More... | |
Multi Jagged coordinate partitioning algorithm.
Definition at line 368 of file Zoltan2_AlgMultiJagged.hpp.
| Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::AlgMJ | ( | ) | 
Multi Jagged coordinate partitioning algorithm default constructor.
Definition at line 1410 of file Zoltan2_AlgMultiJagged.hpp.
| void Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::multi_jagged_part | ( | const RCP< const Environment > & | env, | 
| RCP< const Comm< int > > & | problemComm, | ||
| double | imbalance_tolerance_, | ||
| int | num_teams_, | ||
| size_t | num_global_parts_, | ||
| Kokkos::View< mj_part_t *, Kokkos::HostSpace > & | part_no_array_, | ||
| int | recursion_depth_, | ||
| int | coord_dim_, | ||
| mj_lno_t | num_local_coords_, | ||
| mj_gno_t | num_global_coords_, | ||
| Kokkos::View< const mj_gno_t *, device_t > & | initial_mj_gnos_, | ||
| Kokkos::View< mj_scalar_t **, Kokkos::LayoutLeft, device_t > & | mj_coordinates_, | ||
| int | num_weights_per_coord_, | ||
| Kokkos::View< bool *, Kokkos::HostSpace > & | mj_uniform_weights_, | ||
| Kokkos::View< mj_scalar_t **, device_t > & | mj_weights_, | ||
| Kokkos::View< bool *, Kokkos::HostSpace > & | mj_uniform_parts_, | ||
| Kokkos::View< mj_part_t *, device_t > & | result_assigned_part_ids_, | ||
| Kokkos::View< mj_gno_t *, device_t > & | result_mj_gnos_ | ||
| ) | 
Multi Jagged coordinate partitioning algorithm.
| env,: | library configuration and problem parameters | 
| problemComm,: | the communicator for the problem | 
| imbalance_tolerance,: | the input provided imbalance tolerance. | 
| num_teams,: | number of teams for CUDA kernels. | 
| num_global_parts,: | number of target global parts. | 
| part_no_array,: | part no array, if provided this will be used for partitioning. | 
| recursion_depth,: | if part no array is provided, it is the length of part no array, if part no is not provided than it is the number of steps that algorithm will divide into num_global_parts parts. | 
| coord_dim,: | coordinate dimension | 
| num_local_coords,: | number of local coordinates | 
| num_global_coords,: | number of global coordinates | 
| initial_mj_gnos,: | the list of initial global id's | 
| mj_coordinates,: | the two dimensional coordinate array. | 
| num_weights_per_coord,: | number of weights per coordinate | 
| mj_uniform_weights,: | if weight index [i] has uniform weight or not. | 
| mj_weights,: | the two dimensional array for weights | 
| mj_uniform_parts,: | if the target partitioning aims uniform parts | 
| result_assigned_part_ids,: | Output - the result partids corresponding to the coordinates given im result_mj_gnos. | 
| result_mj_gnos,: | Output - the result coordinate global id's corresponding to the part_ids array. | 
| env | library configuration and problem parameters | 
| problemComm | the communicator for the problem | 
| imbalance_tolerance | : the input provided imbalance tolerance. | 
| num_teams | : number of teams for CUDA kernels | 
| num_global_parts,: | number of target global parts. | 
| part_no_array,: | part no array, if provided this will be used for partitioning. | 
| recursion_depth,: | if part no array is provided, it is the length of part no array, if part no is not provided than it is the number of steps that algorithm will divide into num_global_parts parts. | 
| coord_dim,: | coordinate dimension | 
| num_local_coords,: | number of local coordinates | 
| num_global_coords,: | number of global coordinates | 
| initial_mj_gnos,: | the list of initial global id's | 
| mj_coordinates,: | the two dimensional coordinate array. | 
| num_weights_per_coord,: | number of weights per coordinate | 
| mj_uniform_weights,: | if weight index [i] has uniform weight or not. | 
| mj_weights,: | the two dimensional array for weights | 
| mj_uniform_parts,: | if the target partitioning aims uniform parts | 
| result_assigned_part_ids,: | Output - 1D pointer, should be provided as null. Memory is given in the function. the result partids corresponding to the coordinates given in result_mj_gnos. | 
| result_mj_gnos,: | Output - 1D pointer, should be provided as null. Memory is given in the function. the result coordinate global id's corresponding to the part_ids array. | 
Definition at line 7701 of file Zoltan2_AlgMultiJagged.hpp.
| void Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::set_partitioning_parameters | ( | bool | distribute_points_on_cut_lines_, | 
| int | max_concurrent_part_calculation_, | ||
| int | check_migrate_avoid_migration_option_, | ||
| double | minimum_migration_imbalance_, | ||
| int | migration_type_ = 0 | ||
| ) | 
Multi Jagged coordinate partitioning algorithm.
| distribute_points_on_cut_lines_ | : if partitioning can distribute points on same coordinate to different parts. | 
| max_concurrent_part_calculation_ | : how many parts we can calculate concurrently. | 
| check_migrate_avoid_migration_option_ | : whether to migrate=1, avoid migrate=2, or leave decision to MJ=0 | 
| minimum_migration_imbalance_ | : when MJ decides whether to migrate, the minimum imbalance for migration. | 
| migration_type_ | : when MJ migration whether to migrate for perfect load-imbalance or less messages | 
| distribute_points_on_cut_lines_ | : if partitioning can distribute points on same coordinate to different parts. | 
| max_concurrent_part_calculation_ | : how many parts we can calculate concurrently. | 
| check_migrate_avoid_migration_option_ | : whether to migrate=1, avoid migrate=2, or leave decision to MJ=0 | 
| minimum_migration_imbalance_ | : when MJ decides whether to migrate, the minimum imbalance for migration. | 
| migration_type | : whether to migrate for perfect load imbalance (0) or less messages. | 
Definition at line 7656 of file Zoltan2_AlgMultiJagged.hpp.
| void Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::set_to_keep_part_boxes | ( | ) | 
Function call, if the part boxes are intended to be kept.
Definition at line 2072 of file Zoltan2_AlgMultiJagged.hpp.
| RCP< typename AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::mj_partBox_t > Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::get_global_box | ( | ) | const | 
DOCWORK: Documentation.
Function returns the part boxes stored returns null if boxes are not stored, and prints warning mesage.
Definition at line 2062 of file Zoltan2_AlgMultiJagged.hpp.
| RCP< typename AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::mj_partBoxVector_t > Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::get_kept_boxes | ( | ) | const | 
DOCWORK: Documentation.
Definition at line 8300 of file Zoltan2_AlgMultiJagged.hpp.
| RCP< typename AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::mj_partBoxVector_t > Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::compute_global_box_boundaries | ( | RCP< mj_partBoxVector_t > & | localPartBoxes | ) | const | 
DOCWORK: Documentation.
Definition at line 8315 of file Zoltan2_AlgMultiJagged.hpp.
| void Zoltan2::AlgMJ< mj_scalar_t, mj_lno_t, mj_gno_t, mj_part_t, mj_node_t >::sequential_task_partitioning | ( | const RCP< const Environment > & | env, | 
| mj_lno_t | num_total_coords, | ||
| mj_lno_t | num_selected_coords, | ||
| size_t | num_target_part, | ||
| int | coord_dim_, | ||
| Kokkos::View< mj_scalar_t **, Kokkos::LayoutLeft, device_t > & | mj_coordinates_, | ||
| Kokkos::View< mj_lno_t *, device_t > & | initial_adjList_output_adjlist, | ||
| mj_lno_t * | output_xadj, | ||
| int | recursion_depth_, | ||
| const Kokkos::View< mj_part_t *, Kokkos::HostSpace > & | part_no_array_, | ||
| bool | partition_along_longest_dim, | ||
| int | num_ranks_per_node, | ||
| bool | divide_to_prime_first_, | ||
| mj_part_t | num_first_level_parts_ = 1, | ||
| const Kokkos::View< mj_part_t *, Kokkos::HostSpace > & | first_level_distribution_ = Kokkos::View<mj_part_t *, Kokkos::HostSpace>() | ||
| ) | 
Special function for partitioning for task mapping. Runs sequential, and performs deterministic partitioning for the partitioning the points along a cutline.
| env | library configuration and problem parameters | 
| num_total_coords | number of total coordinates | 
| num_selected_coords | : the number of selected coordinates. This is to set, if there are n processors, but only m<n processors are selected for mapping. | 
| num_target_part,: | number of target global parts. | 
| coord_dim_,: | coordinate dimension for coordinates | 
| mj_coordinates_,: | the coordinates | 
| initial_selected_coords_output_permutation,: | Array allocated by caller, in the size of num_total_coords, first num_selected_coords elements should list the indices of the selected processors. This is output for output permutation array. | 
| output_xadj,: | The output part xadj array, pointing beginning and end of each part on output permutation array (inital_adjList_output_adjlist). Returned in CSR format: part i's info in output_xadj[i] : output_xadj[i+1] | 
| recursion_depth,: | recursion depth | 
| part_no_array_,: | possibly null part_no_array, specifying how many parts each should be divided during partitioning. | 
| partition_along_longest_dim | DOCWORK: Documentation | 
| divide_to_prime_first_ | DOCWORK: Documentation Nonuniform first level partitioning (Currently available only for sequential_task_partitioning): Currently used for Dragonfly task mapping by partitioning Dragonfly RCA machine coordinates and application coordinates. An optimization that completely partitions the most important machine dimension first (i.e. the Dragonfly group coordinate, or RCA's x coordinate). The standard MJ alg follows after the nonuniform first level partitioning. | 
| num_first_level_parts_,: | the number of parts after the first level of partitioning (may be nonuniform) | 
| first_level_distribution_,: | a view containing the distribution of elements in each part after the first cut (used for nonuniform first cuts) Ex. (first level partitioning): If we have 120 elements, num_first_level_parts = 3, first_level_distribution = [4, 10, 6], then part sizes after first level will be [24, 60, 36]. Standard uniform MJ continues for all subsequent levels. | 
| env | library configuration and problem parameters | 
| num_total_coords | number of total coordinates | 
| num_selected_coords | : the number of selected coordinates. This is to set, if there are n processors, but only m<n processors are selected for mapping. | 
| num_target_part,: | number of target global parts. | 
| coord_dim_,: | coordinate dimension for coordinates | 
| mj_coordinates_,: | the coordinates | 
| inital_adjList_output_adjlist,: | Array allocated by caller, in the size of num_total_coords, first num_selected_coords elements should list the indices of the selected processors. This is output for output permutation array. | 
| output_xadj,: | The output part xadj array, pointing beginning and end of each part on output permutation array (inital_adjList_output_adjlist). Returned in CSR format: part i's info in output_xadj[i] : output_xadj[i+1] | 
| rd,: | recursion depth | 
| part_no_array_,: | possibly null part_no_array, specifying how many parts each should be divided during partitioning. | 
| partition_along_longest_dim | DOCWORK: Documentation | 
| num_ranks_per_node | DOCWORK: Documentation | 
| divide_to_prime_first_ | DOCWORK: Documentation | 
Nonuniform first level partitioning: Currently used for Dragonfly task mapping by partitioning Dragonfly RCA machine coordinates and application coordinates. An optimization that completely partitions the most important machine dimension first (i.e. the Dragonfly group coordinate, or RCA's x coordinate). The standard MJ alg follows after the nonuniform first level partitioning.
| num_target_first_level_parts,: | the number of parts requested after the first level of partitioning (resulting parts may be imbalanced) | 
| first_level_dist,: | an array requesting the distribution of elements in each part after the first cut (used for nonuniform first cuts) | 
Ex. (first level partitioning): If we have 120 elements, num_first_level_parts = 3, first_level_distribution = [4, 10, 6], then part sizes after first level will be [24, 60, 36]. Standard uniform MJ continues for all subsequent levels.
Definition at line 1483 of file Zoltan2_AlgMultiJagged.hpp.
 1.8.5
 1.8.5