From 7b850b8ee43fb7f8a0b2a1565ed01102d40b0a14 Mon Sep 17 00:00:00 2001 From: Arnur Nigmetov Date: Mon, 5 Sep 2016 13:32:05 +0200 Subject: Technical changes for R integration Avoid including iostream (R complains about that). All output protected by preprocessor directive (R checker should not see an instance of std::cout << in your code). Also added getWassersteinCost to be in line with the Dionysus implementation used in TDA. --- geom_matching/wasserstein/src/auction_oracle.cpp | 57 ++++--------------- .../wasserstein/src/auction_runner_gs.cpp | 61 +++++++++++++++++--- .../wasserstein/src/auction_runner_jac.cpp | 52 ++++++++++++----- geom_matching/wasserstein/src/basic_defs.cpp | 10 +++- geom_matching/wasserstein/src/wasserstein.cpp | 65 ++++++++++++++++++++-- 5 files changed, 170 insertions(+), 75 deletions(-) (limited to 'geom_matching/wasserstein/src') diff --git a/geom_matching/wasserstein/src/auction_oracle.cpp b/geom_matching/wasserstein/src/auction_oracle.cpp index 088936b..9c764da 100644 --- a/geom_matching/wasserstein/src/auction_oracle.cpp +++ b/geom_matching/wasserstein/src/auction_oracle.cpp @@ -68,8 +68,6 @@ AuctionOracleLazyHeap::AuctionOracleLazyHeap(const std::vector& b, weightMatrix.reserve(b.size()); //const double maxDistUpperBound = 3 * getFurthestDistance3Approx(b, g); //weightAdjConst = pow(maxDistUpperBound, wassersteinPower); - //std::cout << "3getFurthestDistance3Approx = " << getFurthestDistance3Approx(b, g) << std::endl; - //std::cout << "in AuctionOracleLazyHeap weightAdjConst = " << weightAdjConst << std::endl; // init weight matrix for(const auto& pointA : bidders) { std::vector weightVec; @@ -134,7 +132,9 @@ void AuctionOracleLazyHeap::setPrice(IdxType itemIdx, double newPrice) { assert( prices.at(itemIdx) < newPrice ); #ifdef DEBUG_AUCTION +#ifndef FOR_R_TDA std::cout << "price incremented by " << prices.at(itemIdx) - newPrice << std::endl; +#endif #endif prices[itemIdx] = newPrice; // lazy: record the moment we updated the price of the items, @@ -169,8 +169,10 @@ DebugOptimalBid AuctionOracleLazyHeap::getOptimalBidDebug(IdxType bidderIdx) result.secondBestItemValue = topIter->second; result.secondBestItemIdx = topIter->first; +#ifndef FOR_R_TDA //std::cout << "getOptimalBid: bidderIdx = " << bidderIdx << "; bestItemIdx = " << bestItemIdx << "; bestItemValue = " << bestItemValue << "; bestItemsPrice = " << prices[bestItemIdx] << "; secondBestItemIdx = " << topIter->first << "; secondBestValue = " << secondBestItemValue << "; secondBestPrice = " << prices[topIter->first] << "; bid = " << prices[bestItemIdx] + ( bestItemValue - secondBestItemValue ) + epsilon << "; epsilon = " << epsilon << std::endl; //std::cout << "getOptimalBid: bidderIdx = " << bidderIdx << "; bestItemIdx = " << bestItemIdx << "; bestItemsDist= " << (weightAdjConst - bestItemValue) << "; bestItemsPrice = " << prices[bestItemIdx] << "; secondBestItemIdx = " << topIter->first << "; secondBestDist= " << (weightAdjConst - secondBestItemValue) << "; secondBestPrice = " << prices[topIter->first] << "; bid = " << prices[bestItemIdx] + ( bestItemValue - secondBestItemValue ) + epsilon << "; epsilon = " << epsilon << std::endl; +#endif return result; } @@ -190,8 +192,10 @@ IdxValPair AuctionOracleLazyHeap::getOptimalBid(const IdxType bidderIdx) ++topIter; // now points to the second-best items double secondBestItemValue = topIter->second; +#ifndef FOR_R_TDA //std::cout << "getOptimalBid: bidderIdx = " << bidderIdx << "; bestItemIdx = " << bestItemIdx << "; bestItemValue = " << bestItemValue << "; bestItemsPrice = " << prices[bestItemIdx] << "; secondBestItemIdx = " << topIter->first << "; secondBestValue = " << secondBestItemValue << "; secondBestPrice = " << prices[topIter->first] << "; bid = " << prices[bestItemIdx] + ( bestItemValue - secondBestItemValue ) + epsilon << "; epsilon = " << epsilon << std::endl; //std::cout << "getOptimalBid: bidderIdx = " << bidderIdx << "; bestItemIdx = " << bestItemIdx << "; bestItemsDist= " << (weightAdjConst - bestItemValue) << "; bestItemsPrice = " << prices[bestItemIdx] << "; secondBestItemIdx = " << topIter->first << "; secondBestDist= " << (weightAdjConst - secondBestItemValue) << "; secondBestPrice = " << prices[topIter->first] << "; bid = " << prices[bestItemIdx] + ( bestItemValue - secondBestItemValue ) + epsilon << "; epsilon = " << epsilon << std::endl; +#endif // bid value: price + value difference + epsilon return std::make_pair(bestItemIdx, @@ -221,8 +225,6 @@ AuctionOracleLazyHeapRestricted::AuctionOracleLazyHeapRestricted(const std::vect weightMatrix.reserve(b.size()); //const double maxDistUpperBound = 3 * getFurthestDistance3Approx(b, g); //weightAdjConst = pow(maxDistUpperBound, wassersteinPower); - //std::cout << "3getFurthestDistance3Approx = " << getFurthestDistance3Approx(b, g) << std::endl; - //std::cout << "in AuctionOracleLazyHeapRestricted weightAdjConst = " << weightAdjConst << std::endl; // init weight matrix for(const auto& pointA : bidders) { std::vector weightVec; @@ -319,7 +321,9 @@ void AuctionOracleLazyHeapRestricted::setPrice(IdxType itemIdx, double newPrice) { assert( prices.at(itemIdx) < newPrice ); #ifdef DEBUG_AUCTION +#ifndef FOR_R_TDA std::cout << "price incremented by " << prices.at(itemIdx) - newPrice << std::endl; +#endif #endif prices[itemIdx] = newPrice; if (items[itemIdx].isNormal() ) { @@ -426,7 +430,6 @@ DebugOptimalBid AuctionOracleLazyHeapRestricted::getOptimalBidDebug(IdxType bidd //debugNaiveResult.secondBestItemIdx = itemIdx; //} //} - ////std::cout << "got naive result" << std::endl; //if ( fabs( debugMyResult.bestItemValue - debugNaiveResult.bestItemValue ) > 1e-6 or //fabs( debugNaiveResult.secondBestItemValue - debugMyResult.secondBestItemValue) > 1e-6 ) { @@ -648,7 +651,6 @@ AuctionOracleKDTree::AuctionOracleKDTree(const std::vector& _bidde } DnnTraits traits; traits.internal_p = internal_p; - //std::cout << "kdtree: " << dnnPointHandles.size() << " points" << std::endl; kdtree = new dnn::KDTree(traits, dnnPointHandles, wassersteinPower); size_t dnnItemIdxAll { 0 }; @@ -658,7 +660,6 @@ AuctionOracleKDTree::AuctionOracleKDTree(const std::vector& _bidde DnnPoint p(dnnItemIdxAll++); p[0] = g.getRealX(); p[1] = g.getRealY(); - //std::cout << "to all tree: " << p[0] << ", " << p[1] << std::endl; dnnPointsAll.push_back(p); assert(dnnItemIdxAll == dnnPointsAll.size()); } @@ -666,7 +667,6 @@ AuctionOracleKDTree::AuctionOracleKDTree(const std::vector& _bidde for(size_t i = 0; i < dnnPointsAll.size(); ++i) { dnnPointHandlesAll.push_back(&dnnPointsAll[i]); } - //std::cout << "kdtreeAll: " << dnnPointHandlesAll.size() << " points" << std::endl; kdtreeAll = new dnn::KDTree(traits, dnnPointHandlesAll, wassersteinPower); size_t handleIdx {0}; @@ -677,13 +677,9 @@ AuctionOracleKDTree::AuctionOracleKDTree(const std::vector& _bidde } } //to-do: remove maxVal from - //std::cout << "3getFurthestDistance3Approx = " << getFurthestDistance3Approx(_bidders, _items) << std::endl; maxVal = 3*getFurthestDistance3Approx(_bidders, _items); maxVal = pow(maxVal, wassersteinPower); weightAdjConst = maxVal; - //std::cout << "AuctionOracleKDTree: weightAdjConst = " << weightAdjConst << std::endl; - //std::cout << "AuctionOracleKDTree constructor done" << std::endl; - // for debug } DebugOptimalBid AuctionOracleKDTree::getOptimalBidDebug(IdxType bidderIdx) @@ -694,16 +690,12 @@ DebugOptimalBid AuctionOracleKDTree::getOptimalBidDebug(IdxType bidderIdx) bidderDnn[0] = bidder.getRealX(); bidderDnn[1] = bidder.getRealY(); - //std::cout << "bidder.x = " << bidderDnn[0] << std::endl; - //std::cout << "bidder.y = " << bidderDnn[1] << std::endl; - std::vector candItems; if ( bidder.isDiagonal() ) { // auto twoBestItems = kdtree->findK(bidderDnn, 2); - //std::cout << "twoBestItems for non-diag: " << twoBestItems[0].d << " " << twoBestItems[1].d << std::endl; candItems.push_back( std::make_pair(twoBestItems[0].p->id(), twoBestItems[0].d) ); candItems.push_back( std::make_pair(twoBestItems[1].p->id(), twoBestItems[1].d) ); assert(diagItemsHeap.size() > 1); @@ -719,7 +711,6 @@ DebugOptimalBid AuctionOracleKDTree::getOptimalBidDebug(IdxType bidderIdx) assert(candItems[1].second >= candItems[0].second); } else { auto twoBestItems = kdtreeAll->findK(bidderDnn, 2); - //std::cout << "twoBestItems for all: " << twoBestItems[0].d << " " << twoBestItems[1].d << std::endl; candItems.push_back( std::make_pair(twoBestItems[0].p->id(), twoBestItems[0].d) ); candItems.push_back( std::make_pair(twoBestItems[1].p->id(), twoBestItems[1].d) ); //size_t projItemIdx { biddersProjIndices.at(bidderIdx) }; @@ -734,7 +725,6 @@ DebugOptimalBid AuctionOracleKDTree::getOptimalBidDebug(IdxType bidderIdx) result.secondBestItemIdx = candItems[1].first; result.bestItemValue = candItems[0].second; result.secondBestItemValue = candItems[1].second; - //std::cout << "got result: " << result << std::endl; //double bestItemsPrice = prices[bestItemIdx]; //if (items[result.bestItemIdx].type == DiagramPoint::DIAG) { //double bestItemValue1 = pow( distLp(bidder, items[result.bestItemIdx]), q) + prices[result.bestItemIdx]; @@ -820,18 +810,12 @@ void AuctionOracleKDTree::setPrice(IdxType itemIdx, double newPrice) assert(newPrice > prices.at(itemIdx)); prices[itemIdx] = newPrice; if ( items[itemIdx].isNormal() ) { - //std::cout << "before increasing weight in kdtree " << std::endl; - //std::cout << kdtreeItems.at(itemIdx) << std::endl; assert(0 <= itemIdx and itemIdx < kdtreeItems.size()); assert(0 <= kdtreeItems[itemIdx] and kdtreeItems[itemIdx] < dnnPointHandles.size()); kdtree->increase_weight( dnnPointHandles[kdtreeItems[itemIdx]], newPrice); kdtreeAll->increase_weight( dnnPointHandlesAll[itemIdx], newPrice); - //std::cout << "after increasing weight in kdtree" << std::endl; } else { - //std::cout << "before decreasing weight in heap" << std::endl; - //std::cout << "diagHeapHandles.size = " << diagHeapHandles.size() << std::endl; kdtreeAll->increase_weight( dnnPointHandlesAll[itemIdx], newPrice); - //std::cout << "after increasing weight in kdtree" << std::endl; assert(diagHeapHandles.size() > heapHandlesIndices.at(itemIdx)); diagItemsHeap.decrease(diagHeapHandles[heapHandlesIndices[itemIdx]], std::make_pair(itemIdx, newPrice)); } @@ -986,7 +970,6 @@ AuctionOracleKDTreeRestricted::AuctionOracleKDTreeRestricted(const std::vector(traits, dnnPointHandles, wassersteinPower); size_t handleIdx {0}; @@ -997,12 +980,9 @@ AuctionOracleKDTreeRestricted::AuctionOracleKDTreeRestricted(const std::vectorid() }; double bestNormalItemValue { twoBestItems[0].d }; size_t secondBestNormalItemIdx { twoBestItems[1].p->id() }; @@ -1162,9 +1138,6 @@ IdxValPair AuctionOracleKDTreeRestricted::getOptimalBid(IdxType bidderIdx) DiagramPoint bidder = bidders[bidderIdx]; - //std::cout << "bidder.x = " << bidderDnn[0] << std::endl; - //std::cout << "bidder.y = " << bidderDnn[1] << std::endl; - // corresponding point is always considered as a candidate // if bidder is a diagonal point, projItem is a normal point, // and vice versa. @@ -1226,7 +1199,6 @@ IdxValPair AuctionOracleKDTreeRestricted::getOptimalBid(IdxType bidderIdx) bidderDnn[0] = bidder.getRealX(); bidderDnn[1] = bidder.getRealY(); auto twoBestItems = kdtree->findK(bidderDnn, 2); - //std::cout << "twoBestItems for all: " << twoBestItems[0].d << " " << twoBestItems[1].d << std::endl; size_t bestNormalItemIdx { twoBestItems[0].p->id() }; double bestNormalItemValue { twoBestItems[0].d }; // if there is only one off-diagonal point in the second diagram, @@ -1268,15 +1240,10 @@ void AuctionOracleKDTreeRestricted::setPrice(IdxType itemIdx, double newPrice) assert(newPrice > prices.at(itemIdx)); prices[itemIdx] = newPrice; if ( items[itemIdx].isNormal() ) { - //std::cout << "before increasing weight in kdtree " << std::endl; - //std::cout << kdtreeItems.at(itemIdx) << std::endl; assert(0 <= itemIdx and itemIdx < kdtreeItems.size()); assert(0 <= kdtreeItems[itemIdx] and kdtreeItems[itemIdx] < dnnPointHandles.size()); kdtree->increase_weight( dnnPointHandles[kdtreeItems[itemIdx]], newPrice); - //std::cout << "after increasing weight in kdtree" << std::endl; } else { - //std::cout << "before decreasing weight in heap" << std::endl; - //std::cout << "diagHeapHandles.size = " << diagHeapHandles.size() << std::endl; assert(diagHeapHandles.size() > heapHandlesIndices.at(itemIdx)); diagItemsHeap.decrease(diagHeapHandles[heapHandlesIndices[itemIdx]], std::make_pair(itemIdx, newPrice)); bestDiagonalItemsComputed = false; @@ -1300,10 +1267,10 @@ void AuctionOracleKDTreeRestricted::setEpsilon(double newVal) std::ostream& operator<< (std::ostream& output, const DebugOptimalBid& db) { - std::cout << "bestItemValue = " << db.bestItemValue; - std::cout << "; bestItemIdx = " << db.bestItemIdx; - std::cout << "; secondBestItemValue = " << db.secondBestItemValue; - std::cout << "; secondBestItemIdx = " << db.secondBestItemIdx; + output << "bestItemValue = " << db.bestItemValue; + output << "; bestItemIdx = " << db.bestItemIdx; + output << "; secondBestItemValue = " << db.secondBestItemValue; + output << "; secondBestItemIdx = " << db.secondBestItemIdx; return output; } diff --git a/geom_matching/wasserstein/src/auction_runner_gs.cpp b/geom_matching/wasserstein/src/auction_runner_gs.cpp index 4ff495f..bd25442 100644 --- a/geom_matching/wasserstein/src/auction_runner_gs.cpp +++ b/geom_matching/wasserstein/src/auction_runner_gs.cpp @@ -27,6 +27,7 @@ derivative works thereof, in binary and source code form. #include +#include #include #include #include @@ -36,6 +37,10 @@ derivative works thereof, in binary and source code form. #include "auction_runner_gs.h" #include "wasserstein.h" +#ifdef FOR_R_TDA +#include "Rcpp.h" +#endif + //#define PRINT_DETAILED_TIMING namespace geom_ws { @@ -144,11 +149,13 @@ void AuctionRunnerGS::runAuction(void) double denominator = currentResult - numBidders * oracle->getEpsilon(); currentResult = pow(currentResult, 1.0 / wassersteinPower); #ifdef PRINT_DETAILED_TIMING +#ifndef FOR_R_TDA iterResults.push_back(currentResult); iterTimes.push_back(hrClock.now()); std::cout << "Iteration " << iterNum << " finished. "; std::cout << "Current result is " << currentResult << ", epsilon = " << oracle->getEpsilon() << std::endl; std::cout << "Number of rounds (cumulative): " << numRounds << std::endl; +#endif #endif if ( denominator <= 0 ) { //std::cout << "Epsilon is too big." << std::endl; @@ -157,8 +164,10 @@ void AuctionRunnerGS::runAuction(void) denominator = pow(denominator, 1.0 / wassersteinPower); double numerator = currentResult - denominator; #ifdef PRINT_DETAILED_TIMING +#ifndef FOR_R_TDA std::cout << " numerator: " << numerator << " denominator: " << denominator; std::cout << "; error bound: " << numerator / denominator << std::endl; +#endif #endif // if relative error is greater than delta, continue notDone = ( numerator / denominator > delta ); @@ -166,13 +175,16 @@ void AuctionRunnerGS::runAuction(void) // decrease epsilon for the next iteration oracle->setEpsilon( oracle->getEpsilon() / epsilonCommonRatio ); if (iterNum > maxIterNum) { +#ifndef FOR_R_TDA std::cerr << "Maximum iteration number exceeded, exiting. Current result is:"; std::cerr << wassersteinDistance << std::endl; - std::exit(1); +#endif + throw std::runtime_error("Maximum iteration number exceeded"); } } while ( notDone ); //printMatching(); #ifdef PRINT_DETAILED_TIMING +#ifndef FOR_R_TDA for(size_t iterIdx = 0; iterIdx < iterResults.size(); ++iterIdx) { double trueRelError = ( iterResults.at(iterIdx) - currentResult ) / currentResult; auto iterCumulativeTime = iterTimes.at(iterIdx) - startMoment; @@ -183,6 +195,7 @@ void AuctionRunnerGS::runAuction(void) ", iteration time " << iterTime.count() << std::endl; } #endif +#endif } void AuctionRunnerGS::runAuctionPhase(void) @@ -200,15 +213,22 @@ void AuctionRunnerGS::runAuctionPhase(void) assignItemToBidder(optimalBid.first, bidderIdx); oracle->setPrice(optimalItemIdx, bidValue); //printDebug(); +#ifdef FOR_R_TDA + if ( numRounds % 10000 == 0 ) { + Rcpp::checkUserInterrupt(); + } +#endif } while (not unassignedBidders.empty()); //std::cout << "runAuctionPhase finished" << std::endl; #ifdef DEBUG_AUCTION for(size_t bidderIdx = 0; bidderIdx < numBidders; ++bidderIdx) { if ( biddersToItems[bidderIdx] < 0) { +#ifndef FOR_R_TDA std::cerr << "After auction terminated bidder " << bidderIdx; std::cerr << " has no items assigned" << std::endl; - throw "Auction did not give a perfect matching"; +#endif + throw std::runtime_error("Auction did not give a perfect matching"); } } #endif @@ -225,6 +245,7 @@ double AuctionRunnerGS::getDistanceToQthPowerInternal(void) auto pB = items[biddersToItems[bIdx]]; result += pow(distLp(pA, pB, internal_p), wassersteinPower); } + wassersteinCost = result; wassersteinDistance = pow(result, 1.0 / wassersteinPower); return result; } @@ -235,11 +256,20 @@ double AuctionRunnerGS::getWassersteinDistance(void) return wassersteinDistance; } +double AuctionRunnerGS::getWassersteinCost(void) +{ + runAuction(); + return wassersteinCost; +} + + + // Debug routines void AuctionRunnerGS::printDebug(void) { #ifdef DEBUG_AUCTION +#ifndef FOR_R_TDA sanityCheck(); std::cout << "**********************" << std::endl; std::cout << "Current assignment:" << std::endl; @@ -259,6 +289,7 @@ void AuctionRunnerGS::printDebug(void) } std::cout << "**********************" << std::endl; #endif +#endif } @@ -266,13 +297,17 @@ void AuctionRunnerGS::sanityCheck(void) { #ifdef DEBUG_AUCTION if (biddersToItems.size() != numBidders) { +#ifndef FOR_R_TDA std::cerr << "Wrong size of biddersToItems, must be " << numBidders << ", is " << biddersToItems.size() << std::endl; - throw "Wrong size of biddersToItems"; +#endif + throw std::runtime_error("Wrong size of biddersToItems"); } if (itemsToBidders.size() != numBidders) { +#ifndef FOR_R_TDA std::cerr << "Wrong size of itemsToBidders, must be " << numBidders << ", is " << itemsToBidders.size() << std::endl; - throw "Wrong size of itemsToBidders"; +#endif + throw std::runtime_error("Wrong size of itemsToBidders"); } for(size_t bidderIdx = 0; bidderIdx < numBidders; ++bidderIdx) { @@ -281,18 +316,22 @@ void AuctionRunnerGS::sanityCheck(void) if ( std::count(biddersToItems.begin(), biddersToItems.end(), biddersToItems[bidderIdx]) > 1 ) { +#ifndef FOR_R_TDA std::cerr << "Item " << biddersToItems[bidderIdx]; std::cerr << " appears in biddersToItems more than once" << std::endl; - throw "Duplicate in biddersToItems"; +#endif + throw std::runtime_error("Duplicate in biddersToItems"); } if (itemsToBidders.at(biddersToItems[bidderIdx]) != static_cast(bidderIdx)) { +#ifndef FOR_R_TDA std::cerr << "Inconsitency: bidderIdx = " << bidderIdx; std::cerr << ", itemIdx in biddersToItems = "; std::cerr << biddersToItems[bidderIdx]; std::cerr << ", bidderIdx in itemsToBidders = "; std::cerr << itemsToBidders[biddersToItems[bidderIdx]] << std::endl; - throw "inconsistent mapping"; +#endif + throw std::runtime_error("inconsistent mapping"); } } } @@ -304,18 +343,22 @@ void AuctionRunnerGS::sanityCheck(void) if ( std::count(itemsToBidders.begin(), itemsToBidders.end(), itemsToBidders[itemIdx]) > 1 ) { +#ifndef FOR_R_TDA std::cerr << "Bidder " << itemsToBidders[itemIdx]; std::cerr << " appears in itemsToBidders more than once" << std::endl; - throw "Duplicate in itemsToBidders"; +#endif + throw std::runtime_error("Duplicate in itemsToBidders"); } // check for consistency if (biddersToItems.at(itemsToBidders[itemIdx]) != static_cast(itemIdx)) { +#ifndef FOR_R_TDA std::cerr << "Inconsitency: itemIdx = " << itemIdx; std::cerr << ", bidderIdx in itemsToBidders = "; std::cerr << itemsToBidders[itemIdx]; std::cerr << ", itemIdx in biddersToItems= "; std::cerr << biddersToItems[itemsToBidders[itemIdx]] << std::endl; - throw "inconsistent mapping"; +#endif + throw std::runtime_error("inconsistent mapping"); } } } @@ -325,6 +368,7 @@ void AuctionRunnerGS::sanityCheck(void) void AuctionRunnerGS::printMatching(void) { //#ifdef DEBUG_AUCTION +#ifndef FOR_R_TDA sanityCheck(); for(size_t bIdx = 0; bIdx < biddersToItems.size(); ++bIdx) { if (biddersToItems[bIdx] >= 0) { @@ -335,6 +379,7 @@ void AuctionRunnerGS::printMatching(void) assert(false); } } +#endif //#endif } diff --git a/geom_matching/wasserstein/src/auction_runner_jac.cpp b/geom_matching/wasserstein/src/auction_runner_jac.cpp index 08ed202..b892643 100644 --- a/geom_matching/wasserstein/src/auction_runner_jac.cpp +++ b/geom_matching/wasserstein/src/auction_runner_jac.cpp @@ -37,10 +37,10 @@ derivative works thereof, in binary and source code form. namespace geom_ws { // ***************************** -// AuctionRunnerJak +// AuctionRunnerJac // ***************************** -AuctionRunnerJak::AuctionRunnerJak(const std::vector& A, const std::vector& B, const double q, const double _delta, const double _internal_p) : +AuctionRunnerJac::AuctionRunnerJac(const std::vector& A, const std::vector& B, const double q, const double _delta, const double _internal_p) : bidders(A), items(B), numBidders(A.size()), @@ -57,7 +57,7 @@ AuctionRunnerJak::AuctionRunnerJak(const std::vector& A, const std oracle = std::unique_ptr(new AuctionOracle(bidders, items, wassersteinPower, internal_p)); } -void AuctionRunnerJak::assignGoodToBidder(IdxType itemIdx, IdxType bidderIdx) +void AuctionRunnerJac::assignGoodToBidder(IdxType itemIdx, IdxType bidderIdx) { //sanityCheck(); IdxType myOldItem = biddersToItems[bidderIdx]; @@ -79,7 +79,9 @@ void AuctionRunnerJak::assignGoodToBidder(IdxType itemIdx, IdxType bidderIdx) // RE: this cannot be necessary. I submitted the best bid, hence I was // an unassigned bidder. if (myOldItem != -1) { +#ifndef FOR_R_TDA std::cout << "This is not happening" << std::endl; +#endif assert(false); itemsToBidders[myOldItem] = -1; } @@ -88,7 +90,9 @@ void AuctionRunnerJak::assignGoodToBidder(IdxType itemIdx, IdxType bidderIdx) biddersToItems[currItemOwner] = myOldItem; // add the old owner of bids to the list of if ( -1 != myOldItem ) { +#ifndef FOR_R_TDA std::cout << "This is not happening" << std::endl; +#endif assert(false); // if I had something, update itemsToBidders, too // RE: nonsense: if I had something, I am not unassigned and did not @@ -103,18 +107,17 @@ void AuctionRunnerJak::assignGoodToBidder(IdxType itemIdx, IdxType bidderIdx) } -void AuctionRunnerJak::assignToBestBidder(IdxType itemIdx) +void AuctionRunnerJac::assignToBestBidder(IdxType itemIdx) { assert( itemIdx >= 0 and itemIdx < static_cast(numItems) ); assert( bidTable[itemIdx].first != -1); IdxValPair bestBid { bidTable[itemIdx] }; assignGoodToBidder(itemIdx, bestBid.first); - //std::cout << "About to call setPrice" << std::endl; oracle->setPrice(itemIdx, bestBid.second); //dynamic_cast(oracle)->setNai } -void AuctionRunnerJak::clearBidTable(void) +void AuctionRunnerJac::clearBidTable(void) { for(auto& itemWithBidIdx : itemsWithBids) { itemReceivedBidVec[itemWithBidIdx] = 0; @@ -124,7 +127,7 @@ void AuctionRunnerJak::clearBidTable(void) itemsWithBids.clear(); } -void AuctionRunnerJak::submitBid(IdxType bidderIdx, const IdxValPair& itemsBidValuePair) +void AuctionRunnerJac::submitBid(IdxType bidderIdx, const IdxValPair& itemsBidValuePair) { IdxType itemIdx = itemsBidValuePair.first; double bidValue = itemsBidValuePair.second; @@ -139,9 +142,10 @@ void AuctionRunnerJak::submitBid(IdxType bidderIdx, const IdxValPair& itemsBidVa } } -void AuctionRunnerJak::printDebug(void) +void AuctionRunnerJac::printDebug(void) { #ifdef DEBUG_AUCTION +#ifndef FOR_R_TDA sanityCheck(); std::cout << "**********************" << std::endl; std::cout << "Current assignment:" << std::endl; @@ -168,9 +172,10 @@ void AuctionRunnerJak::printDebug(void) //} std::cout << "**********************" << std::endl; #endif +#endif } -void AuctionRunnerJak::flushAssignment(void) +void AuctionRunnerJac::flushAssignment(void) { for(auto& b2g : biddersToItems) { b2g = -1; @@ -181,7 +186,7 @@ void AuctionRunnerJak::flushAssignment(void) //oracle->adjustPrices(); } -void AuctionRunnerJak::runAuction(void) +void AuctionRunnerJac::runAuction(void) { // relative error // choose some initial epsilon @@ -213,15 +218,19 @@ void AuctionRunnerJak::runAuction(void) // decrease epsilon for the next iteration oracle->setEpsilon( oracle->getEpsilon() / epsilonCommonRatio ); if (iterNum > maxIterNum) { +#ifndef FOR_R_TDA std::cerr << "Maximum iteration number exceeded, exiting. Current result is:"; std::cerr << wassersteinDistance << std::endl; std::exit(1); +#else + throw std::runtime_error("Max. iter. exceeded"); +#endif } } while ( notDone ); //printMatching(); } -void AuctionRunnerJak::runAuctionPhase(void) +void AuctionRunnerJac::runAuctionPhase(void) { //std::cout << "Entered runAuctionPhase" << std::endl; //int numUnassignedBidders { 0 }; @@ -266,7 +275,7 @@ void AuctionRunnerJak::runAuctionPhase(void) } // assertion: the matching must be perfect -double AuctionRunnerJak::getDistanceToQthPowerInternal(void) +double AuctionRunnerJac::getDistanceToQthPowerInternal(void) { sanityCheck(); double result = 0.0; @@ -276,19 +285,29 @@ double AuctionRunnerJak::getDistanceToQthPowerInternal(void) auto pB = items[biddersToItems[bIdx]]; result += pow(distLp(pA, pB, internal_p), wassersteinPower); } + wassersteinCost = result; wassersteinDistance = pow(result, 1.0 / wassersteinPower); return result; } -double AuctionRunnerJak::getWassersteinDistance(void) +double AuctionRunnerJac::getWassersteinDistance(void) { runAuction(); return wassersteinDistance; } -void AuctionRunnerJak::sanityCheck(void) +double AuctionRunnerJac::getWassersteinCost(void) +{ + runAuction(); + return wassersteinCost; +} + + + +void AuctionRunnerJac::sanityCheck(void) { #ifdef DEBUG_AUCTION +#ifndef FOR_R_TDA if (biddersToItems.size() != numBidders) { std::cerr << "Wrong size of biddersToItems, must be " << numBidders << ", is " << biddersToItems.size() << std::endl; throw "Wrong size of biddersToItems"; @@ -344,11 +363,13 @@ void AuctionRunnerJak::sanityCheck(void) } } #endif +#endif } -void AuctionRunnerJak::printMatching(void) +void AuctionRunnerJac::printMatching(void) { //#ifdef DEBUG_AUCTION +#ifndef FOR_R_TDA sanityCheck(); for(size_t bIdx = 0; bIdx < biddersToItems.size(); ++bIdx) { if (biddersToItems[bIdx] >= 0) { @@ -359,6 +380,7 @@ void AuctionRunnerJak::printMatching(void) assert(false); } } +#endif //#endif } diff --git a/geom_matching/wasserstein/src/basic_defs.cpp b/geom_matching/wasserstein/src/basic_defs.cpp index d228123..a46e6aa 100644 --- a/geom_matching/wasserstein/src/basic_defs.cpp +++ b/geom_matching/wasserstein/src/basic_defs.cpp @@ -30,6 +30,10 @@ derivative works thereof, in binary and source code form. #include #include "basic_defs_ws.h" +#ifndef FOR_R_TDA +#include +#endif + namespace geom_ws { // Point @@ -43,11 +47,14 @@ bool Point::operator!=(const Point& other) const return !(*this == other); } + +#ifndef FOR_R_TDA std::ostream& operator<<(std::ostream& output, const Point p) { output << "(" << p.x << ", " << p.y << ")"; return output; } +#endif double sqrDist(const Point& a, const Point& b) { @@ -97,6 +104,7 @@ double distLp(const DiagramPoint& a, const DiagramPoint& b, const double p) } +#ifndef FOR_R_TDA std::ostream& operator<<(std::ostream& output, const DiagramPoint p) { if ( p.type == DiagramPoint::DIAG ) { @@ -106,7 +114,7 @@ std::ostream& operator<<(std::ostream& output, const DiagramPoint p) } return output; } - +#endif DiagramPoint::DiagramPoint(double xx, double yy, Type ttype) : x(xx), diff --git a/geom_matching/wasserstein/src/wasserstein.cpp b/geom_matching/wasserstein/src/wasserstein.cpp index cac811b..fc1b662 100644 --- a/geom_matching/wasserstein/src/wasserstein.cpp +++ b/geom_matching/wasserstein/src/wasserstein.cpp @@ -36,7 +36,7 @@ derivative works thereof, in binary and source code form. #ifdef GAUSS_SEIDEL_AUCTION #include "auction_runner_gs.h" #else -#include "auction_runner_jak.h" +#include "auction_runner_jac.h" #endif namespace geom_ws { @@ -50,29 +50,78 @@ double wassersteinDistVec(const std::vector& A, const double _epsFactor) { if (q < 1) { +#ifndef FOR_R_TDA std::cerr << "Wasserstein distance not defined for q = " << q << ", must be >= 1" << std::endl; - throw "Bad q in Wasserstein"; +#endif + throw std::runtime_error("Bad q in Wasserstein " + std::to_string(q)); } if (delta < 0.0) { +#ifndef FOR_R_TDA std::cerr << "Relative error " << delta << ", must be > 0" << std::endl; - throw "Bad delta in Wasserstein"; +#endif + throw std::runtime_error("Bad delta in Wasserstein " + std::to_string(delta)); } if (_initialEpsilon < 0.0) { +#ifndef FOR_R_TDA std::cerr << "Initial epsilon = " << _initialEpsilon << ", must be non-negative" << std::endl; - throw "Bad delta in Wasserstein"; +#endif + throw std::runtime_error("Bad initial epsilon in Wasserstein" + std::to_string(_initialEpsilon)); } if (_epsFactor < 0.0) { +#ifndef FOR_R_TDA std::cerr << "Epsilon factor = " << _epsFactor << ", must be non-negative" << std::endl; - throw "Bad delta in Wasserstein"; +#endif + throw std::runtime_error("Bad epsilon factor in Wasserstein " + std::to_string(_epsFactor)); } + #ifdef GAUSS_SEIDEL_AUCTION AuctionRunnerGS auction(A, B, q, delta, _internal_p, _initialEpsilon, _epsFactor); #else - AuctionRunnerJak auction(A, B, q, delta, _internal_p); + AuctionRunnerJac auction(A, B, q, delta, _internal_p); #endif return auction.getWassersteinDistance(); } +double wassersteinCostVec(const std::vector& A, + const std::vector& B, + const double q, + const double delta, + const double _internal_p, + const double _initialEpsilon, + const double _epsFactor) +{ + if (q < 1) { +#ifndef FOR_R_TDA + std::cerr << "Wasserstein distance not defined for q = " << q << ", must be >= 1" << std::endl; +#endif + throw std::runtime_error("Bad q in Wasserstein " + std::to_string(q)); + } + if (delta < 0.0) { +#ifndef FOR_R_TDA + std::cerr << "Relative error " << delta << ", must be > 0" << std::endl; +#endif + throw std::runtime_error("Bad delta in Wasserstein " + std::to_string(delta)); + } + if (_initialEpsilon < 0.0) { +#ifndef FOR_R_TDA + std::cerr << "Initial epsilon = " << _initialEpsilon << ", must be non-negative" << std::endl; +#endif + throw std::runtime_error("Bad initial epsilon in Wasserstein" + std::to_string(_initialEpsilon)); + } + if (_epsFactor < 0.0) { +#ifndef FOR_R_TDA + std::cerr << "Epsilon factor = " << _epsFactor << ", must be non-negative" << std::endl; +#endif + throw std::runtime_error("Bad epsilon factor in Wasserstein " + std::to_string(_epsFactor)); + } +#ifdef GAUSS_SEIDEL_AUCTION + AuctionRunnerGS auction(A, B, q, delta, _internal_p, _initialEpsilon, _epsFactor); +#else + AuctionRunnerJac auction(A, B, q, delta, _internal_p); +#endif + return auction.getWassersteinCost(); +} + bool readDiagramPointSet(const std::string& fname, std::vector>& result) { return readDiagramPointSet(fname.c_str(), result); @@ -84,7 +133,9 @@ bool readDiagramPointSet(const char* fname, std::vector> x >> y)) { +#ifndef FOR_R_TDA std::cerr << "Error in file " << fname << ", line number " << lineNumber << ": cannot parse \"" << line << "\"" << std::endl; +#endif return false; } result.push_back(std::make_pair(x,y)); -- cgit v1.2.3