]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
Intersect2DMeshWith1DLine : Bug correction concerning cells in mesh1D colinear to...
authorabn <adrien.bruneton@cea.fr>
Thu, 11 Dec 2014 13:24:37 +0000 (14:24 +0100)
committervsr <vsr@opencascade.com>
Thu, 18 Dec 2014 14:24:05 +0000 (17:24 +0300)
src/INTERP_KERNEL/CellModel.cxx
src/INTERP_KERNEL/CellModel.hxx
src/MEDCoupling/MEDCouplingMemArray.cxx
src/MEDCoupling/MEDCouplingMemArray.hxx
src/MEDCoupling/MEDCouplingUMesh.cxx
src/MEDCoupling/MEDCouplingUMesh.hxx
src/MEDCoupling_Swig/MEDCouplingBasicsTest.py
src/MEDCoupling_Swig/MEDCouplingCommon.i
src/MEDCoupling_Swig/MEDCouplingMemArray.i

index 5a22ab1c409f83e1959165ea5d244775863a6f19..2bdfb826bb2a2c47f2edbec4f29c516605262778 100644 (file)
@@ -591,6 +591,52 @@ namespace INTERP_KERNEL
       throw INTERP_KERNEL::Exception("CellModel::fillSonEdgesNodalConnectivity3D : not implemented yet for NORM_POLYHED !");   
   }
 
