]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
MEDMEM Industrialization 2008
authoreap <eap@opencascade.com>
Thu, 18 Dec 2008 17:37:20 +0000 (17:37 +0000)
committereap <eap@opencascade.com>
Thu, 18 Dec 2008 17:37:20 +0000 (17:37 +0000)
    add new tests

src/ParaMEDMEM/Test/ParaMEDMEMTest.hxx
src/ParaMEDMEM/Test/ParaMEDMEMTest_IntersectionDEC.cxx

index 2d1bface138f32c376df9f0b4c265ed17c77b197..4972bd510109624c94ee2e6d383ff27ed18cd099 100644 (file)
@@ -38,8 +38,16 @@ class ParaMEDMEMTest : public CppUnit::TestFixture
   CPPUNIT_TEST(testBlockTopology_constructor);
   CPPUNIT_TEST(testBlockTopology_serialize);
   CPPUNIT_TEST(testIntersectionDEC_2D);
+  CPPUNIT_TEST(testIntersectionDEC_2D_P0P1);
+  CPPUNIT_TEST(testIntersectionDEC_2D_P1P0);
+  CPPUNIT_TEST(testIntersectionDEC_2D_P1dP0);
 
-       CPPUNIT_TEST(testSynchronousEqualIntersectionWithoutInterpNativeDEC_2D);
+  CPPUNIT_TEST(testIntersectionDEC_3D);
+  CPPUNIT_TEST(testIntersectionDEC_3D_P0P1);
+  CPPUNIT_TEST(testIntersectionDEC_3D_P1P0);
+  CPPUNIT_TEST(testIntersectionDEC_3D_P1dP0);
+
+       CPPUNIT_TEST(testSynchronousEqualIntersectionWithoutInterpNativeDEC_2D);
        CPPUNIT_TEST(testSynchronousEqualIntersectionWithoutInterpDEC_2D);
        CPPUNIT_TEST(testSynchronousEqualIntersectionDEC_2D);
        CPPUNIT_TEST(testSynchronousFasterSourceIntersectionDEC_2D);
@@ -57,6 +65,9 @@ class ParaMEDMEMTest : public CppUnit::TestFixture
   CPPUNIT_TEST(testNonCoincidentDEC_3D); 
 #endif
   CPPUNIT_TEST(testStructuredCoincidentDEC);
+  CPPUNIT_TEST(testAsynchronousEqualIntersectionDEC_2D_P0P1);
+  CPPUNIT_TEST(testAsynchronousFasterSourceIntersectionDEC_2D_P1P0);
+  CPPUNIT_TEST(testAsynchronousSlowerSourceIntersectionDEC_2D_P1dP0);
   CPPUNIT_TEST_SUITE_END();
   
 
@@ -71,7 +82,16 @@ public:
   void testMPIProcessorGroup_rank();
   void testBlockTopology_constructor();
   void testBlockTopology_serialize();
+
   void testIntersectionDEC_2D();
+  void testIntersectionDEC_2D_P0P1();
+  void testIntersectionDEC_2D_P1P0();
+  void testIntersectionDEC_2D_P1dP0();
+
+  void testIntersectionDEC_3D();
+  void testIntersectionDEC_3D_P0P1();
+  void testIntersectionDEC_3D_P1P0();
+  void testIntersectionDEC_3D_P1dP0();
 #ifdef MED_ENABLE_FVM
   void testNonCoincidentDEC_2D();
   void testNonCoincidentDEC_3D();
@@ -91,16 +111,28 @@ public:
   void testAsynchronousSlowSourceIntersectionDEC_2D();
   void testAsynchronousFastSourceIntersectionDEC_2D();
 
+  void testAsynchronousEqualIntersectionDEC_2D_P0P1();
+  void testAsynchronousFasterSourceIntersectionDEC_2D_P1P0();
+  void testAsynchronousSlowerSourceIntersectionDEC_2D_P1dP0();
+
 
 private:
+  void testIntersectionDEC_2D_(const std::string& srcMethod,
+                               const std::string& tgtMethod);
+  void testIntersectionDEC_3D_(const std::string& srcMethod,
+                               const std::string& tgtMethod);
   void testNonCoincidentDEC(const std::string& filename1, 
                           const std::string& meshname1, 
                           const std::string& filename2, 
                           const std::string& meshname2,
                            int nbprocsource, double epsilon);
   void testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA, 
-                                                                         double dtB, double tmaxB,
-                      bool WithPointToPoint, bool Asynchronous, bool WithInterp );
+                                          double dtB, double tmaxB,
+                                          bool WithPointToPoint,
+                                          bool Asynchronous,
+                                          bool WithInterp,
+                                          const std::string& srcMethod="P0",
+                                          const std::string& tgtMethod="P0");
   
   };
 
