]> SALOME platform Git repositories - modules/med.git/commitdiff
Salome HOME
Modification of C++ API of MEDCouplingMesh::getCellsContainingPoints to avoid copy.
authorageay <ageay>
Thu, 8 Aug 2013 13:16:37 +0000 (13:16 +0000)
committerageay <ageay>
Thu, 8 Aug 2013 13:16:37 +0000 (13:16 +0000)
src/MEDCoupling/MEDCouplingFieldDiscretization.cxx
src/MEDCoupling/MEDCouplingMesh.cxx
src/MEDCoupling/MEDCouplingMesh.hxx
src/MEDCoupling/MEDCouplingRemapper.cxx
src/MEDCoupling/MEDCouplingUMesh.cxx
src/MEDCoupling/MEDCouplingUMesh.hxx
src/MEDCoupling/Test/MEDCouplingBasicsTest1.cxx
src/MEDCoupling/Test/MEDCouplingExamplesTest.cxx
src/MEDCoupling_Swig/MEDCouplingCommon.i

index a927a731c068e65c8c6dafb02d2c7831e300dc77..ff36cb4628c0f6fce1edd4ae4dd45b984c265596 100644 (file)
@@ -621,8 +621,9 @@ DataArrayDouble *MEDCouplingFieldDiscretizationP0::getValueOnMulti(const DataArr
 {
   if(!mesh)
     throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP0::getValueOnMulti : NULL input mesh !");
-  std::vector<int> elts,eltsIndex;
-  mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,elts,eltsIndex);
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsArr,eltsIndexArr;
+  mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,eltsArr,eltsIndexArr);
+  const int *elts(eltsArr->begin()),*eltsIndex(eltsIndexArr->begin());
   int spaceDim=mesh->getSpaceDimension();
   int nbOfComponents=arr->getNumberOfComponents();
   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
