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

A PartitioningSolution is a solution to a partitioning problem. More...

#include <Zoltan2_PartitioningSolution.hpp>

Inheritance diagram for Zoltan2::PartitioningSolution< Adapter >:
Inheritance graph
[legend]
Collaboration diagram for Zoltan2::PartitioningSolution< Adapter >:
Collaboration graph
[legend]

Public Member Functions

 PartitioningSolution (const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, int nUserWeights, const RCP< Algorithm< Adapter > > &algorithm=Teuchos::null)
 Constructor when part sizes are not supplied. More...
 
 PartitioningSolution (const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, int nUserWeights, ArrayView< ArrayRCP< part_t > > reqPartIds, ArrayView< ArrayRCP< scalar_t > > reqPartSizes, const RCP< Algorithm< Adapter > > &algorithm=Teuchos::null)
 Constructor when part sizes are supplied. More...
 
size_t getTargetGlobalNumberOfParts () const
 Returns the global number of parts desired in the solution. More...
 
size_t getActualGlobalNumberOfParts () const
 Returns the actual global number of parts provided in setParts(). More...
 
size_t getLocalNumberOfParts () const
 Returns the number of parts to be assigned to this process. More...
 
scalar_t getLocalFractionOfPart () const
 If parts are divided across processes, return the fraction of a part on this process. More...
 
bool oneToOnePartDistribution () const
 Is the part-to-process distribution is one-to-one. More...
 
const int * getPartDistribution () const
 Return a distribution by part. More...
 
const part_tgetProcDistribution () const
 Return a distribution by process. More...
 
int getNumberOfCriteria () const
 Get the number of criteria (object weights) More...
 
bool criteriaHasUniformPartSizes (int idx) const
 Determine if balancing criteria has uniform part sizes. (User can specify differing part sizes.) More...
 
scalar_t getCriteriaPartSize (int idx, part_t part) const
 Get the size for a given weight index and a given part. More...
 
bool criteriaHaveSamePartSizes (int c1, int c2) const
 Return true if the two weight indices have the same part size information. More...
 
void setParts (ArrayRCP< part_t > &partList)
 The algorithm uses setParts to set the solution. More...
 
void RemapParts ()
 Remap a new partition for maximum overlap with an input partition. More...
 
long measure_stays (part_t *remap, int *idx, part_t *adj, long *wgt, part_t nrhs, part_t nlhs)
 
const RCP< const Comm< int > > & getCommunicator () const
 Return the communicator associated with the solution. More...
 
const RCP< const Environment > & getEnvironment () const
 Return the environment associated with the solution. More...
 
const part_tgetPartListView () const
 Returns the part list corresponding to the global ID list. More...
 
const int * getProcListView () const
 Returns the process list corresponding to the global ID list. More...
 
virtual bool isPartitioningTreeBinary () const
 calculate if partition tree is binary. More...
 
void getPartitionTree (part_t &numTreeVerts, std::vector< part_t > &permPartNums, std::vector< part_t > &splitRangeBeg, std::vector< part_t > &splitRangeEnd, std::vector< part_t > &treeVertParents) const
 get the partition tree - fill the relevant arrays More...
 
std::vector
< Zoltan2::coordinateModelPartBox > & 
getPartBoxesView () const
 returns the part box boundary list. More...
 
part_t pointAssign (int dim, scalar_t *point) const
 Return the part overlapping a given point in space;. More...
 
void boxAssign (int dim, scalar_t *lower, scalar_t *upper, size_t &nPartsFound, part_t **partsFound) const
 Return an array of all parts overlapping a given box in space. More...
 
void getCommunicationGraph (ArrayRCP< part_t > &comXAdj, ArrayRCP< part_t > &comAdj) const
 returns communication graph resulting from geometric partitioning. More...
 
void getPartsForProc (int procId, double &numParts, part_t &partMin, part_t &partMax) const
 Get the parts belonging to a process. More...
 
void getProcsForPart (part_t partId, part_t &procMin, part_t &procMax) const
 Get the processes containing a part. More...
 
- Public Member Functions inherited from Zoltan2::Solution
virtual ~Solution ()
 