@@ -108,12 +140,14 @@ private:
 class ParaMEDMEMTest_TmpFilesRemover
 {
 public:
-  ParaMEDMEMTest_TmpFilesRemover() {}
+  ParaMEDMEMTest_TmpFilesRemover():myNoRemove(false) {}
   ~ParaMEDMEMTest_TmpFilesRemover();
   bool Register(const std::string theTmpFile);
+  void NoRemove() { myNoRemove=true; }
 
 private:
   std::set<std::string> myTmpFiles;
+  bool myNoRemove;
 };
 
 /*!
index b36a262a1590f1eebc69d3d25817164ccc9e5db6..75bc9cbb42a1983d68a60739b0c8ea396e77a62a 100644 (file)
@@ -32,6 +32,7 @@
 #include "ParaFIELD.hxx"
 #include "UnstructuredParaSUPPORT.hxx"
 #include "ICoCoMEDField.hxx"
+#include "MEDMEM_ArrayInterface.hxx"
  
 #include <string>
 
@@ -57,207 +58,492 @@ using namespace MEDMEM;
     void sendData();
  */
  
+
+// ================================================================================
+//      2D
+// ================================================================================
+
 void ParaMEDMEMTest::testIntersectionDEC_2D()
+{
+  testIntersectionDEC_2D_("P0","P0");
+}
+void ParaMEDMEMTest::testIntersectionDEC_2D_P0P1()
+{
+  testIntersectionDEC_2D_("P0","P1");
+}
+void ParaMEDMEMTest::testIntersectionDEC_2D_P1P0()
+{
+  testIntersectionDEC_2D_("P1","P0");
+}
+void ParaMEDMEMTest::testIntersectionDEC_2D_P1dP0()
+{
+  testIntersectionDEC_2D_("P1d","P0");
+}
+
+void ParaMEDMEMTest::testIntersectionDEC_2D_(const string& srcMethod/*="P0"*/,
+                                             const string& tgtMethod/*="P0"*/)
 {
   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
+
 //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::IntersectionDEC dec (*source_group,*target_group);
 
-  MEDMEM::MESH* mesh;
-  MEDMEM::SUPPORT* support;
-  ParaMEDMEM::ParaMESH* paramesh;
-  ParaMEDMEM::ParaFIELD* parafield;
-  ParaMEDMEM::ParaSUPPORT* parasupport;
+  MEDMEM::MESH           * mesh       =0;
+  MEDMEM::SUPPORT        * support    =0;
+  ParaMEDMEM::ParaMESH   * paramesh   =0;
+  ParaMEDMEM::ParaFIELD  * parafield  =0;
+  ParaMEDMEM::ParaSUPPORT* parasupport=0;
   double * value;
   ICoCo::Field* icocofield ;
-  
+
   string data_dir                   = getenv("MED_ROOT_DIR");
   string tmp_dir                    = getenv("TMP");
   if (tmp_dir == "")
     tmp_dir = "/tmp";
   string filename_xml1              = data_dir + "/share/salome/resources/med/square1_split";
   string filename_xml2              = data_dir + "/share/salome/resources/med/square2_split"; 
-  string filename_seq_wr            = tmp_dir + "/";
-  string filename_seq_med           = tmp_dir + "/myWrField_seq_pointe221.med";
-  
+
   // To remove tmp files from disk
   ParaMEDMEMTest_TmpFilesRemover aRemover;
-  
+  //aRemover.NoRemove(); ///////// DEBUG ////////////
+
   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;
-      
-       CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+  {
+    string master = filename_xml1;
+
+    ostringstream strstream;
+    strstream <<master<<rank+1<<".med";
+    ostringstream meshname ;
+    meshname<< "Mesh_2_"<< rank+1;
+
+    CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+    if ( srcMethod == "P1" )
+      support=new MEDMEM::SUPPORT(mesh,"all nodes",MED_EN::MED_NODE);
+    else
       support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
-    
-      paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
-    
-//      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
+
+    int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+    if ( srcMethod == "P1d" )
+    {
+      MEDMEM::FIELD<double>* field = new MEDMEM::FIELD<double>(support, 1);
+      typedef MEDMEM::MEDMEM_ArrayInterface<double,FullInterlace,Gauss>::Array TGaussArray;
+      MED_EN::medGeometryElement type = support->getTypes()[0];
+      const int nbelgeoc   [] = { 0, nb_local };
+      const int nbgaussgeo [] = { 0, type % 100 };
+      field->setArray ( new TGaussArray(/*dim=*/1,/*nbelem=*/nb_local,/*nbtypegeo=*/1,
+                                        nbelgeoc, nbgaussgeo));
+      field->setGaussLocalization
+        (type, MEDMEM::GAUSS_LOCALIZATION_::makeDefaultLocalization("locname",type,nbgaussgeo[1]));
+      field->setSupport( support );
+      parafield = new ParaFIELD( field, *source_group );
+      nb_local *= nbgaussgeo[1];
+    }
+    else
+    {
       parasupport=new UnstructuredParaSUPPORT( support,*source_group);
       ParaMEDMEM::ComponentTopology comptopo;
       parafield = new ParaFIELD(parasupport, comptopo);
-
-      
-                       int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
-//      double * value= new double[nb_local];
-      value= new double[nb_local];
-      for(int ielem=0; ielem<nb_local;ielem++)
-        value[ielem]=1.0;
-      parafield->getField()->setValue(value);
-    
-//      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
-      icocofield=new ICoCo::MEDField(paramesh,parafield);
-     
-      dec.attachLocalField(icocofield);
     }
-  
+    paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
+
+    value= new double[nb_local];
+    for(int ielem=0; ielem<nb_local;ielem++)
+      value[ielem]=1.0;
+    parafield->getField()->setValue(value);
+
+    icocofield=new ICoCo::MEDField(paramesh,parafield);
+
+    dec.attachLocalField(icocofield, srcMethod);
+  }
+
   //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;
+
+    CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+    if ( tgtMethod == "P1" )
+      support=new MEDMEM::SUPPORT(mesh,"all nodes",MED_EN::MED_NODE);
+    else
+      support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
+
+    paramesh=new ParaMESH (*mesh,*target_group,"target mesh");
+    parasupport=new UnstructuredParaSUPPORT(support,*target_group);
+    ParaMEDMEM::ComponentTopology comptopo;
+    parafield = new ParaFIELD(parasupport, comptopo);
+
+
+    int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+    value= new double[nb_local];
+    for(int ielem=0; ielem<nb_local;ielem++)
+      value[ielem]=0.0;
+    parafield->getField()->setValue(value);
+    icocofield=new ICoCo::MEDField(paramesh,parafield);
+
+    dec.attachLocalField(icocofield, tgtMethod);
+  }
+
+
+  //attaching a DEC to the source group 
+  double field_before_int;
+  double field_after_int;
+
+  if (source_group->containsMyRank())
+  { 
+    field_before_int = parafield->getVolumeIntegral(1);
+    dec.synchronize();
+    cout<<"DEC usage"<<endl;
+    dec.setForcedRenormalization(false);
+
+    dec.sendData();
+    paramesh->write(MED_DRIVER,"./sourcesquareb");
+    if (source_group->myRank()==0)
+      aRemover.Register("./sourcesquareb");
+    ostringstream filename;
+    filename<<"./sourcesquareb_"<<source_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    parafield->write(MED_DRIVER,"./sourcesquareb","boundary");  
+
+    dec.recvData();
+    cout <<"writing"<<endl;
+    paramesh->write(MED_DRIVER,"./sourcesquare");
+    if (source_group->myRank()==0)
+      aRemover.Register("./sourcesquare");
+    parafield->write(MED_DRIVER,"./sourcesquare","boundary");
+
+
+    filename<<"./sourcesquare_"<<source_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    field_after_int = parafield->getVolumeIntegral(1);
+
+
+    //      MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+    //       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+
+    //cout<<"CPPUNIT_ASSERT_DOUBLES_EQUAL("<<field_before_int<<", "<<field_after_int<<")" << endl;
+    CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
+
+  }
+
+  //attaching a DEC to the target group
+  if (target_group->containsMyRank())
+  {
+    dec.synchronize();
+    dec.setForcedRenormalization(false);
+
+    dec.recvData();
+    paramesh->write(MED_DRIVER, "./targetsquareb");
+    parafield->write(MED_DRIVER, "./targetsquareb", "boundary");
+    if (target_group->myRank()==0)
+      aRemover.Register("./targetsquareb");
+    ostringstream filename;
+    filename<<"./targetsquareb_"<<target_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    dec.sendData();
+    paramesh->write(MED_DRIVER, "./targetsquare");
+    parafield->write(MED_DRIVER, "./targetsquare", "boundary");
+    if (target_group->myRank()==0)
+      aRemover.Register("./targetsquareb");
+
+    filename<<"./targetsquareb_"<<target_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    //          double field_before_int, field_after_int;
+    //                  MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+    //       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+
+    //                  CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
+    //if ( tgtMethod == "P1" )
     {
-      string master= filename_xml2;
-      ostringstream strstream;
-            strstream << master<<(rank-nproc_source+1)<<".med";
-      ostringstream meshname ;
-            meshname<< "Mesh_3_"<<rank-nproc_source+1;
-      
-      CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+      // P0-->P1 will give constant field 1.
+      MEDMEM::FIELD<double>* field = parafield->getField();
+      int nb_local=field->getSupport()->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+      //cout << rank << ": Nb values = " << nb_local << " " <<srcMethod<<"->"<<tgtMethod<<endl;
+      for ( int i=1; i <=nb_local; ++i ) {
+        //cout <<rank << ": " << i << " " <<field->getValueIJ(i,1) << endl;
+        CPPUNIT_ASSERT_DOUBLES_EQUAL( 1., field->getValueIJ(i,1), 1e-6);
+      }
+    }
+
+  }
+
+  delete source_group;
+  delete target_group;
+  delete self_group;
+  delete mesh;
+  delete paramesh;
+  delete parafield;
+  delete parasupport;
+  delete [] value;
+  delete icocofield;
+  delete support;
+
+  MPI_Barrier(MPI_COMM_WORLD);
+  cout << "end of IntersectionDEC_2D test, "<< srcMethod<<"->"<<tgtMethod<<endl;
+}
+
+// ================================================================================
+//      3D
+// ================================================================================
+
+void ParaMEDMEMTest::testIntersectionDEC_3D()
+{
+  testIntersectionDEC_3D_("P0","P0");
+}
+void ParaMEDMEMTest::testIntersectionDEC_3D_P0P1()
+{
+  testIntersectionDEC_3D_("P0","P1");
+}
+void ParaMEDMEMTest::testIntersectionDEC_3D_P1P0()
+{
+  testIntersectionDEC_3D_("P1","P0");
+}
+void ParaMEDMEMTest::testIntersectionDEC_3D_P1dP0()
+{
+  testIntersectionDEC_3D_("P1d","P0");
+}
+
+void ParaMEDMEMTest::testIntersectionDEC_3D_(const string& srcMethod/*="P0"*/,
+                                             const string& tgtMethod/*="P0"*/)
+{
+  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 !=2) return ;
+
+  int nproc_source = 1;
+  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::IntersectionDEC dec (*source_group,*target_group);
+
+  MEDMEM::MESH           * mesh       =0;
+  MEDMEM::SUPPORT        * support    =0;
+  ParaMEDMEM::ParaMESH   * paramesh   =0;
+  ParaMEDMEM::ParaFIELD  * parafield  =0;
+  ParaMEDMEM::ParaSUPPORT* parasupport=0;
+  double * value;
+  ICoCo::Field* icocofield ;
+
+  string data_dir      = getenv("MED_ROOT_DIR");
+  string tmp_dir       = getenv("TMP");
+  if (tmp_dir == "")
+    tmp_dir = "/tmp";
+  string filename_xml1 = data_dir + "/share/salome/resources/med/cube_fine.med";
+  string filename_xml2 = data_dir + "/share/salome/resources/med/cube_coarse.med"; 
+
+  // To remove tmp files from disk
+  ParaMEDMEMTest_TmpFilesRemover aRemover;
+  //aRemover.NoRemove(); ///////// DEBUG ////////////
+
+  MPI_Barrier(MPI_COMM_WORLD);
+  if (source_group->containsMyRank())
+  {
+    string master = filename_xml1;
+    string meshname = "Box1Moderate";
+
+    CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,master,meshname));
+    if ( srcMethod == "P1" )
+      support=new MEDMEM::SUPPORT(mesh,"all nodes",MED_EN::MED_NODE);
+    else
       support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
