]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
Implementation of mergefields on Gauss Points + unification of mergefields(f1,f2...
authorAnthony Geay <anthony.geay@edf.fr>
Fri, 4 Jan 2019 09:41:53 +0000 (10:41 +0100)
committerAnthony Geay <anthony.geay@edf.fr>
Fri, 4 Jan 2019 10:30:08 +0000 (11:30 +0100)
src/MEDCoupling/MCAuto.hxx
src/MEDCoupling/MEDCouplingField.cxx
src/MEDCoupling/MEDCouplingFieldDiscretization.cxx
src/MEDCoupling/MEDCouplingFieldDiscretization.hxx
src/MEDCoupling/MEDCouplingFieldDouble.cxx
src/MEDCoupling/Test/MEDCouplingBasicsTest1.cxx
src/MEDCoupling_Swig/MEDCouplingBasicsTest1.py
src/MEDCoupling_Swig/MEDCouplingBasicsTest6.py

index b680659fa521a152a962a8f2546bb16d5533fd2f..23081e94dc51fcdab7864377c83d67a534497e5c 100644 (file)
@@ -24,6 +24,8 @@
 #include "MEDCouplingRefCountObject.hxx"
 #include "InterpKernelException.hxx"
 
+#include <vector>
+
 namespace MEDCoupling
 {
   template<class T>
@@ -80,6 +82,16 @@ namespace MEDCoupling
     return ret;
   }
 
+  template<class T>
+  typename std::vector<const T *> ToConstVect(const typename std::vector< MCAuto<T> >& vec)
+  {
+    std::size_t sz(vec.size());
+    std::vector<const T *> ret(sz);
+    for(std::size_t i=0;i<sz;i++)
+      ret[i]=(const T *)vec[i];
+    return ret;
+  }
+  
   template<class T>
   class MCConstAuto
   {
index 94c0a12fad091179c86d4e2c1fc3586e300646dd..d3c037c799d748076d3799a9342ac064757f79de 100644 (file)
@@ -114,7 +114,9 @@ bool MEDCouplingField::areCompatibleForMerge(const MEDCouplingField *other) cons
 {
   if(!other)
     throw INTERP_KERNEL::Exception("MEDCouplingField::areCompatibleForMerge : input field is NULL !");
-  if(!_type->isEqual(other->_type,1.))
+  if(!_type || !other->_type)
+    throw INTERP_KERNEL::Exception("MEDCouplingField::areCompatibleForMerge : this or other has a nullptr spatial discretization !");
+  if(_type->getEnum()!=other->_type->getEnum())
     return false;
   if(_nature!=other->_nature)
     return false;
index ebec0b1967c40c643440a11f2b23eb6cecf54c82..0c93445c52dc0cb5b05d58b3d8c154c8186810f5 100644 (file)
@@ -476,6 +476,20 @@ void MEDCouplingFieldDiscretization::RenumberEntitiesFromN2OArr(const int *new2O
     }
 }
 
+template<class FIELD_DISC>
+MCAuto<MEDCouplingFieldDiscretization> MEDCouplingFieldDiscretization::EasyAggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds)
+{
+  if(fds.empty())
+    throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::aggregate : input array is empty");
+  for(const MEDCouplingFieldDiscretization * it : fds)
+    {
+      const FIELD_DISC *itc(dynamic_cast<const FIELD_DISC *>(it));
+      if(!itc)
+        throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretization::aggregate : same field discretization expected for all input discretizations !");
+    }
+  return fds[0]->clone();
+}
+
 MEDCouplingFieldDiscretization::~MEDCouplingFieldDiscretization()
 {
 }
@@ -763,6 +777,11 @@ MEDCouplingMesh *MEDCouplingFieldDiscretizationP0::buildSubMeshDataRange(const M
   return ret.retn();
 }
 
+MCAuto<MEDCouplingFieldDiscretization> MEDCouplingFieldDiscretizationP0::aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const
+{
+  return EasyAggregate<MEDCouplingFieldDiscretizationP0>(fds);
+}
+
 int MEDCouplingFieldDiscretizationOnNodes::getNumberOfTuples(const MEDCouplingMesh *mesh) const
 {
   if(!mesh)
@@ -1082,6 +1101,11 @@ void MEDCouplingFieldDiscretizationP1::reprQuickOverview(std::ostream& stream) c
   stream << "P1 spatial discretization.";
 }
 
+MCAuto<MEDCouplingFieldDiscretization> MEDCouplingFieldDiscretizationP1::aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const
+{
+  return EasyAggregate<MEDCouplingFieldDiscretizationP1>(fds);
+}
+
 MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell():_discr_per_cell(0)
 {
 }
@@ -1116,6 +1140,12 @@ MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(con
     }
 }
 
