Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Zoltan2::CoordinateTaskMapper< Adapter, part_t > Class Template Reference

#include <Zoltan2_TaskMapping.hpp>

Inheritance diagram for Zoltan2::CoordinateTaskMapper< Adapter, part_t >:
Inheritance graph
[legend]
Collaboration diagram for Zoltan2::CoordinateTaskMapper< Adapter, part_t >:
Collaboration graph
[legend]

Public Member Functions

void getProcTask (part_t *&proc_to_task_xadj_, part_t *&proc_to_task_adj_)
 
virtual void map (const RCP< MappingSolution< Adapter > > &mappingsoln)
 Mapping method. More...
 
virtual ~CoordinateTaskMapper ()
 
void create_local_task_to_rank (const lno_t num_local_coords, const part_t *local_coord_parts, const ArrayRCP< part_t > task_to_proc_)
 
 CoordinateTaskMapper (const Teuchos::RCP< const Teuchos::Comm< int > > comm_, const Teuchos::RCP< const MachineRepresentation< pcoord_t, part_t > > machine_, const Teuchos::RCP< const Adapter > input_adapter_, const Teuchos::RCP< const Zoltan2::PartitioningSolution< Adapter > > soln_, const Teuchos::RCP< const Environment > envConst, bool is_input_adapter_distributed=true, int num_ranks_per_node=1, bool divide_to_prime_first=false, bool reduce_best_mapping=true)
 Constructor. When this constructor is called, in order to calculate the communication metric, the task adjacency graph is created based on the coordinate model input and partitioning of it. If the communication graph is already calculated, use the other constructors. More...
 
 CoordinateTaskMapper (const Teuchos::RCP< const Teuchos::Comm< int > > comm_, const Teuchos::RCP< const MachineRepresentation< pcoord_t, part_t > > machine_, const Teuchos::RCP< const Adapter > input_adapter_, const part_t num_parts_, const part_t *result_parts, const Teuchos::RCP< const Environment > envConst, bool is_input_adapter_distributed=true, int num_ranks_per_node=1, bool divide_to_prime_first=false, bool reduce_best_mapping=true)
 Constructor. Instead of Solution we have two parameters, numparts. More...
 
 CoordinateTaskMapper (const Environment *env_const_, const Teuchos::Comm< int > *problemComm, int proc_dim, int num_processors, pcoord_t **machine_coords, int task_dim, part_t num_tasks, tcoord_t **task_coords, ArrayRCP< part_t >task_comm_xadj, ArrayRCP< part_t >task_comm_adj, pcoord_t *task_communication_edge_weight_, int recursion_depth, Kokkos::View< part_t *, Kokkos::HostSpace > part_no_array, const part_t *machine_dimensions, int num_ranks_per_node=1, bool divide_to_prime_first=false, bool reduce_best_mapping=true)
 Constructor The mapping constructor which will also perform the mapping operation. The result mapping can be obtained by –getAssignedProcForTask function: which returns the assigned processor id for the given task –getPartsForProc: which returns the assigned tasks with the number of tasks. More...
 
virtual size_t getLocalNumberOfParts () const
 Returns the number of parts to be assigned to this process. More...
 
pcoord_t ** shiftMachineCoordinates (int machine_dim, const part_t *machine_dimensions, bool *machine_extent_wrap_around, part_t numProcs, pcoord_t **mCoords)
 Using the machine dimensions provided, create virtual machine coordinates by assigning the largest gap to be as the wrap around link. More...
 
virtual void getProcsForPart (part_t taskId, part_t &numProcs, part_t *&procs) const
 getAssignedProcForTask function, returns the assigned tasks with the number of tasks. More...
 
part_t getAssignedProcForTask (part_t taskId)
 getAssignedProcForTask function, returns the assigned processor id for the given task More...
 
virtual void getPartsForProc (int procId, part_t &numParts, part_t *&parts) const
 getAssignedProcForTask function, returns the assigned tasks with the number of tasks. More...
 
ArrayView< part_tgetAssignedTasksForProc (part_t procId)
 
