src/INTERPOLATION/MEDMEM_MappingTools.hxx

Go to the documentation of this file.
00001 // Copyright (C) 2005  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
00002 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
00003 // 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either 
00007 // version 2.1 of the License.
00008 // 
00009 // This library is distributed in the hope that it will be useful 
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
00012 // Lesser General Public License for more details.
00013 //
00014 // You should have received a copy of the GNU Lesser General Public  
00015 // License along with this library; if not, write to the Free Software 
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00017 //
00018 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00019 //
00020 #ifndef COORDONNEES_BARYCENTRIQUES_HPP
00021 #define COORDONNEES_BARYCENTRIQUES_HPP
00022 
00023 #define _TEMPLATE_SPE_ template <class NUAGEMAILLE, class NUAGENOEUD, class NOEUD>
00024 #define _COORDBARYC_ Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,DIMENSION>
00025 #define _COORDBARY_2D_ Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,2>
00026 #define _COORDBARY_3D_ Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,3>
00027 
00033 
00034 /*********************************************************/
00035 /*                                                       */
00036 /*           Classe Coordonnees_Barycentriques           */
00037 /*                                                       */
00038 /*********************************************************/
00039 
00040 // C'est la définition de la classe générique template qui n'est utilisée qu'à travers ses spécialisations
00041 // vu que le nombre de spécialisations a faire est petit (nombre de dimensions utilisée, 3 pour MED)
00042 // et vu que l'acces a ces classes doit etre rapide, car ce sont des classes de calcul
00043 // la technique de spécialisation, plus lourde, mais plus rapide, a été préférée aux techniques d'héritage
00044 
00045 template <class NUAGEMAILLE, class NUAGENOEUD, class NOEUD, int DIMENSION> class Coordonnees_Barycentriques
00046 {
00047 // TEMPLATE GENERIQUE VIDE OBLIGE DE PASSER PAR UNE SPECIALISATION
00048 };
00049 
00050 /*********************************************************/
00051 /*                                                       */
00052 /*                   Spécialisation 2D                   */
00053 /*                                                       */
00054 /*********************************************************/
00055 
00056 _TEMPLATE_SPE_ class _COORDBARY_2D_
00057 {
00058 protected :
00059      NUAGEMAILLE * mailles;
00060      NUAGENOEUD * sommets;
00061      
00062      vector<int> etat_coord_baryc;
00063      vector< vector< vector<double> > > coord_baryc;
00064      
00065 public :
00066      
00067      Coordonnees_Barycentriques():mailles(NULL),sommets(NULL) {}
00068      Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n);
00069      ~Coordonnees_Barycentriques() {}
00070      // donne les pseudos coordonnées barycentriques de M dans ma maille de numéro global num_maille dans mailles
00071      // la pseudo coordonnées barycentrique par rapport a une face est la distance normalisée a cette face, 
00072      // dans le cas d'une face triangulaire, c'est la coordonnées ba
00073      vector<double> Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M);
00074      vector<double> Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base); 
00075      vector<double> Calcule_Coord_Baryc(int num_maille, const NOEUD & M);
00076 };
00077 
00078 /*********************************************************/
00079 /*                                                       */
00080 /*                   Spécialisation 3D                   */
00081 /*                                                       */
00082 /*********************************************************/
00083 
00084 
00085 _TEMPLATE_SPE_ class _COORDBARY_3D_
00086 {
00087 protected :
00088      NUAGEMAILLE * mailles;
00089      NUAGENOEUD * sommets;
00090      
00091      vector<int> etat_coord_baryc;
00092      vector< vector< vector<double> > > coord_baryc;
00093      
00094 public :
00095      
00096      Coordonnees_Barycentriques():mailles(NULL),sommets(NULL) {}
00097      Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n);
00098      ~Coordonnees_Barycentriques() {}
00099      vector<double> Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M);
00100      vector<double> Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base); 
00101      vector<double> Calcule_Coord_Baryc(int num_maille, const NOEUD & M);
00102 };
00103 
00109 
00110 _TEMPLATE_SPE_ _COORDBARY_2D_::Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n):mailles(m),sommets(n)
00111           {
00112           cout<<"Creation des Coordonnées Barycentriques : "<<flush;
00113           int nbr_mailles=mailles->SIZE();
00114           etat_coord_baryc=vector<int>(nbr_mailles,FAUX);
00115           coord_baryc=vector< vector< vector<double> > >(nbr_mailles);
00116           cout<<"OK ! "<<endl;
00117           }
00118 
00119 _TEMPLATE_SPE_ vector<double> _COORDBARY_2D_::Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M)
00120      {
00121      int i,nbr_faces;
00122      if (etat_coord_baryc[num_maille]==FAUX) 
00123           {
00124           nbr_faces=(*mailles)[num_maille].DONNE_NBR_FACES();
00125           
00126           coord_baryc[num_maille]=vector< vector<double> >(nbr_faces);
00127           
00128           type_retour simplexe_base;
00129           
00130           for (i=0;i<nbr_faces;i++)
00131                {
00132                (*mailles)[num_maille].DONNE_SIMPLEXE_BASE(i,simplexe_base);
00133                coord_baryc[num_maille][i]=Calcule_Base_Coord_Baryc(vector<int>(&simplexe_base.quoi[0],&simplexe_base.quoi[simplexe_base.combien]));
00134                etat_coord_baryc[num_maille]=VRAI;
00135                }
00136           }    
00137      return Calcule_Coord_Baryc(num_maille,M);
00138      }
00139 
00140 _TEMPLATE_SPE_ vector<double> _COORDBARY_2D_::Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base)
00141      {
00142      const vector<int> &ref=simplexe_base;
00143      vector<double> retour(3);
00144           
00145      double x0=(*sommets)[ref[0]][0];
00146      double y0=(*sommets)[ref[0]][1];
00147      double x1=(*sommets)[ref[1]][0];
00148      double y1=(*sommets)[ref[1]][1];
00149      double x2=(*sommets)[ref[2]][0];
00150      double y2=(*sommets)[ref[2]][1];
00151           
00152      double delta=(x1*y2-x2*y1)+(x2*y0-x0*y2)+(x0*y1-x1*y0);
00153 
00154      retour[0]=(y1-y2)/delta;
00155      retour[1]=(x2-x1)/delta;
00156      retour[2]=(x1*y2-x2*y1)/delta;
00157      
00158      return retour;
00159      }
00160 
00161 _TEMPLATE_SPE_ vector<double> _COORDBARY_2D_::Calcule_Coord_Baryc(int num_maille, const NOEUD & M)
00162      {
00163      int i,j;
00164         // for PAL11458
00165         double zero = 0. ;
00166      //vector<double> coord_baryc_M(3,0);
00167      int nbr_faces=coord_baryc[num_maille].size();
00168      vector<double> coord_baryc_M(nbr_faces,zero);
00169      for (i=0;i</*3*/nbr_faces;i++) 
00170           {
00171           for (j=0;j<2;j++) coord_baryc_M[i]+=coord_baryc[num_maille][i][j]*M[j];
00172           coord_baryc_M[i]+=coord_baryc[num_maille][i][2];
00173           }
00174      return coord_baryc_M;
00175      }
00176 
00177 _TEMPLATE_SPE_ _COORDBARY_3D_::Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n):mailles(m),sommets(n)
00178           {
00179           cout<<"Creation des Coordonnées Barycentriques : "<<flush;
00180           int nbr_mailles=mailles->SIZE();
00181           etat_coord_baryc=vector<int>(nbr_mailles,FAUX);
00182           coord_baryc=vector< vector< vector<double> > >(nbr_mailles);
00183           cout<<"OK ! "<<endl;
00184           }
00185      
00186 _TEMPLATE_SPE_ vector<double> _COORDBARY_3D_::Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M)
00187      {
00188      int i,nbr_faces;
00189      if (etat_coord_baryc[num_maille]==FAUX) 
00190           {
00191           nbr_faces=(*mailles)[num_maille].DONNE_NBR_FACES();
00192           
00193           coord_baryc[num_maille]=vector< vector<double> >(nbr_faces);
00194           
00195           type_retour simplexe_base;
00196           
00197           for (i=0;i<nbr_faces;i++)
00198                {
00199                (*mailles)[num_maille].DONNE_SIMPLEXE_BASE(i,simplexe_base);
00200                coord_baryc[num_maille][i]=Calcule_Base_Coord_Baryc(vector<int>(&simplexe_base.quoi[0],&simplexe_base.quoi[simplexe_base.combien]));
00201                etat_coord_baryc[num_maille]=VRAI;
00202                }
00203           }    
00204      return Calcule_Coord_Baryc(num_maille,M);
00205      }
00206 
00207 
00208 _TEMPLATE_SPE_ vector<double> _COORDBARY_3D_::Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base)
00209      {
00210      const vector<int> &ref=simplexe_base;
00211      vector<double> retour(4);
00212           
00213      double x0=(*sommets)[ref[0]][0];
00214      double y0=(*sommets)[ref[0]][1];
00215      double z0=(*sommets)[ref[0]][2];
00216      double x1=(*sommets)[ref[1]][0];
00217      double y1=(*sommets)[ref[1]][1];
00218      double z1=(*sommets)[ref[1]][2];
00219      double x2=(*sommets)[ref[2]][0];
00220      double y2=(*sommets)[ref[2]][1];
00221      double z2=(*sommets)[ref[2]][2];
00222      double x3=(*sommets)[ref[3]][0];
00223      double y3=(*sommets)[ref[3]][1];
00224      double z3=(*sommets)[ref[3]][2];
00225           
00226      double delta1=((y2-y1)*(z3-z1)-(z2-z1)*(y3-y1));
00227      double delta2=((x3-x1)*(z2-z1)-(x2-x1)*(z3-z1));
00228      double delta3=((x2-x1)*(y3-y1)-(x3-x1)*(y2-y1));
00229      
00230      double delta=delta1*(x0-x1)+delta2*(y0-y1)+delta3*(z0-z1);
00231 
00232      retour[0]=delta1/delta;
00233      retour[1]=delta2/delta;
00234      retour[2]=delta3/delta;
00235      retour[3]=-(delta1*x1+delta2*y1+delta3*z1)/delta;
00236      
00237      return retour;
00238      }
00239 
00240 _TEMPLATE_SPE_ vector<double> _COORDBARY_3D_::Calcule_Coord_Baryc(int num_maille, const NOEUD & M)
00241      {
00242      int i,j;
00243      int nbr_faces=coord_baryc[num_maille].size();
00244      vector<double> coord_baryc_M(nbr_faces);
00245      for (i=0;i<nbr_faces;i++) 
00246           {//CCRT Porting : constructor of vector<T>(int,const T&) not supported on CCRT
00247           coord_baryc_M[i]=0.;
00248           for (j=0;j<3;j++) coord_baryc_M[i]+=coord_baryc[num_maille][i][j]*M[j];
00249           coord_baryc_M[i]+=coord_baryc[num_maille][i][3];
00250           }
00251      return coord_baryc_M;
00252      }
00253 
00254 //*/
00255 
00256 #undef _TEMPLATE_SPE_
00257 // template <class NUAGEMAILLE, class NUAGENOEUD, class NOEUD, int DIMENSION>
00258 #undef _COORDBARYC_
00259 // Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,DIMENSION>
00260 #undef _COORDBARY_2D_
00261 // Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,2>
00262 #undef _COORDBARY_3D_
00263 // Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,3>
00264 
00265 #endif