50 #include <Teuchos_DefaultComm.hpp> 
   51 #include <Teuchos_CommandLineProcessor.hpp> 
   52 #include <Teuchos_ParameterList.hpp> 
   55 int main(
int narg, 
char **arg)
 
   57   Tpetra::ScopeGuard tscope(&narg, &arg);
 
   58   Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
 
   62   int rank = comm->getRank();
 
   63   int nprocs = comm->getSize();
 
   66   int numGlobalIdentifiers = 100;
 
   67   int numMyIdentifiers = numGlobalIdentifiers / nprocs;
 
   68   if (rank < numGlobalIdentifiers % nprocs)
 
   69     numMyIdentifiers += 1;
 
   76   if (!myIds || !myWeights){
 
   84       std::cout << 
"Memory allocation failure" << std::endl;
 
   85       std::cout << 
"FAIL" << std::endl;
 
   91   for (
int i=0; i < numMyIdentifiers; i++){
 
   92     myIds[i] = myBaseId+i;
 
   93     myWeights[i] = rank%3 + 1;
 
   94     origsumwgts += myWeights[i];
 
   98   int *origcnt = 
new int[nprocs];
 
  100   Teuchos::gather<int, int>(&numMyIdentifiers, 1, origcnt, 1, 0, *comm);
 
  101   Teuchos::gather<int, zscalar_t>(&origsumwgts, 1, origwgts, 1, 0, *comm);
 
  103     std::cout << 
"BEFORE PART CNTS: ";
 
  104     for (
int i = 0; i < nprocs; i++)
 
  105       std::cout << origcnt[i] << 
" ";
 
  106     std::cout << std::endl;
 
  107     std::cout << 
"BEFORE PART WGTS: ";
 
  108     for (
int i = 0; i < nprocs; i++)
 
  109       std::cout << origwgts[i] << 
" ";
 
  110     std::cout << std::endl;
 
  116   std::vector<const zscalar_t *> weightValues;
 
  117   std::vector<int> weightStrides;   
 
  118   weightValues.push_back(const_cast<const zscalar_t *>(myWeights));
 
  126     new adapter_t(
zlno_t(numMyIdentifiers),myIds,weightValues,weightStrides);
 
  129   bool useWeights = 
true;
 
  130   Teuchos::CommandLineProcessor cmdp (
false, 
false);
 
  131   cmdp.setOption(
"weights", 
"no-weights", &useWeights,
 
  132                 "Indicated whether to use identifier weights in partitioning");
 
  133   cmdp.parse(narg, arg);
 
  135   Teuchos::ParameterList params(
"test parameters");
 
  137   params.set(
"num_global_parts", nprocs);
 
  138   params.set(
"algorithm", 
"block");
 
  139   params.set(
"partitioning_approach", 
"partition");
 
  140   if (!useWeights) params.set(
"partitioning_objective", 
"balance_object_count");
 
  150   quality_t *metricObject = 
new quality_t(adapter, ¶ms, comm, &solution);
 
  155   memset(totalWeight, 0, nprocs * 
sizeof(
zscalar_t));
 
  156   int *totalCnt = 
new int [nprocs];
 
  157   int *sumCnt = 
new int [nprocs];
 
  158   memset(totalCnt, 0, nprocs * 
sizeof(
int));
 
  161   zscalar_t libImbalance = metricObject->getWeightImbalance(0);
 
  165   for (
int i=0; i < numMyIdentifiers; i++){
 
  166     totalCnt[partList[i]]++;
 
  167     totalWeight[partList[i]] += myWeights[i];
 
  170   Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_SUM, nprocs, 
 
  172   Teuchos::reduceAll<int, zscalar_t>(*comm, Teuchos::REDUCE_SUM, nprocs, 
 
  173                                     totalWeight, sumWeight);
 
  178     std::cout << 
"AFTER PART CNTS: ";
 
  179     for (
int i=0; i < nprocs; i++)
 
  180       std::cout << sumCnt[i] << 
" ";
 
  181     std::cout << std::endl;
 
  184     std::cout << 
"AFTER PART WGTS: ";
 
  185     for (
int i=0; i < nprocs; i++){
 
  186       std::cout << sumWeight[i] << 
" ";
 
  187       total += sumWeight[i];
 
  189     std::cout << std::endl;
 
  195     for (
int i=0; i < nprocs; i++){
 
  197       if (sumWeight[i] > avg)
 
  198         imb = (sumWeight[i] - avg) / avg;
 
  200         imb = (avg - sumWeight[i]) / avg;
 
  207     std::cout << 
"Computed imbalance: " << imbalance << std::endl;
 
  208     std::cout << 
"Library's imbalance: " << libImbalance << std::endl;
 
  211     if (imbalance > libImbalance)
 
  212       err = imbalance - libImbalance;
 
  214       err = libImbalance - imbalance;
 
  227       std::cout << 
"failure in solution's imbalance data" << std::endl;
 
  228       std::cout << 
"FAIL" << std::endl;
 
  234     std::cout << 
"PASS" << std::endl;
 
  241   delete [] totalWeight;
 
int main(int narg, char *arg[])
 
Defines the PartitioningSolution class. 
 
common code used by tests 
 
SparseMatrixAdapter_t::part_t part_t
 
This class represents a collection of global Identifiers and their associated weights, if any. 
 
A PartitioningSolution is a solution to a partitioning problem. 
 
const part_t * getPartListView() const 
Returns the part list corresponding to the global ID list. 
 
static const std::string fail
 
Defines the BasicIdentifierAdapter class. 
 
int globalFail(const Comm< int > &comm, int fail)
 
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem. 
 
PartitioningProblem sets up partitioning problems for the user. 
 
Defines the PartitioningProblem class. 
 
A class that computes and returns quality metrics. 
 
void solve(bool updateInputData=true)
Direct the problem to create a solution.