Zoltan2
Zoltan2_Metric.hpp
Go to the documentation of this file.
00001 // @HEADER
00002 //
00003 // ***********************************************************************
00004 //
00005 //   Zoltan2: A package of combinatorial algorithms for scientific computing
00006 //                  Copyright 2012 Sandia Corporation
00007 //
00008 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00009 // the U.S. Government retains certain rights in this software.
00010 //
00011 // Redistribution and use in source and binary forms, with or without
00012 // modification, are permitted provided that the following conditions are
00013 // met:
00014 //
00015 // 1. Redistributions of source code must retain the above copyright
00016 // notice, this list of conditions and the following disclaimer.
00017 //
00018 // 2. Redistributions in binary form must reproduce the above copyright
00019 // notice, this list of conditions and the following disclaimer in the
00020 // documentation and/or other materials provided with the distribution.
00021 //
00022 // 3. Neither the name of the Corporation nor the names of the
00023 // contributors may be used to endorse or promote products derived from
00024 // this software without specific prior written permission.
00025 //
00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00037 //
00038 // Questions? Contact Karen Devine      (kddevin@sandia.gov)
00039 //                    Erik Boman        (egboman@sandia.gov)
00040 //                    Siva Rajamanickam (srajama@sandia.gov)
00041 //
00042 // ***********************************************************************
00043 //
00044 // @HEADER
00045 
00051 #ifndef ZOLTAN2_METRIC_HPP
00052 #define ZOLTAN2_METRIC_HPP
00053 
00054 #include <Zoltan2_StridedData.hpp>
00055 #include <Zoltan2_PartitioningSolution.hpp>
00056 
00057 #include <Epetra_SerialDenseVector.h>
00058 #include <Epetra_DataAccess.h>
00059 
00060 #include <cmath>
00061 #include <iomanip>
00062 #include <iostream>
00063 
00064 namespace Zoltan2{
00065 
00067 // Classes
00069 
00073 template <typename scalar_t>
00074   class MetricValues{
00075 
00076 private:
00077   void resetValues(){
00078     scalar_t *tmp = new scalar_t [evalNumMetrics];
00079     memset(tmp, 0, sizeof(scalar_t) * evalNumMetrics);
00080     values_ = arcp(tmp, 0, evalNumMetrics, true);
00081   }
00082   ArrayRCP<scalar_t> values_;
00083   std::string metricName_;
00084   multiCriteriaNorm mcnorm_;   // store "actualNorm + 1"
00085 
00086 public:
00087 
00095 enum metricOffset{
00096   evalLocalSum,    
00097   evalGlobalSum,   
00098   evalGlobalMin,   
00099   evalGlobalMax,   
00100   evalGlobalAvg,   
00101   evalMinImbalance, 
00102   evalAvgImbalance, 
00103   evalMaxImbalance, 
00104   evalNumMetrics    
00105 };
00106 
00114 static void printHeader(std::ostream &os);
00115 
00117 void printLine(std::ostream &os) const;
00118 
00120 MetricValues(std::string mname) : 
00121   values_(), metricName_(mname), mcnorm_(multiCriteriaNorm(0)) { 
00122     resetValues();}
00123 
00125 MetricValues() : 
00126   values_(), metricName_("unset"), mcnorm_(multiCriteriaNorm(0)) { 
00127     resetValues();}
00128 
00130 void setName(std::string name) { metricName_ = name;}
00131 
00133 void setNorm(multiCriteriaNorm normVal) { 
00134   mcnorm_ = multiCriteriaNorm(normVal+1);}
00135 
00137 void setLocalSum(scalar_t x) { values_[evalLocalSum] = x;}
00138 
00140 void setGlobalSum(scalar_t x) { values_[evalGlobalSum] = x;}
00141 
00143 void setGlobalMin(scalar_t x) { values_[evalGlobalMin] = x;}
00144 
00146 void setGlobalMax(scalar_t x) { values_[evalGlobalMax] = x;}
00147 
00149 void setGlobalAvg(scalar_t x) { values_[evalGlobalAvg] = x;}
00150 
00152 void setMinImbalance(scalar_t x) { values_[evalMinImbalance] = x;}
00153 
00157 void setMaxImbalance(scalar_t x) { values_[evalMaxImbalance] = x;}
00158 
00160 void setAvgImbalance(scalar_t x) { values_[evalAvgImbalance] = x;}
00161 
00163 const std::string &getName() const { return metricName_; }
00164 
00166 multiCriteriaNorm getNorm() { return multiCriteriaNorm(mcnorm_-1);}
00167 
00169 scalar_t getLocalSum() const { return values_[evalLocalSum];}
00170 
00172 scalar_t getGlobalSum() const { return values_[evalGlobalSum];}
00173 
00175 scalar_t getGlobalMin() const { return values_[evalGlobalMin];}
00176 
00178 scalar_t getGlobalMax() const { return values_[evalGlobalMax];}
00179 
00181 scalar_t getGlobalAvg() const { return values_[evalGlobalAvg];}
00182 
00184 scalar_t getMinImbalance() const { return values_[evalMinImbalance];}
00185 
00189 scalar_t getMaxImbalance() const { return values_[evalMaxImbalance];}
00190 
00192 scalar_t getAvgImbalance() const { return values_[evalAvgImbalance];}
00193 
00194 };  // end class
00195 
00196 
00197 template <typename scalar_t>
00198   void MetricValues<scalar_t>::printLine(std::ostream &os) const
00199 {
00200   std::string label(metricName_);
00201   if (mcnorm_ > 0){
00202     multiCriteriaNorm realNorm = multiCriteriaNorm(mcnorm_ - 1);
00203     std::ostringstream oss;
00204     switch (realNorm){
00205       case normMinimizeTotalWeight:   // 1-norm = Manhattan norm 
00206         oss << metricName_ << " (1)";
00207         break;
00208       case normBalanceTotalMaximum:   // 2-norm = sqrt of sum of squares
00209         oss << metricName_ << " (2)";
00210         break;
00211       case normMinimizeMaximumWeight: // inf-norm = maximum norm 
00212         oss << metricName_ << " (inf)";
00213         break;
00214       default:
00215         oss << metricName_ << " (?)";
00216         break;
00217     }
00218 
00219     label = oss.str();
00220   }
00221 
00222   os << std::setw(20) << label;
00223   os << std::setw(12) << std::setprecision(4) << values_[evalGlobalMin];
00224   os << std::setw(12) << std::setprecision(4) << values_[evalGlobalMax];
00225   os << std::setw(12) << std::setprecision(4) << values_[evalGlobalAvg];
00226   os << std::setw(2) << " ";
00227   os << std::setw(6) << std::setprecision(4) << values_[evalMinImbalance];
00228   os << std::setw(6) << std::setprecision(4) << values_[evalMaxImbalance];
00229   os << std::setw(6) << std::setprecision(4) << values_[evalAvgImbalance];
00230   os << std::endl;
00231 }
00232 
00233 template <typename scalar_t>
00234   void MetricValues<scalar_t>::printHeader(std::ostream &os)
00235 {
00236   os << std::setw(20) << " ";
00237   os << std::setw(36) << "------------SUM PER PART-----------";
00238   os << std::setw(2) << " ";
00239   os << std::setw(18) << "IMBALANCE PER PART";
00240   os << std::endl;
00241 
00242   os << std::setw(20) << " ";
00243   os << std::setw(12) << "min" << std::setw(12) << "max" << std::setw(12) << "avg";
00244   os << std::setw(2) << " ";
00245   os << std::setw(6) << "min" << std::setw(6) << "max" << std::setw(6) << "avg";
00246   os << std::endl;
00247 }
00248 
00250 // Namespace methods to compute metric values
00252 
00262 template <typename scalar_t>
00263  void getStridedStats(const ArrayView<scalar_t> &v, int stride, 
00264    int offset, scalar_t &min, scalar_t &max, scalar_t &sum)
00265 {
00266   if (v.size() < 1) return;
00267   min = max = sum = v[offset];
00268 
00269   for (int i=offset+stride; i < v.size(); i += stride){
00270     if (v[i] < min) min = v[i];
00271     else if (v[i] > max) max = v[i];
00272     sum += v[i];
00273   }
00274 }
00275 
00294 template <typename scalar_t, typename pnum_t, typename lno_t, typename part_t>
00295   void normedPartWeights(
00296     const RCP<const Environment> &env,
00297     part_t numberOfParts,
00298     const ArrayView<const pnum_t> &parts,
00299     const ArrayView<StridedData<lno_t, scalar_t> > &vwgts,
00300     multiCriteriaNorm mcNorm,
00301     scalar_t *weights)
00302 {
00303   env->localInputAssertion(__FILE__, __LINE__, "parts or weights", 
00304     numberOfParts > 0 && vwgts.size() > 0, BASIC_ASSERTION);
00305 
00306   int numObjects = parts.size();
00307   int vwgtDim = vwgts.size();
00308 
00309   memset(weights, 0, sizeof(scalar_t) * numberOfParts);
00310 
00311   if (numObjects == 0)
00312     return;
00313 
00314   if (vwgtDim == 0) {
00315     for (lno_t i=0; i < numObjects; i++){
00316       weights[parts[i]]++;
00317     }
00318   }
00319   else if (vwgtDim == 1){
00320     for (lno_t i=0; i < numObjects; i++){
00321       weights[parts[i]] += vwgts[0][i];
00322     }
00323   }
00324   else{
00325     switch (mcNorm){
00326       case normMinimizeTotalWeight:   
00327         for (int wdim=0; wdim < vwgtDim; wdim++){
00328           for (lno_t i=0; i < numObjects; i++){
00329             weights[parts[i]] += vwgts[wdim][i];
00330           }
00331         }  // next weight index
00332         break;
00333        
00334       case normBalanceTotalMaximum:   
00335         for (lno_t i=0; i < numObjects; i++){
00336           scalar_t ssum = 0;
00337           for (int wdim=0; wdim < vwgtDim; wdim++)
00338             ssum += (vwgts[wdim][i] * vwgts[wdim][i]);
00339           weights[parts[i]] += sqrt(ssum);
00340         }
00341         break;
00342 
00343       case normMinimizeMaximumWeight: 
00344         for (lno_t i=0; i < numObjects; i++){
00345           scalar_t max = 0;
00346           for (int wdim=0; wdim < vwgtDim; wdim++)
00347             if (vwgts[wdim][i] > max)
00348               max = vwgts[wdim][i];
00349           weights[parts[i]] += max;
00350         }
00351         break;
00352 
00353       default:
00354         env->localBugAssertion(__FILE__, __LINE__, "invalid norm", false,
00355           BASIC_ASSERTION);
00356         break;
00357     } 
00358   } 
00359 }
00360 
00402 template <typename scalar_t, typename pnum_t, typename lno_t, typename part_t>
00403   void globalSumsByPart( 
00404     const RCP<const Environment> &env,
00405     const RCP<const Comm<int> > &comm, 
00406     const ArrayView<const pnum_t> &part, 
00407     int vwgtDim,
00408     const ArrayView<StridedData<lno_t, scalar_t> > &vwgts,
00409     multiCriteriaNorm mcNorm,
00410     part_t &numParts, 
00411     part_t &numNonemptyParts,
00412     ArrayRCP<MetricValues<scalar_t> > &metrics,
00413     ArrayRCP<scalar_t> &globalSums)
00414 {
00415   env->debug(DETAILED_STATUS, "Entering globalSumsByPart");
00417   // Initialize return values
00418 
00419   numParts = numNonemptyParts = 0;
00420 
00421   int numMetrics = 1;                       // "object count"
00422   if (vwgtDim) numMetrics++;                // "normed weight" or "weight 1"
00423   if (vwgtDim > 1) numMetrics += vwgtDim;   // "weight n"
00424 
00425   typedef MetricValues<scalar_t> mv_t;
00426   mv_t *newMetrics = new mv_t [numMetrics];
00427   env->localMemoryAssertion(__FILE__, __LINE__, numMetrics, newMetrics); 
00428   ArrayRCP<mv_t> metricArray(newMetrics, 0, numMetrics, true);
00429 
00430   metrics = metricArray;
00431 
00433   // Figure out the global number of parts in use.
00434   // Verify number of vertex weights is the same everywhere.
00435 
00436   lno_t localNumObj = part.size();
00437   part_t localNum[2], globalNum[2];
00438   localNum[0] = static_cast<part_t>(vwgtDim);  
00439   localNum[1] = 0;
00440 
00441   for (lno_t i=0; i < localNumObj; i++)
00442     if (part[i] > localNum[1]) localNum[1] = part[i];
00443 
00444   try{
00445     reduceAll<int, part_t>(*comm, Teuchos::REDUCE_MAX, 2, 
00446       localNum, globalNum);
00447   }
00448   Z2_THROW_OUTSIDE_ERROR(*env)
00449 
00450   env->globalBugAssertion(__FILE__,__LINE__,
00451     "inconsistent number of vertex weights",
00452     globalNum[0] == localNum[0], DEBUG_MODE_ASSERTION, comm);
00453 
00454   part_t nparts = globalNum[1] + 1;
00455 
00456   part_t globalSumSize = nparts * numMetrics;
00457   scalar_t * sumBuf = new scalar_t [globalSumSize];
00458   env->localMemoryAssertion(__FILE__, __LINE__, globalSumSize, sumBuf);
00459   globalSums = arcp(sumBuf, 0, globalSumSize);
00460 
00462   // Calculate the local totals by part.
00463 
00464   scalar_t *localBuf = new scalar_t [globalSumSize];
00465   env->localMemoryAssertion(__FILE__, __LINE__, globalSumSize, localBuf);
00466   memset(localBuf, 0, sizeof(scalar_t) * globalSumSize);
00467 
00468   scalar_t *obj = localBuf;              // # of objects
00469 
00470   for (lno_t i=0; i < localNumObj; i++)
00471     obj[part[i]]++;
00472 
00473   if (numMetrics > 1){
00474 
00475     scalar_t *wgt = localBuf + nparts; // single normed weight
00476     try{
00477       normedPartWeights<scalar_t, pnum_t, lno_t, part_t>(env, nparts, 
00478         part, vwgts, mcNorm, wgt);
00479     }
00480     Z2_FORWARD_EXCEPTIONS
00481   
00482 //KDDKDD TODO  This code assumes the solution has the part ordered the
00483 //KDDKDD TODO  same way as the user input.  That assumption is not
00484 //KDDKDD TODO  currently true, although we plan to make it true.
00485 //KDDKDD TODO  As a results, currently the weight metrics may be wrong.
00486 //KDDKDD TODO  See bug 5891.  April 5, 2013
00487     if (vwgtDim > 1){
00488       wgt += nparts;         // individual weights
00489       for (int vdim = 0; vdim < vwgtDim; vdim++){
00490         for (lno_t i=0; i < localNumObj; i++)
00491           wgt[part[i]] += vwgts[vdim][i];
00492         wgt += nparts;
00493       }
00494     }
00495   }
00496 
00497   // Metric: local sums on process
00498 
00499   int next = 0;
00500 
00501   metrics[next].setName("object count");
00502   metrics[next].setLocalSum(localNumObj);
00503   next++;
00504 
00505   if (numMetrics > 1){
00506     scalar_t *wgt = localBuf + nparts; // single normed weight
00507     scalar_t total = 0.0;
00508   
00509     for (int p=0; p < nparts; p++){
00510       total += wgt[p];
00511     }
00512 
00513     if (vwgtDim == 1)
00514       metrics[next].setName("weight 1");
00515     else
00516       metrics[next].setName("normed weight");
00517 
00518     metrics[next].setLocalSum(total);
00519     next++;
00520   
00521     if (vwgtDim > 1){
00522       for (int vdim = 0; vdim < vwgtDim; vdim++){
00523         wgt += nparts;
00524         total = 0.0;
00525         for (int p=0; p < nparts; p++){
00526           total += wgt[p];
00527         }
00528 
00529         std::ostringstream oss;
00530         oss << "weight " << vdim+1;
00531 
00532         metrics[next].setName(oss.str());
00533         metrics[next].setLocalSum(total);
00534         next++;
00535       }
00536     }
00537   }
00538 
00540   // Obtain global totals by part.
00541 
00542   try{
00543     reduceAll<int, scalar_t>(*comm, Teuchos::REDUCE_SUM, globalSumSize,
00544       localBuf, sumBuf);
00545   }
00546   Z2_THROW_OUTSIDE_ERROR(*env);
00547 
00548   delete [] localBuf;
00549 
00551   // Global sum, min, max, and average over all parts
00552 
00553   obj = sumBuf;                     // # of objects
00554   scalar_t min=0, max=0, sum=0;
00555   next = 0;
00556 
00557   ArrayView<scalar_t> objVec(obj, nparts);
00558   getStridedStats<scalar_t>(objVec, 1, 0, min, max, sum);
00559 
00560   metrics[next].setGlobalMin(min);
00561   metrics[next].setGlobalMax(max);
00562   metrics[next].setGlobalSum(sum);
00563   metrics[next].setGlobalAvg(sum / nparts);
00564   next++;
00565 
00566   if (numMetrics > 1){
00567     scalar_t *wgt = sumBuf + nparts;        // single normed weight
00568   
00569     ArrayView<scalar_t> normedWVec(wgt, nparts);
00570     getStridedStats<scalar_t>(normedWVec, 1, 0, min, max, sum);
00571 
00572     if (vwgtDim > 1)
00573       metrics[next].setNorm(multiCriteriaNorm(mcNorm));
00574 
00575     metrics[next].setGlobalMin(min);
00576     metrics[next].setGlobalMax(max);
00577     metrics[next].setGlobalSum(sum);
00578     metrics[next].setGlobalAvg(sum / nparts);
00579     next++;
00580   
00581     if (vwgtDim > 1){
00582       for (int vdim=0; vdim < vwgtDim; vdim++){
00583         wgt += nparts;       // individual weights
00584         ArrayView<scalar_t> fromVec(wgt, nparts);
00585         getStridedStats<scalar_t>(fromVec, 1, 0, min, max, sum);
00586 
00587         metrics[next].setGlobalMin(min);
00588         metrics[next].setGlobalMax(max);
00589         metrics[next].setGlobalSum(sum);
00590         metrics[next].setGlobalAvg(sum / nparts);
00591         next++;
00592       }
00593     }
00594   }
00595 
00597   // How many parts do we actually have.
00598 
00599   numParts = nparts;
00600   obj = sumBuf;               // # of objects
00601 
00602   for (part_t p=nparts-1; p > 0; p--){
00603     if (obj[p] > 0) break;
00604     numParts--;
00605   }
00606 
00607   numNonemptyParts = numParts; 
00608 
00609   for (part_t p=0; p < numParts; p++)
00610     if (obj[p] == 0) numNonemptyParts--;
00611 
00612   env->debug(DETAILED_STATUS, "Exiting globalSumsByPart");
00613 }
00614 
00645 template <typename scalar_t, typename part_t>
00646   void computeImbalances(part_t numParts, part_t targetNumParts,
00647     const scalar_t *psizes, scalar_t sumVals , const scalar_t *vals, 
00648     scalar_t &min, scalar_t &max, scalar_t &avg)
00649 {
00650   min = sumVals;
00651   max = avg = 0;
00652 
00653   if (sumVals <= 0 || targetNumParts < 1 || numParts < 1)
00654     return;
00655 
00656   if (targetNumParts==1 || numParts==1){
00657     min = max = avg = 0;  // 0 imbalance
00658     return;
00659   }
00660 
00661   if (!psizes){
00662     scalar_t target = sumVals / targetNumParts;
00663     for (part_t p=0; p < numParts; p++){
00664       scalar_t diff = abs(vals[p] - target);
00665       scalar_t tmp = diff / target;
00666       avg += tmp;
00667       if (tmp > max) max = tmp;
00668       if (tmp < min) min = tmp;
00669     }
00670     part_t emptyParts = targetNumParts - numParts;  
00671     if (emptyParts > 0){
00672       if (max < 1.0)
00673         max = 1.0;       // target divided by target
00674       avg += emptyParts;
00675     }
00676   }
00677   else{
00678     for (part_t p=0; p < targetNumParts; p++){
00679       if (psizes[p] > 0){
00680         if (p < numParts){
00681           scalar_t target = sumVals * psizes[p];
00682           scalar_t diff = abs(vals[p] - target);
00683           scalar_t tmp = diff / target;
00684           avg += tmp;
00685           if (tmp > max) max = tmp;
00686           if (tmp < min) min = tmp;
00687         }
00688         else{
00689           if (max < 1.0)
00690             max = 1.0;  // target divided by target
00691           avg += 1.0;
00692         }
00693       }
00694     }
00695   }
00696 
00697   avg /= targetNumParts;
00698 }
00699 
00733 template <typename scalar_t, typename part_t>
00734  void computeImbalances(part_t numParts, part_t targetNumParts,
00735    int numSizes, ArrayView<ArrayRCP<scalar_t> > psizes,
00736    scalar_t sumVals , const scalar_t *vals, 
00737    scalar_t &min, scalar_t &max, scalar_t &avg)
00738 {
00739   min = sumVals;
00740   max = avg = 0;
00741 
00742   if (sumVals <= 0 || targetNumParts < 1 || numParts < 1)
00743     return;
00744 
00745   if (targetNumParts==1 && numParts==1){
00746     min = max = avg = 0;  // 0 imbalance
00747     return;
00748   }
00749 
00750   bool allUniformParts = true;
00751   for (int i=0; i < numSizes; i++){
00752     if (psizes[i].size() != 0){
00753       allUniformParts = false;
00754       break;
00755     }
00756   }
00757 
00758   if (allUniformParts){
00759     computeImbalances<scalar_t, part_t>(numParts, targetNumParts, NULL,
00760       sumVals, vals, min, max, avg);
00761     return;
00762   }
00763 
00764   double uniformSize = 1.0 / targetNumParts;
00765   ArrayRCP<double> sizeVec(new double [numSizes], 0, numSizes, true);
00766   for (int i=0; i < numSizes; i++){
00767     sizeVec[i] = uniformSize;
00768   }
00769 
00770   for (part_t p=0; p < numParts; p++){
00771 
00772     // If we have objects in parts that should have 0 objects,
00773     // we don't compute an imbalance.  It means that other
00774     // parts have too few objects, and the imbalance will be
00775     // picked up there.
00776 
00777     if (p >= targetNumParts)
00778       break;
00779 
00780     // Vector of target amounts: T
00781 
00782     for (int i=0; i < numSizes; i++)
00783       if (psizes[i].size() > 0)
00784         sizeVec[i] = psizes[i][p];
00785 
00786     Epetra_SerialDenseVector target(::View, sizeVec.getRawPtr(), numSizes);
00787     target.Scale(sumVals);
00788     double targetNorm = target.Norm2();
00789 
00790     // If part is supposed to be empty, we don't compute an
00791     // imbalance.  Same argument as above.
00792 
00793     if (targetNorm > 0){
00794 
00795       // Vector of actual amounts: A
00796 
00797       Epetra_SerialDenseVector actual(numSizes);
00798       for (int i=0; i < numSizes; i++)
00799         actual[i] = vals[p] * -1.0;
00800       
00801       actual += target;
00802 
00803       //  |A - T| / |T|
00804 
00805       scalar_t imbalance = actual.Norm2() / targetNorm;
00806 
00807       if (imbalance < min)
00808         min = imbalance;
00809       else if (imbalance > max)
00810         max = imbalance;
00811       avg += imbalance; 
00812     }
00813   }
00814 
00815   part_t numEmptyParts = 0;
00816 
00817   for (part_t p=numParts; p < targetNumParts; p++){
00818     bool nonEmptyPart = false;
00819     for (int i=0; !nonEmptyPart && i < numSizes; i++)
00820       if (psizes[i].size() > 0 && psizes[i][p] > 0.0)
00821         nonEmptyPart = true;
00822 
00823     if (nonEmptyPart){
00824       // The partitioning has no objects for this part, which
00825       // is supposed to be non-empty.
00826       numEmptyParts++;
00827     }
00828   }
00829 
00830   if (numEmptyParts > 0){
00831     avg += numEmptyParts;
00832     if (max < 1.0)
00833       max = 1.0;       // target divided by target
00834   }
00835 
00836   avg /= targetNumParts;
00837 }
00838 
00868 template <typename Adapter>
00869   void objectMetrics(
00870     const RCP<const Environment> &env,
00871     const RCP<const Comm<int> > &comm,
00872     multiCriteriaNorm mcNorm,
00873     const RCP<const Adapter> &ia,
00874     const RCP<const PartitioningSolution<Adapter> > &solution,
00875     typename Adapter::part_t &numParts,
00876     typename Adapter::part_t &numNonemptyParts,
00877     ArrayRCP<MetricValues<typename Adapter::scalar_t> > &metrics)
00878 {
00879   env->debug(DETAILED_STATUS, "Entering objectMetrics");
00880 
00881   typedef typename Adapter::scalar_t scalar_t;
00882   typedef typename Adapter::lno_t lno_t;
00883   typedef typename Adapter::part_t part_t;
00884   typedef StridedData<lno_t, scalar_t> sdata_t;
00885 
00886   // Local number of objects.
00887 
00888   size_t numLocalObjects = ia->getLocalNumIDs();
00889 
00890   // Parts to which objects are assigned.
00891 
00892   const part_t *parts = solution->getPartList();
00893   env->localInputAssertion(__FILE__, __LINE__, "parts not set", 
00894     parts, BASIC_ASSERTION);
00895   ArrayView<const part_t> partArray(parts, numLocalObjects);
00896 
00897   // Weights, if any, for each object.
00898 
00899   int nWeights = ia->getNumWeightsPerID();
00900   int numCriteria = (nWeights > 0 ? nWeights : 1);
00901   Array<sdata_t> weights(numCriteria);
00902 
00903   if (nWeights == 0){
00904     // One set of uniform weights is implied.
00905     // StridedData default constructor creates length 0 strided array.
00906     weights[0] = sdata_t();
00907   }
00908   else{
00909     for (int i=0; i < nWeights; i++){
00910       int stride;
00911       const scalar_t *wgt;
00912       ia->getWeightsView(wgt, stride, i); 
00913       ArrayRCP<const scalar_t> wgtArray(wgt, 0, stride*numLocalObjects, false);
00914       weights[i] = sdata_t(wgtArray, stride);
00915     }
00916   }
00917 
00918   // Relative part sizes, if any, assigned to the parts.
00919 
00920   part_t targetNumParts = solution->getTargetGlobalNumberOfParts();
00921   scalar_t *psizes = NULL;
00922 
00923   ArrayRCP<ArrayRCP<scalar_t> > partSizes(numCriteria);
00924   for (int dim=0; dim < numCriteria; dim++){
00925     if (solution->criteriaHasUniformPartSizes(dim) != true){
00926       psizes = new scalar_t [targetNumParts];
00927       env->localMemoryAssertion(__FILE__, __LINE__, numParts, psizes);
00928       for (part_t i=0; i < targetNumParts; i++){
00929         psizes[i] = solution->getCriteriaPartSize(dim, i);
00930       }
00931       partSizes[dim] = arcp(psizes, 0, targetNumParts, true);
00932     }
00933   }
00934 
00936   // Get number of parts, and the number that are non-empty.
00937   // Get sums per part of objects, individual weights, and normed weight sums.
00938 
00939   ArrayRCP<scalar_t> globalSums;
00940 
00941   try{
00942     globalSumsByPart<scalar_t, part_t, lno_t, part_t>(env, comm, 
00943       partArray, nWeights, weights.view(0, numCriteria), mcNorm,
00944       numParts, numNonemptyParts, metrics, globalSums);
00945   }
00946   Z2_FORWARD_EXCEPTIONS
00947 
00949   // Compute imbalances for the object count.  
00950   // (Use first index of part sizes.) 
00951 
00952   scalar_t *objCount  = globalSums.getRawPtr();
00953   scalar_t min, max, avg;
00954   psizes=NULL;
00955 
00956   if (partSizes[0].size() > 0)
00957     psizes = partSizes[0].getRawPtr();
00958 
00959   computeImbalances<scalar_t, part_t>(numParts, targetNumParts, psizes,
00960       metrics[0].getGlobalSum(), objCount, 
00961       min, max, avg);
00962 
00963   metrics[0].setMinImbalance(1.0 + min);
00964   metrics[0].setMaxImbalance(1.0 + max);
00965   metrics[0].setAvgImbalance(1.0 + avg);
00966 
00968   // Compute imbalances for the normed weight sum.
00969 
00970   scalar_t *wgts = globalSums.getRawPtr() + numParts;
00971 
00972   if (metrics.size() > 1){
00973   
00974     computeImbalances<scalar_t, part_t>(numParts, targetNumParts, 
00975       numCriteria, partSizes.view(0, numCriteria),
00976       metrics[1].getGlobalSum(), wgts,
00977       min, max, avg);
00978   
00979     metrics[1].setMinImbalance(1.0 + min);
00980     metrics[1].setMaxImbalance(1.0 + max);
00981     metrics[1].setAvgImbalance(1.0 + avg);
00982 
00983     if (metrics.size() > 2){
00984 
00986     // Compute imbalances for each individual weight.
00987 
00988       int next = 2;
00989   
00990       for (int vdim=0; vdim < numCriteria; vdim++){
00991         wgts += numParts;
00992         psizes = NULL;
00993   
00994         if (partSizes[vdim].size() > 0)
00995            psizes = partSizes[vdim].getRawPtr();
00996          
00997         computeImbalances<scalar_t, part_t>(numParts, targetNumParts, psizes,
00998           metrics[next].getGlobalSum(), wgts, min, max, avg);
00999   
01000         metrics[next].setMinImbalance(1.0 + min);
01001         metrics[next].setMaxImbalance(1.0 + max);
01002         metrics[next].setAvgImbalance(1.0 + avg);
01003         next++;
01004       }
01005     }
01006 
01007   }
01008   env->debug(DETAILED_STATUS, "Exiting objectMetrics");
01009 }
01010 
01014 template <typename scalar_t, typename part_t>
01015   void printMetrics( std::ostream &os,
01016     part_t targetNumParts, part_t numParts, part_t numNonemptyParts, 
01017     const ArrayView<MetricValues<scalar_t> > &infoList)
01018 {
01019   os << "NUMBER OF PARTS IS " << numParts;
01020   if (numNonemptyParts < numParts)
01021     os << " (" << numNonemptyParts << " of which are non-empty)";
01022   os << std::endl;
01023   if (targetNumParts != numParts)
01024     os << "TARGET NUMBER OF PARTS IS " << targetNumParts << std::endl;
01025 
01026   std::string unset("unset");
01027 
01028   MetricValues<scalar_t>::printHeader(os);
01029 
01030   for (int i=0; i < infoList.size(); i++)
01031     if (infoList[i].getName() != unset)
01032       infoList[i].printLine(os);
01033 
01034   os << std::endl;
01035 }
01036 
01040 template <typename scalar_t, typename part_t>
01041   void printMetrics( std::ostream &os,
01042     part_t targetNumParts, part_t numParts, part_t numNonemptyParts, 
01043     const MetricValues<scalar_t> &info)
01044 {
01045   ArrayView<MetricValues<scalar_t> > infoList(&info, 1);
01046   printMetrics( os, targetNumParts, numParts, numNonemptyParts, infoList);
01047 }
01048 
01049 
01052 template <typename scalar_t>
01053   scalar_t normedWeight(ArrayView <scalar_t> weights, 
01054     multiCriteriaNorm norm)
01055 {
01056   size_t dim = weights.size();
01057   if (dim == 0)
01058     return 0.0;
01059   else if (dim == 1)
01060     return weights[0];
01061  
01062   scalar_t nweight = 0;
01063 
01064   switch (norm){
01065     case normMinimizeTotalWeight:   
01067       for (size_t i=0; i <dim; i++)
01068         nweight += weights[i];
01069 
01070       break;
01071      
01072     case normBalanceTotalMaximum:   
01073       for (size_t i=0; i <dim; i++)
01074         nweight += (weights[i] * weights[i]);
01075 
01076       nweight = sqrt(nweight);
01077 
01078       break;
01079 
01080     case normMinimizeMaximumWeight: 
01082       nweight = weights[0];
01083       for (size_t i=1; i <dim; i++)
01084         if (weights[i] > nweight)
01085           nweight = weights[i];
01086 
01087       break;
01088 
01089     default:
01090       Environment env;  // default environment
01091       env.localBugAssertion(__FILE__, __LINE__, "invalid norm", false,
01092         BASIC_ASSERTION);
01093       break;
01094   } 
01095 
01096   return nweight;
01097 }
01098 
01101 template<typename lno_t, typename scalar_t>
01102   scalar_t normedWeight(ArrayView<StridedData<lno_t, scalar_t> > weights, 
01103      lno_t idx, multiCriteriaNorm norm)
01104 {
01105   size_t dim = weights.size();
01106   if (dim < 1)
01107     return 0;
01108 
01109   Array<scalar_t> vec(dim, 1.0);
01110   for (size_t i=0; i < dim; i++)
01111     if (weights[i].size() > 0)
01112        vec[dim] = weights[i][idx];
01113 
01114   return normedWeight(vec, norm);
01115 }
01116 
01117 } //namespace Zoltan2
01118 #endif