- Public Member Functions inherited from Zoltan2::PartitionMapping< Adapter >
 PartitionMapping (const Teuchos::RCP< const Teuchos::Comm< int > >comm_, const Teuchos::RCP< const Zoltan2::MachineRepresentation< pcoord_t, part_t > >machine_, const Teuchos::RCP< const Adapter > input_adapter_, const Teuchos::RCP< const Zoltan2::PartitioningSolution< Adapter > >soln_, const Teuchos::RCP< const Environment > envConst_)
 Constructor Constructor builds the map from parts to ranks. KDDKDD WILL NEED THE SOLUTION FOR INTELLIGENT MAPPING KDDKDD BUT MAY WANT TO SET PART SIZES BASED ON CAPABILITY OF A RANK. KDDKDD SO WHEN SHOULD THE MAP BE CREATED? More...
 
 PartitionMapping (const Teuchos::RCP< const Teuchos::Comm< int > >comm_, const Teuchos::RCP< const Zoltan2::MachineRepresentation< pcoord_t, part_t > >machine_, const Teuchos::RCP< const Adapter > input_adapter_, const part_t num_parts_, const part_t *result_parts, const Teuchos::RCP< const Environment > envConst_)
 
 PartitionMapping (const Teuchos::RCP< const Teuchos::Comm< int > >comm_, const Teuchos::RCP< const Environment > envConst_)
 
 PartitionMapping ()
 
 PartitionMapping (const Teuchos::RCP< const Environment >envConst_)
 
 PartitionMapping (const Teuchos::RCP< const Environment > envConst_, const Teuchos::RCP< const Teuchos::Comm< int > >comm_, const Teuchos::RCP< const MachineRepresentation< pcoord_t, part_t > >machine_)
 
virtual ~PartitionMapping ()
 
- Public Member Functions inherited from Zoltan2::Algorithm< Adapter >
virtual ~Algorithm ()
 
virtual int localOrder (const RCP< LocalOrderingSolution< lno_t > > &)
 Ordering method. More...
 
virtual int globalOrder (const RCP< GlobalOrderingSolution< gno_t > > &)
 Ordering method. More...
 
virtual void color (const RCP< ColoringSolution< Adapter > > &)
 Coloring method. More...
 
virtual void match ()
 Matching method. More...
 
virtual void partition (const RCP< PartitioningSolution< Adapter > > &)
 Partitioning method. More...
 
virtual void partitionMatrix (const RCP< MatrixPartitioningSolution< Adapter > > &)
 Matrix Partitioning method. More...
 
virtual bool isPartitioningTreeBinary () const
 return if algorithm determins tree to be binary More...
 
virtual void getPartitionTree (part_t, part_t &, std::vector< part_t > &, std::vector< part_t > &, std::vector< part_t > &, std::vector< part_t > &) const
 for partitioning methods, fill arrays with partition tree info More...
 
virtual std::vector
< coordinateModelPartBox > & 
getPartBoxesView () const
 for partitioning methods, return bounding boxes of the More...
 
virtual part_t pointAssign (int, scalar_t *) const
 pointAssign method: Available only for some partitioning algorithms More...
 
virtual void boxAssign (int, scalar_t *, scalar_t *, size_t &, part_t **) const
 boxAssign method: Available only for some partitioning algorithms More...
 
virtual void getCommunicationGraph (const PartitioningSolution< Adapter > *, ArrayRCP< part_t > &, ArrayRCP< part_t > &)
 returns serial communication graph of a computed partition More...
 
virtual int getRankForPart (part_t)
 In mapping, returns the rank to which a part is assigned. More...
 
virtual void getMyPartsView (part_t &, part_t *&)
 In mapping, returns a view of parts assigned to the current rank. More...
 

Protected Member Functions

void doMapping (int myRank, const Teuchos::RCP< const Teuchos::Comm< int > > comm_)
 doMapping function, calls getMapping function of communicationModel object. More...
 
RCP< Comm< int > > create_subCommunicator ()
 creates and returns the subcommunicator for the processor group. More...
 
void getBestMapping ()
 finds the lowest cost mapping and broadcasts solution to everyone. More...
 
void writeMapping ()
 
void writeMapping2 (int myRank)
 

Protected Attributes

ArrayRCP< part_tproc_to_task_xadj
 
ArrayRCP< part_tproc_to_task_adj
 
