Salome HOME
18fd7e8735e26f564b6b39c5c279af44b9f00a7f
[tools/medcoupling.git] / src / INTERP_KERNELTest / MEDMeshMaker.cxx
1 //  Copyright (C) 2007-2008  CEA/DEN, EDF R&D
2 //
3 //  This library is free software; you can redistribute it and/or
4 //  modify it under the terms of the GNU Lesser General Public
5 //  License as published by the Free Software Foundation; either
6 //  version 2.1 of the License.
7 //
8 //  This library is distributed in the hope that it will be useful,
9 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
10 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 //  Lesser General Public License for more details.
12 //
13 //  You should have received a copy of the GNU Lesser General Public
14 //  License along with this library; if not, write to the Free Software
15 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 //  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 #include "MEDMEM_Mesh.hxx"
20 #include "MEDMEM_Meshing.hxx"
21
22 MEDMEM::MESH* MEDMeshMaker(int dim, int nbedge, MED_EN::medGeometryElement type)
23 {
24   MEDMEM::MESHING* mesh=new MEDMEM::MESHING();
25   mesh->setSpaceDimension(dim);
26   int nbnodes;
27   int nbelems;
28   switch (dim)
29     {
30     case 2: 
31       nbnodes=(nbedge+1)*(nbedge+1);
32       if(type==MED_EN::MED_QUAD4)
33         nbelems=(nbedge*nbedge);
34       else
35         throw MEDMEM::MEDEXCEPTION("MEDMeshMaker: type not impletmented");
36       break;
37     case 3:
38       nbnodes=(nbedge+1)*(nbedge+1)*(nbedge+1);
39       if (type==MED_EN::MED_HEXA8)
40         nbelems= nbedge*nbedge*nbedge;
41       else
42         throw MEDMEM::MEDEXCEPTION("MEDMeshMaker: type not impletmented");
43       break;
44     }
45   double* coords = new double[dim*nbnodes];
46   int nz;
47   if (dim==2) nz =1; else nz=nbedge+1;
48   {
49     for (int ix=0; ix < nbedge+1; ix++)
50       for (int iy=0; iy<nbedge+1; iy++)
51         for (int iz=0; iz<nz;iz++)
52           {
53             int inode=(ix*(nbedge+1)*nz+iy*nz+iz);
54             coords[inode*dim]=double(ix)/double(nbedge);
55             coords[inode*dim+1]=double(iy)/double(nbedge);
56             if (dim==3)
57               coords[inode*dim+2]=double(iz)/double(nbedge);
58           }
59   }
60   mesh->setCoordinates(dim, nbnodes,coords,"CARTESIAN",MED_EN::MED_FULL_INTERLACE);
61   delete [] coords;
62   mesh->setNumberOfTypes(1,MED_EN::MED_CELL);
63   mesh->setTypes(&type,MED_EN::MED_CELL);
64   mesh->setNumberOfElements(&nbelems,MED_EN::MED_CELL);
65   
66   int* conn = new int [nbelems*(type%100)];
67   if (dim==2)
68     {
69       for (int ix=0; ix<nbedge; ix++)
70         for (int iy=0; iy<nbedge; iy++)
71           {
72             int ielem=(ix*nbedge+iy);
73             conn [ielem*4]=ix*(nbedge+1)+iy+1;
74             conn [ielem*4+1]=ix*(nbedge+1)+iy+1+1;
75             conn [ielem*4+2]=(ix+1)*(nbedge+1)+iy+1+1;
76             conn [ielem*4+3]=(ix+1)*(nbedge+1)+iy+1;                               
77           }
78     }
79   if (dim==3)
80     {
81       for (int ix=0; ix<nbedge; ix++)
82         for (int iy=0; iy<nbedge; iy++)
83           for (int iz=0; iz<nbedge; iz++)
84             {
85               int ielem=(ix*nbedge*nbedge+iy*nbedge+iz);
86               conn [ielem*8]=ix*(nbedge+1)*(nbedge+1)+iy*(nbedge+1)+iz+1;
87               conn [ielem*8+1]=(ix+1)*(nbedge+1)*(nbedge+1)+iy*(nbedge+1)+iz+1;
88               conn [ielem*8+2]=(ix+1)*(nbedge+1)*(nbedge+1)+(iy+1)*(nbedge+1)+iz+1;
89               conn [ielem*8+3]=ix*(nbedge+1)*(nbedge+1)+(iy+1)*(nbedge+1)+iz+1;
90               conn [ielem*8+4]=ix*(nbedge+1)*(nbedge+1)+iy*(nbedge+1)+iz+1+1;
91               conn [ielem*8+5]=(ix+1)*(nbedge+1)*(nbedge+1)+iy*(nbedge+1)+iz+1+1;
92               conn [ielem*8+6]=(ix+1)*(nbedge+1)*(nbedge+1)+(iy+1)*(nbedge+1)+iz+1+1;
93               conn [ielem*8+7]=ix*(nbedge+1)*(nbedge+1)+(iy+1)*(nbedge+1)+iz+1+1;
94             }
95     }
96   mesh->setConnectivity(conn, MED_EN::MED_CELL,type);
97   delete [] conn;
98   mesh->setMeshDimension(dim);
99   return mesh;
100 }