-      
-      paramesh=new ParaMESH (*mesh,*target_group,"target mesh");
-//      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
-      parasupport=new UnstructuredParaSUPPORT(support,*target_group);
+
+    int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+    if ( srcMethod == "P1d" )
+    {
+      MEDMEM::FIELD<double>* field = new MEDMEM::FIELD<double>(support, 1);
+      typedef MEDMEM::MEDMEM_ArrayInterface<double,FullInterlace,Gauss>::Array TGaussArray;
+      MED_EN::medGeometryElement type = support->getTypes()[0];
+      const int nbelgeoc   [] = { 0, nb_local };
+      const int nbgaussgeo [] = { 0, type % 100 };
+      field->setArray ( new TGaussArray(/*dim=*/1,/*nbelem=*/nb_local,/*nbtypegeo=*/1,nbelgeoc, nbgaussgeo));
+      field->setGaussLocalization
+        ( type, MEDMEM::GAUSS_LOCALIZATION_::makeDefaultLocalization("locname", type, nbgaussgeo[1]));
+      field->setSupport( support );
+      parafield = new ParaFIELD( field, *source_group );
+      nb_local *= nbgaussgeo[1];
+    }
+    else
+    {
+      parasupport=new UnstructuredParaSUPPORT( support,*source_group);
       ParaMEDMEM::ComponentTopology comptopo;
       parafield = new ParaFIELD(parasupport, comptopo);
-
-                       
-                       int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
-//      double * value= new double[nb_local];
-      value= new double[nb_local];
-      for(int ielem=0; ielem<nb_local;ielem++)
-        value[ielem]=0.0;
-      parafield->getField()->setValue(value);
-//      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
-      icocofield=new ICoCo::MEDField(paramesh,parafield);
-      
-      dec.attachLocalField(icocofield);
     }