ArrayRCP< part_ttask_to_proc
 
ArrayRCP< part_tlocal_task_to_rank
 
bool isOwnerofModel
 
CoordinateCommunicationModel
< pcoord_t, tcoord_t, part_t,
node_t > * 
proc_task_comm
 
part_t nprocs
 
part_t ntasks
 
ArrayRCP< part_ttask_communication_xadj
 
ArrayRCP< part_ttask_communication_adj
 
ArrayRCP< scalar_ttask_communication_edge_weight
 

Additional Inherited Members

- Public Types inherited from Zoltan2::Algorithm< Adapter >
typedef Adapter::lno_t lno_t
 
typedef Adapter::gno_t gno_t
 
typedef Adapter::scalar_t scalar_t
 
typedef Adapter::part_t part_t
 
- Public Attributes inherited from Zoltan2::PartitionMapping< Adapter >
const Teuchos::RCP< const
Teuchos::Comm< int > > 
comm
 
const Teuchos::RCP< const
Zoltan2::MachineRepresentation
< pcoord_t, part_t > > 
machine
 
const Teuchos::RCP< const Adapter > input_adapter
 
const Teuchos::RCP< const
Zoltan2::PartitioningSolution
< Adapter > > 
soln
 
const Teuchos::RCP< const
Environment
env
 
const part_t num_parts
 
const part_tsolution_parts
 

Detailed Description

template<typename Adapter, typename part_t>
class Zoltan2::CoordinateTaskMapper< Adapter, part_t >

Definition at line 1750 of file Zoltan2_TaskMapping.hpp.

Constructor & Destructor Documentation

template<typename Adapter, typename part_t>
virtual Zoltan2::CoordinateTaskMapper< Adapter, part_t >::~CoordinateTaskMapper ( )
inlinevirtual

Definition at line 2218 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
Zoltan2::CoordinateTaskMapper< Adapter, part_t >::CoordinateTaskMapper ( const Teuchos::RCP< const Teuchos::Comm< int > >  comm_,
const Teuchos::RCP< const MachineRepresentation< pcoord_t, part_t > >  machine_,
const Teuchos::RCP< const Adapter >  input_adapter_,
const Teuchos::RCP< const Zoltan2::PartitioningSolution< Adapter > >  soln_,
const Teuchos::RCP< const Environment envConst,
bool  is_input_adapter_distributed = true,
int  num_ranks_per_node = 1,
bool  divide_to_prime_first = false,
bool  reduce_best_mapping = true 
)
inline

Constructor. When this constructor is called, in order to calculate the communication metric, the task adjacency graph is created based on the coordinate model input and partitioning of it. If the communication graph is already calculated, use the other constructors.

Parameters
comm_is the communication object.
machine_is the machineRepresentation object. Stores the coordinates of machines.
model_is the input adapter.
soln_is the solution object. Holds the assignment of points.
envConst_is the environment object.

Definition at line 2255 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
Zoltan2::CoordinateTaskMapper< Adapter, part_t >::CoordinateTaskMapper ( const Teuchos::RCP< const Teuchos::Comm< int > >  comm_,
const Teuchos::RCP< const MachineRepresentation< pcoord_t, part_t > >  machine_,
const Teuchos::RCP< const Adapter >  input_adapter_,
const part_t  num_parts_,
const part_t result_parts,
const Teuchos::RCP< const Environment envConst,
bool  is_input_adapter_distributed = true,
int  num_ranks_per_node = 1,
bool  divide_to_prime_first = false,
bool  reduce_best_mapping = true 
)
inline

Constructor. Instead of Solution we have two parameters, numparts.

When this constructor is called, in order to calculate the communication metric, the task adjacency graph is created based on the coordinate model input and partitioning of it. If the communication graph is already calculated, use the other constructors.

Parameters
comm_is the communication object.
machine_is the machineRepresentation object. Stores the coordinates of machines.
model_is the input adapter.
soln_is the solution object. Holds the assignment of points.
envConst_is the environment object.

