]> SALOME platform Git repositories - modules/med.git/blob - src/INTERPOLATION/MEDMEM_MappingTools.hxx
Salome HOME
updating the main trunk with the CEA debug devellopment from the branch
[modules/med.git] / src / INTERPOLATION / MEDMEM_MappingTools.hxx
1 #ifndef COORDONNEES_BARYCENTRIQUES_HPP
2 #define COORDONNEES_BARYCENTRIQUES_HPP
3
4 #define _TEMPLATE_SPE_ template <class NUAGEMAILLE, class NUAGENOEUD, class NOEUD>
5 #define _COORDBARYC_ Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,DIMENSION>
6 #define _COORDBARY_2D_ Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,2>
7 #define _COORDBARY_3D_ Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,3>
8
9 //////////////////////////////////////////////////////////////////
10 ///                                                            ///
11 ///                        DECLARATIONS                        ///
12 ///                                                            ///
13 //////////////////////////////////////////////////////////////////
14
15 /*********************************************************/
16 /*                                                       */
17 /*           Classe Coordonnees_Barycentriques           */
18 /*                                                       */
19 /*********************************************************/
20
21 // C'est la définition de la classe générique template qui n'est utilisée qu'à travers ses spécialisations
22 // vu que le nombre de spécialisations a faire est petit (nombre de dimensions utilisée, 3 pour MED)
23 // et vu que l'acces a ces classes doit etre rapide, car ce sont des classes de calcul
24 // la technique de spécialisation, plus lourde, mais plus rapide, a été préférée aux techniques d'héritage
25
26 template <class NUAGEMAILLE, class NUAGENOEUD, class NOEUD, int DIMENSION> class Coordonnees_Barycentriques
27 {
28 // TEMPLATE GENERIQUE VIDE OBLIGE DE PASSER PAR UNE SPECIALISATION
29 };
30
31 /*********************************************************/
32 /*                                                       */
33 /*                   Spécialisation 2D                   */
34 /*                                                       */
35 /*********************************************************/
36
37 _TEMPLATE_SPE_ class _COORDBARY_2D_
38 {
39 protected :
40         NUAGEMAILLE * mailles;
41         NUAGENOEUD * sommets;
42         
43         vector<int> etat_coord_baryc;
44         vector< vector< vector<double> > > coord_baryc;
45         
46 public :
47         
48         Coordonnees_Barycentriques():mailles(NULL),sommets(NULL) {}
49         Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n);
50         ~Coordonnees_Barycentriques() {}
51         // donne les pseudos coordonnées barycentriques de M dans ma maille de numéro global num_maille dans mailles
52         // la pseudo coordonnées barycentrique par rapport a une face est la distance normalisée a cette face, 
53         // dans le cas d'une face triangulaire, c'est la coordonnées ba
54         vector<double> Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M);
55         vector<double> Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base);      
56         vector<double> Calcule_Coord_Baryc(int num_maille, const NOEUD & M);
57 };
58
59 /*********************************************************/
60 /*                                                       */
61 /*                   Spécialisation 3D                   */
62 /*                                                       */
63 /*********************************************************/
64
65
66 _TEMPLATE_SPE_ class _COORDBARY_3D_
67 {
68 protected :
69         NUAGEMAILLE * mailles;
70         NUAGENOEUD * sommets;
71         
72         vector<int> etat_coord_baryc;
73         vector< vector< vector<double> > > coord_baryc;
74         
75 public :
76         
77         Coordonnees_Barycentriques():mailles(NULL),sommets(NULL) {}
78         Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n);
79         ~Coordonnees_Barycentriques() {}
80         vector<double> Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M);
81         vector<double> Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base);      
82         vector<double> Calcule_Coord_Baryc(int num_maille, const NOEUD & M);
83 };
84
85 //////////////////////////////////////////////////////////////////
86 ///                                                            ///
87 ///                            CODE                            ///
88 ///                                                            ///
89 //////////////////////////////////////////////////////////////////
90
91 _TEMPLATE_SPE_ _COORDBARY_2D_::Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n):mailles(m),sommets(n)
92                 {
93                 cout<<"Creation des Coordonnées Barycentriques : "<<flush;
94                 int nbr_mailles=mailles->SIZE();
95                 etat_coord_baryc=vector<int>(nbr_mailles,FAUX);
96                 coord_baryc=vector< vector< vector<double> > >(nbr_mailles);
97                 cout<<"OK ! "<<endl;
98                 }
99
100 _TEMPLATE_SPE_ vector<double> _COORDBARY_2D_::Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M)
101         {
102         int i,nbr_faces;
103         if (etat_coord_baryc[num_maille]==FAUX) 
104                 {
105                 nbr_faces=(*mailles)[num_maille].DONNE_NBR_FACES();
106                 
107                 coord_baryc[num_maille]=vector< vector<double> >(nbr_faces);
108                 
109                 type_retour simplexe_base;
110                 
111                 for (i=0;i<nbr_faces;i++)
112                         {
113                         (*mailles)[num_maille].DONNE_SIMPLEXE_BASE(i,simplexe_base);
114                         coord_baryc[num_maille][i]=Calcule_Base_Coord_Baryc(vector<int>(&simplexe_base.quoi[0],&simplexe_base.quoi[simplexe_base.combien]));
115                         etat_coord_baryc[num_maille]=VRAI;
116                         }
117                 }       
118         return Calcule_Coord_Baryc(num_maille,M);
119         }
120
121 _TEMPLATE_SPE_ vector<double> _COORDBARY_2D_::Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base)
122         {
123         const vector<int> &ref=simplexe_base;
124         vector<double> retour(3);
125                 
126         double x0=(*sommets)[ref[0]][0];
127         double y0=(*sommets)[ref[0]][1];
128         double x1=(*sommets)[ref[1]][0];
129         double y1=(*sommets)[ref[1]][1];
130         double x2=(*sommets)[ref[2]][0];
131         double y2=(*sommets)[ref[2]][1];
132                 
133         double delta=(x1*y2-x2*y1)+(x2*y0-x0*y2)+(x0*y1-x1*y0);
134
135         retour[0]=(y1-y2)/delta;
136         retour[1]=(x2-x1)/delta;
137         retour[2]=(x1*y2-x2*y1)/delta;
138         
139         return retour;
140         }
141
142 _TEMPLATE_SPE_ vector<double> _COORDBARY_2D_::Calcule_Coord_Baryc(int num_maille, const NOEUD & M)
143         {
144         int i,j;
145         vector<double> coord_baryc_M(3,0);
146         for (i=0;i<3;i++) 
147                 {
148                 for (j=0;j<2;j++) coord_baryc_M[i]+=coord_baryc[num_maille][i][j]*M[j];
149                 coord_baryc_M[i]+=coord_baryc[num_maille][i][2];
150                 }
151         return coord_baryc_M;
152         }
153
154 _TEMPLATE_SPE_ _COORDBARY_3D_::Coordonnees_Barycentriques(NUAGEMAILLE * m, NUAGENOEUD *n):mailles(m),sommets(n)
155                 {
156                 cout<<"Creation des Coordonnées Barycentriques : "<<flush;
157                 int nbr_mailles=mailles->SIZE();
158                 etat_coord_baryc=vector<int>(nbr_mailles,FAUX);
159                 coord_baryc=vector< vector< vector<double> > >(nbr_mailles);
160                 cout<<"OK ! "<<endl;
161                 }
162         
163 _TEMPLATE_SPE_ vector<double> _COORDBARY_3D_::Donne_Pseudo_Coord_Baryc(int num_maille,const NOEUD &M)
164         {
165         int i,nbr_faces;
166         if (etat_coord_baryc[num_maille]==FAUX) 
167                 {
168                 nbr_faces=(*mailles)[num_maille].DONNE_NBR_FACES();
169                 
170                 coord_baryc[num_maille]=vector< vector<double> >(nbr_faces);
171                 
172                 type_retour simplexe_base;
173                 
174                 for (i=0;i<nbr_faces;i++)
175                         {
176                         (*mailles)[num_maille].DONNE_SIMPLEXE_BASE(i,simplexe_base);
177                         coord_baryc[num_maille][i]=Calcule_Base_Coord_Baryc(vector<int>(&simplexe_base.quoi[0],&simplexe_base.quoi[simplexe_base.combien]));
178                         etat_coord_baryc[num_maille]=VRAI;
179                         }
180                 }       
181         return Calcule_Coord_Baryc(num_maille,M);
182         }
183
184
185 _TEMPLATE_SPE_ vector<double> _COORDBARY_3D_::Calcule_Base_Coord_Baryc(const vector<int> &simplexe_base)
186         {
187         const vector<int> &ref=simplexe_base;
188         vector<double> retour(4);
189                 
190         double x0=(*sommets)[ref[0]][0];
191         double y0=(*sommets)[ref[0]][1];
192         double z0=(*sommets)[ref[0]][2];
193         double x1=(*sommets)[ref[1]][0];
194         double y1=(*sommets)[ref[1]][1];
195         double z1=(*sommets)[ref[1]][2];
196         double x2=(*sommets)[ref[2]][0];
197         double y2=(*sommets)[ref[2]][1];
198         double z2=(*sommets)[ref[2]][2];
199         double x3=(*sommets)[ref[3]][0];
200         double y3=(*sommets)[ref[3]][1];
201         double z3=(*sommets)[ref[3]][2];
202                 
203         double delta1=((y2-y1)*(z3-z1)-(z2-z1)*(y3-y1));
204         double delta2=((x3-x1)*(z2-z1)-(x2-x1)*(z3-z1));
205         double delta3=((x2-x1)*(y3-y1)-(x3-x1)*(y2-y1));
206         
207         double delta=delta1*(x0-x1)+delta2*(y0-y1)+delta3*(z0-z1);
208
209         retour[0]=delta1/delta;
210         retour[1]=delta2/delta;
211         retour[2]=delta3/delta;
212         retour[3]=-(delta1*x1+delta2*y1+delta3*z1)/delta;
213         
214         return retour;
215         }
216
217 _TEMPLATE_SPE_ vector<double> _COORDBARY_3D_::Calcule_Coord_Baryc(int num_maille, const NOEUD & M)
218         {
219         int i,j;
220         int nbr_faces=coord_baryc[num_maille].size();
221         vector<double> coord_baryc_M(nbr_faces);
222         for (i=0;i<nbr_faces;i++) 
223                 {//CCRT Porting : constructor of vector<T>(int,const T&) not supported on CCRT
224                 coord_baryc_M[i]=0.;
225                 for (j=0;j<3;j++) coord_baryc_M[i]+=coord_baryc[num_maille][i][j]*M[j];
226                 coord_baryc_M[i]+=coord_baryc[num_maille][i][3];
227                 }
228         return coord_baryc_M;
229         }
230
231 //*/
232
233 #undef _TEMPLATE_SPE_
234 // template <class NUAGEMAILLE, class NUAGENOEUD, class NOEUD, int DIMENSION>
235 #undef _COORDBARYC_
236 // Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,DIMENSION>
237 #undef _COORDBARY_2D_
238 // Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,2>
239 #undef _COORDBARY_3D_
240 // Coordonnees_Barycentriques<NUAGEMAILLE,NUAGENOEUD,NOEUD,3>
241
242 #endif