-    
-  
+    paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
+
+    value= new double[nb_local];
+    for(int ielem=0; ielem<nb_local;ielem++)
+      value[ielem]=1.0;
+    parafield->getField()->setValue(value);
+
+    icocofield=new ICoCo::MEDField(paramesh,parafield);
+
+    dec.attachLocalField(icocofield, srcMethod);
+  }
+
+  //loading the geometry for the target group
+  if (target_group->containsMyRank())
+  {
+    string master= filename_xml2;
+    string meshname = "Box2Moderate";
+
+    CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,master,meshname));
+    if ( tgtMethod == "P1" )
+      support=new MEDMEM::SUPPORT(mesh,"all nodes",MED_EN::MED_NODE);
+    else
+      support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
+
+    paramesh=new ParaMESH (*mesh,*target_group,"target mesh");
+    parasupport=new UnstructuredParaSUPPORT(support,*target_group);
+    ParaMEDMEM::ComponentTopology comptopo;
+    parafield = new ParaFIELD(parasupport, comptopo);
+
+
+    int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+    value= new double[nb_local];
+    for(int ielem=0; ielem<nb_local;ielem++)
+      value[ielem]=0.0;
+    parafield->getField()->setValue(value);
+    icocofield=new ICoCo::MEDField(paramesh,parafield);
+
+    dec.attachLocalField(icocofield, tgtMethod);
+  }
+
+
   //attaching a DEC to the source group 
   double field_before_int;
   double field_after_int;