+MEDCouplingFieldDiscretizationPerCell::MEDCouplingFieldDiscretizationPerCell(DataArrayInt *dpc):_discr_per_cell(dpc)
+{
+  if(_discr_per_cell)
+    _discr_per_cell->incrRef();
+}
+
 void MEDCouplingFieldDiscretizationPerCell::updateTime() const
 {
   if(_discr_per_cell)
@@ -1843,6 +1873,50 @@ void MEDCouplingFieldDiscretizationGauss::renumberValuesOnCellsR(const MEDCoupli
   throw INTERP_KERNEL::Exception("Number of cells has changed and becomes higher with some cells that have been split ! Unable to conserve the Gauss field !");
 }
 
+MCAuto<MEDCouplingFieldDiscretization> MEDCouplingFieldDiscretizationGauss::aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const
+{
+  if(fds.empty())
+    throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::aggregate : input array is empty");
+  std::vector<MEDCouplingGaussLocalization> loc;//store the localizations for the output GaussDiscretization object
+  std::vector< MCAuto<DataArrayInt> > discPerCells(fds.size());
+  std::size_t i(0);
+  for(auto it=fds.begin();it!=fds.end();++it,++i)
+    {
+      const MEDCouplingFieldDiscretizationGauss *itc(dynamic_cast<const MEDCouplingFieldDiscretizationGauss *>(*it));
+      if(!itc)
+        throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::aggregate : same field discretization expected for all input discretizations !");
+      //
+      std::vector<MEDCouplingGaussLocalization> loc2(itc->_loc);
+      std::vector<int> newLocId(loc2.size());
+      for(std::size_t j=0;j<loc2.size();++j)
+        {
+          std::size_t k(0);
+          for(;k<loc.size();++k)
+            {
+              if(loc2[j].isEqual(loc[k],1e-10))
+                {
+                  newLocId[j]=(int)k;
+                  break;
+                }
+            }
+          if(k==loc.size())// current loc2[j]
+            {
+              newLocId[j]=(int)loc.size();
+              loc.push_back(loc2[j]);
+            }
+        }
+      const DataArrayInt *dpc(itc->_discr_per_cell);
+      if(!dpc)
+        throw INTERP_KERNEL::Exception("MEDCouplingFieldDiscretizationGauss::aggregate : Presence of nullptr array of disc per cell !");
+      MCAuto<DataArrayInt> dpc2(dpc->deepCopy());
+      dpc2->transformWithIndArr(newLocId.data(),newLocId.data()+newLocId.size());
+      discPerCells[i]=dpc2;
+    }
+  MCAuto<DataArrayInt> dpc3(DataArrayInt::Aggregate(ToConstVect(discPerCells)));
+  MCAuto<MEDCouplingFieldDiscretizationGauss> ret(new MEDCouplingFieldDiscretizationGauss(dpc3,loc));
+  return DynamicCast<MEDCouplingFieldDiscretizationGauss,MEDCouplingFieldDiscretization>(ret);
+}
+
 void MEDCouplingFieldDiscretizationGauss::setGaussLocalizationOnType(const MEDCouplingMesh *mesh, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
                                                                      const std::vector<double>& gsCoo, const std::vector<double>& wg)
 {
@@ -2718,6 +2792,11 @@ void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCells(double epsOnVa
   throw INTERP_KERNEL::Exception("Not implemented yet !");
 }
 
+MCAuto<MEDCouplingFieldDiscretization> MEDCouplingFieldDiscretizationGaussNE::aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const
+{
+  return EasyAggregate<MEDCouplingFieldDiscretizationGaussNE>(fds);
+}
+
 void MEDCouplingFieldDiscretizationGaussNE::renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const
 {
   throw INTERP_KERNEL::Exception("Not implemented yet !");
@@ -2813,6 +2892,11 @@ void MEDCouplingFieldDiscretizationKriging::reprQuickOverview(std::ostream& stre
   stream << "Kriging spatial discretization.";
 }
 
+MCAuto<MEDCouplingFieldDiscretization> MEDCouplingFieldDiscretizationKriging::aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const
+{
+  return EasyAggregate<MEDCouplingFieldDiscretizationKriging>(fds);
+}
+
 /*!
  * Returns the matrix of size nbRows = \a nbOfTargetPoints and \a nbCols = \a nbCols. This matrix is useful if 
  * 
index 37362a8a0662c42e0052923e3cd917314bb32f62..ece607e77124938d7a1b89fc80729e7b49d89a29 100644 (file)
@@ -87,6 +87,7 @@ namespace MEDCoupling
     MEDCOUPLING_EXPORT virtual void renumberValuesOnNodes(double epsOnVals, const int *old2New, int newNbOfNodes, DataArrayDouble *arr) const = 0;
     MEDCOUPLING_EXPORT virtual void renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const = 0;
     MEDCOUPLING_EXPORT virtual void renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const = 0;
+    MEDCOUPLING_EXPORT virtual MCAuto<MEDCouplingFieldDiscretization> aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const = 0;
     MEDCOUPLING_EXPORT virtual void getSerializationIntArray(DataArrayInt *& arr) const;
     MEDCOUPLING_EXPORT virtual void getTinySerializationIntInformation(std::vector<int>& tinyInfo) const;
     MEDCOUPLING_EXPORT virtual void getTinySerializationDbleInformation(std::vector<double>& tinyInfo) const;
@@ -111,6 +112,8 @@ namespace MEDCoupling
     MEDCOUPLING_EXPORT MEDCouplingFieldDiscretization();
     MEDCOUPLING_EXPORT static void RenumberEntitiesFromO2NArr(double epsOnVals, const int *old2NewPtr, int newNbOfEntity, DataArrayDouble *arr, const std::string& msg);
     MEDCOUPLING_EXPORT static void RenumberEntitiesFromN2OArr(const int *new2OldPtr, int new2OldSz, DataArrayDouble *arr, const std::string& msg);
+    template<class FIELD_DISC>
+    static MCAuto<MEDCouplingFieldDiscretization> EasyAggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds);
   protected:
     double _precision;
     static const double DFLT_PRECISION;
@@ -142,6 +145,7 @@ namespace MEDCoupling
     MEDCOUPLING_EXPORT void renumberValuesOnNodes(double epsOnVals, const int *old2New, int newNbOfNodes, DataArrayDouble *arr) const;
     MEDCOUPLING_EXPORT void renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const;
     MEDCOUPLING_EXPORT void renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const;
+    MEDCOUPLING_EXPORT MCAuto<MEDCouplingFieldDiscretization> aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const override;
     MEDCOUPLING_EXPORT MEDCouplingMesh *buildSubMeshData(const MEDCouplingMesh *mesh, const int *start, const int *end, DataArrayInt *&di) const;
     MEDCOUPLING_EXPORT MEDCouplingMesh *buildSubMeshDataRange(const MEDCouplingMesh *mesh, int beginCellIds, int endCellIds, int stepCellIds, int& beginOut, int& endOut, int& stepOut, DataArrayInt *&di) const;
     MEDCOUPLING_EXPORT DataArrayInt *computeTupleIdsToSelectFromCellIds(const MEDCouplingMesh *mesh, const int *startCellIds, const int *endCellIds) const;
@@ -187,6 +191,7 @@ namespace MEDCoupling
     MEDCOUPLING_EXPORT void getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const;
     MEDCOUPLING_EXPORT DataArrayDouble *getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const;
     MEDCOUPLING_EXPORT void reprQuickOverview(std::ostream& stream) const;
+    MEDCOUPLING_EXPORT MCAuto<MEDCouplingFieldDiscretization> aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const override;
   public:
     static const char REPR[];
     static const TypeOfField TYPE;
@@ -209,6 +214,7 @@ namespace MEDCoupling
     MEDCouplingFieldDiscretizationPerCell();
     MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, const int *startCellIds, const int *endCellIds);
     MEDCouplingFieldDiscretizationPerCell(const MEDCouplingFieldDiscretizationPerCell& other, int beginCellIds, int endCellIds, int stepCellIds);
+    MEDCouplingFieldDiscretizationPerCell(DataArrayInt *dpc);
     ~MEDCouplingFieldDiscretizationPerCell();
     void updateTime() const;
     std::size_t getHeapMemorySizeWithoutChildren() const;
@@ -265,6 +271,7 @@ namespace MEDCoupling
     MEDCOUPLING_EXPORT void renumberValuesOnNodes(double epsOnVals, const int *old2New, int newNbOfNodes, DataArrayDouble *arr) const;
     MEDCOUPLING_EXPORT void renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const;
     MEDCOUPLING_EXPORT void renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const;
+    MEDCOUPLING_EXPORT MCAuto<MEDCouplingFieldDiscretization> aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const override;
     MEDCOUPLING_EXPORT void setGaussLocalizationOnType(const MEDCouplingMesh *mesh, INTERP_KERNEL::NormalizedCellType type, const std::vector<double>& refCoo,
                                                        const std::vector<double>& gsCoo, const std::vector<double>& wg);
     MEDCOUPLING_EXPORT void setGaussLocalizationOnCells(const MEDCouplingMesh *mesh, const int *begin, const int *end, const std::vector<double>& refCoo,
@@ -284,6 +291,7 @@ namespace MEDCoupling
   protected:
     MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, const int *startCellIds=0, const int *endCellIds=0);
     MEDCouplingFieldDiscretizationGauss(const MEDCouplingFieldDiscretizationGauss& other, int beginCellIds, int endCellIds, int stepCellIds);
+    MEDCouplingFieldDiscretizationGauss(DataArrayInt *dpc, const std::vector<MEDCouplingGaussLocalization>& loc):MEDCouplingFieldDiscretizationPerCell(dpc),_loc(loc) { }
     void zipGaussLocalizations();
     int getOffsetOfCell(int cellId) const;
     void checkLocalizationId(int locId) const;
@@ -330,6 +338,7 @@ namespace MEDCoupling
     MEDCOUPLING_EXPORT void renumberValuesOnNodes(double epsOnVals, const int *old2New, int newNbOfNodes, DataArrayDouble *arr) const;
     MEDCOUPLING_EXPORT void renumberValuesOnCells(double epsOnVals, const MEDCouplingMesh *mesh, const int *old2New, int newSz, DataArrayDouble *arr) const;
     MEDCOUPLING_EXPORT void renumberValuesOnCellsR(const MEDCouplingMesh *mesh, const int *new2old, int newSz, DataArrayDouble *arr) const;
+    MEDCOUPLING_EXPORT MCAuto<MEDCouplingFieldDiscretization> aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const override;
     MEDCOUPLING_EXPORT void reprQuickOverview(std::ostream& stream) const;
     MEDCOUPLING_EXPORT static const double *GetWeightArrayFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth);
     MEDCOUPLING_EXPORT static const double *GetRefCoordsFromGeometricType(INTERP_KERNEL::NormalizedCellType geoType, std::size_t& lgth);
@@ -409,6 +418,7 @@ namespace MEDCoupling
     MEDCOUPLING_EXPORT void getValueOn(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, double *res) const;
     MEDCOUPLING_EXPORT DataArrayDouble *getValueOnMulti(const DataArrayDouble *arr, const MEDCouplingMesh *mesh, const double *loc, int nbOfPoints) const;
     MEDCOUPLING_EXPORT void reprQuickOverview(std::ostream& stream) const;
+    MEDCOUPLING_EXPORT MCAuto<MEDCouplingFieldDiscretization> aggregate(std::vector<const MEDCouplingFieldDiscretization *>& fds) const override;
   public://specific part
     MEDCOUPLING_EXPORT DataArrayDouble *computeEvaluationMatrixOnGivenPts(const MEDCouplingMesh *mesh, const double *loc, int nbOfTargetPoints, int& nbCols) const;
     MEDCOUPLING_EXPORT DataArrayDouble *computeInverseMatrix(const MEDCouplingMesh *mesh, int& isDrift, int& matSz) const;
index 6e6d5dc21c5c55c1b0e7802e9846071fa86940a4..e3e9d8d02a8709c56d11f110af1c0ad921e67619 100644 (file)
@@ -2361,24 +2361,9 @@ void MEDCouplingFieldDouble::sortPerTuple(bool asc)
  */
 MEDCouplingFieldDouble *MEDCouplingFieldDouble::MergeFields(const MEDCouplingFieldDouble *f1, const MEDCouplingFieldDouble *f2)
 {
-  if(!f1->areCompatibleForMerge(f2))
-    throw INTERP_KERNEL::Exception("Fields are not compatible. Unable to apply MergeFields on them ! Check support mesh, field nature, and spatial and time discretisation.");
-  const MEDCouplingMesh *m1(f1->getMesh()),*m2(f2->getMesh());
-  if(!f1->timeDiscr())
-    throw INTERP_KERNEL::Exception("MEDCouplingFieldDouble::MergeFields : no time discr of f1 !");
-  if(!f1->_type)
-    throw INTERP_KERNEL::Exception("MEDCouplingFieldDouble::MergeFields : no spatial discr of f1 !");
-  MEDCouplingTimeDiscretization *td(f1->timeDiscr()->aggregate(f2->timeDiscr()));
-  td->copyTinyAttrFrom(*f1->timeDiscr());
-  MCAuto<MEDCouplingFieldDouble> ret(new MEDCouplingFieldDouble(f1->getNature(),td,f1->_type->clone()));
-  ret->setName(f1->getName());
-  ret->setDescription(f1->getDescription());
-  if(m1)
-    {
-      MCAuto<MEDCouplingMesh> m=m1->mergeMyselfWith(m2);
-      ret->setMesh(m);
-    }
-  return ret.retn();
+  std::vector<const MEDCouplingFieldDouble *> a(2);
+  a[0]=f1; a[1]=f2;
+  return MergeFields(a);
 }
 
 /*!
@@ -2402,29 +2387,36 @@ MEDCouplingFieldDouble *MEDCouplingFieldDouble::MergeFields(const MEDCouplingFie
  */
 MEDCouplingFieldDouble *MEDCouplingFieldDouble::MergeFields(const std::vector<const MEDCouplingFieldDouble *>& a)
 {
-  if(a.size()<1)
-    throw INTERP_KERNEL::Exception("FieldDouble::MergeFields : size of array must be >= 1 !");
+  if(a.empty())
+    throw INTERP_KERNEL::Exception("FieldDouble::MergeFields : input array is empty !");
   std::vector< MCAuto<MEDCouplingUMesh> > ms(a.size());
   std::vector< const MEDCouplingUMesh *> ms2(a.size());
   std::vector< const MEDCouplingTimeDiscretization *> tds(a.size());
-  std::vector<const MEDCouplingFieldDouble *>::const_iterator it=a.begin();
-  const MEDCouplingFieldDouble *ref=(*it++);
+  std::vector< const MEDCouplingFieldDouble *>::const_iterator it=a.begin();
+  std::vector<const MEDCouplingFieldDiscretization *> fds(a.size());
+  const MEDCouplingFieldDouble *ref((*it++));
   if(!ref)
-    throw INTERP_KERNEL::Exception("MEDCouplingFieldDouble::MergeFields : presence of NULL instance in first place of input vector !");
+    throw INTERP_KERNEL::Exception("MEDCouplingFieldDouble::MergeFields : presence of nullptr instance in first place of input vector !");
+  if(!ref->getDiscretization())
+    throw INTERP_KERNEL::Exception("MEDCouplingFieldDouble::MergeFields : nullptr spatial discretization !");
   for(;it!=a.end();it++)
     if(!ref->areCompatibleForMerge(*it))
       throw INTERP_KERNEL::Exception("Fields are not compatible. Unable to apply MergeFields on them! Check support mesh, field nature, and spatial and time discretisation.");
   for(int i=0;i<(int)a.size();i++)
     {
+      if(!a[i])
+        throw INTERP_KERNEL::Exception("MEDCouplingFieldDouble::MergeFields : presence of nullptr instance in input vector !");
       if(a[i]->getMesh())
         { ms[i]=a[i]->getMesh()->buildUnstructured(); ms2[i]=ms[i]; }
       else
         { ms[i]=0; ms2[i]=0; }
       tds[i]=a[i]->timeDiscr();
+      fds[i]=a[i]->getDiscretization();
     }
   MEDCouplingTimeDiscretization *td(tds[0]->aggregate(tds));
+  MCAuto<MEDCouplingFieldDiscretization> fda(fds[0]->aggregate(fds));
   td->copyTinyAttrFrom(*(a[0]->timeDiscr()));
-  MCAuto<MEDCouplingFieldDouble> ret(new MEDCouplingFieldDouble(a[0]->getNature(),td,a[0]->_type->clone()));
+  MCAuto<MEDCouplingFieldDouble> ret(new MEDCouplingFieldDouble(a[0]->getNature(),td,fda.retn()));
   ret->setName(a[0]->getName());
   ret->setDescription(a[0]->getDescription());
   if(ms2[0])
index 99abe4394994778fa51f475a6dffebe686d2898c..4e859d3225f69168f53e482843643b9e08640f58 100644 (file)
@@ -1373,6 +1373,7 @@ void MEDCouplingBasicsTest1::testMergeField1()
   MEDCouplingFieldDouble *f3=MEDCouplingFieldDouble::MergeFields(f1,f2);
   f3->checkConsistencyLight();
   MEDCouplingUMesh *m4=build2DTargetMeshMerged_1();
+  m4->setName(f1->getMesh()->getName());
   CPPUNIT_ASSERT(f3->getMesh()->isEqual(m4,1.e-12));
   std::string name=f3->getName();
   CPPUNIT_ASSERT(name=="MeasureOfMesh_");
index 2d6a4c0cb3b4294c5ad28f1d6d7e93ad0f5e5dea..ae07b6d34464c5997133b83b26d69f62b672d2dc 100644 (file)
@@ -1030,6 +1030,7 @@ class MEDCouplingBasicsTest1(unittest.TestCase):
         f3=MEDCouplingFieldDouble.MergeFields(f1,f2);
         f3.checkConsistencyLight();
         m4=MEDCouplingDataForTest.build2DTargetMeshMerged_1();
+        m4.setName(f1.getMesh().getName())
         self.assertTrue(f3.getMesh().isEqual(m4,1.e-12));
         name=f3.getName();
         self.assertEqual(name,"MeasureOfMesh_");
index 6c48d9b45bea850a551a13e180e51a20078b4443..0bf6ee2455af109c800668a37386ee5c1dee67fc 100644 (file)
@@ -325,6 +325,77 @@ If true geometry (with curve as edges) is considered the result of getCellsConta
         self.assertTrue(m.getNodalConnectivityIndex().isEqual(expConnI))
         pass
 
+    def testMergeFieldsOnGauss1(self):
+        mName="mesh"
+        fieldName="field"
+        #
+        _a=0.446948490915965;
+        _b=0.091576213509771;
+        _p1=0.11169079483905;
+        _p2=0.0549758718227661;
+        refCoo1=[ 0.,0., 1.,0., 0.,1. ]
+        gsCoo1=[ 2*_b-1, 1-4*_b, 2*_b-1, 2.07*_b-1, 1-4*_b,
+                 2*_b-1, 1-4*_a, 2*_a-1, 2*_a-1, 1-4*_a, 2*_a-1, 2*_a-1 ]
+        wg1=[ 4*_p2, 4*_p2, 4*_p2, 4*_p1, 4*_p1, 4*_p1 ]
+        #
+        refCoo2=[ -1.,-1., 1.,-1., 1.,1., -1.,1. ]
+        gsCoo2=[0.1,0.1, 0.2,0.2, 0.5,0.5, 0.6,0.6, 0.7,0.7]
+        wg2=[0.1,0.2,0.3,0.4,0.5]
+        #
+        coo=DataArrayDouble([0,0,1,0,2,0,0,1,1,1,2,1,0,2,1,2,2,2],9,2)
+        m1=MEDCouplingUMesh(mName,2)
+        m1.allocateCells() ; m1.setCoords(coo)
+        m1.insertNextCell(NORM_TRI3,[1,4,2])
+        m1.insertNextCell(NORM_TRI3,[4,5,2])
+        m1.insertNextCell(NORM_QUAD4,[4,7,8,5])
+        f1=MEDCouplingFieldDouble(ON_GAUSS_PT) ; f1.setName(fieldName)
+        f1.setMesh(m1)
+        f1.setGaussLocalizationOnType(NORM_TRI3,refCoo1,gsCoo1,wg1)
+        f1.setGaussLocalizationOnType(NORM_QUAD4,refCoo2,gsCoo2,wg2)
+        arr=DataArrayDouble(f1.getNumberOfTuplesExpected())
+        arr.iota()
+        f1.setArray(arr)
+        f1.checkConsistencyLight()
+        #
+        m2=MEDCouplingUMesh(mName,2)
+        m2.allocateCells() ; m2.setCoords(coo)
+        m2.insertNextCell(NORM_QUAD4,[0,3,4,1])
+        m2.insertNextCell(NORM_QUAD4,[3,6,7,4])
+        ###################
+        self.assertTrue(f1.getMesh().getCoords().isEqual(m2.getCoords(),1e-12))
+        f1.getMesh().setCoords(m2.getCoords())
+        #
+        f2=MEDCouplingFieldDouble(ON_GAUSS_PT)
+        f2.setMesh(m2)
+        for gt in m2.getAllGeoTypes(): # on recopie les localisation en utilisant f1
+            glt=f1.getGaussLocalizationIdOfOneType(gt)
+            gloc=f1.getGaussLocalization(glt)
+            f2.setGaussLocalizationOnType(gt,gloc.getRefCoords(),gloc.getGaussCoords(),gloc.getWeights())
+        arr2=DataArrayDouble(f2.getNumberOfTuplesExpected())
+        arr2[:]=0
+        f2.setArray(arr2)
+        f2.checkConsistencyLight()
+        #
+        fout1=MEDCouplingFieldDouble.MergeFields([f1,f2])
+        fout2=MEDCouplingFieldDouble.MergeFields(f1,f2)
+        #
+        fOut=MEDCouplingFieldDouble(ON_GAUSS_PT)
+        mOut=MEDCouplingUMesh.MergeUMeshes([f1.getMesh(),m2])
+        mOut.setName(f1.getMesh().getName())
+        fOut.setMesh(mOut)
+        for gt in f1.getMesh().getAllGeoTypes(): # on recopie les localisation en utilisant f1
+            glt=f1.getGaussLocalizationIdOfOneType(gt)
+            gloc=f1.getGaussLocalization(glt)
+            fOut.setGaussLocalizationOnType(gt,gloc.getRefCoords(),gloc.getGaussCoords(),gloc.getWeights())
+        fOut.setArray(DataArrayDouble.Aggregate([f1.getArray(),arr2]))
+        fOut.checkConsistencyLight()
+        fOut.setName(f1.getName())
+        fOut.getMesh().setName(f1.getMesh().getName())
+        #
+        self.assertTrue(fout1.isEqual(fOut,1e-12,1e-12))
+        self.assertTrue(fout2.isEqual(fOut,1e-12,1e-12))
+        pass
+
     pass
 
 if __name__ == '__main__':