+  void CellModel::changeOrientationOf2D(int *nodalConn, unsigned int sz) const
+  {
+    if(sz<1)
+      return ;
+    if(!isQuadratic())
+      {
+        std::vector<int> tmp(sz-1);
+        std::copy(nodalConn+1,nodalConn+sz,tmp.rbegin());
+        std::copy(tmp.begin(),tmp.end(),nodalConn+1);
+      }
+    else
+      {
+        unsigned int sz2(sz/2);
+        std::vector<int> tmp0(sz2-1),tmp1(sz2);
+        std::copy(nodalConn+1,nodalConn+sz2,tmp0.rbegin());
+        std::copy(nodalConn+sz2,nodalConn+sz,tmp1.rbegin());
+        std::copy(tmp0.begin(),tmp0.end(),nodalConn+1);
+        std::copy(tmp1.begin(),tmp1.end(),nodalConn+sz2);
+      }
+  }
+
+  void CellModel::changeOrientationOf1D(int *nodalConn, unsigned int sz) const
+  {
+    if(!isDynamic())
+      {
+        if(sz==2 || sz==3)
+          {
+            std::swap(nodalConn[0],nodalConn[1]);
+            return ;
+          }
+        else if(sz==4)
+          {
+            std::swap(nodalConn[0],nodalConn[1]);
+            std::swap(nodalConn[2],nodalConn[3]);
+          }
+        else
+          throw Exception("CellModel::changeOrientationOf1D : unrecognized 1D cell type !");
+      }
+    else
+      {
+        std::vector<int> tmp(sz-1);
+        std::copy(nodalConn+1,nodalConn+sz,tmp.rbegin());
+        std::copy(tmp.begin(),tmp.end(),nodalConn+1);
+      }
+  }
+
   //================================================================================
   /*!
    * \brief Return number of nodes in sonId-th son of a Dynamic() cell
index 32974d57cd6dba21b04b4413ed7fd531080802a3..ab49dcf9978eb1f344d204479d8f6fd9818a0c70 100644 (file)
@@ -72,6 +72,8 @@ namespace INTERP_KERNEL
     INTERPKERNEL_EXPORT unsigned fillSonCellNodalConnectivity2(int sonId, const int *nodalConn, int lgth, int *sonNodalConn, NormalizedCellType& typeOfSon) const;
     INTERPKERNEL_EXPORT unsigned fillSonCellNodalConnectivity4(int sonId, const int *nodalConn, int lgth, int *sonNodalConn, NormalizedCellType& typeOfSon) const;
     INTERPKERNEL_EXPORT unsigned fillSonEdgesNodalConnectivity3D(int sonId, const int *nodalConn, int lgth, int *sonNodalConn, NormalizedCellType& typeOfSon) const;
+    INTERPKERNEL_EXPORT void changeOrientationOf2D(int *nodalConn, unsigned int sz) const;
+    INTERPKERNEL_EXPORT void changeOrientationOf1D(int *nodalConn, unsigned int sz) const;
   private:
     bool _dyn;
     bool _quadratic;
index d3f460529a0a2edb09ea81fa3377fd89b5b2c3c6..e9058f908f1dfc4f688794108a3b0cfdbffe8fa6 100644 (file)
@@ -6345,7 +6345,7 @@ void DataArrayInt::reprQuickOverviewData(std::ostream& stream, std::size_t maxNb
 }
 
 /*!
- * Modifies \a this one-dimensional array so that each value \a v = \a indArrBg[ \a v ],
+ * Modifies in place \a this one-dimensional array so that each value \a v = \a indArrBg[ \a v ],
  * i.e. a current value is used as in index to get a new value from \a indArrBg.
  *  \param [in] indArrBg - pointer to the first element of array of new values to assign
  *         to \a this array.
@@ -6354,15 +6354,15 @@ void DataArrayInt::reprQuickOverviewData(std::ostream& stream, std::size_t maxNb
  *  \throw If \a this->getNumberOfComponents() != 1
  *  \throw If any value of \a this can't be used as a valid index for 
  *         [\a indArrBg, \a indArrEnd).
+ *
+ *  \sa replaceOneValByInThis
  */
 void DataArrayInt::transformWithIndArr(const int *indArrBg, const int *indArrEnd)
 {
   checkAllocated();
   if(getNumberOfComponents()!=1)
     throw INTERP_KERNEL::Exception("Call transformWithIndArr method on DataArrayInt with only one component, you can call 'rearrange' method before !");
-  int nbElemsIn=(int)std::distance(indArrBg,indArrEnd);
-  int nbOfTuples=getNumberOfTuples();
-  int *pt=getPointer();
+  int nbElemsIn((int)std::distance(indArrBg,indArrEnd)),nbOfTuples(getNumberOfTuples()),*pt(getPointer());
   for(int i=0;i<nbOfTuples;i++,pt++)
     {
       if(*pt>=0 && *pt<nbElemsIn)
@@ -6376,6 +6376,29 @@ void DataArrayInt::transformWithIndArr(const int *indArrBg, const int *indArrEnd
   declareAsNew();
 }
 
+/*!
+ * Modifies in place \a this one-dimensional array like this : each id in \a this so that this[id] equal to \a valToBeReplaced will be replaced at the same place by \a replacedBy.
+ *
+ * \param [in] valToBeReplaced - the value in \a this to be replaced.
+ * \param [in] replacedBy - the value taken by each tuple previously equal to \a valToBeReplaced.
+ *
+ * \sa DataArrayInt::transformWithIndArr
+ */
+void DataArrayInt::replaceOneValByInThis(int valToBeReplaced, int replacedBy)
+{
+  checkAllocated();
+  if(getNumberOfComponents()!=1)
+    throw INTERP_KERNEL::Exception("Call replaceOneValByInThis method on DataArrayInt with only one component, you can call 'rearrange' method before !");
+  if(valToBeReplaced==replacedBy)
+    return ;
+  int nbOfTuples(getNumberOfTuples()),*pt(getPointer());
+  for(int i=0;i<nbOfTuples;i++,pt++)
+    {
+      if(*pt==valToBeReplaced)
+        *pt=replacedBy;
+    }
+}
+
 /*!
  * Computes distribution of values of \a this one-dimensional array between given value
  * ranges (casts). This method is typically useful for entity number spliting by types,
@@ -9249,6 +9272,32 @@ int DataArrayInt::getMinValueInArray() const
   return *loc;
 }
 
+/*!
+ * Returns in a single walk in \a this the min value and the max value in \a this.
+ * \a this is expected to be single component array.
+ *
+ * \param [out] minValue - the min value in \a this.
+ * \param [out] maxValue - the max value in \a this.
+ *
+ * \sa getMinValueInArray, getMinValue, getMaxValueInArray, getMaxValue
+ */
+void DataArrayInt::getMinMaxValues(int& minValue, int& maxValue) const
+{
+  checkAllocated();
+  if(getNumberOfComponents()!=1)
+    throw INTERP_KERNEL::Exception("DataArrayInt::getMinMaxValues : must be applied on DataArrayInt with only one component !");
+  int nbTuples(getNumberOfTuples());
+  const int *pt(begin());
+  minValue=std::numeric_limits<int>::max(); maxValue=-std::numeric_limits<int>::max();
+  for(int i=0;i<nbTuples;i++,pt++)
+    {
+      if(*pt<minValue)
+        minValue=*pt;
+      if(*pt>maxValue)
+        maxValue=*pt;
+    }
+}
+
 /*!
  * Converts every value of \a this array to its absolute value.
  * \b WARNING this method is non const. If a new DataArrayInt instance should be built containing the result of abs DataArrayInt::computeAbs
@@ -9422,7 +9471,7 @@ void DataArrayInt::applyModulus(int val)
  * \param [in] vmax end of range. This value is \b not included in range (excluded).
  * \return a newly allocated data array that the caller should deal with.
  *
- * \sa DataArrayInt::getIdsNotInRange
+ * \sa DataArrayInt::getIdsNotInRange , DataArrayInt::getIdsStrictlyNegative
  */
 DataArrayInt *DataArrayInt::getIdsInRange(int vmin, int vmax) const
 {
@@ -9447,7 +9496,7 @@ DataArrayInt *DataArrayInt::getIdsInRange(int vmin, int vmax) const
  * \param [in] vmax end of range. This value is included in range (included).
  * \return a newly allocated data array that the caller should deal with.
  * 
- * \sa DataArrayInt::getIdsInRange
+ * \sa DataArrayInt::getIdsInRange , DataArrayInt::getIdsStrictlyNegative
  */
 DataArrayInt *DataArrayInt::getIdsNotInRange(int vmin, int vmax) const
 {
@@ -9463,6 +9512,26 @@ DataArrayInt *DataArrayInt::getIdsNotInRange(int vmin, int vmax) const
   return ret.retn();
 }
 
+/*!
+ * This method works only on data array with one component. This method returns a newly allocated array storing stored ascendantly of tuple ids in \a this so that this[id]<0.
+ *
+ * \return a newly allocated data array that the caller should deal with.
+ * \sa DataArrayInt::getIdsInRange
+ */
+DataArrayInt *DataArrayInt::getIdsStrictlyNegative() const
+{
+  checkAllocated();
+  if(getNumberOfComponents()!=1)
+    throw INTERP_KERNEL::Exception("DataArrayInt::getIdsStrictlyNegative : this must have exactly one component !");
+  const int *cptr(getConstPointer());
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret(DataArrayInt::New()); ret->alloc(0,1);
+  int nbOfTuples(getNumberOfTuples());
+  for(int i=0;i<nbOfTuples;i++,cptr++)
+    if(*cptr<0)
+      ret->pushBackSilent(i);
+  return ret.retn();
+}
+
 /*!
  * This method works only on data array with one component.
  * This method checks that all ids in \b this are in [ \b vmin, \b vmax ). If there is at least one element in \a this not in [ \b vmin, \b vmax ) an exception will be thrown.
@@ -10039,6 +10108,7 @@ DataArrayInt *DataArrayInt::buildIntersection(const DataArrayInt *other) const
  * 
  * \return a newly allocated array that contain the result of the unique operation applied on \a this.
  * \throw if \a this is not allocated or if \a this has not exactly one component.
+ * \sa DataArrayInt::buildUniqueNotSorted
  */
 DataArrayInt *DataArrayInt::buildUnique() const
 {
@@ -10055,6 +10125,38 @@ DataArrayInt *DataArrayInt::buildUnique() const
   return ret.retn();
 }
 
+/*!
+ * This method can be applied on allocated with one component DataArrayInt instance.
+ * This method keep elements only once by keeping the same order in \a this that is not expected to be sorted.
+ *
+ * \return a newly allocated array that contain the result of the unique operation applied on \a this.
+ *
+ * \throw if \a this is not allocated or if \a this has not exactly one component.
+ *
+ * \sa DataArrayInt::buildUnique
+ */
+DataArrayInt *DataArrayInt::buildUniqueNotSorted() const
+{
+  checkAllocated();
+    if(getNumberOfComponents()!=1)
+      throw INTERP_KERNEL::Exception("DataArrayInt::buildUniqueNotSorted : only single component allowed !");
+  int minVal,maxVal;
+  getMinMaxValues(minVal,maxVal);
+  std::vector<bool> b(maxVal-minVal+1,false);
+  const int *ptBg(begin()),*endBg(end());
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret(DataArrayInt::New()); ret->alloc(0,1);
+  for(const int *pt=ptBg;pt!=endBg;pt++)
+    {
+      if(!b[*pt-minVal])
+        {
+          ret->pushBackSilent(*pt);
+          b[*pt-minVal]=true;
+        }
+    }
+  ret->copyStringInfoFrom(*this);
+  return ret.retn();
+}
+
 /*!
  * Returns a new DataArrayInt which contains size of every of groups described by \a this
  * "index" array. Such "index" array is returned for example by 
index 8c25de2fb6e94cf3c859a79924322ff7aab17c0a..b47a14ebc95741ab39134bafd08986c317535be0 100644 (file)
@@ -476,6 +476,7 @@ namespace ParaMEDMEM
     MEDCOUPLING_EXPORT void reprQuickOverview(std::ostream& stream) const;
     MEDCOUPLING_EXPORT void reprQuickOverviewData(std::ostream& stream, std::size_t maxNbOfByteInRepr) const;
     MEDCOUPLING_EXPORT void transformWithIndArr(const int *indArrBg, const int *indArrEnd);
+    MEDCOUPLING_EXPORT void replaceOneValByInThis(int valToBeReplaced, int replacedBy);
     MEDCOUPLING_EXPORT DataArrayInt *transformWithIndArrR(const int *indArrBg, const int *indArrEnd) const;
     MEDCOUPLING_EXPORT void splitByValueRange(const int *arrBg, const int *arrEnd,
                                               DataArrayInt *& castArr, DataArrayInt *& rankInsideCast, DataArrayInt *& castsPresent) const;
@@ -555,6 +556,7 @@ namespace ParaMEDMEM
     MEDCOUPLING_EXPORT int getMaxValueInArray() const;
     MEDCOUPLING_EXPORT int getMinValue(int& tupleId) const;
     MEDCOUPLING_EXPORT int getMinValueInArray() const;
+    MEDCOUPLING_EXPORT void getMinMaxValues(int& minValue, int& maxValue) const;
     MEDCOUPLING_EXPORT void abs();
     MEDCOUPLING_EXPORT DataArrayInt *computeAbs() const;
     MEDCOUPLING_EXPORT void applyLin(int a, int b, int compoId);
@@ -568,6 +570,7 @@ namespace ParaMEDMEM
     MEDCOUPLING_EXPORT void applyRPow(int val);
     MEDCOUPLING_EXPORT DataArrayInt *getIdsInRange(int vmin, int vmax) const;
     MEDCOUPLING_EXPORT DataArrayInt *getIdsNotInRange(int vmin, int vmax) const;
+    MEDCOUPLING_EXPORT DataArrayInt *getIdsStrictlyNegative() const;
     MEDCOUPLING_EXPORT bool checkAllIdsInRange(int vmin, int vmax) const;
     MEDCOUPLING_EXPORT static DataArrayInt *Aggregate(const DataArrayInt *a1, const DataArrayInt *a2, int offsetA2);
     MEDCOUPLING_EXPORT static DataArrayInt *Aggregate(const std::vector<const DataArrayInt *>& arr);
@@ -586,6 +589,7 @@ namespace ParaMEDMEM
     MEDCOUPLING_EXPORT DataArrayInt *buildUnion(const DataArrayInt *other) const;
     MEDCOUPLING_EXPORT DataArrayInt *buildIntersection(const DataArrayInt *other) const;
     MEDCOUPLING_EXPORT DataArrayInt *buildUnique() const;
+    MEDCOUPLING_EXPORT DataArrayInt *buildUniqueNotSorted() const;
     MEDCOUPLING_EXPORT DataArrayInt *deltaShiftIndex() const;
     MEDCOUPLING_EXPORT void computeOffsets();
     MEDCOUPLING_EXPORT void computeOffsets2();
index f9c16dcbfd8fa186c86f9b10e6dfb830d2b122ca..04636183b67e6759a106f37bd3dcc26753a86e77 100644 (file)
@@ -4243,11 +4243,6 @@ namespace ParaMEDMEM
     // end
   };
 
-
-
-  /*!
-   * Warning the nodes in \a m should be decrRefed ! To avoid that Node * pointer be replaced by another instance.
-   */
   INTERP_KERNEL::Edge *MEDCouplingUMeshBuildQPFromEdge2(INTERP_KERNEL::NormalizedCellType typ, const int *bg, const double *coords2D, std::map< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Node>,int>& m)
   {
     INTERP_KERNEL::Edge *ret(0);
@@ -5978,40 +5973,28 @@ void MEDCouplingUMesh::are2DCellsNotCorrectlyOriented(const double *vec, bool po
  *  \ref cpp_mcumesh_are2DCellsNotCorrectlyOriented "Here is a C++ example".<br>
  *  \ref  py_mcumesh_are2DCellsNotCorrectlyOriented "Here is a Python example".
  *  \endif
+ *
+ *  \sa changeOrientationOfCells
  */
 void MEDCouplingUMesh::orientCorrectly2DCells(const double *vec, bool polyOnly)
 {
   if(getMeshDimension()!=2 || getSpaceDimension()!=3)
     throw INTERP_KERNEL::Exception("Invalid mesh to apply orientCorrectly2DCells on it : must be meshDim==2 and spaceDim==3 !");
-  int nbOfCells=getNumberOfCells();
-  int *conn=_nodal_connec->getPointer();
-  const int *connI=_nodal_connec_index->getConstPointer();
-  const double *coordsPtr=_coords->getConstPointer();
-  bool isModified=false;
+  int nbOfCells(getNumberOfCells()),*conn(_nodal_connec->getPointer());
+  const int *connI(_nodal_connec_index->getConstPointer());
+  const double *coordsPtr(_coords->getConstPointer());
+  bool isModified(false);
   for(int i=0;i<nbOfCells;i++)
     {
-      INTERP_KERNEL::NormalizedCellType type=(INTERP_KERNEL::NormalizedCellType)conn[connI[i]];
+      INTERP_KERNEL::NormalizedCellType type((INTERP_KERNEL::NormalizedCellType)conn[connI[i]]);
       if(!polyOnly || (type==INTERP_KERNEL::NORM_POLYGON || type==INTERP_KERNEL::NORM_QPOLYG))
         {
-          bool isQuadratic(INTERP_KERNEL::CellModel::GetCellModel(type).isQuadratic());
+          const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(type));
+          bool isQuadratic(cm.isQuadratic());
           if(!IsPolygonWellOriented(isQuadratic,vec,conn+connI[i]+1,conn+connI[i+1],coordsPtr))
             {
               isModified=true;
-              if(!isQuadratic)
-                {
-                  std::vector<int> tmp(connI[i+1]-connI[i]-2);
-                  std::copy(conn+connI[i]+2,conn+connI[i+1],tmp.rbegin());
-                  std::copy(tmp.begin(),tmp.end(),conn+connI[i]+2);
-                }
-              else
-                {
-                  int sz(((int)(connI[i+1]-connI[i]-1))/2);
-                  std::vector<int> tmp0(sz-1),tmp1(sz);
-                  std::copy(conn+connI[i]+2,conn+connI[i]+1+sz,tmp0.rbegin());
-                  std::copy(conn+connI[i]+1+sz,conn+connI[i+1],tmp1.rbegin());
-                  std::copy(tmp0.begin(),tmp0.end(),conn+connI[i]+2);
-                  std::copy(tmp1.begin(),tmp1.end(),conn+connI[i]+1+sz);
-                }
+              cm.changeOrientationOf2D(conn+connI[i]+1,(unsigned int)(connI[i+1]-connI[i]-1));
             }
         }
     }
@@ -6020,6 +6003,38 @@ void MEDCouplingUMesh::orientCorrectly2DCells(const double *vec, bool polyOnly)
   updateTime();
 }
 
+/*!
+ * This method change the orientation of cells in \a this without any consideration of coordinates. Only connectivity is impacted.
+ *
+ * \sa orientCorrectly2DCells
+ */
+void MEDCouplingUMesh::changeOrientationOfCells()
+{
+  int mdim(getMeshDimension());
+  if(mdim!=2 && mdim!=1)
+    throw INTERP_KERNEL::Exception("Invalid mesh to apply changeOrientationOfCells on it : must be meshDim==2 or meshDim==1 !");
+  int nbOfCells(getNumberOfCells()),*conn(_nodal_connec->getPointer());
+  const int *connI(_nodal_connec_index->getConstPointer());
+  if(mdim==2)
+    {//2D
+      for(int i=0;i<nbOfCells;i++)
+        {
+          INTERP_KERNEL::NormalizedCellType type((INTERP_KERNEL::NormalizedCellType)conn[connI[i]]);
+          const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(type));
+          cm.changeOrientationOf2D(conn+connI[i]+1,(unsigned int)(connI[i+1]-connI[i]-1));
+        }
+    }
+  else
+    {//1D
+      for(int i=0;i<nbOfCells;i++)
+        {
+          INTERP_KERNEL::NormalizedCellType type((INTERP_KERNEL::NormalizedCellType)conn[connI[i]]);
+          const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(type));
+          cm.changeOrientationOf1D(conn+connI[i]+1,(unsigned int)(connI[i+1]-connI[i]-1));
+        }
+    }
+}
+
 /*!
  * Finds incorrectly oriented polyhedral cells, i.e. polyhedrons having correctly
  * oriented facets. The normal vector of the facet should point out of the cell.
@@ -8954,7 +8969,7 @@ MEDCouplingUMesh *BuildMesh1DCutFrom(const MEDCouplingUMesh *mesh1D, const std::
   return ret.retn();
 }
 
-MEDCouplingUMesh *BuildRefined2DCell(const DataArrayDouble *coords, const int *descBg, const int *descEnd, const std::vector< std::vector<int> >& intersectEdge1)
+MEDCouplingUMesh *BuildRefined2DCellLinear(const DataArrayDouble *coords, const int *descBg, const int *descEnd, const std::vector< std::vector<int> >& intersectEdge1)
 {
   std::vector<int> allEdges;
   for(const int *it2(descBg);it2!=descEnd;it2++)
@@ -8967,7 +8982,7 @@ MEDCouplingUMesh *BuildRefined2DCell(const DataArrayDouble *coords, const int *d
     }
   std::size_t nb(allEdges.size());
   if(nb%2!=0)
-    throw INTERP_KERNEL::Exception("BuildRefined2DCell : internal error 1 !");
+    throw INTERP_KERNEL::Exception("BuildRefined2DCellLinear : internal error 1 !");
   std::size_t nbOfEdgesOf2DCellSplit(nb/2);
   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> ret(MEDCouplingUMesh::New("",2));
   ret->setCoords(coords);
@@ -8979,6 +8994,84 @@ MEDCouplingUMesh *BuildRefined2DCell(const DataArrayDouble *coords, const int *d
   return ret.retn();
 }
 
+MEDCouplingUMesh *BuildRefined2DCellQuadratic(const DataArrayDouble *coords, const MEDCouplingUMesh *mesh2D, int cellIdInMesh2D, const int *descBg, const int *descEnd, const std::vector< std::vector<int> >& intersectEdge1)
+{
+  const int *c(mesh2D->getNodalConnectivity()->begin()),*ci(mesh2D->getNodalConnectivityIndex()->begin());
+  const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)c[ci[cellIdInMesh2D]]));
+  std::size_t ii(0);
+  unsigned sz(cm.getNumberOfSons2(c+ci[cellIdInMesh2D]+1,ci[cellIdInMesh2D+1]-ci[cellIdInMesh2D]-1));
+  if(sz!=std::distance(descBg,descEnd))
+    throw INTERP_KERNEL::Exception("BuildRefined2DCellQuadratic : internal error 1 !");
+  INTERP_KERNEL::AutoPtr<int> tmpPtr(new int[ci[cellIdInMesh2D+1]-ci[cellIdInMesh2D]]);
+  std::vector<int> allEdges,centers;
+  const double *coordsPtr(coords->begin());
+  MEDCouplingAutoRefCountObjectPtr<DataArrayDouble> addCoo(DataArrayDouble::New()); addCoo->alloc(0,1);
+  int offset(coords->getNumberOfTuples());
+  for(const int *it2(descBg);it2!=descEnd;it2++,ii++)
+    {
+      INTERP_KERNEL::NormalizedCellType typeOfSon;
+      cm.fillSonCellNodalConnectivity2(ii,c+ci[cellIdInMesh2D]+1,ci[cellIdInMesh2D+1]-ci[cellIdInMesh2D]-1,tmpPtr,typeOfSon);
+      const std::vector<int>& edge1(intersectEdge1[std::abs(*it2)-1]);
+      if(*it2>0)
+        allEdges.insert(allEdges.end(),edge1.begin(),edge1.end());
+      else
+        allEdges.insert(allEdges.end(),edge1.rbegin(),edge1.rend());
+      if(edge1.size()==2)
+        centers.push_back(tmpPtr[2]);//special case where no subsplit of edge -> reuse the original center.
+      else
+        {//the current edge has been subsplit -> create corresponding centers.
+          std::size_t nbOfCentersToAppend(edge1.size()/2);
+          std::map< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Node>,int> m;
+          MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> ee(MEDCouplingUMeshBuildQPFromEdge2(typeOfSon,tmpPtr,coordsPtr,m));
+          std::vector<int>::const_iterator it3(allEdges.end()-edge1.size());
+          for(std::size_t k=0;k<nbOfCentersToAppend;k++)
+            {
+              double tmpp[2];
+              const double *aa(coordsPtr+2*(*it3++));
+              const double *bb(coordsPtr+2*(*it3++));
+              ee->getMiddleOfPoints(aa,bb,tmpp);
+              addCoo->insertAtTheEnd(tmpp,tmpp+2);
+              centers.push_back(offset+k);
+            }
+        }
+    }
+  std::size_t nb(allEdges.size());
+  if(nb%2!=0)
+    throw INTERP_KERNEL::Exception("BuildRefined2DCellQuadratic : internal error 2 !");
+  std::size_t nbOfEdgesOf2DCellSplit(nb/2);
+  MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> ret(MEDCouplingUMesh::New("",2));
+  if(addCoo->empty())
+    ret->setCoords(coords);
+  else
+    {
+      addCoo->rearrange(2);
+      addCoo=DataArrayDouble::Aggregate(coords,addCoo);
+      ret->setCoords(addCoo);
+    }
+  ret->allocateCells(1);
+  std::vector<int> connOut(nbOfEdgesOf2DCellSplit);
+  for(std::size_t kk=0;kk<nbOfEdgesOf2DCellSplit;kk++)
+    connOut[kk]=allEdges[2*kk];
+  connOut.insert(connOut.end(),centers.begin(),centers.end());
+  ret->insertNextCell(INTERP_KERNEL::NORM_QPOLYG,connOut.size(),&connOut[0]);
+  return ret.retn();
+}
+
+/*!
+ * This method creates a refinement of a cell in \a mesh2D. Those cell is defined by descending connectivity and the sorted subdivided nodal connectivity
+ * of those edges.
+ *
+ * \param [in] mesh2D - The origin 2D mesh. \b Warning \b coords are not those of \a mesh2D. But mesh2D->getCoords()==coords[:mesh2D->getNumberOfNodes()]
+ */
+MEDCouplingUMesh *BuildRefined2DCell(const DataArrayDouble *coords, const MEDCouplingUMesh *mesh2D, int cellIdInMesh2D, const int *descBg, const int *descEnd, const std::vector< std::vector<int> >& intersectEdge1)
+{
+  const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(mesh2D->getTypeOfCell(cellIdInMesh2D)));
+  if(!cm.isQuadratic())
+    return BuildRefined2DCellLinear(coords,descBg,descEnd,intersectEdge1);
+  else
+    return BuildRefined2DCellQuadratic(coords,mesh2D,cellIdInMesh2D,descBg,descEnd,intersectEdge1);
+}
+
 void AddCellInMesh2D(MEDCouplingUMesh *mesh2D, const std::vector<int>& conn, const std::vector< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> >& edges)
 {
   bool isQuad(false);
@@ -9000,7 +9093,7 @@ void AddCellInMesh2D(MEDCouplingUMesh *mesh2D, const std::vector<int>& conn, con
       for(std::size_t i=0;i<sz;i++)
         {
           double tmp[2];
-          edges[(i+1)%sz]->getMiddleOfPoints(coo+2*conn[i],coo+2*conn[(i+1)%sz],tmp);
+          edges[(i+1)%sz]->getMiddleOfPoints(coo+2*conn[i],coo+2*conn[(i+1)%sz],tmp);// tony a chier i+1 -> i
           addCoo.insert(addCoo.end(),tmp,tmp+2);
           conn2.push_back(offset+(int)i);
         }
@@ -9013,6 +9106,9 @@ void AddCellInMesh2D(MEDCouplingUMesh *mesh2D, const std::vector<int>& conn, con
 
 /*!
  * \b WARNING edges in out1 coming from \a splitMesh1D are \b NOT oriented because only used for equation of curve.
+ *
+ * This method cuts in 2 parts the input 2D cell given using boundaries description (\a edge1Bis and \a edge1BisPtr) using
+ * a set of edges defined in \a splitMesh1D.
  */
 void BuildMesh2DCutInternal2(const MEDCouplingUMesh *splitMesh1D, const std::vector<int>& edge1Bis, const std::vector< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> >& edge1BisPtr,
                              std::vector< std::vector<int> >& out0, std::vector< std::vector< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> > >& out1)
@@ -9051,20 +9147,20 @@ void BuildMesh2DCutInternal2(const MEDCouplingUMesh *splitMesh1D, const std::vec
       for(std::size_t k=ii;k<jj+1;k++)
         { connOutLeft.push_back(edge1Bis[2*k+1]); eleft.push_back(edge1BisPtr[2*k+1]); }
       std::vector< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> > ees(iEnd);
-      for(int ik=iEnd-1;ik>=0;ik--)
+      for(int ik=0;ik<iEnd;ik++)
         {
           std::map< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Node>,int> m;
           MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> ee(MEDCouplingUMeshBuildQPFromEdge2((INTERP_KERNEL::NormalizedCellType)cSplitPtr[ciSplitPtr[ik]],cSplitPtr+ciSplitPtr[ik]+1,splitMesh1D->getCoords()->begin(),m));
-          ees[iEnd-1-ik]=ee;
+          ees[ik]=ee;
         }
       for(int ik=iEnd-1;ik>=0;ik--)
         connOutLeft.push_back(cSplitPtr[ciSplitPtr[ik]+1]);
       for(std::size_t k=jj+1;k<nbOfEdgesOf2DCellSplit+ii;k++)
         { connOutRight.push_back(edge1Bis[2*k+1]); eright.push_back(edge1BisPtr[2*k+1]); }
-      eleft.insert(eleft.end(),ees.begin(),ees.end());
+      eleft.insert(eleft.end(),ees.rbegin(),ees.rend());
       for(int ik=0;ik<iEnd;ik++)
         connOutRight.push_back(cSplitPtr[ciSplitPtr[ik]+2]);
-      eright.insert(eright.end(),ees.rbegin(),ees.rend());
+      eright.insert(eright.end(),ees.begin(),ees.end());
     }
 }
 
@@ -9089,7 +9185,7 @@ CellInfo::CellInfo(const std::vector<int>& edges, const std::vector< MEDCoupling
   for(std::size_t i=0;i<nbe;i++)
     {
       edges2[2*i]=edges[i]; edges2[2*i+1]=edges[(i+1)%nbe];
-      edgesPtr2[2*i]=edgesPtr[i]; edgesPtr2[2*i+1]=edgesPtr[i];
+      edgesPtr2[2*i]=edgesPtr[(i+1)%nbe]; edgesPtr2[2*i+1]=edgesPtr[(i+1)%nbe];//tony a chier
     }
   _edges.resize(4*nbe); _edges_ptr.resize(4*nbe);
   std::copy(edges2.begin(),edges2.end(),_edges.begin()); std::copy(edges2.begin(),edges2.end(),_edges.begin()+2*nbe);
@@ -9237,7 +9333,7 @@ void VectorOfCellInfo::setMeshAt(int pos, const MEDCouplingAutoRefCountObjectPtr
   for(std::size_t j=0;j<sz;j++)
     pool[pos+j]=CellInfo(edges[j],edgePtrs[j]);
   for(int i=pos+1;i<(int)_pool.size();i++)
-    pool[pos+sz-1]=_pool[i];
+    pool[i+sz-1]=_pool[i];
   _pool=pool;
   //
   if(sz==2)
@@ -9309,6 +9405,17 @@ CellInfo& VectorOfCellInfo::get(int pos)
   return _pool[pos];
 }
 
+/*!
+ * Given :
+ * - a \b closed set of edges ( \a allEdges and \a allEdgesPtr ) that defines the split descending 2D cell.
+ * - \a splitMesh1D a split 2D curve mesh contained into 2D cell defined above.
+ *
+ * This method returns the 2D mesh and feeds \a idsLeftRight using offset.
+ *
+ * Algorithm : \a splitMesh1D is cut into contiguous parts. Each contiguous parts will build incrementally the output 2D cells.
+ *
+ * \param [in] allEdges a list of pairs (beginNode, endNode). Linked with \a allEdgesPtr to get the equation of edge.
+ */
 MEDCouplingUMesh *BuildMesh2DCutInternal(double eps, const MEDCouplingUMesh *splitMesh1D, const std::vector<int>& allEdges, const std::vector< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> >& allEdgesPtr, int offset,
                                          MEDCouplingAutoRefCountObjectPtr<DataArrayInt>& idsLeftRight)
 {
@@ -9372,8 +9479,8 @@ MEDCouplingUMesh *BuildMesh2DCutFrom(double eps, int cellIdInMesh2D, const MEDCo
   const int *cdescPtr(mesh2DDesc->getNodalConnectivity()->begin()),*cidescPtr(mesh2DDesc->getNodalConnectivityIndex()->begin());
   //
   std::vector<int> allEdges;
-  std::vector< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> > allEdgesPtr;
-  for(const int *it(descBg);it!=descEnd;it++)
+  std::vector< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Edge> > allEdgesPtr; // for each sub edge in splitMesh2D the uncut Edge object of the original mesh2D
+  for(const int *it(descBg);it!=descEnd;it++) // for all edges in the descending connectivity of the 2D mesh in relative Fortran mode
     {
       int edgeId(std::abs(*it)-1);
       std::map< MEDCouplingAutoRefCountObjectPtr<INTERP_KERNEL::Node>,int> m;
@@ -9391,6 +9498,76 @@ MEDCouplingUMesh *BuildMesh2DCutFrom(double eps, int cellIdInMesh2D, const MEDCo
   return BuildMesh2DCutInternal(eps,splitMesh1D,allEdges,allEdgesPtr,offset,idsLeftRight);
 }
 
+bool AreEdgeEqual(const double *coo2D, const INTERP_KERNEL::CellModel& typ1, const int *conn1, const INTERP_KERNEL::CellModel& typ2, const int *conn2, double eps)
+{
+  if(!typ1.isQuadratic() && !typ2.isQuadratic())
+    {//easy case comparison not
+      return conn1[0]==conn2[0] && conn1[1]==conn2[1];
+    }
+  else if(typ1.isQuadratic() && typ2.isQuadratic())
+    {
+      bool status0(conn1[0]==conn2[0] && conn1[1]==conn2[1]);
+      if(!status0)
+        return false;
+      if(conn1[2]==conn2[2])
+        return true;
+      const double *a(coo2D+2*conn1[2]),*b(coo2D+2*conn2[2]);
+      double dist(sqrt((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])));
+      return dist<eps;
+    }
+  else
+    {//only one is quadratic
+      bool status0(conn1[0]==conn2[0] && conn1[1]==conn2[1]);
+      if(!status0)
+        return false;
+      if(typ1.isQuadratic())
+        {
+          const double *a(coo2D+2*conn1[2]),*bb(coo2D+2*conn2[0]),*be(coo2D+2*conn2[1]);
+          double b[2]; b[0]=(be[0]+bb[0])/2.; b[1]=(be[1]+bb[1])/2.;
+          double dist(sqrt((a[0]-b[0])*(a[0]-b[0])+(a[1]-b[1])*(a[1]-b[1])));
+          return dist<eps;
+        }
+    }
+}
+
+/*!
+ * This method returns among the cellIds [ \a candidatesIn2DBg , \a candidatesIn2DEnd ) in \a mesh2DSplit those exactly sharing \a cellIdInMesh1DSplitRelative in \a mesh1DSplit.
+ * \a mesh2DSplit and \a mesh1DSplit are expected to share the coordinates array.
+ *
+ * \param [in] cellIdInMesh1DSplitRelative is in Fortran mode using sign to specify direction.
+ */
+int FindRightCandidateAmong(const MEDCouplingUMesh *mesh2DSplit, const int *candidatesIn2DBg, const int *candidatesIn2DEnd, const MEDCouplingUMesh *mesh1DSplit, int cellIdInMesh1DSplitRelative, double eps)
+{
+  if(candidatesIn2DEnd==candidatesIn2DBg)
+    throw INTERP_KERNEL::Exception("FindRightCandidateAmong : internal error 1 !");
+  const double *coo(mesh2DSplit->getCoords()->begin());
+  if(std::distance(candidatesIn2DBg,candidatesIn2DEnd)==1)
+    return *candidatesIn2DBg;
+  int edgeId(std::abs(cellIdInMesh1DSplitRelative)-1);
+  MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> cur1D(static_cast<MEDCouplingUMesh *>(mesh1DSplit->buildPartOfMySelf(&edgeId,&edgeId+1,true)));
+  if(cellIdInMesh1DSplitRelative<0)
+    cur1D->changeOrientationOfCells();
+  const int *c1D(cur1D->getNodalConnectivity()->begin());
+  const INTERP_KERNEL::CellModel& ref1DType(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)c1D[0]));
+  for(const int *it=candidatesIn2DBg;it!=candidatesIn2DEnd;it++)
+    {
+      MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> cur2D(static_cast<MEDCouplingUMesh *>(mesh2DSplit->buildPartOfMySelf(it,it+1,true)));
+      const int *c(cur2D->getNodalConnectivity()->begin()),*ci(cur2D->getNodalConnectivityIndex()->begin());
+      const INTERP_KERNEL::CellModel &cm(INTERP_KERNEL::CellModel::GetCellModel((INTERP_KERNEL::NormalizedCellType)c[ci[0]]));
+      unsigned sz(cm.getNumberOfSons2(c+ci[0]+1,ci[1]-ci[0]-1));
+      INTERP_KERNEL::AutoPtr<int> tmpPtr(new int[ci[1]-ci[0]]);
+      for(unsigned it2=0;it2<sz;it2++)
+        {
+          INTERP_KERNEL::NormalizedCellType typeOfSon;
+          cm.fillSonCellNodalConnectivity2(it2,c+ci[0]+1,ci[1]-ci[0]-1,tmpPtr,typeOfSon);
+          const INTERP_KERNEL::CellModel &curCM(INTERP_KERNEL::CellModel::GetCellModel(typeOfSon));
+          if(AreEdgeEqual(coo,ref1DType,c1D+1,curCM,tmpPtr,eps))
+            return *it;
+        }
+    }
+  throw INTERP_KERNEL::Exception("FindRightCandidateAmong : internal error 2 ! Unable to find the edge among split cell !");
+}
+
 /// @endcond
 
 /*!
@@ -9460,7 +9637,7 @@ void MEDCouplingUMesh::Intersect2DMeshWith1DLine(const MEDCouplingUMesh *mesh2D,
   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret2(DataArrayInt::New()); ret2->alloc(0,1);
   MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> ret1(BuildMesh1DCutFrom(mesh1D,intersectEdge2,mesh2D->getCoords(),addCoo,mergedNodes,colinear2,intersectEdge1,
       idsInRet1Colinear,idsInDescMesh2DForIdsInRetColinear));
-  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret3(DataArrayInt::New()); ret3->alloc(ret1->getNumberOfCells()*2,1); ret3->fillWithValue(-1); ret3->rearrange(2);
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> ret3(DataArrayInt::New()); ret3->alloc(ret1->getNumberOfCells()*2,1); ret3->fillWithValue(std::numeric_limits<int>::max()); ret3->rearrange(2);
   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> idsInRet1NotColinear(idsInRet1Colinear->buildComplement(ret1->getNumberOfCells()));
   // deal with cells in mesh2D that are not cut but only some of their edges are
   MEDCouplingAutoRefCountObjectPtr<DataArrayInt> idsInDesc2DToBeRefined(idsInDescMesh2DForIdsInRetColinear->deepCpy());
@@ -9503,12 +9680,13 @@ void MEDCouplingUMesh::Intersect2DMeshWith1DLine(const MEDCouplingUMesh *mesh2D,
       const int *rdptr(dd3->begin()),*rdiptr(dd4->begin()),*dptr(dd1->begin()),*diptr(dd2->begin());
       for(const int *it=fewModifiedCells->begin();it!=fewModifiedCells->end();it++)
         {
-          outMesh2DSplit.push_back(BuildRefined2DCell(ret1->getCoords(),dptr+diptr[*it],dptr+diptr[*it+1],intersectEdge1));
+          outMesh2DSplit.push_back(BuildRefined2DCell(ret1->getCoords(),mesh2D,*it,dptr+diptr[*it],dptr+diptr[*it+1],intersectEdge1));
+          ret1->setCoords(outMesh2DSplit.back()->getCoords());
         }
       int offset(ret2->getNumberOfTuples());
       ret2->pushBackValsSilent(fewModifiedCells->begin(),fewModifiedCells->end());
       MEDCouplingAutoRefCountObjectPtr<DataArrayInt> partOfRet3(DataArrayInt::New()); partOfRet3->alloc(2*idsInRet1Colinear->getNumberOfTuples(),1);
-      partOfRet3->fillWithValue(-1); partOfRet3->rearrange(2);
+      partOfRet3->fillWithValue(std::numeric_limits<int>::max()); partOfRet3->rearrange(2);
       int kk(0),*ret3ptr(partOfRet3->getPointer());
       for(const int *it=idsInDescMesh2DForIdsInRetColinear->begin();it!=idsInDescMesh2DForIdsInRetColinear->end();it++,kk++)
         {
@@ -9524,11 +9702,25 @@ void MEDCouplingUMesh::Intersect2DMeshWith1DLine(const MEDCouplingUMesh *mesh2D,
                     ret3ptr[2*kk+1]=tmp+offset;
                 }
               else
-                throw INTERP_KERNEL::Exception("MEDCouplingUMesh::Intersect2DMeshWith1DLine : internal error 1 !");
+                {//the current edge is shared by a 2D cell that will be split just after
+                  if(std::find(dptr+diptr[*it2],dptr+diptr[*it2+1],-(*it))!=dptr+diptr[*it2+1])
+                    ret3ptr[2*kk]=-(*it2+1);
+                  if(std::find(dptr+diptr[*it2],dptr+diptr[*it2+1],(*it))!=dptr+diptr[*it2+1])
+                    ret3ptr[2*kk+1]=-(*it2+1);
+                }
             }
         }
+      m1Desc->setCoords(ret1->getCoords());
+      ret1NonCol->setCoords(ret1->getCoords());
       ret3->setPartOfValues3(partOfRet3,idsInRet1Colinear->begin(),idsInRet1Colinear->end(),0,2,1,true);
+      if(!outMesh2DSplit.empty())
+        {
+          DataArrayDouble *da(outMesh2DSplit.back()->getCoords());
+          for(std::vector< MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> >::iterator itt=outMesh2DSplit.begin();itt!=outMesh2DSplit.end();itt++)
+            (*itt)->setCoords(da);
+        }
     }
+  cellsToBeModified=cellsToBeModified->buildUniqueNotSorted();
   for(const int *it=cellsToBeModified->begin();it!=cellsToBeModified->end();it++)
     {
       MEDCouplingAutoRefCountObjectPtr<DataArrayInt> idsNonColPerCell(elts->getIdsEqual(*it));
@@ -9549,9 +9741,21 @@ void MEDCouplingUMesh::Intersect2DMeshWith1DLine(const MEDCouplingUMesh *mesh2D,
     tmp[i]=outMesh2DSplit[i];
   //
   ret1->getCoords()->setInfoOnComponents(compNames);
+  MEDCouplingAutoRefCountObjectPtr<MEDCouplingUMesh> ret2D(MEDCouplingUMesh::MergeUMeshesOnSameCoords(tmp));
+  // To finish - filter ret3 - std::numeric_limits<int>::max() -> -1 - negate values must be resolved.
+  ret3->rearrange(1);
+  MEDCouplingAutoRefCountObjectPtr<DataArrayInt> edgesToDealWith(ret3->getIdsStrictlyNegative());
+  for(const int *it=edgesToDealWith->begin();it!=edgesToDealWith->end();it++)
+    {
+      int old2DCellId(-ret3->getIJ(*it,0)-1);
+      MEDCouplingAutoRefCountObjectPtr<DataArrayInt> candidates(ret2->getIdsEqual(old2DCellId));
+      ret3->setIJ(*it,0,FindRightCandidateAmong(ret2D,candidates->begin(),candidates->end(),ret1,*it%2==0?-((*it)/2+1):(*it)/2+1,eps));// div by 2 because 2 components natively in ret3
+    }
+  ret3->replaceOneValByInThis(std::numeric_limits<int>::max(),-1);
+  ret3->rearrange(2);
   //
   splitMesh1D=ret1.retn();
-  splitMesh2D=MEDCouplingUMesh::MergeUMeshesOnSameCoords(tmp);
+  splitMesh2D=ret2D.retn();
   cellIdInMesh2D=ret2.retn();
   cellIdInMesh1D=ret3.retn();
 }
index e1b2c93d742da6149c3d72066243b5af65ebf97f..202f0f33481847d61247a0faa029929bf32c20bc 100644 (file)
@@ -185,6 +185,7 @@ namespace ParaMEDMEM
     MEDCOUPLING_EXPORT void convertDegeneratedCells();
     MEDCOUPLING_EXPORT void are2DCellsNotCorrectlyOriented(const double *vec, bool polyOnly, std::vector<int>& cells) const;
     MEDCOUPLING_EXPORT void orientCorrectly2DCells(const double *vec, bool polyOnly);
+    MEDCOUPLING_EXPORT void changeOrientationOfCells();
     MEDCOUPLING_EXPORT void arePolyhedronsNotCorrectlyOriented(std::vector<int>& cells) const;
     MEDCOUPLING_EXPORT void orientCorrectlyPolyhedrons();
     MEDCOUPLING_EXPORT void getFastAveragePlaneOfThis(double *vec, double *pos) const;
index 1178b5846fa33ef7b05759e185a1aae3355a548c..8a5c0e8f263e5e9c1a4183cccfa5371b887028f2 100644 (file)
@@ -15667,24 +15667,24 @@ class MEDCouplingBasicsTest(unittest.TestCase):
         pass
 
     def testSwig2Intersect2DMeshWith1DLine2(self):
-         """A basic test with colinearity between m1 and m2 and the last cell of m2 outside m1."""
-         i=MEDCouplingIMesh("mesh",2,[5,5],[0.,0.],[1.,1.])
-         m1=i.buildUnstructured()
-         m2=MEDCouplingUMesh("mesh",1) ; m2.setCoords(DataArrayDouble([0.5,2.,2.25,2.,2.5,2.,2.75,2.,3.,2.,4.,2.,5.,2.],7,2)) ; m2.allocateCells()
-         for i in xrange(6):
-             m2.insertNextCell(NORM_SEG2,[i,i+1])
-             pass
-         a,b,c,d=MEDCouplingUMesh.Intersect2DMeshWith1DLine(m1,m2,1e-12)
-         self.assertTrue(a.getNodalConnectivity().isEqual(DataArrayInt([4,1,0,5,6,4,2,1,6,7,4,3,2,7,8,4,4,3,8,9,4,16,15,20,21,4,17,16,21,22,4,18,17,22,23,4,19,18,23,24,5,6,5,10,25,11,5,7,6,11,12,5,8,7,12,26,27,28,13,5,9,8,13,14,5,11,25,10,15,16,5,12,11,16,17,5,13,28,27,26,12,17,18,5,14,13,18,19])))
-         self.assertTrue(a.getNodalConnectivityIndex().isEqual(DataArrayInt([0,5,10,15,20,25,30,35,40,46,51,59,64,70,75,83,88])))
-         self.assertTrue(b.getNodalConnectivity().isEqual(DataArrayInt([1,25,11,1,11,12,1,12,26,1,26,27,1,27,28,1,28,13,1,13,14,1,14,31])))
-         self.assertTrue(b.getNodalConnectivityIndex().isEqual(DataArrayInt([0,3,6,9,12,15,18,21,24])))
-         self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
-         self.assertTrue(a.getCoords()[:25].isEqual(m1.getCoords(),1e-12))
-         self.assertTrue(a.getCoords()[25:].isEqualWithoutConsideringStr(m2.getCoords(),1e-12))
-         self.assertTrue(c.isEqual(DataArrayInt([0,1,2,3,12,13,14,15,4,5,6,7,8,9,10,11])))
-         self.assertTrue(d.isEqual(DataArrayInt([(12,8),(13,9),(14,10),(14,10),(14,10),(14,10),(15,11),(-1,-1)])))
-         pass
+        """A basic test with colinearity between m1 and m2 and the last cell of m2 outside m1."""
+        i=MEDCouplingIMesh("mesh",2,[5,5],[0.,0.],[1.,1.])
+        m1=i.buildUnstructured()
+        m2=MEDCouplingUMesh("mesh",1) ; m2.setCoords(DataArrayDouble([0.5,2.,2.25,2.,2.5,2.,2.75,2.,3.,2.,4.,2.,5.,2.],7,2)) ; m2.allocateCells()
+        for i in xrange(6):
+            m2.insertNextCell(NORM_SEG2,[i,i+1])
+            pass
+        a,b,c,d=MEDCouplingUMesh.Intersect2DMeshWith1DLine(m1,m2,1e-12)
+        self.assertTrue(a.getNodalConnectivity().isEqual(DataArrayInt([4,1,0,5,6,4,2,1,6,7,4,3,2,7,8,4,4,3,8,9,4,16,15,20,21,4,17,16,21,22,4,18,17,22,23,4,19,18,23,24,5,6,5,10,25,11,5,7,6,11,12,5,8,7,12,26,27,28,13,5,9,8,13,14,5,11,25,10,15,16,5,12,11,16,17,5,13,28,27,26,12,17,18,5,14,13,18,19])))
+        self.assertTrue(a.getNodalConnectivityIndex().isEqual(DataArrayInt([0,5,10,15,20,25,30,35,40,46,51,59,64,70,75,83,88])))
+        self.assertTrue(b.getNodalConnectivity().isEqual(DataArrayInt([1,25,11,1,11,12,1,12,26,1,26,27,1,27,28,1,28,13,1,13,14,1,14,31])))
+        self.assertTrue(b.getNodalConnectivityIndex().isEqual(DataArrayInt([0,3,6,9,12,15,18,21,24])))
+        self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
+        self.assertTrue(a.getCoords()[:25].isEqual(m1.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[25:].isEqualWithoutConsideringStr(m2.getCoords(),1e-12))
+        self.assertTrue(c.isEqual(DataArrayInt([0,1,2,3,12,13,14,15,4,5,6,7,8,9,10,11])))
+        self.assertTrue(d.isEqual(DataArrayInt([(12,8),(13,9),(14,10),(14,10),(14,10),(14,10),(15,11),(-1,-1)])))
+        pass
 
     def testSwig2Intersect2DMeshWith1DLine3(self):
         """m2 fully included in cell #12. of m1"""
@@ -15864,23 +15864,181 @@ class MEDCouplingBasicsTest(unittest.TestCase):
         m_circ = MEDCouplingDataForTest.buildCircle2(0.0, 0.0, 2.0)
         coords = [0.0,3.0,0.0,-3.0]
         connec = [0,1]
-        m_line = MEDCouplingUMesh.New("seg", 1)  
+        m_line = MEDCouplingUMesh("seg", 1)  
         m_line.allocateCells(1)
         meshCoords = DataArrayDouble.New(coords, len(coords)/2, 2)
         m_line.setCoords(meshCoords)
         m_line.insertNextCell(NORM_SEG2, connec)
-        a, b, c, d = MEDCouplingUMesh.Intersect2DMeshWith1DLine(m_circ, m_line, eps)    
-        self.assertEqual([32, 1, 7, 10, 11, 12, 13, 14, 15, 32, 5, 3, 11, 10, 16, 17, 18, 19], a.getNodalConnectivity().getValues())
-        self.assertEqual([0, 9, 18],   a.getNodalConnectivityIndex().getValues())
-        self.assertEqual([1, 8, 11, 1, 11, 10, 1, 10, 9], b.getNodalConnectivity().getValues())
-        self.assertEqual([0, 3, 6, 9], b.getNodalConnectivityIndex().getValues())
+        a, b, c, d = MEDCouplingUMesh.Intersect2DMeshWith1DLine(m_circ, m_line, eps)
+        self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
+        self.assertTrue(a.getCoords()[:m_circ.getNumberOfNodes()].isEqual(m_circ.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[m_circ.getNumberOfNodes():m_circ.getNumberOfNodes()+m_line.getNumberOfNodes()].isEqual(m_line.getCoords(),1e-12))
+        self.assertTrue(a.getCoords().isEqual(DataArrayDouble([(2.,0.),(1.4142135623730951,1.414213562373095),(0.,2.),(-1.414213562373095,1.4142135623730951),(-2.,0.),(-1.4142135623730954,-1.414213562373095),(0.,-2.),(1.4142135623730947,-1.4142135623730954),(0.,3.),(0.,-3.),(0.,-2.),(0.,2.),(2.,0.),(0.7653668647301797,-1.8477590650225735),(0.,0.),(0.7653668647301797,1.8477590650225735),(-2,0.),(-0.7653668647301795,1.8477590650225735),(0.,0.),(-0.7653668647301795,-1.8477590650225735)]),1e-12))
+        self.assertEqual([32,1,7,10,11,12,13,14,15,32,5,3,11,10,16,17,18,19],a.getNodalConnectivity().getValues())
+        self.assertEqual([0,9,18],  a.getNodalConnectivityIndex().getValues())
+        self.assertEqual([1,8,11,1,11,10,1,10,9],b.getNodalConnectivity().getValues())
+        self.assertEqual([0,3,6,9],b.getNodalConnectivityIndex().getValues())
         self.assertTrue(a.getCoords()[:8].isEqual(m_circ.getCoords(),1e-12))
         self.assertTrue(a.getCoords()[8:10].isEqual(m_line.getCoords(),1e-12))
-        coo_tgt = DataArrayDouble([2.0, 0.0, 1.4142135623730951, 1.414213562373095, 1.2246467991473532e-16, 2.0, -1.414213562373095, 1.4142135623730951, -2.0, 2.4492935982947064e-16, -1.4142135623730954, -1.414213562373095, -3.6739403974420594e-16, -2.0, 1.4142135623730947, -1.4142135623730954, 0.0, 3.0, 0.0, -3.0, 0.0, -2.0, 0.0, 2.0, 2.0, -2.220446049250313e-16, 0.7653668647301797, -1.8477590650225735, 0.0, 0.0, 0.7653668647301797, 1.8477590650225735, -1.9999999999999998, -3.2343398276365944e-16, -0.7653668647301795, 1.8477590650225735, 0.0, 0.0, -0.7653668647301795, -1.8477590650225735])
-        self.assertTrue(a.getCoords().isEqualWithoutConsideringStr(coo_tgt, 1.0e-12))
+        coo_tgt = DataArrayDouble([2.,0.,1.4142135623730951,1.414213562373095,1.2246467991473532e-16,2.,-1.414213562373095,1.4142135623730951,-2.,0.,-1.4142135623730954,-1.414213562373095,-3.6739403974420594e-16,-2.,1.4142135623730947,-1.4142135623730954,0.,3.,0.,-3.,0.,-2.,0.,2.,2.,0.,0.7653668647301797,-1.8477590650225735,0.,0.,0.7653668647301797,1.8477590650225735,-2.,0.,-0.7653668647301795,1.8477590650225735,0.,0.,-0.7653668647301795,-1.8477590650225735])
+        self.assertTrue(a.getCoords().isEqualWithoutConsideringStr(coo_tgt,1.0e-12))
         self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
-        self.assertEqual([0, 0], c.getValues())
-        self.assertEqual([-1, -1, 0, 1, -1, -1], d.getValues())
+        self.assertEqual([0,0],c.getValues())
+        self.assertEqual([-1,-1,0,1,-1,-1],d.getValues())
+
+    def testSwig2Intersect2DMeshWith1DLine11(self):
+        """ Quad line re-entering a square cell """
+        eps = 1.0e-8
+        m = MEDCouplingUMesh("box", 2)
+        m.setCoords(DataArrayDouble([-1., -1., -1., 1., 1., 1., 1., -1.0],4,2))
+        c, cI = [NORM_POLYGON, 0, 1, 2, 3], [0, 5]
+        m.setConnectivity(DataArrayInt(c), DataArrayInt(cI))
+        m.checkCoherency()
+        coords2 = [0., 1.3, -1.3, 0., -0.6, 0.6, 0., -1.3, -0.5, -0.5]
+        connec2, cI2 = [NORM_SEG3, 0, 1, 2, NORM_SEG3, 1, 3, 4], [0,4,8]
+        m_line = MEDCouplingUMesh("seg", 1)  
+        m_line.setCoords(DataArrayDouble(coords2, len(coords2)/2, 2))
+        m_line.setConnectivity(DataArrayInt(connec2), DataArrayInt(cI2))
+        a, b, c, d = MEDCouplingUMesh.Intersect2DMeshWith1DLine(m, m_line, eps)
+        self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
+        self.assertTrue(a.getCoords()[:m.getNumberOfNodes()].isEqual(m.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[m.getNumberOfNodes():m.getNumberOfNodes()+m_line.getNumberOfNodes()].isEqual(m_line.getCoords(),1e-12))
+        self.assertTrue(a.getCoords().isEqual(DataArrayDouble([(-1.,-1.),(-1.,1.),(1.,1.),(1.,-1.),(0.,1.3),(-1.3,0.),(-0.6,0.6),(0.,-1.3),(-0.5,-0.5),(-1.,0.23453685964236054),(-1.,-0.13033276368660177),(-0.2345368596423598,1.),(-0.1303327636866019,-1.),(-0.11489196370692323,1.1481421036683868),(-0.6,0.6),(-1.1481421036683859,0.11489196370692323),(-1.147455889106615,-0.0593103465193594),(-0.5,-0.5),(-0.0593103465193594,-1.147455889106615),(1.,0.),(0.4348336181566991,-1.),(-0.5651663818433009,-1.),(-1.,-0.5651663818433009),(-1.,0.05210204797787939),(-0.6,0.6),(0.3827315701788201,1.),(-0.6172684298211799,1.),(-0.6,0.6),(-1.,0.6172684298211802),(-0.6,0.6),(0.3827315701788201,1.),(1.,0.),(0.4348336181566991,-1.),(-0.5,-0.5),(-1.,0.05210204797787939),(-1.,-0.5651663818433009),(-0.5,-0.5),(-0.5651663818433009,-1.)]),1e-12))
+        self.assertEqual([32,9,11,2,3,12,10,29,30,31,32,33,34,32,0,10,12,35,36,37,32,1,11,9,26,27,28],a.getNodalConnectivity().getValues())
+        self.assertEqual([0,13,20,27],a.getNodalConnectivityIndex().getValues())
+        self.assertEqual([2,4,11,13,2,11,9,14,2,9,5,15,2,5,10,16,2,10,12,17,2,12,7,18],b.getNodalConnectivity().getValues())
+        self.assertEqual([0,4,8,12,16,20,24],b.getNodalConnectivityIndex().getValues())
+        self.assertTrue(a.getCoords()[:4].isEqual(m.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[4:9].isEqual(m_line.getCoords(),1e-12))
+        self.assertTrue(DataArrayInt([0,0,0]).isEqual(c))
+        self.assertTrue(DataArrayInt([(-1,-1),(0,2),(-1,-1),(-1,-1),(0,1),(-1,-1)]).isEqual(d))
+        pass
+
+    def testSwig2Intersect2DMeshWith1DLine12(self):
+        """ Two squares one in the other intersected by an horizontal line """
+        eps = 1.0e-8
+        m = MEDCouplingUMesh("boxbox", 2)
+        m.setCoords(DataArrayDouble([-0.5,-0.5,-0.5,0.5,0.5,0.5,0.5,-0.5,-0.25,-0.25,-0.25,0.25,0.25,0.25,0.25,-0.25],8,2))
+        c = [NORM_POLYGON, 4, 5, 6, 7, NORM_POLYGON, 0, 1, 5, 4, NORM_POLYGON, 1, 2, 3, 0, 4, 7, 6, 5]
+        cI = [0, 5, 10, 19]
+        m.setConnectivity(DataArrayInt(c), DataArrayInt(cI))
+        m.checkCoherency()
+        coords2 = [-1., 0.25, 1., 0.25]
+        connec2, cI2 = [NORM_SEG2, 0, 1], [0,3]
+        m_line = MEDCouplingUMesh.New("seg", 1)  
+        m_line.setCoords(DataArrayDouble(coords2, len(coords2)/2, 2))
+        m_line.setConnectivity(DataArrayInt(connec2), DataArrayInt(cI2))
+        m_line2 = m_line.deepCpy()
+        m2 = m.deepCpy()
+        a, b, c, d = MEDCouplingUMesh.Intersect2DMeshWith1DLine(m, m_line, eps)
+        self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
+        self.assertTrue(a.getCoords()[:m.getNumberOfNodes()].isEqual(m.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[m.getNumberOfNodes():m.getNumberOfNodes()+m_line.getNumberOfNodes()].isEqual(m_line.getCoords(),1e-12))
+        self.assertTrue(a.getCoords().isEqual(DataArrayDouble([(-0.5,-0.5),(-0.5,0.5),(0.5,0.5),(0.5,-0.5),(-0.25,-0.25),(-0.25,0.25),(0.25,0.25),(0.25,-0.25),(-1.,0.25),(1.,0.25),(-0.5,0.25),(0.5,0.25)]),1e-12))
+        self.assertEqual([5,4,5,6,7,5,1,5,10,5,4,0,10,5,5,5,1,2,11,6,5,3,0,4,7,6,11],a.getNodalConnectivity().getValues())
+        self.assertEqual([0,5,9,14,20,27],a.getNodalConnectivityIndex().getValues())
+        self.assertEqual([1,8,10,1,10,5,1,5,6,1,6,11,1,11,9],b.getNodalConnectivity().getValues())
+        self.assertEqual([0,3,6,9,12,15],b.getNodalConnectivityIndex().getValues())
+        self.assertTrue(c.isEqual(DataArrayInt([0,1,1,2,2])))
+        self.assertTrue(d.isEqual(DataArrayInt([(-1,-1),(1,2),(3,0),(3,4),(-1,-1)])))
+        pass
+
+    def testSwig2Intersect2DMeshWith1DLine13(self):
+        """ A square (side length) in a circle intersected by a simple horizontal line """
+        import math
+        eps = 1.0e-8
+        m = MEDCouplingUMesh("boxcircle", 2)
+        sq2 = math.sqrt(2.0)
+        soth = (sq2+1.0)/2.0
+        coo = [2., 0., sq2, sq2, 0., 2., -sq2, sq2, -2., 0., -sq2, -sq2, 0., -2., sq2, -sq2, -1., -1., -1., 1., 1., 
+         1., 1., -1., -1., 0., 0., 1., 1., 0., 0., -1., -soth, soth, soth,soth]
+        coo = DataArrayDouble(coo); coo.rearrange(2) 
+        m.setCoords(coo)
+        c = [NORM_QPOLYG, 8, 9, 10, 11, 12, 13, 14, 15, NORM_QPOLYG, 3, 1, 10, 9, 2, 17, 13, 16, NORM_QPOLYG, 1, 7, 5, 3, 9, 8, 11, 10, 0, 6, 4, 16, 12, 15, 14, 17]
+        cI = [0, 9, 18, 35]
+        m.setConnectivity(DataArrayInt(c), DataArrayInt(cI))
+        m.checkCoherency()
+        coords2 = [-2., 1., 2., 1.0]
+        connec2, cI2 = [NORM_SEG2, 0, 1], [0,3]
+        m_line = MEDCouplingUMesh("seg", 1)  
+        m_line.setCoords(DataArrayDouble(coords2, len(coords2)/2, 2))
+        m_line.setConnectivity(DataArrayInt(connec2), DataArrayInt(cI2))
+        a, b, c, d = MEDCouplingUMesh.Intersect2DMeshWith1DLine(m, m_line, eps)
+        self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
+        self.assertTrue(a.getCoords()[:m.getNumberOfNodes()].isEqual(m.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[m.getNumberOfNodes():m.getNumberOfNodes()+m_line.getNumberOfNodes()].isEqual(m_line.getCoords(),1e-12))
+        self.assertTrue(a.getCoords().isEqual(DataArrayDouble([(2.,0.),(1.4142135623730951,1.4142135623730951),(0.,2.),(-1.4142135623730951,1.4142135623730951),(-2.,0.),(-1.4142135623730951,-1.4142135623730951),(0.,-2.),(1.4142135623730951,-1.4142135623730951),(-1.,-1.),(-1.,1.),(1.,1.),(1.,-1.),(-1.,0.),(0.,1.),(1.,0.),(0.,-1.),(-1.2071067811865475,1.2071067811865475),(1.2071067811865475,1.2071067811865475),(-2.,1.),(2.,1.),(1.7320508075688772,1.),(-1.7320508075688772,1.),(-1.2071067811865475,1.2071067811865475),(-1.3660254037844386,1.),(-1.58670668058247,1.2175228580174415),(0.,-1.),(1.,0.),(1.2071067811865475,1.2071067811865475),(1.5867066805824703,1.2175228580174413),(1.9828897227476205,-0.26105238444010315),(0.,-2.),(-1.9828897227476205,-0.2610523844401032),(-1.3660254037844386,1.),(-1.,0.),(1.5867066805824703,1.2175228580174413),(1.3660254037844386,1.),(1.2071067811865475,1.2071067811865475),(0.,-2.),(-1.9828897227476205,-0.2610523844401032),(-1.3660254037844386,1.),(-1.,0.),(0.,-1.),(1.,0.),(1.3660254037844386,1.),(1.9828897227476205,-0.26105238444010315)]),1e-12))
+        self.assertEqual([32,8,9,10,11,12,13,14,15,32,3,1,10,9,2,17,13,16,32,3,9,21,22,23,24,32,1,20,10,34,35,36,32,7,5,21,9,8,11,10,20,37,38,39,40,41,42,43,44],a.getNodalConnectivity().getValues())
+        self.assertEqual([0,9,18,25,32,49],a.getNodalConnectivityIndex().getValues())
+        self.assertEqual([1,18,21,1,21,9,1,9,10,1,10,20,1,20,19],b.getNodalConnectivity().getValues())
+        self.assertEqual([0,3,6,9,12,15],b.getNodalConnectivityIndex().getValues())
+        self.assertTrue(c.isEqual(DataArrayInt([0,1,2,2,2])))
+        self.assertTrue(d.isEqual(DataArrayInt([(-1,-1),(2,4),(1,0),(3,4),(-1,-1)])))
+        pass
+
+    def testSwig2Intersect2DMeshWith1DLine14(self):
+        """ A circle in a circle intersected by a simple horizontal line, not tangent to the circles """
+        eps = 1.0e-8
+        m = MEDCouplingUMesh("boxcircle", 2)
+        coo = [2.,0.,1.4142135623730951,1.414213562373095,0.,2.,-1.414213562373095,1.4142135623730951,-2.,0.,-1.4142135623730954,-1.414213562373095,0.,-2.,
+               1.4142135623730947,-1.4142135623730954,1.,0.,0.7071067811865476,0.7071067811865475,0.,1.,-0.7071067811865475,0.7071067811865476,-1.,0.,-0.7071067811865477,-0.7071067811865475,
+               0.,-1.,0.7071067811865474,-0.7071067811865477,1.060660171779821,-1.0606601717798214,-1.0606601717798214,-1.0606601717798212]
+        coo = DataArrayDouble(coo); coo.rearrange(2) 
+        m.setCoords(coo)
+        c = [NORM_QPOLYG, 15, 13, 11, 9, 14, 12, 10, 8, NORM_QPOLYG, 7, 5, 13, 15, 6, 17, 14, 16, NORM_QPOLYG, 5, 3, 1, 7, 15, 9, 11, 13, 4, 2, 0, 16, 8, 10, 12, 17]
+        cI = [0, 9, 18, 35] 
+        m.setConnectivity(DataArrayInt(c), DataArrayInt(cI))
+        m.checkCoherency()
+        coords2 = [-2., 0., 2., 0.]
+        connec2, cI2 = [NORM_SEG2, 0, 1], [0,3]
+        m_line = MEDCouplingUMesh.New("seg", 1)  
+        m_line.setCoords(DataArrayDouble(coords2, len(coords2)/2, 2))
+        m_line.setConnectivity(DataArrayInt(connec2), DataArrayInt(cI2))
+        a, b, c, d = MEDCouplingUMesh.Intersect2DMeshWith1DLine(m, m_line, eps)
+        self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
+        self.assertTrue(a.getCoords()[:m.getNumberOfNodes()].isEqual(m.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[m.getNumberOfNodes():m.getNumberOfNodes()+m_line.getNumberOfNodes()].isEqual(m_line.getCoords(),1e-12))
+        self.assertTrue(a.getCoords().isEqual(DataArrayDouble([(2.,0.),(1.4142135623730951,1.414213562373095),(0.,2.),(-1.414213562373095,1.4142135623730951),(-2.,0.),(-1.4142135623730954,-1.414213562373095),(0.,-2.),(1.4142135623730947,-1.4142135623730954),(1.,0.),(0.7071067811865476,0.7071067811865475),(0.,1.),(-0.7071067811865475,0.7071067811865476),(-1.,0.),(-0.7071067811865477,-0.7071067811865475),(0.,-1.),(0.7071067811865474,-0.7071067811865477),(1.060660171779821,-1.0606601717798214),(-1.0606601717798214,-1.0606601717798212),(-2.,0.),(2.,0.),(-1.,0.),(1.,0.),(0.,2.),(1.8477590650225735,0.7653668647301795),(1.8477590650225735,-0.7653668647301797),(1.060660171779821,-1.0606601717798214),(0.9238795325112867,-0.38268343236508984),(0.9238795325112867,0.3826834323650897),(0.,1.),(-0.9238795325112867,0.3826834323650896),(-1.5,0.),(-1.8477590650225735,0.7653668647301792),(-1.0606601717798214,-1.0606601717798212),(-1.8477590650225733,-0.7653668647301799),(-1.5,0.),(-0.9238795325112866,-0.38268343236508995),(0.,1.),(-0.9238795325112867,0.3826834323650896),(-1.5,0.),(-1.8477590650225735,0.7653668647301792),(0.,2.),(1.8477590650225735,0.7653668647301795),(1.5,0.),(0.9238795325112867,0.3826834323650897),(1.060660171779821,-1.0606601717798214),(0.9238795325112867,-0.38268343236508984),(1.5,0.),(1.8477590650225735,-0.7653668647301797),(0.,1.),(0.9238795325112867,0.3826834323650897),(0.,0.),(-0.9238795325112867,0.3826834323650896),(0.,-1.),(-0.9238795325112866,-0.38268343236508995),(0.,0.),(0.9238795325112867,-0.38268343236508984)]),1e-12))
+        self.assertEqual([32,7,5,13,15,6,17,14,16,32,9,11,20,18,3,1,19,21,36,37,38,39,40,41,42,43,32,7,15,21,19,44,45,46,47,32,13,5,18,20,32,33,34,35,32,11,9,21,20,48,49,50,51,32,15,13,20,21,52,53,54,55],a.getNodalConnectivity().getValues())
+        self.assertEqual([0,9,26,35,44,53,62],a.getNodalConnectivityIndex().getValues())
+        self.assertEqual([1,18,20,1,20,21,1,21,19],b.getNodalConnectivity().getValues())
+        self.assertEqual([0,3,6,9],b.getNodalConnectivityIndex().getValues())
+        self.assertTrue(c.isEqual(DataArrayInt([1,2,2,2,0,0])))
+        self.assertTrue(d.isEqual(DataArrayInt([(1,3),(4,5),(1,2)])))
+        pass
+
+    def testSwig2Intersect2DMeshWith1DLine15(self):
+        """ Same as testSwig2Intersect2DMeshWith1DLine13 except that the line is colinear AND splits on of the common edge of 2D mesh."""
+        import math
+        eps = 1.0e-8
+        m = MEDCouplingUMesh("boxcircle", 2)
+        sq2 = math.sqrt(2.0)
+        soth = (sq2+1.0)/2.0
+        coo = [2., 0., sq2, sq2, 0., 2., -sq2, sq2, -2., 0., -sq2, -sq2, 0., -2., sq2, -sq2, -1., -1., -1., 1., 1., 
+         1., 1., -1., -1., 0., 0., 1., 1., 0., 0., -1., -soth, soth, soth,soth]
+        coo = DataArrayDouble(coo); coo.rearrange(2) 
+        m.setCoords(coo)
+        c = [NORM_QPOLYG, 8, 9, 10, 11, 12, 13, 14, 15, NORM_QPOLYG, 3, 1, 10, 9, 2, 17, 13, 16, NORM_QPOLYG, 1, 7, 5, 3, 9, 8, 11, 10, 0, 6, 4, 16, 12, 15, 14, 17]
+        cI = [0, 9, 18, 35]
+        m.setConnectivity(DataArrayInt(c), DataArrayInt(cI))
+        m.checkCoherency()
+        coords2 = [(-2., 1.),(2.,1.),(0.,1)]
+        connec2, cI2 = [NORM_SEG2, 0, 2, NORM_SEG2, 2, 1], [0,3,6]
+        m_line = MEDCouplingUMesh("seg", 1)  
+        m_line.setCoords(DataArrayDouble(coords2))
+        m_line.setConnectivity(DataArrayInt(connec2), DataArrayInt(cI2))
+        a, b, c, d = MEDCouplingUMesh.Intersect2DMeshWith1DLine(m, m_line, eps)
+        self.assertTrue(a.getCoords().getHiddenCppPointer()==b.getCoords().getHiddenCppPointer())
+        self.assertTrue(a.getCoords()[:m.getNumberOfNodes()].isEqual(m.getCoords(),1e-12))
+        self.assertTrue(a.getCoords()[m.getNumberOfNodes():m.getNumberOfNodes()+m_line.getNumberOfNodes()].isEqual(m_line.getCoords(),1e-12))
+        self.assertTrue(a.getCoords().isEqual(DataArrayDouble([(2.,0.),(1.4142135623730951,1.4142135623730951),(0.,2.),(-1.4142135623730951,1.4142135623730951),(-2.,0.),(-1.4142135623730951,-1.4142135623730951),(0.,-2.),(1.4142135623730951,-1.4142135623730951),(-1.,-1.),(-1.,1.),(1.,1.),(1.,-1.),(-1.,0.),(0.,1.),(1.,0.),(0.,-1.),(-1.2071067811865475,1.2071067811865475),(1.2071067811865475,1.2071067811865475),(-2.,1.),(2.,1.),(0.,1.),(1.7320508075688776,1.),(-1.7320508075688776,1.),(-0.5,1.),(0.5,1.),(0.5,1.),(-0.5,1.),(-1.2071067811865475,1.2071067811865475),(-1.3660254037844388,1.),(-1.58670668058247,1.2175228580174415),(0.,-1.),(1.,0.),(1.2071067811865475,1.2071067811865475),(1.5867066805824703,1.2175228580174413),(1.9828897227476205,-0.26105238444010315),(0.,-2.),(-1.9828897227476205,-0.2610523844401032),(-1.3660254037844388,1.),(-1.,0.),(1.5867066805824703,1.2175228580174413),(1.3660254037844388,1.),(1.2071067811865475,1.2071067811865475),(0.,-2.),(-1.9828897227476205,-0.2610523844401032),(-1.3660254037844388,1.),(-1.,0.),(0.,-1.),(1.,0.),(1.3660254037844388,1.),(1.9828897227476205,-0.26105238444010315)]),1e-12))
+        self.assertEqual([32,8,9,20,10,11,12,23,24,14,15,32,3,1,10,20,9,2,17,25,26,16,32,3,9,22,27,28,29,32,1,21,10,39,40,41,32,7,5,22,9,8,11,10,21,42,43,44,45,46,47,48,49],a.getNodalConnectivity().getValues())
+        self.assertEqual([0,11,22,29,36,53],a.getNodalConnectivityIndex().getValues())
+        self.assertEqual([1,18,22,1,22,9,1,9,20,1,20,10,1,10,21,1,21,19],b.getNodalConnectivity().getValues())
+        self.assertEqual([0,3,6,9,12,15,18],b.getNodalConnectivityIndex().getValues())
+        self.assertTrue(c.isEqual(DataArrayInt([0,1,2,2,2])))
+        self.assertTrue(d.isEqual(DataArrayInt([(-1,-1),(2,4),(1,0),(1,0),(3,4),(-1,-1)])))
+        pass
 
     def testOrderConsecutiveCells1D1(self):
         """A line in several unconnected pieces:"""
@@ -16046,6 +16204,47 @@ class MEDCouplingBasicsTest(unittest.TestCase):
         self.assertEqual(p2.getSlice(),slice(2,11,4))
         pass
 
+    def testSwig2DAIGetIdsStrictlyNegative1(self):
+        d=DataArrayInt([4,-5,-1,0,3,99,-7])
+        self.assertTrue(d.getIdsStrictlyNegative().isEqual(DataArrayInt([1,2,6])))
+        pass
+
+    def testSwig2DAIReplaceOneValByInThis1(self):
+        d=DataArrayInt([4,-5,-1,0,-5,99,-7,5])
+        d.replaceOneValByInThis(-5,900)
+        self.assertTrue(d.isEqual(DataArrayInt([4,900,-1,0,900,99,-7,5])))
+        pass
+
+    def testSwig2DAIGetMinMaxValues1(self):
+        d=DataArrayInt([4,-5,-1,0,3,99,-7])
+        a,b=d.getMinMaxValues()
+        self.assertEqual(a,-7)
+        self.assertEqual(b,99)
+        pass
+
+    def testSwig2DAIBuildUniqueNotSorted1(self):
+        d=DataArrayInt([-5,3,2,-1,2,3,-6,4,2,-5,3,7])
+        self.assertTrue(d.buildUniqueNotSorted().isEqual(DataArrayInt([-5,3,2,-1,-6,4,7])))
+        pass
+
+    def testSwig2UMeshChangeOrientationOfCells1(self):
+        """ Here testing changeOrientationOfCell method on unstructured meshes lying on no coords."""
+        m=MEDCouplingUMesh("mesh",1)
+        c=DataArrayInt([NORM_SEG2,4,5,NORM_SEG2,10,8,NORM_SEG3,20,7,33,NORM_SEG3,13,15,12,NORM_SEG2,3,2,NORM_SEG4,5,6,8,10,NORM_SEG4,34,33,3,2])
+        cI=DataArrayInt([0,3,6,10,14,17,22,27])
+        m.setConnectivity(c,cI)
+        m.changeOrientationOfCells()
+        self.assertTrue(m.getNodalConnectivity().isEqual(DataArrayInt([NORM_SEG2,5,4,NORM_SEG2,8,10,NORM_SEG3,7,20,33,NORM_SEG3,15,13,12,NORM_SEG2,2,3,NORM_SEG4,6,5,10,8,NORM_SEG4,33,34,2,3])))
+        self.assertTrue(m.getNodalConnectivityIndex().isEqual(cI))
+        # testing 2D cells
+        m=MEDCouplingUMesh("mesh",2)
+        c=DataArrayInt([NORM_TRI3,0,1,2,NORM_QUAD4,3,4,5,6,NORM_POLYGON,7,8,9,10,11,NORM_TRI6,12,13,14,15,16,17,NORM_QUAD8,18,19,20,21,22,23,24,25,NORM_QPOLYG,26,27,28,29,30,31,32,33,34,35])
+        cI=DataArrayInt([0,4,9,15,22,31,42])
+        m.setConnectivity(c,cI)
+        m.changeOrientationOfCells()
+        self.assertTrue(m.getNodalConnectivity().isEqual(DataArrayInt([NORM_TRI3,0,2,1,NORM_QUAD4,3,6,5,4,NORM_POLYGON,7,11,10,9,8,NORM_TRI6,12,14,13,17,16,15,NORM_QUAD8,18,21,20,19,25,24,23,22,NORM_QPOLYG,26,30,29,28,27,35,34,33,32,31])))
+        self.assertTrue(m.getNodalConnectivityIndex().isEqual(cI))
+        pass
     pass
 
 if __name__ == '__main__':
index b78056bfb231b87ddf145d4a5f946c0b1dcee346..3866bfe1b405c2d86062db448ebe905b323a58a4 100644 (file)
@@ -1687,6 +1687,7 @@ namespace ParaMEDMEM
     DataArrayDouble *getBoundingBoxForBBTreeFast() const throw(INTERP_KERNEL::Exception);
     DataArrayDouble *getBoundingBoxForBBTree2DQuadratic(double arcDetEps=1e-12) const throw(INTERP_KERNEL::Exception);
     DataArrayDouble *getBoundingBoxForBBTree1DQuadratic(double arcDetEps=1e-12) const throw(INTERP_KERNEL::Exception);
+    void changeOrientationOfCells() throw(INTERP_KERNEL::Exception);
     int split2DCells(const DataArrayInt *desc, const DataArrayInt *descI, const DataArrayInt *subNodesInSeg, const DataArrayInt *subNodesInSegI, const DataArrayInt *midOpt=0, const DataArrayInt *midOptI=0) throw(INTERP_KERNEL::Exception);
     static MEDCouplingUMesh *Build0DMeshFromCoords(DataArrayDouble *da) throw(INTERP_KERNEL::Exception);
     static MEDCouplingUMesh *MergeUMeshes(const MEDCouplingUMesh *mesh1, const MEDCouplingUMesh *mesh2) throw(INTERP_KERNEL::Exception);
index 911cb7a39a666c74833891b9b993b26057ba5b3a..544ba336793293f23bfd15d02459562e5bb1857a 100644 (file)
@@ -73,6 +73,7 @@
 %newobject ParaMEDMEM::DataArrayInt::computeAbs;
 %newobject ParaMEDMEM::DataArrayInt::getIdsInRange;
 %newobject ParaMEDMEM::DataArrayInt::getIdsNotInRange;
+%newobject ParaMEDMEM::DataArrayInt::getIdsStrictlyNegative;
 %newobject ParaMEDMEM::DataArrayInt::Aggregate;
 %newobject ParaMEDMEM::DataArrayInt::AggregateIndexes;
 %newobject ParaMEDMEM::DataArrayInt::Meld;
@@ -92,6 +93,7 @@
 %newobject ParaMEDMEM::DataArrayInt::buildSubstractionOptimized;
 %newobject ParaMEDMEM::DataArrayInt::buildIntersection;
 %newobject ParaMEDMEM::DataArrayInt::buildUnique;
+%newobject ParaMEDMEM::DataArrayInt::buildUniqueNotSorted;
 %newobject ParaMEDMEM::DataArrayInt::deltaShiftIndex;
 %newobject ParaMEDMEM::DataArrayInt::buildExplicitArrByRanges;
 %newobject ParaMEDMEM::DataArrayInt::buildExplicitArrOfSliceOnScaledArr;
@@ -2507,6 +2509,7 @@ namespace ParaMEDMEM
     void fillWithZero() throw(INTERP_KERNEL::Exception);
     void fillWithValue(int val) throw(INTERP_KERNEL::Exception);
     void iota(int init=0) throw(INTERP_KERNEL::Exception);
+    void replaceOneValByInThis(int valToBeReplaced, int replacedBy) throw(INTERP_KERNEL::Exception);
     std::string repr() const throw(INTERP_KERNEL::Exception);
     std::string reprZip() const throw(INTERP_KERNEL::Exception);
     DataArrayInt *invertArrayO2N2N2O(int newNbOfElem) const throw(INTERP_KERNEL::Exception);
@@ -2568,6 +2571,7 @@ namespace ParaMEDMEM
     void applyRPow(int val) throw(INTERP_KERNEL::Exception);
     DataArrayInt *getIdsInRange(int vmin, int vmax) const throw(INTERP_KERNEL::Exception);
     DataArrayInt *getIdsNotInRange(int vmin, int vmax) const throw(INTERP_KERNEL::Exception);
+    DataArrayInt *getIdsStrictlyNegative() const throw(INTERP_KERNEL::Exception);
     bool checkAllIdsInRange(int vmin, int vmax) const throw(INTERP_KERNEL::Exception);
     static DataArrayInt *Aggregate(const DataArrayInt *a1, const DataArrayInt *a2, int offsetA2) throw(INTERP_KERNEL::Exception);
     static DataArrayInt *Meld(const DataArrayInt *a1, const DataArrayInt *a2) throw(INTERP_KERNEL::Exception);
@@ -2581,6 +2585,7 @@ namespace ParaMEDMEM
     DataArrayInt *buildUnion(const DataArrayInt *other) const throw(INTERP_KERNEL::Exception);
     DataArrayInt *buildIntersection(const DataArrayInt *other) const throw(INTERP_KERNEL::Exception);
     DataArrayInt *buildUnique() const throw(INTERP_KERNEL::Exception);
+    DataArrayInt *buildUniqueNotSorted() const throw(INTERP_KERNEL::Exception);
     DataArrayInt *deltaShiftIndex() const throw(INTERP_KERNEL::Exception);
     void computeOffsets() throw(INTERP_KERNEL::Exception);
     void computeOffsets2() throw(INTERP_KERNEL::Exception);
@@ -2782,6 +2787,16 @@ namespace ParaMEDMEM
           throw INTERP_KERNEL::Exception("DataArrayInt::buildExplicitArrOfSliceOnScaledArr (wrap) : the input slice contains some unknowns that can't be determined in static method ! Call DataArray::getSlice (non static) instead !");
         return self->buildExplicitArrOfSliceOnScaledArr(strt,stp,step);
       }
+
+      PyObject *getMinMaxValues() const throw(INTERP_KERNEL::Exception)
+      {
+        int a,b;
+        self->getMinMaxValues(a,b);
+        PyObject *ret=PyTuple_New(2);
+        PyTuple_SetItem(ret,0,PyInt_FromLong(a));
+        PyTuple_SetItem(ret,1,PyInt_FromLong(b));
+        return ret;
+      }
    
       static PyObject *BuildOld2NewArrayFromSurjectiveFormat2(int nbOfOldTuples, PyObject *arr, PyObject *arrI) throw(INTERP_KERNEL::Exception)
       {