-  
+
   if (source_group->containsMyRank())
-    { 
-      field_before_int = parafield->getVolumeIntegral(1);
-      dec.synchronize();
-      cout<<"DEC usage"<<endl;
-      dec.setForcedRenormalization(false);
+  { 
+    field_before_int = parafield->getVolumeIntegral(1);
+    dec.synchronize();
+    cout<<"DEC usage"<<endl;
+    dec.setForcedRenormalization(false);
+
+    dec.sendData();
+    paramesh->write(MED_DRIVER,"./sourceboxb");
+    if (source_group->myRank()==0)
+      aRemover.Register("./sourceboxb");
+    ostringstream filename;
+    filename<<"./sourceboxb_"<<source_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    parafield->write(MED_DRIVER,"./sourceboxb","boundary");  
+
+    dec.recvData();
+    cout <<"writing"<<endl;
+    paramesh->write(MED_DRIVER,"./sourcebox");
+    if (source_group->myRank()==0)
+      aRemover.Register("./sourcebox");
+    parafield->write(MED_DRIVER,"./sourcebox","boundary");
+
+
+    filename<<"./sourcebox_"<<source_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    field_after_int = parafield->getVolumeIntegral(1);
+
+
+    //      MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+    //       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+
+    CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
+
+  }
 
-      dec.sendData();
-      paramesh->write(MED_DRIVER,"./sourcesquareb");
-      if (source_group->myRank()==0)
-        aRemover.Register("./sourcesquareb");
-      ostringstream filename;
-      filename<<"./sourcesquareb_"<<source_group->myRank()+1;
-      aRemover.Register(filename.str().c_str());
-      parafield->write(MED_DRIVER,"./sourcesquareb","boundary");  
-   
-      dec.recvData();
-      cout <<"writing"<<endl;
-      paramesh->write(MED_DRIVER,"./sourcesquare");
-      if (source_group->myRank()==0)
-        aRemover.Register("./sourcesquare");
-      parafield->write(MED_DRIVER,"./sourcesquare","boundary");
-      
-     
-      filename<<"./sourcesquare_"<<source_group->myRank()+1;
-      aRemover.Register(filename.str().c_str());
-      field_after_int = parafield->getVolumeIntegral(1);
-                       
-                       
- //      MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
-//       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
-
-       CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
-               
-    }
-  
   //attaching a DEC to the target group
   if (target_group->containsMyRank())