Detailed Description

template<typename Adapter>
class Zoltan2::PartitioningSolution< Adapter >

A PartitioningSolution is a solution to a partitioning problem.

It is initialized by a PartitioningProblem, written to by an algorithm, and may be read by the user or by a data migration routine in an input adapter.

Todo:

Problem computes metrics using the Solution. Should Solution have a pointer to the metrics, since it may persist after the Problem is gone?

save an RCB tree, so it can be used in repartitioning, and supplied to the caller.

doxyfy the comments in this file.

Examples:
rcb_C.cpp.

Definition at line 55 of file Zoltan2_PartitioningSolution.hpp.

Constructor & Destructor Documentation

template<typename Adapter >
Zoltan2::PartitioningSolution< Adapter >::PartitioningSolution ( const RCP< const Environment > &  env,
const RCP< const Comm< int > > &  comm,
int  nUserWeights,
const RCP< Algorithm< Adapter > > &  algorithm = Teuchos::null 
)

Constructor when part sizes are not supplied.

The Solution constructor may require global communication. The rest of the Solution methods do not.

Parameters
envthe environment for the application
commthe communicator for the problem associated with this solution
nUserWeightsthe number of weights supplied by the application for each object.
algorithmAlgorithm, if any, used to compute the solution.

It is possible that part sizes were supplied on other processes, so this constructor does do a check to see if part sizes need to be globally calculated.

Definition at line 660 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter >
Zoltan2::PartitioningSolution< Adapter >::PartitioningSolution ( const RCP< const Environment > &  env,
const RCP< const Comm< int > > &  comm,
int  nUserWeights,
ArrayView< ArrayRCP< part_t > >  reqPartIds,
ArrayView< ArrayRCP< scalar_t > >  reqPartSizes,
const RCP< Algorithm< Adapter > > &  algorithm = Teuchos::null 
)

Constructor when part sizes are supplied.

The Solution constructor may require global communication. The rest of the Solution methods do not.

Parameters
envthe environment for the application
commthe communicator for the problem associated with this solution
nUserWeightsthe number of weights supplied by the application
reqPartIdsreqPartIds[i] is a list of of part numbers for weight index i.
reqPartSizesreqPartSizes[i] is the list of part sizes for weight i corresponding to parts in reqPartIds[i]
algorithmAlgorithm, if any, used to compute the solution.

If reqPartIds[i].size() and reqPartSizes[i].size() are zero for all processes, it is assumed that part sizes for weight dimension "i" are uniform.

If across the application there are some part numbers that are not included in the reqPartIds lists, then those part sizes are assumed to be 1.0.

Todo:
handle errors that may arise - like duplicate part numbers

Definition at line 693 of file Zoltan2_PartitioningSolution.hpp.

Member Function Documentation

template<typename Adapter>
size_t Zoltan2::PartitioningSolution< Adapter >::getTargetGlobalNumberOfParts ( ) const
inline

Returns the global number of parts desired in the solution.

Definition at line 165 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
size_t Zoltan2::PartitioningSolution< Adapter >::getActualGlobalNumberOfParts ( ) const
inline

Returns the actual global number of parts provided in setParts().

Definition at line 169 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
size_t Zoltan2::PartitioningSolution< Adapter >::getLocalNumberOfParts ( ) const
inline

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

Definition at line 173 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
scalar_t Zoltan2::PartitioningSolution< Adapter >::getLocalFractionOfPart ( ) const
inline

If parts are divided across processes, return the fraction of a part on this process.

Returns
zero if parts are not split across processes, approximate fraction of the part otherwise.
Todo:
More useful to get number of processes owning part? Or not useful at all - remove this?

Definition at line 182 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
bool Zoltan2::PartitioningSolution< Adapter >::oneToOnePartDistribution ( ) const
inline

Is the part-to-process distribution is one-to-one.

Returns
true if Process p owns part p for all p, and false if the part to process distribution is more complex.

If this is true, then getPartDistribution() and getProcDistribution() return NULL pointers. If either of the latter two methods is non-NULL, then this method returns false.

Definition at line 193 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
const int* Zoltan2::PartitioningSolution< Adapter >::getPartDistribution ( ) const
inline

