Salome HOME
[bos #38048] [EDF] (2023-T3) PARAMEDMEM Ergonomy.
[tools/medcoupling.git] / src / ParaMEDMEMTest / ParaMEDMEMTest_ByStringMPIProcessorGroup.cxx
1 // Copyright (C) 2007-2023  CEA, EDF
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "ParaMEDMEMTest.hxx"
21 #include <cppunit/TestAssert.h>
22 #include "CommInterface.hxx"
23 #include "ProcessorGroup.hxx"
24 #include "ByStringMPIProcessorGroup.hxx"
25
26 #include <string>
27
28 // use this define to enable lines, execution of which leads to Segmentation Fault
29 #define ENABLE_FAULTS
30
31 // use this define to enable CPPUNIT asserts and fails, showing bugs
32 #define ENABLE_FORCED_FAILURES
33
34
35 using namespace std;
36 using namespace MEDCoupling;
37  
38 /*
39  * Check methods defined in MPPIProcessorGroup.hxx
40  *
41  (+) ByStringMPIProcessorGroup(const CommInterface& interface);
42  (+) ByStringMPIProcessorGroup(const CommInterface& interface, std::string& codeTag, const MPI_Comm& world_comm );
43  (+) ByStringMPIProcessorGroup(ByStringMPIProcessorGroup& other );
44 */
45  
46 void ParaMEDMEMTest::testByStringMPIProcessorGroup_constructor()
47 {
48   CommInterface comm_interface;
49   ByStringMPIProcessorGroup* group = new ByStringMPIProcessorGroup(comm_interface);
50   int size;
51   MPI_Comm_size(MPI_COMM_WORLD, &size);
52   CPPUNIT_ASSERT_EQUAL(size,group->size());
53   int size2;
54   const MPI_Comm* communicator=group->getComm();
55   MPI_Comm_size(*communicator, &size2);
56   CPPUNIT_ASSERT_EQUAL(size,size2);
57   delete group;  
58 }
59  
60 void ParaMEDMEMTest::testByStringMPIProcessorGroup_stringconstructor()
61 {
62   int size, rankId;
63   MPI_Comm_size(MPI_COMM_WORLD, &size);
64   MPI_Comm_rank(MPI_COMM_WORLD, &rankId);
65
66   if (size != 3)
67     return;
68   
69   std::string myTag;
70   if ( rankId == 0 || rankId == 2 ) 
71     myTag = "group0";
72   else
73     myTag = "gr1";
74
75   CommInterface comm_interface;
76   ByStringMPIProcessorGroup * group = new ByStringMPIProcessorGroup(comm_interface,myTag,MPI_COMM_WORLD);
77   ByStringMPIProcessorGroup * copygroup = new ByStringMPIProcessorGroup(*group);
78   CPPUNIT_ASSERT(group);
79   CPPUNIT_ASSERT(copygroup);
80
81   std::set<int> ranksInGroup = group->getProcIDs();
82   std::set<int> ranksInCopiedGroup = group->getProcIDs();
83   if ( rankId == 0 || rankId == 2 )  
84   {
85     CPPUNIT_ASSERT_EQUAL( (int)ranksInGroup.size(), 2 );
86     CPPUNIT_ASSERT_EQUAL( (int)ranksInCopiedGroup.size(), 2 );
87   }
88   else
89   {
90     CPPUNIT_ASSERT_EQUAL( (int)ranksInGroup.size(), 1 );
91     CPPUNIT_ASSERT_EQUAL( (int)ranksInCopiedGroup.size(), 1 );
92   }  
93   CPPUNIT_ASSERT( group->contains(rankId) );
94   CPPUNIT_ASSERT( copygroup->contains(rankId) );
95   delete group;
96   delete copygroup;
97 }