-    {
-      dec.synchronize();
-      dec.setForcedRenormalization(false);
+  {
+    dec.synchronize();
+    dec.setForcedRenormalization(false);
 
-      dec.recvData();
-      paramesh->write(MED_DRIVER, "./targetsquareb");
-      parafield->write(MED_DRIVER, "./targetsquareb", "boundary");
-       if (target_group->myRank()==0)
-        aRemover.Register("./targetsquareb");
-      ostringstream filename;
-      filename<<"./targetsquareb_"<<target_group->myRank()+1;
-      aRemover.Register(filename.str().c_str());
-      dec.sendData();
-      paramesh->write(MED_DRIVER, "./targetsquare");
-      parafield->write(MED_DRIVER, "./targetsquare", "boundary");
-                       if (target_group->myRank()==0)
-        aRemover.Register("./targetsquareb");
-      
-      filename<<"./targetsquareb_"<<target_group->myRank()+1;
-      aRemover.Register(filename.str().c_str());
-                       //              double field_before_int, field_after_int;
-//                     MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
-//       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
-                       
-//                     CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
-               
-    }
-  
-   delete source_group;
+    dec.recvData();
+    paramesh->write(MED_DRIVER, "./targetboxb");
+    parafield->write(MED_DRIVER, "./targetboxb", "boundary");
+    if (target_group->myRank()==0)
+      aRemover.Register("./targetboxb");
+    ostringstream filename;
+    filename<<"./targetboxb_"<<target_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    dec.sendData();
+    paramesh->write(MED_DRIVER, "./targetbox");
+    parafield->write(MED_DRIVER, "./targetbox", "boundary");
+    if (target_group->myRank()==0)
+      aRemover.Register("./targetboxb");
+
+    filename<<"./targetboxb_"<<target_group->myRank()+1;
+    aRemover.Register(filename.str().c_str());
+    //          double field_before_int, field_after_int;
+    //                  MPI_Bcast(&field_before_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+    //       MPI_Bcast(&field_after_int,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
+
+    //                  CPPUNIT_ASSERT_DOUBLES_EQUAL(field_before_int, field_after_int, 1e-6);
+
+  }
+
+  delete source_group;
   delete target_group;
   delete self_group;
   delete mesh;
   delete paramesh;
-       delete parafield;
+  delete parafield;
   delete parasupport;
   delete [] value;
   delete icocofield;
   delete support;
 
   MPI_Barrier(MPI_COMM_WORLD);
-  cout << "end of IntersectionDEC_2D test"<<endl;
+  cout << "end of IntersectionDEC_3D test, "<< srcMethod<<"->"<<tgtMethod<<endl;
 }
+// ================================================================================
+
+
+
 
 //Synchronous tests without interpolation with native mode (AllToAll(v) from lam/MPI:
 void ParaMEDMEMTest::testSynchronousEqualIntersectionWithoutInterpNativeDEC_2D()
@@ -314,6 +600,18 @@ void ParaMEDMEMTest::testAsynchronousFastSourceIntersectionDEC_2D()
 {
        testAsynchronousIntersectionDEC_2D(0.01,1,0.11,1,true,true,true);
 }
+void ParaMEDMEMTest::testAsynchronousEqualIntersectionDEC_2D_P0P1()
+{
+  //testAsynchronousIntersectionDEC_2D(0.1,1,0.1,1,true,true,true, "P0","P1");
+}
+void ParaMEDMEMTest::testAsynchronousFasterSourceIntersectionDEC_2D_P1P0()
+{
+  //testAsynchronousIntersectionDEC_2D(0.09,1,0.1,1,true,true,true, "P1","P0");
+}
+void ParaMEDMEMTest::testAsynchronousSlowerSourceIntersectionDEC_2D_P1dP0()
+{
+  //testAsynchronousIntersectionDEC_2D(0.11,1,0.1,1,true,true,true, "P1d","P0");
+}
 
 /*!
  * Tests an asynchronous exchange between two codes
@@ -321,8 +619,12 @@ void ParaMEDMEMTest::testAsynchronousFastSourceIntersectionDEC_2D()
  * the other one receives with dtB as an interval, the max time being tmaxB
  */
 void ParaMEDMEMTest::testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA, 