Return a distribution by part.

Returns
If any parts are divided across processes, then a mapping A is returned. A such that A[i] is the lowest numbered process owning part i. The length of the array is one greater than the global number of parts. The value of the last element is the global number of processes.

Parts are divided across processes only if there are fewer parts than processes and the caller did not define "num_local_parts" for each process. In this case, parts are divided somewhat evenly across the processes. This situation is more likely to arise in Zoltan2 algorithms than in user applications.

If either oneToOnePartDistribution() is true or getProcDistribution() is non-NULL, then this method returns a NULL pointer. The three are mutually exclusive and collective exhaustive.

Definition at line 214 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
const part_t* Zoltan2::PartitioningSolution< Adapter >::getProcDistribution ( ) const
inline

Return a distribution by process.

Returns
If the mapping of parts to processes is not one-to-one, and if parts are not divided across processes, then the mapping A is returned. A such that A[i] is the first part assigned to process i. (Parts are assigned sequentially to processes.) However if A[i+1] is equal to A[i], there are no parts assigned to process i. The length of the array is one greater than the number of processes. The last element of the array is the global number of parts.

If the mapping is one-to-one, or if parts are divided across processes, then this method returns NULL pointer, and either oneToOnePartDistribution() or getPartDistribution() describes the mapping.

Definition at line 235 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
int Zoltan2::PartitioningSolution< Adapter >::getNumberOfCriteria ( ) const
inline

Get the number of criteria (object weights)

Returns
the number of criteria for which the solution has part sizes.

Definition at line 243 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
bool Zoltan2::PartitioningSolution< Adapter >::criteriaHasUniformPartSizes ( int  idx) const
inline

Determine if balancing criteria has uniform part sizes. (User can specify differing part sizes.)

Parameters
idxA value from 0 to one less than the number of weights per object.
Returns
true if part sizes are uniform for this criteria.

Definition at line 252 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
scalar_t Zoltan2::PartitioningSolution< Adapter >::getCriteriaPartSize ( int  idx,
part_t  part 
) const
inline

Get the size for a given weight index and a given part.

Parameters
idxA value from 0 to one less than the number of weights per object.
partA value from 0 to one less than the global number of parts to be computed
Returns
The size for that part. Part sizes for a given weight dimension sum to 1.0.
Todo:
It would be useful to algorithms to get the sum of part sizes from a to b, or the sum or a list of parts.

Definition at line 266 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter >
bool Zoltan2::PartitioningSolution< Adapter >::criteriaHaveSamePartSizes ( int  c1,
int  c2 
) const

Return true if the two weight indices have the same part size information.

Parameters
c1A value from 0 through one less than the number of weights.
c2A value from 0 through one less than the number of weights.
Returns
If weight index c1 and weight index c2 have the same part size information, the true is returned, otherwise false is returned.

It may be a problem for some algorithms if there are multiple weight dimensions with differing part size constraints to be satisfied.

Definition at line 1456 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter >
void Zoltan2::PartitioningSolution< Adapter >::setParts ( ArrayRCP< part_t > &  partList)

The algorithm uses setParts to set the solution.

Parameters
gnoListA list of global numbers.
partListThe part assigned to gnoList[i] by the algorithm should be in partList[i]. The partList is allocated and written by the algorithm.
dataDidNotMoveThe algorithm did not change the order of the data provided by the model; that is, it did not move the data to other processes or reorganize within the process. Thus, the gnoList and partList are ordered in the same way as the array provided by the model. Setting this flag to true avoids processing to remap the data to the original process that provided the data.

The global numbers supplied by the algorithm do not need to be those representing the global Ids of that process. But all global numbers should be assigned a part by exactly one process.

setParts() must be called by all processes in the problem, as the part for each global identifier supplied by each process in its InputAdapter is found and saved in this PartitioningSolution.

Definition at line 1255 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter >
void Zoltan2::PartitioningSolution< Adapter >::RemapParts ( )

Remap a new partition for maximum overlap with an input partition.

