]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
*** empty log message *** PARAVIS_29012010
authorageay <ageay>
Wed, 27 Jan 2010 14:56:34 +0000 (14:56 +0000)
committerageay <ageay>
Wed, 27 Jan 2010 14:56:34 +0000 (14:56 +0000)
src/ParaMEDMEM/DEC.cxx
src/ParaMEDMEM/ParaFIELD.cxx
src/ParaMEDMEM/ParaFIELD.hxx
src/ParaMEDMEM/Test/ParaMEDMEMTest.hxx
src/ParaMEDMEM/Test/ParaMEDMEMTest_InterpKernelDEC.cxx

index 8bb1f2d27fe3639d3d6f8058e9b82564c7a0397b..fe0d2f87442831b6556d28f969d000a800e4074f 100644 (file)
@@ -21,6 +21,7 @@
 #include "BlockTopology.hxx"
 #include "ComponentTopology.hxx"
 #include "ParaFIELD.hxx"
+#include "ParaMESH.hxx"
 #include "DEC.hxx"
 #include "ICoCoField.hxx"
 #include "ICoCoMEDField.hxx"
@@ -128,10 +129,12 @@ namespace ParaMEDMEM
       local_group=_target_group;
     else
       throw INTERP_KERNEL::Exception("Invalid procgroup for field attachment to DEC");
-        
-    _local_field= new ParaFIELD(field, *local_group);
+    ParaMESH *paramesh=new ParaMESH((MEDCouplingPointSet *)field->getMesh(),*local_group,field->getMesh()->getName());
+    _local_field = new ParaFIELD(field, paramesh, *local_group);
     _owns_field=true;
-    _comm_interface=&(local_group->getCommInterface());
+    _local_field->setOwnSupport(true);
+    attachLocalField(_local_field);
+    //_comm_interface=&(local_group->getCommInterface());
   }
 
   /*! 
index 10fc964aaaaa0ec922b7f20d5fbfcab3e72bdb9f..a0fb5695bc3f201c5e671009a9fdd027a5e29b65 100644 (file)
@@ -64,7 +64,7 @@ namespace ParaMEDMEM
   */
   ParaFIELD::ParaFIELD(TypeOfField type, TypeOfTimeDiscretization td, ParaMESH* para_support, const ComponentTopology& component_topology)
     :_field(0),
