From: abn Date: Tue, 1 Dec 2015 15:31:44 +0000 (+0100) Subject: OverlapDEC: further imp of the load balancing algo to try to keep X-Git-Tag: simple_cmake~4^2 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=de86a608553fc1b1e3d9e4e065a130344671e658;p=tools%2Fmedcoupling.git OverlapDEC: further imp of the load balancing algo to try to keep locality. --- diff --git a/src/ParaMEDMEM/OverlapDEC.cxx b/src/ParaMEDMEM/OverlapDEC.cxx index 9e5221dea..af842fb9d 100644 --- a/src/ParaMEDMEM/OverlapDEC.cxx +++ b/src/ParaMEDMEM/OverlapDEC.cxx @@ -354,4 +354,9 @@ namespace ParaMEDMEM return false; return _group->containsMyRank(); } + + void OverlapDEC::debugPrintWorkSharing(std::ostream & ostr) const + { + _locator->debugPrintWorkSharing(ostr); + } } diff --git a/src/ParaMEDMEM/OverlapDEC.hxx b/src/ParaMEDMEM/OverlapDEC.hxx index a1f29cf22..2d63789b7 100644 --- a/src/ParaMEDMEM/OverlapDEC.hxx +++ b/src/ParaMEDMEM/OverlapDEC.hxx @@ -25,6 +25,7 @@ #include "InterpolationOptions.hxx" #include +#include namespace ParaMEDMEM { @@ -50,9 +51,12 @@ namespace ParaMEDMEM bool isInGroup() const; void setDefaultValue(double val) {_default_field_value = val;} + //! 0 means initial algo from Antho, 1 or 2 means Adrien's algo (2 should be better). Make your choice :-)) void setWorkSharingAlgo(int method) { _load_balancing_algo = method; } + + void debugPrintWorkSharing(std::ostream & ostr) const; private: - int _load_balancing_algo; // 0 means initial algo from Antho, 1 means Adrien's algo. Make your choice :-)) + int _load_balancing_algo; bool _own_group; OverlapInterpolationMatrix* _interpolation_matrix; diff --git a/src/ParaMEDMEM/OverlapElementLocator.cxx b/src/ParaMEDMEM/OverlapElementLocator.cxx index 93984ccce..81cb987db 100644 --- a/src/ParaMEDMEM/OverlapElementLocator.cxx +++ b/src/ParaMEDMEM/OverlapElementLocator.cxx @@ -58,13 +58,17 @@ namespace ParaMEDMEM _comm=getCommunicator(); computeBoundingBoxesAndInteractionList(); - if (workSharingAlgo == 0) - computeTodoList_original(); - else - if(workSharingAlgo == 1) - computeTodoList_new(); - else + switch(workSharingAlgo) + { + case 0: + computeTodoList_original(); break; + case 1: + computeTodoList_new(false); break; + case 2: + computeTodoList_new(true); break; + default: throw INTERP_KERNEL::Exception("OverlapElementLocator::OverlapElementLocator(): invalid algorithm selected!"); + } fillProcToSend(); } @@ -125,8 +129,8 @@ namespace ParaMEDMEM _proc_pairs.resize(_group.size()); for(int i=0;i<_group.size();i++) for(int j=0;j<_group.size();j++) - if(intersectsBoundingBox(i,j)) - _proc_pairs[i].push_back(j); + if(intersectsBoundingBox(i,j)) + _proc_pairs[i].push_back(j); } void OverlapElementLocator::computeTodoList_original() @@ -160,13 +164,14 @@ namespace ParaMEDMEM } /* More efficient (?) work sharing algorithm: a job (i,j) is initially assigned twice: to proc#i and to proc#j. - * Then try to reduce as much as possible the variance of the num of jobs per proc: + * Then try to reduce as much as possible the variance of the num of jobs per proc by selecting the right duplicate + * to remove: * - take the most loaded proc i, * + select the job (i,j) for which proc#j is the less loaded - * + remove this job from proc#i + * + remove this job from proc#i, and mark it as 'unremovable' from proc#j * - repeat until no more duplicates are found */ - void OverlapElementLocator::computeTodoList_new() + void OverlapElementLocator::computeTodoList_new(bool revertIter) { using namespace std; int infinity = std::numeric_limits::max(); @@ -214,7 +219,7 @@ namespace ParaMEDMEM int sz = (*itVector).size(); if (proc_valid[procID] && sz > max_sz) { - max_sz = (*itVector).size(); + max_sz = sz; max_id = procID; } } @@ -226,9 +231,21 @@ namespace ParaMEDMEM int min_sz = infinity; map & max_map = full_set[max_id]; ProcCouple hit_cpl = make_pair(-1,-1); - for(itMap=max_map.begin(); itMap != max_map.end(); itMap++) - if ((*itMap).second < min_sz) - hit_cpl = (*itMap).first; + if(revertIter) + { + // Use a reverse iterator here increases our chances to hit a couple of the form (i, myProcId) + // meaning that the final matrix computed won't have to be sent: save some comm. + map ::const_reverse_iterator ritMap; + for(ritMap=max_map.rbegin(); ritMap != max_map.rend(); ritMap++) + if ((*ritMap).second < min_sz) + hit_cpl = (*ritMap).first; + } + else + { + for(itMap=max_map.begin(); itMap != max_map.end(); itMap++) + if ((*itMap).second < min_sz) + hit_cpl = (*itMap).first; + } if (hit_cpl.first == -1) { // Plouf. Current proc 'max_id' can not be reduced. Invalid it: @@ -275,6 +292,15 @@ namespace ParaMEDMEM #endif } + void OverlapElementLocator::debugPrintWorkSharing(std::ostream & ostr) const + { + std::vector< std::vector< ProcCouple > >::const_iterator it = _all_todo_lists.begin(); + ostr << "TODO list lengths: "; + for(; it != _all_todo_lists.end(); ++it) + ostr << (*it).size() << " "; + ostr << "\n"; + } + void OverlapElementLocator::fillProcToSend() { // Feeding now '_procs_to_send*'. A same id can appears twice. The second parameter in pair means what diff --git a/src/ParaMEDMEM/OverlapElementLocator.hxx b/src/ParaMEDMEM/OverlapElementLocator.hxx index bf5580232..3ffd028db 100644 --- a/src/ParaMEDMEM/OverlapElementLocator.hxx +++ b/src/ParaMEDMEM/OverlapElementLocator.hxx @@ -59,10 +59,11 @@ namespace ParaMEDMEM const MEDCouplingPointSet *getTargetMesh(int procId) const; const DataArrayInt *getTargetIds(int procId) const; bool isInMyTodoList(int i, int j) const; + void debugPrintWorkSharing(std::ostream & ostr) const; private: void computeBoundingBoxesAndInteractionList(); void computeTodoList_original(); - void computeTodoList_new(); + void computeTodoList_new(bool revertIter); void fillProcToSend(); bool intersectsBoundingBox(int i, int j) const; void sendLocalMeshTo(int procId, bool sourceOrTarget, OverlapInterpolationMatrix& matrix) const; diff --git a/src/ParaMEDMEMTest/ParaMEDMEMTest.hxx b/src/ParaMEDMEMTest/ParaMEDMEMTest.hxx index 4272b9985..a8d777962 100644 --- a/src/ParaMEDMEMTest/ParaMEDMEMTest.hxx +++ b/src/ParaMEDMEMTest/ParaMEDMEMTest.hxx @@ -49,8 +49,10 @@ class ParaMEDMEMTest : public CppUnit::TestFixture CPPUNIT_TEST(testInterpKernelDEC3DSurfEmptyBBox); // 3 procs CPPUNIT_TEST(testOverlapDEC1); // 3 procs CPPUNIT_TEST(testOverlapDEC1_bis); // 3 procs + CPPUNIT_TEST(testOverlapDEC1_ter); // 3 procs CPPUNIT_TEST(testOverlapDEC2); // 3 procs CPPUNIT_TEST(testOverlapDEC2_bis); // 3 procs + CPPUNIT_TEST(testOverlapDEC2_ter); // 3 procs CPPUNIT_TEST(testOverlapDEC3); // 2 procs CPPUNIT_TEST(testOverlapDEC4); // 2 procs @@ -112,8 +114,10 @@ public: void testInterpKernelDEC3DSurfEmptyBBox(); void testOverlapDEC1(); void testOverlapDEC1_bis(); + void testOverlapDEC1_ter(); void testOverlapDEC2(); void testOverlapDEC2_bis(); + void testOverlapDEC2_ter(); void testOverlapDEC3(); // void testOverlapDEC3_bis(); void testOverlapDEC4(); diff --git a/src/ParaMEDMEMTest/ParaMEDMEMTest_OverlapDEC.cxx b/src/ParaMEDMEMTest/ParaMEDMEMTest_OverlapDEC.cxx index 995c79c87..f6b14cf49 100644 --- a/src/ParaMEDMEMTest/ParaMEDMEMTest_OverlapDEC.cxx +++ b/src/ParaMEDMEMTest/ParaMEDMEMTest_OverlapDEC.cxx @@ -335,7 +335,6 @@ void prepareData2_buildOneSquare(MEDCouplingUMesh* & meshS_0, MEDCouplingUMesh* meshT_0->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT+3); meshT_0->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT+6); meshT_0->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,connT+9); - } /** @@ -491,6 +490,7 @@ void testOverlapDEC_generic(int workSharingAlgo, double bbAdj) dec.attachSourceLocalField(mcfieldS); dec.attachTargetLocalField(mcfieldT); dec.synchronize(); +// dec.debugPrintWorkSharing(std::cout); dec.sendRecvData(true); // MPI_Barrier(MPI_COMM_WORLD); @@ -527,13 +527,17 @@ void ParaMEDMEMTest::testOverlapDEC1() void ParaMEDMEMTest::testOverlapDEC1_bis() { - /*!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! - * Same BB hack as above - * !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! - */ + // Same BB hack as above testOverlapDEC_generic(1,-1.0e-12); } +void ParaMEDMEMTest::testOverlapDEC1_ter() +{ + // Same BB hack as above + testOverlapDEC_generic(2, -1.0e-12); +} + + /*! * Same as testOverlapDEC1() but with regular bounding boxes. If you're looking for a nice debug case, * testOverlapDEC1() is identical in terms of geometry and field values, and more appropriate. @@ -548,6 +552,11 @@ void ParaMEDMEMTest::testOverlapDEC2_bis() testOverlapDEC_generic(1,1.0e-12); } +void ParaMEDMEMTest::testOverlapDEC2_ter() +{ + testOverlapDEC_generic(2,1.0e-12); +} + /*! Test focused on the mapping of cell IDs. * (i.e. when only part of the source/target mesh is transmitted)