Definition at line 2542 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
Zoltan2::CoordinateTaskMapper< Adapter, part_t >::CoordinateTaskMapper ( const Environment env_const_,
const Teuchos::Comm< int > *  problemComm,
int  proc_dim,
int  num_processors,
pcoord_t **  machine_coords,
int  task_dim,
part_t  num_tasks,
tcoord_t **  task_coords,
ArrayRCP< part_t task_comm_xadj,
ArrayRCP< part_t task_comm_adj,
pcoord_t *  task_communication_edge_weight_,
int  recursion_depth,
Kokkos::View< part_t *, Kokkos::HostSpace >  part_no_array,
const part_t machine_dimensions,
int  num_ranks_per_node = 1,
bool  divide_to_prime_first = false,
bool  reduce_best_mapping = true 
)
inline

Constructor The mapping constructor which will also perform the mapping operation. The result mapping can be obtained by –getAssignedProcForTask function: which returns the assigned processor id for the given task –getPartsForProc: which returns the assigned tasks with the number of tasks.

-task_comm_xadj, task_comm_adj, task_communication_edge_weight_
  can be provided NULL. In this case all processors will
  calculate the same mapping.
-If task_comm_xadj, task_comm_adj and provided, algorithm will
  perform rotations and processors will calculate different
  mappings, and best one will be reduced.
-If task_communication_edge_weight_ is provided with
  task_comm_xadj, task_comm_adj this will be used when cost is
  calculated.
-recursion_depth is a mandatory argument. In the case
  part_no_array is not null, this parameter should represent the
  length of part_no_array.
  If part_no_array is given as NULL, then this will give the
  recursion depth for the algorith,
  Maximum number is ceil(log_2(min(num_processors, num_tasks))),
  and providing a higher number will be equivalant to this.
  Partitioning algorithm will work as RCB when maximum number is
  given, which performs the best mapping results.
-part_no_array: The best results are obtained when this parameter
  is given as NULL. But if this is provided, partitioning will
  use this array for partitioning each dimension to the given
  numbers.
  The multiplication of these numbers should be equal to
  min(num_processors, num_tasks).
-machine_dimensions: This can be NULL, but if provided the
  algorithm will perform shift of the machine coords so that
  the largest gap is treated as wrap-around link.
Parameters
env_const_the environment object.
problemCommis the communication object.
proc_dimdimensions of the processor coordinates.
num_processorsis the number of processors
machine_coordsis the coordinates of the processors.
task_dimis the dimension of the tasks.
num_tasksis the number of tasks.
task_coordsis the coordinates of the tasks.
task_comm_xadjis the task communication graphs xadj array. (task i adjacency is between task_comm_xadj[i] and task_comm_xadj[i + 1])
task_comm_adjis task communication graphs adj array.
task_communication_edge_weight_is the weight of the communication in task graph.
recursion_depthis the recursion depth that will be applied to partitioning. If part_no_array is provided, then it is the length of this array.
part_no_arrayif part_no_array is provided, partitioning algorithm will be forced to use this array for partitioning. However, the multiplication of each entries in this array should be equal to min(num_processors, num_tasks).
*machine_dimensions,:the dimensions of the machine network. For example for hopper 17x8x24 This can be NULL, but if provided the algorithm will perform shift of the machine coords so that the largest gap is treated as wrap-around link.

Definition at line 2882 of file Zoltan2_TaskMapping.hpp.

Member Function Documentation

template<typename Adapter, typename part_t>
void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::doMapping ( int  myRank,
const Teuchos::RCP< const Teuchos::Comm< int > >  comm_ 
)
inlineprotected

doMapping function, calls getMapping function of communicationModel object.

Definition at line 1802 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
RCP<Comm<int> > Zoltan2::CoordinateTaskMapper< Adapter, part_t >::create_subCommunicator ( )
inlineprotected

creates and returns the subcommunicator for the processor group.

Definition at line 1829 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::getBestMapping ( )
inlineprotected

finds the lowest cost mapping and broadcasts solution to everyone.

Definition at line 1896 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::writeMapping ( )
inlineprotected

Definition at line 1939 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::writeMapping2 ( int  myRank)
inlineprotected

Definition at line 2008 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::getProcTask ( part_t *&  proc_to_task_xadj_,
part_t *&  proc_to_task_adj_ 
)
inline

Definition at line 2198 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
virtual void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::map ( const RCP< MappingSolution< Adapter > > &  )
inlinevirtual