Assumptions for this version: input part assignment == processor rank for every local object. assuming nGlobalParts <= num ranks TODO: Write a version that takes the input part number as input; this change requires input parts in adapters to be provided in the Adapter. TODO: For repartitioning, compare to old remapping results; see Zoltan1.

Definition at line 1700 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
long Zoltan2::PartitioningSolution< Adapter >::measure_stays ( part_t remap,
int *  idx,
part_t adj,
long *  wgt,
part_t  nrhs,
part_t  nlhs 
)
inline

Definition at line 340 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
const RCP<const Comm<int> >& Zoltan2::PartitioningSolution< Adapter >::getCommunicator ( ) const
inline

Return the communicator associated with the solution.

Definition at line 364 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
const RCP<const Environment>& Zoltan2::PartitioningSolution< Adapter >::getEnvironment ( ) const
inline

Return the environment associated with the solution.

Definition at line 368 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
const part_t* Zoltan2::PartitioningSolution< Adapter >::getPartListView ( ) const
inline

Returns the part list corresponding to the global ID list.

Definition at line 372 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
const int* Zoltan2::PartitioningSolution< Adapter >::getProcListView ( ) const
inline

Returns the process list corresponding to the global ID list.

Returns
The return value is a NULL pointer if part IDs are synonomous with process IDs.

Definition at line 381 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
virtual bool Zoltan2::PartitioningSolution< Adapter >::isPartitioningTreeBinary ( ) const
inlinevirtual

calculate if partition tree is binary.

Definition at line 388 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
void Zoltan2::PartitioningSolution< Adapter >::getPartitionTree ( part_t numTreeVerts,
std::vector< part_t > &  permPartNums,
std::vector< part_t > &  splitRangeBeg,
std::vector< part_t > &  splitRangeEnd,
std::vector< part_t > &  treeVertParents 
) const
inline

get the partition tree - fill the relevant arrays

Definition at line 397 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
std::vector<Zoltan2::coordinateModelPartBox>& Zoltan2::PartitioningSolution< Adapter >::getPartBoxesView ( ) const
inline

returns the part box boundary list.

Definition at line 419 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
part_t Zoltan2::PartitioningSolution< Adapter >::pointAssign ( int  dim,
scalar_t *  point 
) const
inline

Return the part overlapping a given point in space;.

Definition at line 433 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
void Zoltan2::PartitioningSolution< Adapter >::boxAssign ( int  dim,
scalar_t *  lower,
scalar_t *  upper,
size_t &  nPartsFound,
part_t **  partsFound 
) const
inline

Return an array of all parts overlapping a given box in space.

Definition at line 458 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
void Zoltan2::PartitioningSolution< Adapter >::getCommunicationGraph ( ArrayRCP< part_t > &  comXAdj,
ArrayRCP< part_t > &  comAdj 
) const
inline

returns communication graph resulting from geometric partitioning.

Definition at line 473 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
void Zoltan2::PartitioningSolution< Adapter >::getPartsForProc ( int  procId,
double &  numParts,
part_t partMin,
part_t partMax 
) const
inline

Get the parts belonging to a process.

Parameters
procIda process rank
numPartson return will be set the number of parts belonging to the process.
partMinon return will be set to minimum part number
partMaxon return will be set to maximum part number

Normally numParts is at least one. But if there are more processes than parts, one of two things can happen. Either there are processes with no parts, and so numParts will be zero, or a part may be split across more than one process, in which numParts will be non-zero but less than 1.

In the latter case, numParts is 1.0 divided by the number of processes that share the part.

Definition at line 502 of file Zoltan2_PartitioningSolution.hpp.

template<typename Adapter>
void Zoltan2::PartitioningSolution< Adapter >::getProcsForPart ( part_t  partId,
part_t procMin,
part_t procMax 
) const
inline

Get the processes containing a part.

Parameters
partIda part number from 0 to one less than the global number of parts.
procMinon return will be set to minimum proc number
procMaxon return will be set to maximum proc number

Normally procMin and procMax are the same value and a part is assigned to one process. But if there are more processes than parts, it's possible that a part will be divided across more than one process.

Definition at line 522 of file Zoltan2_PartitioningSolution.hpp.


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