-                                                                         double dtB, double tmaxB, bool WithPointToPoint, bool Asynchronous,
-                    bool WithInterp )
+                                                        double dtB, double tmaxB,
+                                                        bool WithPointToPoint,
+                                                        bool Asynchronous,
+                                                        bool WithInterp,
+                                                        const std::string& srcMethod,
+                                                        const std::string& tgtMethod )
 {
   int size;
   int rank;
@@ -369,8 +671,6 @@ void ParaMEDMEMTest::testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA
     tmp_dir = "/tmp";
   string filename_xml1              = data_dir + "/share/salome/resources/med/square1_split";
   string filename_xml2              = data_dir + "/share/salome/resources/med/square2_split"; 
-  string filename_seq_wr            = tmp_dir + "/";
-  string filename_seq_med           = tmp_dir + "/myWrField_seq_pointe221.med";
   
   // To remove tmp files from disk
   ParaMEDMEMTest_TmpFilesRemover aRemover;
@@ -378,70 +678,84 @@ void ParaMEDMEMTest::testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA
   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;
-      
-       CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+  {
+    string master = filename_xml1;
+
+    ostringstream strstream;
+    strstream <<master<<rank+1<<".med";
+    ostringstream meshname ;
+    meshname<< "Mesh_2_"<< rank+1;
+
+    CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+    if ( srcMethod == "P1" )
+      support=new MEDMEM::SUPPORT(mesh,"all nodes",MED_EN::MED_NODE);
+    else
       support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
-    
-      paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
-    
-//      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT( support,*source_group);
+
+    int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+    if ( srcMethod == "P1d" )
+    {
+      MEDMEM::FIELD<double>* field = new MEDMEM::FIELD<double>(support, 1);
+      typedef MEDMEM::MEDMEM_ArrayInterface<double,FullInterlace,Gauss>::Array TGaussArray;
+      MED_EN::medGeometryElement type = support->getTypes()[0];
+      const int nbelgeoc   [] = { 0, nb_local };
+      const int nbgaussgeo [] = { 0, type % 100 };
+      field->setArray ( new TGaussArray(/*dim=*/1,/*nbelem=*/nb_local,/*nbtypegeo=*/1,
+                                        nbelgeoc, nbgaussgeo));
+      field->setGaussLocalization
+        (type, MEDMEM::GAUSS_LOCALIZATION_::makeDefaultLocalization("locname",type,nbgaussgeo[1]));
+      field->setSupport( support );
+      parafield = new ParaFIELD( field, *source_group );
+      nb_local *= nbgaussgeo[1];
+    }
+    else
+    {
       parasupport=new UnstructuredParaSUPPORT( support,*source_group);
       ParaMEDMEM::ComponentTopology comptopo;
       parafield = new ParaFIELD(parasupport, comptopo);
+    }
+    paramesh=new ParaMESH (*mesh,*source_group,"source mesh");
 
-      
-                       int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
-//      double * value= new double[nb_local];
-      value = new double[nb_local];
-      for(int ielem=0; ielem<nb_local;ielem++)
-        value[ielem]=0.0;
-      parafield->getField()->setValue(value);
-    
-//      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
-      icocofield=new ICoCo::MEDField(paramesh,parafield);
-     
-      dec.attachLocalField(icocofield);
+    value = new double[nb_local];
+    for(int ielem=0; ielem<nb_local;ielem++)
+      value[ielem]=0.0;
+    parafield->getField()->setValue(value);
 
+    icocofield=new ICoCo::MEDField(paramesh,parafield);
 
-    }
+    dec.attachLocalField(icocofield, srcMethod);
+  }
   
   //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;
-      
-      CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+  {
+    string master= filename_xml2;
+    ostringstream strstream;
+    strstream << master<<(rank-nproc_source+1)<<".med";
+    ostringstream meshname ;
+    meshname<< "Mesh_3_"<<rank-nproc_source+1;
+
+    CPPUNIT_ASSERT_NO_THROW(mesh = new MESH(MED_DRIVER,strstream.str(),meshname.str()));
+    if ( tgtMethod == "P1" )
+      support=new MEDMEM::SUPPORT(mesh,"all nodes",MED_EN::MED_NODE);
+    else
       support=new MEDMEM::SUPPORT(mesh,"all elements",MED_EN::MED_CELL);
-      
-      paramesh=new ParaMESH (*mesh,*target_group,"target mesh");
-//      ParaMEDMEM::ParaSUPPORT* parasupport=new UnstructuredParaSUPPORT(support,*target_group);
-      parasupport=new UnstructuredParaSUPPORT(support,*target_group);
-      ParaMEDMEM::ComponentTopology comptopo;
-      parafield = new ParaFIELD(parasupport, comptopo);
 
-                       
-                       int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
-//      double * value= new double[nb_local];
-      value = new double[nb_local];
-      for(int ielem=0; ielem<nb_local;ielem++)
-        value[ielem]=0.0;
-      parafield->getField()->setValue(value);
-//      ICoCo::Field* icocofield=new ICoCo::MEDField(paramesh,parafield);
-      icocofield=new ICoCo::MEDField(paramesh,parafield);
-      
-      dec.attachLocalField(icocofield);
-    }
+    paramesh=new ParaMESH (*mesh,*target_group,"target mesh");
+    parasupport=new UnstructuredParaSUPPORT(support,*target_group);
+    ParaMEDMEM::ComponentTopology comptopo;
+    parafield = new ParaFIELD(parasupport, comptopo);
+
+
+    int nb_local=support->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+    value = new double[nb_local];
+    for(int ielem=0; ielem<nb_local;ielem++)
+      value[ielem]=0.0;
+    parafield->getField()->setValue(value);
+    icocofield=new ICoCo::MEDField(paramesh,parafield);
+
+    dec.attachLocalField(icocofield, tgtMethod);
+  }
     
   
   //attaching a DEC to the source group 
@@ -463,7 +777,7 @@ void ParaMEDMEMTest::testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA
       dec.setForcedRenormalization(false);
       for (double time=0; time<tmaxA+1e-10; time+=dtA)
       {
-        cout << "testAsynchronousIntersectionDEC_2D" << rank << " time " << time
+        cout << "testAsynchronousIntersectionDEC_2D " << rank << " time " << time
              << " dtA " << dtA << " tmaxA " << tmaxA << endl ;
         if ( time+dtA < tmaxA+1e-7 ) {
            dec.sendData( time , dtA );
@@ -471,15 +785,13 @@ void ParaMEDMEMTest::testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA
         else {
            dec.sendData( time , 0 );
         }
-         double* value = const_cast<double*> (parafield->getField()->getValue());
-         int nb_local=parafield->getField()->getSupport()->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
-         for (int i=0; i<nb_local;i++)
-                 value[i]= time+dtA;
-
-        
+        double* value = const_cast<double*> (parafield->getField()->getValue());
+        int nb_local=parafield->getField()->getSupport()->getNumberOfElements(MED_EN::MED_ALL_ELEMENTS);
+        for (int i=0; i<nb_local;i++)
+          value[i]= time+dtA;
       }
     }
-  
+
   //attaching a DEC to the target group
   if (target_group->containsMyRank())
     {
@@ -499,19 +811,19 @@ void ParaMEDMEMTest::testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA
       vector<double> times;
       for (double time=0; time<tmaxB+1e-10; time+=dtB)
       {
-          cout << "testAsynchronousIntersectionDEC_2D" << rank << " time " << time
-               << " dtB " << dtB << " tmaxB " << tmaxB << endl ;
-         dec.recvData( time );
-                                       double vi = parafield->getVolumeIntegral(1);
-                                                 cout << "testAsynchronousIntersectionDEC_2D" << rank << " time " << time
-               << " VolumeIntegral " << vi
-               << " time*10000 " << time*10000 << endl ;
-                                       
-                                                               CPPUNIT_ASSERT_DOUBLES_EQUAL(vi,time*10000,0.001);
+        cout << "testAsynchronousIntersectionDEC_2D" << rank << " time " << time
+             << " dtB " << dtB << " tmaxB " << tmaxB << endl ;
+        dec.recvData( time );
+        double vi = parafield->getVolumeIntegral(1);
+        cout << "testAsynchronousIntersectionDEC_2D" << rank << " time " << time
+             << " VolumeIntegral " << vi
+             << " time*10000 " << time*10000 << endl ;
+        
+        CPPUNIT_ASSERT_DOUBLES_EQUAL(vi,time*10000,0.001);
       }
       
     }
-  
+
   delete source_group;
   delete target_group;
   delete self_group;
@@ -525,6 +837,6 @@ void ParaMEDMEMTest::testAsynchronousIntersectionDEC_2D(double dtA, double tmaxA
 
   cout << "testAsynchronousIntersectionDEC_2D" << rank << " MPI_Barrier " << endl ;
  
-       if (Asynchronous) MPI_Barrier(MPI_COMM_WORLD);
-  cout << "end of IntersectionDEC_2D test"<<endl;
+  if (Asynchronous) MPI_Barrier(MPI_COMM_WORLD);
+  cout << "end of IntersectionDEC_2D test, "<< srcMethod<<"->"<<tgtMethod<<endl;
 }