Mapping method.

Reimplemented from Zoltan2::Algorithm< Adapter >.

Definition at line 2204 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::create_local_task_to_rank ( const lno_t  num_local_coords,
const part_t local_coord_parts,
const ArrayRCP< part_t task_to_proc_ 
)
inline

Definition at line 2227 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
virtual size_t Zoltan2::CoordinateTaskMapper< Adapter, part_t >::getLocalNumberOfParts ( ) const
inlinevirtual

Returns the number of parts to be assigned to this process.

Implements Zoltan2::PartitionMapping< Adapter >.

Definition at line 3010 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
pcoord_t** Zoltan2::CoordinateTaskMapper< Adapter, part_t >::shiftMachineCoordinates ( int  machine_dim,
const part_t machine_dimensions,
bool *  machine_extent_wrap_around,
part_t  numProcs,
pcoord_t **  mCoords 
)
inline

Using the machine dimensions provided, create virtual machine coordinates by assigning the largest gap to be as the wrap around link.

Parameters
machine_dim,:the number of dimensions in the machine network.
machine_dimensions,:the dimension of the machien network. For example for hopper, 17,8,24
numProcs,:the number of allocated processors.
mCoords,:allocated machine coordinates.

Definition at line 3024 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
virtual void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::getProcsForPart ( part_t  taskId,
part_t numProcs,
part_t *&  procs 
) const
inlinevirtual

getAssignedProcForTask function, returns the assigned tasks with the number of tasks.

Parameters
procIdprocId being queried.
numProcs(output), the number of processor the part is assigned to.
procs(output), the list of processors assigned to given part..

Implements Zoltan2::PartitionMapping< Adapter >.

Definition at line 3130 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
part_t Zoltan2::CoordinateTaskMapper< Adapter, part_t >::getAssignedProcForTask ( part_t  taskId)
inline

getAssignedProcForTask function, returns the assigned processor id for the given task

Parameters
taskIdtaskId being queried.

Definition at line 3140 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
virtual void Zoltan2::CoordinateTaskMapper< Adapter, part_t >::getPartsForProc ( int  procId,
part_t numParts,
part_t *&  parts 
) const
inlinevirtual

getAssignedProcForTask function, returns the assigned tasks with the number of tasks.

Parameters
procIdprocId being queried.
numParts(output), the number of parts the processor is assigned to.
parts(output), the list of parts assigned to given processor..

Implements Zoltan2::PartitionMapping< Adapter >.

Definition at line 3151 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
ArrayView<part_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::getAssignedTasksForProc ( part_t  procId)
inline

Definition at line 3161 of file Zoltan2_TaskMapping.hpp.

Member Data Documentation

template<typename Adapter, typename part_t>
ArrayRCP<part_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::proc_to_task_xadj
protected

Definition at line 1775 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
ArrayRCP<part_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::proc_to_task_adj
protected

Definition at line 1779 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
ArrayRCP<part_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::task_to_proc
protected

Definition at line 1783 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
ArrayRCP<part_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::local_task_to_rank
protected

Definition at line 1787 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
bool Zoltan2::CoordinateTaskMapper< Adapter, part_t >::isOwnerofModel
protected

Definition at line 1790 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
CoordinateCommunicationModel<pcoord_t,tcoord_t,part_t,node_t>* Zoltan2::CoordinateTaskMapper< Adapter, part_t >::proc_task_comm
protected

Definition at line 1791 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
part_t Zoltan2::CoordinateTaskMapper< Adapter, part_t >::nprocs
protected

Definition at line 1792 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
part_t Zoltan2::CoordinateTaskMapper< Adapter, part_t >::ntasks
protected

Definition at line 1793 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
ArrayRCP<part_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::task_communication_xadj
protected

Definition at line 1794 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
ArrayRCP<part_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::task_communication_adj
protected

Definition at line 1795 of file Zoltan2_TaskMapping.hpp.

template<typename Adapter, typename part_t>
ArrayRCP<scalar_t> Zoltan2::CoordinateTaskMapper< Adapter, part_t >::task_communication_edge_weight
protected

Definition at line 1796 of file Zoltan2_TaskMapping.hpp.


The documentation for this class was generated from the following file: