]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
Bug fix: InterpKernelDEC was not accepting source mesh with extern coords
authorabn <adrien.bruneton@cea.fr>
Tue, 20 Aug 2024 14:38:39 +0000 (16:38 +0200)
committerabn <adrien.bruneton@cea.fr>
Wed, 28 Aug 2024 12:38:08 +0000 (14:38 +0200)
+ improper use of getPointer() instead of getConstPointer) in
  ElementLocator

src/ParaMEDMEM/CommInterface.hxx
src/ParaMEDMEM/ElementLocator.cxx
src/ParaMEDMEMTest/ParaMEDMEMTest.hxx
src/ParaMEDMEMTest/ParaMEDMEMTest_InterpKernelDEC.cxx

index 01691c61b1310f8a0b84c709aaf2e53b79fbbff8..3d52fe9a5b502e0d66ca45e06c063a89944f2414 100644 (file)
@@ -106,7 +106,7 @@ namespace MEDCoupling
 
     int send(void* buffer, int count, MPI_Datatype datatype, int target, int tag, MPI_Comm comm) const { return MPI_Send(buffer,count, datatype, target, tag, comm); }
     int recv(void* buffer, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status* status) const { return MPI_Recv(buffer,count, datatype, source, tag, comm, status); }
-    int sendRecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, 
+    int sendRecv(const void* sendbuf, int sendcount, MPI_Datatype sendtype,
                  int dest, int sendtag, void* recvbuf, int recvcount, 
                  MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm,
                  MPI_Status* status) { return MPI_Sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount, recvtype, source, recvtag, comm,status); }
index 31d7033cdb775fa1fd41b43ddff6d039247e34f8..cc2a90cd6a1514ba601fe5373a1229cb2dfda32e 100644 (file)
@@ -267,12 +267,12 @@ namespace MEDCoupling
     distant_mesh_tmp->resizeForUnserialization(tinyInfoDistant,v1Distant,v2Distant,unusedTinyDistantSts);
     mcIdType nbLocalElems=0;
     mcIdType nbDistElem=0;
-    mcIdType *ptLocal=0;
+    const mcIdType *ptLocal=0;
     mcIdType *ptDist=0;
     if(v1Local)
       {
         nbLocalElems=v1Local->getNbOfElems();
-        ptLocal=v1Local->getPointer();
+        ptLocal=v1Local->getConstPointer();
       }
     if(v1Distant)
       {
@@ -285,12 +285,12 @@ namespace MEDCoupling
                             iprocdistant_in_union,1111,
                             *comm, &status);
     nbLocalElems=0;
-    double *ptLocal2=0;
+    const double *ptLocal2=0;
     double *ptDist2=0;
     if(v2Local)
       {
         nbLocalElems=v2Local->getNbOfElems();
-        ptLocal2=v2Local->getPointer();
+        ptLocal2=v2Local->getConstPointer();
       }
     nbDistElem=0;
     if(v2Distant)
index fb2b339ec7f2df06aa7a1920fef042178547027a..9f12434ddb4ed682355b74a3502a13cae0181199 100644 (file)
@@ -36,6 +36,7 @@ class ParaMEDMEMTest : public CppUnit::TestFixture
   CPPUNIT_TEST(testMPIProcessorGroup_rank);         // >=2 procs
   CPPUNIT_TEST(testBlockTopology_constructor);      // >=2 procs
   CPPUNIT_TEST(testBlockTopology_serialize);        // 1 proc
+  CPPUNIT_TEST(testInterpKernelDEC_ext_coords);     // 1 procs
   CPPUNIT_TEST(testInterpKernelDEC_1D);             // 5 procs
   CPPUNIT_TEST(testInterpKernelDEC_2DCurve);        // 5 procs
   CPPUNIT_TEST(testInterpKernelDEC_2D);             // 5 procs
@@ -104,6 +105,7 @@ public:
   void testMPIProcessorGroup_rank();
   void testBlockTopology_constructor();
   void testBlockTopology_serialize();
+  void testInterpKernelDEC_ext_coords();
   void testInterpKernelDEC_1D();
   void testInterpKernelDEC_2DCurve();
   void testInterpKernelDEC_2D();
index 781c4ed84c4ef40acde2456fbc5668c916b7b447..9cb3069d54f8e742e40a96bf0cb9b916fc36d8c4 100644 (file)
@@ -70,6 +70,70 @@ void ParaMEDMEMTest::testInterpKernelDEC_2DP0P1()
   //testInterpKernelDEC_2D_("P0","P1");
 }
 
+/*! Dummy test to check that we can pass a mesh with a coordinate array not owned by the mesh
+ * as a source for the InterpKernelDEC.
+ */
+void ParaMEDMEMTest::testInterpKernelDEC_ext_coords()
+{
+  using namespace MEDCoupling;
+
+  int size;
+  int rank;
+  MPI_Comm_size(MPI_COMM_WORLD,&size);
+  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
+
+  if (size != 2) return;
+
+  set<int> src_procs;
+  set<int> dest_procs;
+  src_procs.insert(0);
+  dest_procs.insert(1);
+  //
+  MCAuto<MEDCouplingUMesh> mesh(MEDCouplingUMesh::New("Source mesh Proc0",1));
+  ParaMESH *paramesh=0;
+  ParaFIELD *parafieldP0=0;
+  CommInterface interface;
+  ProcessorGroup* src_group = new MEDCoupling::MPIProcessorGroup(interface,src_procs);
+  ProcessorGroup* dest_group = new MEDCoupling::MPIProcessorGroup(interface,dest_procs);
+
+  double coords[4]={0.3,0.7, 0.9,1.0};
+  mcIdType conn[4]={0,1,2,3};
+  mesh->allocateCells(2);
+  mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn);
+  mesh->insertNextCell(INTERP_KERNEL::NORM_SEG2,2,conn+2);
+  mesh->finishInsertingCells();
+  MCAuto<DataArrayDouble> myCoords(DataArrayDouble::New());
+  mesh->setCoords(myCoords);
+
+  if (rank == 0)  // the one with an external array for coords
+    {
+      myCoords->useArray(coords, false, DeallocType::CPP_DEALLOC, 4, 1);
+      paramesh=new ParaMESH(mesh,*src_group,"source mesh");
+    }
+  else
+    {
+      myCoords->alloc(4,1);
+      std::copy(coords, coords+4, myCoords->rwBegin());
+      paramesh=new ParaMESH(mesh,*dest_group,"source mesh");
+    }
+
+  MEDCoupling::ComponentTopology comptopo;
+  parafieldP0 = new ParaFIELD(ON_CELLS,NO_TIME,paramesh, comptopo);
+  double *valueP0=parafieldP0->getField()->getArray()->getPointer();
+  parafieldP0->getField()->setNature(IntensiveMaximum);
+
+  MEDCoupling::InterpKernelDEC dec(*src_group, *dest_group);
+  dec.setMethod("P0");
+  dec.attachLocalField(parafieldP0);
+  // Was crashing here, because of a non const pointer in ElementLocator:
+  dec.synchronize();
+
+  delete parafieldP0;
+  delete paramesh;
+  delete src_group;
+  delete dest_group;
+}
+
 void ParaMEDMEMTest::testInterpKernelDEC_1D()
 {
   int size;