-     _component_topology(component_topology),_topology(0),
+     _component_topology(component_topology),_topology(0),_own_support(false),
      _support(para_support)
   {
     if (para_support->isStructured() || (para_support->getTopology()->getProcGroup()->size()==1 && component_topology.nbBlocks()!=1))
@@ -105,21 +105,23 @@ namespace ParaMEDMEM
     This constructor supposes that support underlying \a subdomain_field has no ParaSUPPORT 
     attached and it therefore recreates one. It therefore takes ownership over _support. The component topology associated with the field is a basic one (all components on the same processor). 
   */
-  ParaFIELD::ParaFIELD(MEDCouplingFieldDouble* subdomain_field, const ProcessorGroup& proc_group):
+  ParaFIELD::ParaFIELD(MEDCouplingFieldDouble* subdomain_field, ParaMESH *sup, const ProcessorGroup& proc_group):
     _field(subdomain_field),
-    _component_topology(ComponentTopology(_field->getNumberOfComponents())),_topology(0),
-    _support()
+    _component_topology(ComponentTopology(_field->getNumberOfComponents())),_topology(0),_own_support(false),
+    _support(sup)
   {
     if(_field)
       _field->incrRef();
-    const ExplicitTopology* source_topo=dynamic_cast<const ExplicitTopology*> (_support->getTopology());
-    _topology=new ExplicitTopology(*source_topo,_component_topology.nbLocalComponents());
+    const BlockTopology* source_topo=dynamic_cast<const BlockTopology*> (_support->getTopology());
+    _topology=new BlockTopology(*source_topo,_component_topology.nbLocalComponents());
   }
 
   ParaFIELD::~ParaFIELD()
   {
     if(_field)
       _field->decrRef();
+    if(_own_support)
+      delete _support;
     delete _topology;
   }
 
index 3bd12d4767a8cb1eaab012f7da77d8e429237d46..08229a8372b4b8567595745286aaa66813e7cf86 100644 (file)
@@ -20,6 +20,7 @@
 #define __PARAFIELD_HXX__
 
 #include "MEDCouplingRefCountObject.hxx"
+#include "ComponentTopology.hxx"
 
 namespace ParaMEDMEM
 {
@@ -37,12 +38,13 @@ namespace ParaMEDMEM
     ParaFIELD(TypeOfField type, TypeOfTimeDiscretization td, ParaMESH* mesh, const ComponentTopology& component_topology); 
 
 
-    ParaFIELD(MEDCouplingFieldDouble* field, const ProcessorGroup& group);
+    ParaFIELD(MEDCouplingFieldDouble* field, ParaMESH *sup, const ProcessorGroup& group);
   
     virtual ~ParaFIELD();
     void synchronizeTarget( ParaMEDMEM::ParaFIELD* source_field);
     void synchronizeSource( ParaMEDMEM::ParaFIELD* target_field);
     MEDCouplingFieldDouble* getField() const { return _field; }
+    void setOwnSupport(bool v) const { _own_support=v; }
     DataArrayInt* returnCumulativeGlobalNumbering() const;
     DataArrayInt* returnGlobalNumbering() const;
     Topology* getTopology() const { return _topology; }
@@ -52,9 +54,9 @@ namespace ParaMEDMEM
     double getL2Norm()const { return -1; }
   private:
     MEDCouplingFieldDouble* _field;
-    const  ParaMEDMEM::ComponentTopology& _component_topology;
+    ParaMEDMEM::ComponentTopology _component_topology;
     Topology* _topology; 
-
+    mutable bool _own_support;
     ParaMESH* _support;
   };
 
index 3562684f3b17e6416c05aeb0fe904267b6a4ec05..140911a42046ec9ea6f1a772234c9a5fd1201e07 100644 (file)
@@ -36,6 +36,7 @@ class ParaMEDMEMTest : public CppUnit::TestFixture
   CPPUNIT_TEST(testBlockTopology_constructor);
   CPPUNIT_TEST(testBlockTopology_serialize);
   CPPUNIT_TEST(testInterpKernelDEC_2D);
+  CPPUNIT_TEST(testInterpKernelDEC2_2D);
   CPPUNIT_TEST(testInterpKernelDEC_2DP0P1);
   CPPUNIT_TEST(testInterpKernelDEC_3D);
   CPPUNIT_TEST(testInterpKernelDECNonOverlapp_2D_P0P0);
@@ -83,6 +84,7 @@ public:
   void testBlockTopology_constructor();
   void testBlockTopology_serialize();
   void testInterpKernelDEC_2D();
+  void testInterpKernelDEC2_2D();
   void testInterpKernelDEC_2DP0P1();
   void testInterpKernelDEC_3D();
   void testInterpKernelDECNonOverlapp_2D_P0P0();
@@ -130,6 +132,7 @@ private:
                                           double dtB, double tmaxB,
                                           bool WithPointToPoint, bool Asynchronous, bool WithInterp, const char *srcMeth, const char *targetMeth);
   void testInterpKernelDEC_2D_(const char *srcMeth, const char *targetMeth);
+  void testInterpKernelDEC2_2D_(const char *srcMeth, const char *targetMeth);
   void testInterpKernelDEC_3D_(const char *srcMeth, const char *targetMeth);
 };
 
index ed54396b6ec15edc87bb8ea7a761037d37e433a7..284a027871ea085d25a0b6b51dc4ca089a5e218d 100644 (file)
@@ -50,6 +50,11 @@ void ParaMEDMEMTest::testInterpKernelDEC_2D()
   testInterpKernelDEC_2D_("P0","P0");
 }
 