@@ -1002,8 +1003,9 @@ DataArrayDouble *MEDCouplingFieldDiscretizationP1::getValueOnMulti(const DataArr
 {
   if(!mesh)
     throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationP1::getValueOnMulti : NULL input mesh !");
-  std::vector<int> elts,eltsIndex;
-  mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,elts,eltsIndex);
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsArr,eltsIndexArr;
+  mesh->getCellsContainingPoints(loc,nbOfPoints,_precision,eltsArr,eltsIndexArr);
+  const int *elts(eltsArr->begin()),*eltsIndex(eltsIndexArr->begin());
   int spaceDim=mesh->getSpaceDimension();
   int nbOfComponents=arr->getNumberOfComponents();
   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> ret=DataArrayDouble::New();
index 3408e9e22bddd3eb1de1b3f59ac1cf70b8ef407c..83d2764c87c0d8ff3705ef2cdb91c8c2f4949a2a 100644 (file)
@@ -631,8 +631,8 @@ void MEDCouplingMesh::getCellsContainingPoint(const double *pos, double eps, std
  *         this->getSpaceDimension() * \a nbOfPoints 
  *  \param [in] nbOfPoints - number of points to locate within \a this mesh.
  *  \param [in] eps - radius of balls (i.e. the precision).
- *  \param [in,out] elts - vector returning ids of found cells.
- *  \param [in,out] eltsIndex - an array, of length \a nbOfPoints + 1,
+ *  \param [out] elts - vector returning ids of found cells.
+ *  \param [out] eltsIndex - an array, of length \a nbOfPoints + 1,
  *         dividing cell ids in \a elts into groups each referring to one
  *         point. Its every element (except the last one) is an index pointing to the
  *         first id of a group of cells. For example cells in contact with the *i*-th
@@ -645,23 +645,22 @@ void MEDCouplingMesh::getCellsContainingPoint(const double *pos, double eps, std
  *  \ref cpp_mcumesh_getCellsContainingPoints "Here is a C++ example".<br>
  *  \ref  py_mcumesh_getCellsContainingPoints "Here is a Python example".
  */
-void MEDCouplingMesh::getCellsContainingPoints(const double *pos, int nbOfPoints, double eps, std::vector<int>& elts, std::vector<int>& eltsIndex) const
+void MEDCouplingMesh::getCellsContainingPoints(const double *pos, int nbOfPoints, double eps, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& elts, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& eltsIndex) const
 {
-  eltsIndex.resize(nbOfPoints+1);
-  eltsIndex[0]=0;
-  elts.clear();
-  int spaceDim=getSpaceDimension();
-  const double *work=pos;
+  eltsIndex=DataArrayInt::New(); elts=DataArrayInt::New(); eltsIndex->alloc(nbOfPoints+1,1); eltsIndex->setIJ(0,0,0); elts->alloc(0,1);
+  int *eltsIndexPtr(eltsIndex->getPointer());
+  int spaceDim(getSpaceDimension());
+  const double *work(pos);
   for(int i=0;i<nbOfPoints;i++,work+=spaceDim)
     {
       int ret=getCellContainingPoint(work,eps);
       if(ret>=0)
         {
-          elts.push_back(ret);
-          eltsIndex[i+1]=eltsIndex[i]+1;
+          elts->pushBackSilent(ret);
+          eltsIndexPtr[i+1]=eltsIndexPtr[i]+1;
         }
       else
-        eltsIndex[i+1]=eltsIndex[i];
+        eltsIndexPtr[i+1]=eltsIndexPtr[i];
     }
 }
 
index 053d60e07c4a95565d1ee7efd237a223b0d766ce..e98443682c3025ba18223282571155fc41569c26 100644 (file)
@@ -25,6 +25,8 @@
 #include "MEDCouplingTimeLabel.hxx"
 #include "MEDCouplingRefCountObject.hxx"
 #include "NormalizedUnstructuredMesh.hxx"
+#include "MEDCouplingAutoRefCountObjectPtr.hxx"
+
 #include "InterpKernelException.hxx"
 
 #include <set>
@@ -108,7 +110,7 @@ namespace ParaMEDMEM
     virtual MEDCouplingFieldDouble *getMeasureFieldOnNode(bool isAbs) const = 0;
     virtual int getCellContainingPoint(const double *pos, double eps) const = 0;
     virtual void getCellsContainingPoint(const double *pos, double eps, std::vector<int>& elts) const;
-    virtual void getCellsContainingPoints(const double *pos, int nbOfPoints, double eps, std::vector<int>& elts, std::vector<int>& eltsIndex) const;
+    virtual void getCellsContainingPoints(const double *pos, int nbOfPoints, double eps, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& elts, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& eltsIndex) const;
     virtual MEDCouplingFieldDouble *fillFromAnalytic(TypeOfField t, int nbOfComp, FunctionToEvaluate func) const;
     virtual MEDCouplingFieldDouble *fillFromAnalytic(TypeOfField t, int nbOfComp, const char *func) const;
     virtual MEDCouplingFieldDouble *fillFromAnalytic2(TypeOfField t, int nbOfComp, const char *func) const;
index 2879ada37d4d0f4972567f9e67b67f03c259bdc7..30356fd72d37a9fb0ba54de0e0a6ce8204b6717d 100644 (file)
@@ -739,12 +739,12 @@ int MEDCouplingRemapper::prepareNotInterpKernelOnlyGaussGauss() throw(INTERP_KER
   const int *srcOffsetArrPtr=srcOffsetArr->begin();
   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> srcLoc=_src_ft->getLocalizationOfDiscr();
   const double *srcLocPtr=srcLoc->begin();
-  std::vector<int> elts,eltsIndex;
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsArr,eltsIndexArr;
   int trgNbOfGaussPts=trgLoc->getNumberOfTuples();
   _matrix.resize(trgNbOfGaussPts);
-  _src_ft->getMesh()->getCellsContainingPoints(trgLoc->begin(),trgNbOfGaussPts,getPrecision(),elts,eltsIndex);
-  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsIndex2=DataArrayInt::New(); eltsIndex2->useArray(&eltsIndex[0],false,CPP_DEALLOC,(int)eltsIndex.size(),1);
-  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfSrcCellsShTrgPts=eltsIndex2->deltaShiftIndex();
+  _src_ft->getMesh()->getCellsContainingPoints(trgLoc->begin(),trgNbOfGaussPts,getPrecision(),eltsArr,eltsIndexArr);
+  const int *elts(eltsArr->begin()),*eltsIndex(eltsIndexArr->begin());
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> nbOfSrcCellsShTrgPts(eltsIndexArr->deltaShiftIndex());
   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ids0=nbOfSrcCellsShTrgPts->getIdsNotEqual(0);
   for(const int *trgId=ids0->begin();trgId!=ids0->end();trgId++)
     {
index 8015d9ed5b531292e7713d44d83b91d777ba5118..efa6c63faeba05ec256e9f1b763982739d74f6a4 100644 (file)
@@ -33,7 +33,6 @@
 #include "InterpKernelMeshQuality.hxx"
 #include "InterpKernelCellSimplify.hxx"
 #include "InterpKernelGeo2DEdgeArcCircle.hxx"
-#include "MEDCouplingAutoRefCountObjectPtr.hxx"
 #include "InterpKernelAutoPtr.hxx"
 #include "InterpKernelGeo2DNode.hxx"
 #include "InterpKernelGeo2DEdgeLin.hxx"
@@ -4014,7 +4013,7 @@ int MEDCouplingUMesh::getCellContainingPoint(const double *pos, double eps) cons
  *          faster. 
  *  \param [in] pos - array of coordinates of the ball central point.
  *  \param [in] eps - ball radius.
- *  \param [in,out] elts - vector returning ids of the found cells. It is cleared
+ *  \param [out] elts - vector returning ids of the found cells. It is cleared
  *         before inserting ids.
  *  \throw If the coordinates array is not set.
  *  \throw If \a this->getMeshDimension() != \a this->getSpaceDimension().
@@ -4024,8 +4023,9 @@ int MEDCouplingUMesh::getCellContainingPoint(const double *pos, double eps) cons
  */
 void MEDCouplingUMesh::getCellsContainingPoint(const double *pos, double eps, std::vector<int>& elts) const
 {
-  std::vector<int> eltsIndex;
-  getCellsContainingPoints(pos,1,eps,elts,eltsIndex);
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> eltsUg,eltsIndexUg;
+  getCellsContainingPoints(pos,1,eps,eltsUg,eltsIndexUg);
+  elts.clear(); elts.insert(elts.end(),eltsUg->begin(),eltsUg->end());
 }
 
 /// @cond INTERNAL
@@ -4154,11 +4154,10 @@ namespace ParaMEDMEM
 
 template<int SPACEDIM>
 void MEDCouplingUMesh::getCellsContainingPointsAlg(const double *coords, const double *pos, int nbOfPoints,
-                                                   double eps, std::vector<int>& elts, std::vector<int>& eltsIndex) const
+                                                   double eps, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& elts, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& eltsIndex) const
 {
-  eltsIndex.resize(nbOfPoints+1);
-  eltsIndex[0]=0;
-  elts.clear();
+  elts=DataArrayInt::New(); eltsIndex=DataArrayInt::New(); eltsIndex->alloc(nbOfPoints+1,1); eltsIndex->setIJ(0,0,0); elts->alloc(0,1);
+  int *eltsIndexPtr(eltsIndex->getPointer());
   MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> bboxArr(getBoundingBoxForBBTree());
   const double *bbox(bboxArr->begin());
   int nbOfCells=getNumberOfCells();
@@ -4168,7 +4167,7 @@ void MEDCouplingUMesh::getCellsContainingPointsAlg(const double *coords, const d
   BBTree<SPACEDIM,int> myTree(&bbox[0],0,0,nbOfCells,-eps);
   for(int i=0;i<nbOfPoints;i++)
     {
-      eltsIndex[i+1]=eltsIndex[i];
+      eltsIndexPtr[i+1]=eltsIndexPtr[i];
       for(int j=0;j<SPACEDIM;j++)
         {
           bb[2*j]=pos[SPACEDIM*i+j];
@@ -4183,8 +4182,8 @@ void MEDCouplingUMesh::getCellsContainingPointsAlg(const double *coords, const d
                                                                                                (INTERP_KERNEL::NormalizedCellType)conn[connI[*iter]],
                                                                                                coords,conn+connI[*iter]+1,sz,eps))
             {
-              eltsIndex[i+1]++;
-              elts.push_back(*iter);
+              eltsIndexPtr[i+1]++;
+              elts->pushBackSilent(*iter);
             }
         }
     }
@@ -4198,8 +4197,8 @@ void MEDCouplingUMesh::getCellsContainingPointsAlg(const double *coords, const d
  *         this->getSpaceDimension() * \a nbOfPoints 
  *  \param [in] nbOfPoints - number of points to locate within \a this mesh.
  *  \param [in] eps - radius of balls (i.e. the precision).
- *  \param [in,out] elts - vector returning ids of found cells.
- *  \param [in,out] eltsIndex - an array, of length \a nbOfPoints + 1,
+ *  \param [out] elts - vector returning ids of found cells.
+ *  \param [out] eltsIndex - an array, of length \a nbOfPoints + 1,
  *         dividing cell ids in \a elts into groups each referring to one
  *         point. Its every element (except the last one) is an index pointing to the
  *         first id of a group of cells. For example cells in contact with the *i*-th
@@ -4215,7 +4214,7 @@ void MEDCouplingUMesh::getCellsContainingPointsAlg(const double *coords, const d
  *  \ref  py_mcumesh_getCellsContainingPoints "Here is a Python example".
  */
 void MEDCouplingUMesh::getCellsContainingPoints(const double *pos, int nbOfPoints, double eps,
-                                                std::vector<int>& elts, std::vector<int>& eltsIndex) const
+                                                MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& elts, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& eltsIndex) const
 {
   int spaceDim=getSpaceDimension();
   int mDim=getMeshDimension();
index 64db7ce5742be7ab46ee73852cbf2b006d5cb87d..49d6c89a2aba36d6bc68603df71abbe6aeeac28c 100644 (file)
@@ -158,7 +158,7 @@ namespace ParaMEDMEM
     MEDCOUPLING_EXPORT DataArrayDouble *distanceToPoints(const DataArrayDouble *pts, DataArrayInt *& cellIds) const throw(INTERP_KERNEL::Exception);
     MEDCOUPLING_EXPORT int getCellContainingPoint(const double *pos, double eps) const;
     MEDCOUPLING_EXPORT void getCellsContainingPoint(const double *pos, double eps, std::vector<int>& elts) const;
-    MEDCOUPLING_EXPORT void getCellsContainingPoints(const double *pos, int nbOfPoints, double eps, std::vector<int>& elts, std::vector<int>& eltsIndex) const;
+    MEDCOUPLING_EXPORT void getCellsContainingPoints(const double *pos, int nbOfPoints, double eps, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& elts, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& eltsIndex) const;
     MEDCOUPLING_EXPORT void checkButterflyCells(std::vector<int>& cells, double eps=1e-12) const;
     MEDCOUPLING_EXPORT DataArrayInt *convexEnvelop2D() throw(INTERP_KERNEL::Exception);
     MEDCOUPLING_EXPORT DataArrayInt *findAndCorrectBadOriented3DExtrudedCells() throw(INTERP_KERNEL::Exception);
@@ -283,7 +283,7 @@ namespace ParaMEDMEM
     DataArrayInt *convertLinearCellsToQuadratic3D1(DataArrayInt *&conn, DataArrayInt *&connI, DataArrayDouble *& coords, std::set<INTERP_KERNEL::NormalizedCellType>& types) const throw(INTERP_KERNEL::Exception);
     template<int SPACEDIM>
     void getCellsContainingPointsAlg(const double *coords, const double *pos, int nbOfPoints,
-                                     double eps, std::vector<int>& elts, std::vector<int>& eltsIndex) const;
+                                     double eps, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& elts, MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& eltsIndex) const;
 /// @cond INTERNAL
     static MEDCouplingUMesh *MergeUMeshesLL(std::vector<const MEDCouplingUMesh *>& a) throw(INTERP_KERNEL::Exception);
     typedef int (*DimM1DescNbrer)(int id, unsigned nb, const INTERP_KERNEL::CellModel& cm, bool compute, const int *conn1, const int *conn2);
index 9808e0205d113042f28bc8745a07d098741c4fc9..fc06ef5a3b80687316bed8a358b63bb8b4ea60c7 100644 (file)
@@ -2162,25 +2162,24 @@ void MEDCouplingBasicsTest1::testGetCellsContainingPoint()
 {
   MEDCouplingUMesh *targetMesh=build2DTargetMesh_1();
   double pos[12]={0.,0.,0.4,0.4,0.,0.4,0.1,0.1,0.25,0.,0.65,0.};
-  std::vector<int> t1,t2;
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> t1,t2;
   //2D basic
   targetMesh->getCellsContainingPoints(pos,6,1e-12,t1,t2);
-  CPPUNIT_ASSERT_EQUAL(6,(int)t1.size());
-  CPPUNIT_ASSERT_EQUAL(7,(int)t2.size());
+  CPPUNIT_ASSERT_EQUAL(6,(int)t1->getNbOfElems());
+  CPPUNIT_ASSERT_EQUAL(7,(int)t2->getNbOfElems());
   const int expectedValues1[6]={0,4,3,0,1,2};
   const int expectedValues2[7]={0,1,2,3,4,5,6};
-  CPPUNIT_ASSERT(std::equal(t1.begin(),t1.end(),expectedValues1));
-  CPPUNIT_ASSERT(std::equal(t2.begin(),t2.end(),expectedValues2));
+  CPPUNIT_ASSERT(std::equal(t1->begin(),t1->end(),expectedValues1));
+  CPPUNIT_ASSERT(std::equal(t2->begin(),t2->end(),expectedValues2));
   //2D with no help of bounding box.
   double center[2]={0.2,0.2};
   MEDCouplingPointSet::Rotate2DAlg(center,0.78539816339744830962,6,pos);
   targetMesh->rotate(center,0,0.78539816339744830962);
-  t1.clear(); t2.clear();
   targetMesh->getCellsContainingPoints(pos,6,1e-12,t1,t2);
-  CPPUNIT_ASSERT_EQUAL(6,(int)t1.size());
-  CPPUNIT_ASSERT_EQUAL(7,(int)t2.size());
-  CPPUNIT_ASSERT(std::equal(t1.begin(),t1.end(),expectedValues1));
-  CPPUNIT_ASSERT(std::equal(t2.begin(),t2.end(),expectedValues2));
+  CPPUNIT_ASSERT_EQUAL(6,(int)t1->getNbOfElems());
+  CPPUNIT_ASSERT_EQUAL(7,(int)t2->getNbOfElems());
+  CPPUNIT_ASSERT(std::equal(t1->begin(),t1->end(),expectedValues1));
+  CPPUNIT_ASSERT(std::equal(t2->begin(),t2->end(),expectedValues2));
   //2D outside
   const double pos1bis[2]={-0.3303300858899107,-0.11819805153394641};
   CPPUNIT_ASSERT_EQUAL(-1,targetMesh->getCellContainingPoint(pos1bis,1e-12));
@@ -2188,17 +2187,18 @@ void MEDCouplingBasicsTest1::testGetCellsContainingPoint()
   //test limits 2D
   targetMesh=build2DTargetMesh_1();
   const double pos2[2]={0.2,-0.05};
-  t1.clear();
-  targetMesh->getCellsContainingPoint(pos2,1e-12,t1);
-  CPPUNIT_ASSERT_EQUAL(2,(int)t1.size());
+  std::vector<int> t11;
+  t11.clear();
+  targetMesh->getCellsContainingPoint(pos2,1e-12,t11);
+  CPPUNIT_ASSERT_EQUAL(2,(int)t11.size());
   const int expectedValues3[2]={0,1};
-  CPPUNIT_ASSERT(std::equal(t1.begin(),t1.end(),expectedValues3));
+  CPPUNIT_ASSERT(std::equal(t11.begin(),t11.end(),expectedValues3));
   const double pos3[2]={0.2,0.2};
-  t1.clear();
-  targetMesh->getCellsContainingPoint(pos3,1e-12,t1);
-  CPPUNIT_ASSERT_EQUAL(5,(int)t1.size());
+  t11.clear();
+  targetMesh->getCellsContainingPoint(pos3,1e-12,t11);
+  CPPUNIT_ASSERT_EQUAL(5,(int)t11.size());
   const int expectedValues4[5]={0,1,2,3,4};
-  CPPUNIT_ASSERT(std::equal(t1.begin(),t1.end(),expectedValues4));
+  CPPUNIT_ASSERT(std::equal(t11.begin(),t11.end(),expectedValues4));
   CPPUNIT_ASSERT_EQUAL(0,targetMesh->getCellContainingPoint(pos3,1e-12));
   targetMesh->decrRef();
   //3D
@@ -2206,17 +2206,17 @@ void MEDCouplingBasicsTest1::testGetCellsContainingPoint()
   const double pos4[3]={25.,25.,25.};
   CPPUNIT_ASSERT_EQUAL(0,targetMesh->getCellContainingPoint(pos4,1e-12));
   const double pos5[3]={50.,50.,50.};
-  t1.clear();
-  targetMesh->getCellsContainingPoint(pos5,1e-12,t1);
-  CPPUNIT_ASSERT_EQUAL(8,(int)t1.size());
+  t11.clear();
+  targetMesh->getCellsContainingPoint(pos5,1e-12,t11);
+  CPPUNIT_ASSERT_EQUAL(8,(int)t11.size());
   const int expectedValues5[8]={0,1,2,3,4,5,6,7};
-  CPPUNIT_ASSERT(std::equal(t1.begin(),t1.end(),expectedValues5));
+  CPPUNIT_ASSERT(std::equal(t11.begin(),t11.end(),expectedValues5));
   const double pos6[3]={0., 50., 0.};
-  t1.clear();
-  targetMesh->getCellsContainingPoint(pos6,1e-12,t1);
-  CPPUNIT_ASSERT_EQUAL(2,(int)t1.size());
+  t11.clear();
+  targetMesh->getCellsContainingPoint(pos6,1e-12,t11);
+  CPPUNIT_ASSERT_EQUAL(2,(int)t11.size());
   const int expectedValues6[2]={0,2};
-  CPPUNIT_ASSERT(std::equal(t1.begin(),t1.end(),expectedValues6));
+  CPPUNIT_ASSERT(std::equal(t11.begin(),t11.end(),expectedValues6));
   //3D outside
   const double pos7[3]={-1.0,-1.0,0.};
   CPPUNIT_ASSERT_EQUAL(-1,targetMesh->getCellContainingPoint(pos7,1e-12));
index be1e0bf730614d5460f3ea88d7eebcaa9d1decae..86a9f7e6f19f789de0e9b3f99e9dd16807d14d0e 100644 (file)
@@ -945,12 +945,12 @@ void CppExample_MEDCouplingUMesh_getCellsContainingPoints()
                             0.3, 0.3,              // point located somewhere inside the mesh
                             coords[2], coords[3]}; // point at the node #1
   const double eps = 1e-4; // ball radius
-  std::vector<int> cells, cellsIndex;
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> cells, cellsIndex;
   mesh->getCellsContainingPoints( pos, 3, eps, cells, cellsIndex );
   const int cellsExpected[3]={4, 0, 1};
   const int cellsIndexExpected[4]={0, 0, 1, 3};
-  CPPUNIT_ASSERT(std::equal( cellsExpected,      cellsExpected+3,      &cells[0]));
-  CPPUNIT_ASSERT(std::equal( cellsIndexExpected, cellsIndexExpected+4, &cellsIndex[0]));
+  CPPUNIT_ASSERT(std::equal( cellsExpected,      cellsExpected+3,      cells->begin()));
+  CPPUNIT_ASSERT(std::equal( cellsIndexExpected, cellsIndexExpected+4, cellsIndex->begin()));
   //! [CppSnippet_MEDCouplingUMesh_getCellsContainingPoints_2]
 }
 
index 92bdc3a2a3a74cb96273a669a313d00bd62849d6..3c3e6a0be68688031025c37edb915799b4f084ca 100644 (file)
@@ -755,23 +755,17 @@ namespace ParaMEDMEM
            int spaceDim=self->getSpaceDimension();
            const char msg[]="Python wrap of MEDCouplingMesh::getCellsContainingPoint : ";
            const double *pos=convertObjToPossibleCpp5_Safe(p,sw,val,a,aa,bb,msg,nbOfPoints,spaceDim,true);
-           std::vector<int> elts,eltsIndex;
+           MEDCouplingAutoRefCountObjectPtr<DataArrayInt> elts,eltsIndex;
            self->getCellsContainingPoints(pos,nbOfPoints,eps,elts,eltsIndex);
-           MEDCouplingAutoRefCountObjectPtr<DataArrayInt> d0=DataArrayInt::New();
-           MEDCouplingAutoRefCountObjectPtr<DataArrayInt> d1=DataArrayInt::New();
-           d0->alloc(elts.size(),1);
-           d1->alloc(eltsIndex.size(),1);
-           std::copy(elts.begin(),elts.end(),d0->getPointer());
-           std::copy(eltsIndex.begin(),eltsIndex.end(),d1->getPointer());
            PyObject *ret=PyTuple_New(2);
-           PyTuple_SetItem(ret,0,SWIG_NewPointerObj(SWIG_as_voidptr(d0.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
-           PyTuple_SetItem(ret,1,SWIG_NewPointerObj(SWIG_as_voidptr(d1.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
+           PyTuple_SetItem(ret,0,SWIG_NewPointerObj(SWIG_as_voidptr(elts.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
+           PyTuple_SetItem(ret,1,SWIG_NewPointerObj(SWIG_as_voidptr(eltsIndex.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
            return ret;
          }
 
          PyObject *getCellsContainingPoints(PyObject *p, double eps) const throw(INTERP_KERNEL::Exception)
          {
-           std::vector<int> elts,eltsIndex;
+           MEDCouplingAutoRefCountObjectPtr<DataArrayInt> elts,eltsIndex;
            int spaceDim=self->getSpaceDimension();
            void *da=0;
            int res1=SWIG_ConvertPtr(p,&da,SWIGTYPE_p_ParaMEDMEM__DataArrayDouble, 0 |  0 );
@@ -800,15 +794,9 @@ namespace ParaMEDMEM
                  }
                self->getCellsContainingPoints(da2->getConstPointer(),size,eps,elts,eltsIndex);
              }
-           MEDCouplingAutoRefCountObjectPtr<DataArrayInt> d0=DataArrayInt::New();
-           MEDCouplingAutoRefCountObjectPtr<DataArrayInt> d1=DataArrayInt::New();
-           d0->alloc(elts.size(),1);
-           d1->alloc(eltsIndex.size(),1);
-           std::copy(elts.begin(),elts.end(),d0->getPointer());
-           std::copy(eltsIndex.begin(),eltsIndex.end(),d1->getPointer());
            PyObject *ret=PyTuple_New(2);
-           PyTuple_SetItem(ret,0,SWIG_NewPointerObj(SWIG_as_voidptr(d0.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
-           PyTuple_SetItem(ret,1,SWIG_NewPointerObj(SWIG_as_voidptr(d1.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
+           PyTuple_SetItem(ret,0,SWIG_NewPointerObj(SWIG_as_voidptr(elts.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
+           PyTuple_SetItem(ret,1,SWIG_NewPointerObj(SWIG_as_voidptr(eltsIndex.retn()),SWIGTYPE_p_ParaMEDMEM__DataArrayInt, SWIG_POINTER_OWN | 0 ));
            return ret;
          }