+void ParaMEDMEMTest::testInterpKernelDEC2_2D()
+{
+  testInterpKernelDEC2_2D_("P0","P0");
+}
+
 void ParaMEDMEMTest::testInterpKernelDEC_3D()
 {
   testInterpKernelDEC_3D_("P0","P0");
@@ -277,6 +282,161 @@ void ParaMEDMEMTest::testInterpKernelDEC_2D_(const char *srcMeth, const char *ta
   cout << "end of InterpKernelDEC_2D test"<<endl;
 }
 
+void ParaMEDMEMTest::testInterpKernelDEC2_2D_(const char *srcMeth, const char *targetMeth)
+{
+  std::string srcM(srcMeth);
+  std::string targetM(targetMeth);
+  int size;
+  int rank;
+  MPI_Comm_size(MPI_COMM_WORLD,&size);
+  MPI_Comm_rank(MPI_COMM_WORLD,&rank);
+
+  //the test is meant to run on five processors
+  if (size !=5) return ;
+   
+  int nproc_source = 3;
+  set<int> self_procs;
+  set<int> procs_source;
+  set<int> procs_target;
+  
+  for (int i=0; i<nproc_source; i++)
+    procs_source.insert(i);
+  for (int i=nproc_source; i<size; i++)
+    procs_target.insert(i);
+  self_procs.insert(rank);
+  
+  ParaMEDMEM::CommInterface interface;
+    
+  ParaMEDMEM::ProcessorGroup* self_group = new ParaMEDMEM::MPIProcessorGroup(interface,self_procs);
+  ParaMEDMEM::ProcessorGroup* target_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_target);
+  ParaMEDMEM::ProcessorGroup* source_group = new ParaMEDMEM::MPIProcessorGroup(interface,procs_source);
+  
+  //loading the geometry for the source group
+
+  ParaMEDMEM::InterpKernelDEC dec (*source_group,*target_group);
+
+  ParaMEDMEM::MEDCouplingUMesh* mesh;
+  ParaMEDMEM::MEDCouplingFieldDouble* mcfield;
+  
+  string filename_xml1              = getResourceFile("square1_split");
+  string filename_xml2              = getResourceFile("square2_split");
+  
+  // To remove tmp files from disk
+  ParaMEDMEMTest_TmpFilesRemover aRemover;
+  
+  MPI_Barrier(MPI_COMM_WORLD);
+  if (source_group->containsMyRank())
+    {
+      string master = filename_xml1;
+      
+      ostringstream strstream;
+      strstream <<master<<rank+1<<".med";
+      ostringstream meshname ;
+      meshname<< "Mesh_2_"<< rank+1;
+      
+      mesh=MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
+      ParaMEDMEM::ComponentTopology comptopo;
+      if(srcM=="P0")
+        {
+          mcfield = MEDCouplingFieldDouble::New(ON_CELLS,NO_TIME);
+          mcfield->setMesh(mesh);
+          DataArrayDouble *array=DataArrayDouble::New();
+          array->alloc(mcfield->getNumberOfTuples(),1);
+          mcfield->setArray(array);
+          array->decrRef();
+          mcfield->setNature(ConservativeVolumic);
+        }
+      else
+        {
+          mcfield = MEDCouplingFieldDouble::New(ON_CELLS,NO_TIME);
+          mcfield->setMesh(mesh);
+          DataArrayDouble *array=DataArrayDouble::New();
+          array->alloc(mcfield->getNumberOfTuples(),1);
+          mcfield->setArray(array);
+          array->decrRef();
+        }
+      int nb_local;
+      if(srcM=="P0")
+        nb_local=mesh->getNumberOfCells();
+      else
+        nb_local=mesh->getNumberOfNodes();
+      double *value=mcfield->getArray()->getPointer();
+      for(int ielem=0; ielem<nb_local;ielem++)
+        value[ielem]=1.0;
+      dec.setMethod(srcMeth);
+      dec.attachLocalField(mcfield);
+    }
+  
+  //loading the geometry for the target group
+  if (target_group->containsMyRank())
+    {
+      string master= filename_xml2;
+      ostringstream strstream;
+      strstream << master<<(rank-nproc_source+1)<<".med";
+      ostringstream meshname ;
+      meshname<< "Mesh_3_"<<rank-nproc_source+1;
+      mesh = MEDLoader::ReadUMeshFromFile(strstream.str().c_str(),meshname.str().c_str(),0);
+      ParaMEDMEM::ComponentTopology comptopo;
+      if(targetM=="P0")
+        {
+          mcfield = MEDCouplingFieldDouble::New(ON_CELLS,NO_TIME);
+          mcfield->setMesh(mesh);
+          DataArrayDouble *array=DataArrayDouble::New();
+          array->alloc(mcfield->getNumberOfTuples(),1);
+          mcfield->setArray(array);
+          array->decrRef();
+          mcfield->setNature(ConservativeVolumic);
+        }
+      else
+        {
+          mcfield = MEDCouplingFieldDouble::New(ON_NODES,NO_TIME);
+          mcfield->setMesh(mesh);
+          DataArrayDouble *array=DataArrayDouble::New();
+          array->alloc(mcfield->getNumberOfTuples(),1);
+          mcfield->setArray(array);
+          array->decrRef();
+        }
+      int nb_local;
+      if(targetM=="P0")
+        nb_local=mesh->getNumberOfCells();
+      else
+        nb_local=mesh->getNumberOfNodes();
+      double *value=mcfield->getArray()->getPointer();
+      for(int ielem=0; ielem<nb_local;ielem++)
+        value[ielem]=0.0;
+      dec.setMethod(targetMeth);
+      dec.attachLocalField(mcfield);
+    }
+    
+  
+  //attaching a DEC to the source group 
+
+  if (source_group->containsMyRank())
+    { 
+      dec.synchronize();
+      dec.setForcedRenormalization(false);
+      dec.sendData();
+      dec.recvData();
+    }
+  
+  //attaching a DEC to the target group
+  if (target_group->containsMyRank())
+    {
+      dec.synchronize();
+      dec.setForcedRenormalization(false);
+      dec.recvData();
+      dec.sendData();
+    }
+  delete source_group;
+  delete target_group;
+  delete self_group;
+  mcfield->decrRef();
+  mesh->decrRef();
+
+  MPI_Barrier(MPI_COMM_WORLD);
+  cout << "end of InterpKernelDEC2_2D test"<<endl;
+}
+
 void ParaMEDMEMTest::testInterpKernelDEC_3D_(const char *srcMeth, const char *targetMeth)
 {
   std::string srcM(srcMeth);
@@ -563,7 +723,6 @@ void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P0()
   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
   ParaMEDMEM::ParaMESH *paramesh=0;
   ParaMEDMEM::ParaFIELD* parafield=0;
-  ICoCo::Field* icocofield=0;
   //
   ParaMEDMEM::CommInterface interface;
   //
@@ -793,7 +952,6 @@ void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
   ParaMEDMEM::MEDCouplingUMesh *mesh=0;
   ParaMEDMEM::ParaMESH *paramesh=0;
   ParaMEDMEM::ParaFIELD *parafieldP0=0,*parafieldP1=0;
-  ICoCo::Field* icocofield=0;
   //
   ParaMEDMEM::CommInterface interface;
   //
@@ -808,7 +966,7 @@ void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
         {
           double coords[6]={-0.3,-0.3, 0.7,0.7, 0.7,-0.3};
           int conn[3]={0,1,2};
-          int globalNode[3]={1,2,0};
+          //int globalNode[3]={1,2,0};
           mesh=MEDCouplingUMesh::New("Source mesh Proc0",2);
           mesh->allocateCells(1);
           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
@@ -823,7 +981,7 @@ void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
         {
           double coords[6]={-0.3,-0.3, -0.3,0.7, 0.7,0.7};
           int conn[3]={0,1,2};
-          int globalNode[3]={1,3,2};
+          //int globalNode[3]={1,3,2};
           mesh=MEDCouplingUMesh::New("Source mesh Proc1",2);
           mesh->allocateCells(1);
           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
@@ -860,7 +1018,7 @@ void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
         {
           double coords[10]={-0.3,-0.3, 0.2,-0.3, 0.7,-0.3, -0.3,0.2, 0.2,0.2 };
           int conn[7]={0,3,4,1, 1,4,2};
-          int globalNode[5]={4,3,0,2,1};
+          //int globalNode[5]={4,3,0,2,1};
           mesh=MEDCouplingUMesh::New("Target mesh Proc2",2);
           mesh->allocateCells(2);
           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
@@ -882,7 +1040,7 @@ void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
         {
           double coords[6]={0.2,0.2, 0.7,-0.3, 0.7,0.2};
           int conn[3]={0,2,1};
-          int globalNode[3]={1,0,5};
+          //int globalNode[3]={1,0,5};
           mesh=MEDCouplingUMesh::New("Target mesh Proc3",2);
           mesh->allocateCells(1);
           mesh->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
@@ -903,7 +1061,7 @@ void ParaMEDMEMTest::testInterpKernelDECNonOverlapp_2D_P0P1P1P0()
         {
           double coords[12]={-0.3,0.2, -0.3,0.7, 0.2,0.7, 0.2,0.2, 0.7,0.7, 0.7,0.2};
           int conn[8]={0,1,2,3, 3,2,4,5};
-          int globalNode[6]={2,6,7,1,8,5};
+          //int globalNode[6]={2,6,7,1,8,5};
           mesh=MEDCouplingUMesh::New("Target mesh Proc4",2);
           mesh->allocateCells(2);
           mesh->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
@@ -1037,7 +1195,6 @@ void ParaMEDMEMTest::testAsynchronousInterpKernelDEC_2D(double dtA, double tmaxA
   ParaMEDMEM::ParaMESH* paramesh;
   ParaMEDMEM::ParaFIELD* parafield;
   
-  double * value ;
   ICoCo::Field* icocofield ;
 
   string tmp_dir                    = getenv("TMP");