ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/REPOS_ERICCA/magic/lib/diamesh/src/m3d_struct.h
Revision: 283
Committed: Tue Sep 13 21:11:20 2011 UTC (13 years, 8 months ago) by francois
Content type: text/plain
File size: 10429 byte(s)
Log Message:
structure de l'écriture

File Contents

# User Rev Content
1 francois 283
2     /* preparation au maillage */
3     /* SEGMENT */
4     struct st_segment
5     {
6     int n1 ;/* numero du noeud 1 */
7     int n2 ;/* numero du noeud 2 */
8     int card ;/* cardinalite */
9     int tabele[10] ;
10     struct st_segment *suivant ;/* chainage */
11     } ;
12     typedef struct st_segment SEGMENT ;
13    
14     /* ELEMENT */
15     struct st_element
16     {
17     int num ;/* numero de l'element */
18     int mark ;
19     struct st_element *suivant ; /* chainage */
20     } ;
21     typedef struct st_element ELEMENT ;
22    
23     /* SURFACE */
24     struct st_surface
25     {
26     struct st_element *ele ;/* tete de liste des elements de la surface fermee */
27     struct st_surface *suivant ;/* chainage */
28     struct st_envel *envel ;
29     int etat ;
30     };
31     typedef struct st_surface SURFACE ;
32     /* ENVEL = ensemble de surfaces = volume enveloppe */
33     struct st_envel
34     {
35     struct st_surface *surf_ext ;/* surface exterieure */
36     struct st_surface *surf_int ;/* surfaces interieures */
37     } ;
38     typedef struct st_envel ENVEL ;
39    
40     /* structures employees */
41    
42     struct st_noeud
43     {
44     int num ;
45     int mark ;
46     int mark2 ;
47     int flag ;
48     struct st_noeud *suivant ;
49     struct st_oct *oct ;
50     struct st_connec *lis_con ;
51     } ;
52     typedef struct st_noeud NOEUD ;
53     struct st_flag
54     {
55     unsigned val0:
56     1 ;
57     unsigned val1:
58     1 ;
59     unsigned val2:
60     1 ;
61     unsigned val3:
62     1 ;
63     /* alignement */
64     unsigned val4:
65     1 ;
66     unsigned val5:
67     1 ;
68     unsigned val6:
69     1 ;
70     unsigned val7:
71     1 ;
72     unsigned val8:
73     1 ;
74     unsigned val9:
75     1 ;
76     unsigned val10:
77     1 ;
78     unsigned val11:
79     1 ;
80     unsigned val12:
81     1 ;
82     unsigned val13:
83     1 ;
84     unsigned val14:
85     1 ;
86     unsigned val15:
87     1 ;
88     } ;
89     typedef struct st_flag FLAG ;
90     struct st_face
91     {
92     int nb_qua ;
93     int mark ;
94     int hist ;/* historique */
95     int essai ;
96     FLAG tab ;
97     struct st_noeud *n1 ;
98     struct st_noeud *n2 ;
99     struct st_noeud *n3 ;
100     struct st_tetra *tetra1 ;
101     struct st_tetra *tetra2 ;
102     struct st_face *prec ;
103     struct st_face *suivant ;
104     } ;
105     typedef struct st_face FACE ;
106    
107     /* structure octree */
108     struct st_oct
109     {
110     float taille ;
111     float vec[3] ;
112     struct st_oct *tab_oct[8] ;
113     struct st_noeud *lis_noe ;/* liste des noeuds */
114     struct st_objet *lis_obj ;/* liste des faces */
115     struct st_connec *lis_con ;/* liste des tetraedres */
116     };
117     typedef struct st_oct OCT ;
118    
119     /* structure objet */
120     struct st_objet
121     {
122     struct st_face *ref ;/* objet de reference */
123     struct st_objet *suivant ;/* chainage dans l'octree */
124     };
125     typedef struct st_objet OBJET ;
126    
127     struct st_tetra /* ensemble de faces creees */
128     {
129     int num ;
130     int mark ;
131     struct st_face *face1 ;
132     struct st_face *face2 ;
133     struct st_face *face3 ;
134     struct st_face *face4 ;
135     };
136     typedef struct st_tetra TETRAEDRE ;
137     struct st_connec
138     {
139     struct st_tetra *tetra ;/* tetraedre de reference */
140     struct st_connec *suivant ;/* chainage des tetraedres coupant l'octree */
141     };
142     typedef struct st_connec CONNEC ;
143     struct st_gest_mem
144     {
145     /* preparation au maillage */
146    
147     int lecture_donnees ;
148     int premiere_couche ;
149     int maillage_termine ;
150     int size ;
151     double cpu ;
152     float crit_2d_min ;
153     float xmin, ymin, zmin ;
154     float xmax, ymax, zmax ;
155     float coeff ;
156     int icode ;
157     /* int nb_tetra ;*/
158     int nb_good ;
159     int nb_medium ;
160     int nb_low ;
161     int nb_bad ;
162     float rap ;
163     /*int nb_noeud ;*/
164     SEGMENT **tab_seg ;
165     float *v_inter ;
166    
167     /* maillage */
168     /* tableau de coordonnees */
169     float *vcorg ;
170     float *coord ;
171    
172     /* tableau de connectivite */
173     int *numele ;
174     int nb_init ;/* nombre initial de noeuds */
175     int nb_info ;/* nombre de noeuds info */
176     int nb_2d ;/* nombre initial de mailles 2d */
177     int *tabele ;/* table des mailles 2d */
178     int *number ;/* numeros mosaic des noeuds de peau */
179    
180     int numax ;
181     int numael ;
182     /* maillage quadratique */
183     float *coord2 ;
184     int *tablin ;
185     int *neww ;
186     int *old ;
187     int is_qua ;
188     int nb_noe_qua ;
189     int *numele2 ;
190     int *connec ;
191     int *tabcor ;
192     int *tabcoul ;
193     int *tab_rec ;
194     int *tab_ref ;
195     int *tab_card ;
196     int nb_rec ;
197     int card ;
198     /* m3d_maidft */
199     float *val_noeud ;
200     int *corresp ;
201     int *tab_voisin ;
202     OCT *root ;
203     char name[256] ;
204     char buffer[256] ;
205     char extension[256] ;
206     char start[256];
207     char end[256];
208     char release[256] ;
209    
210     /* nouvelles structures : alloc par paquets */
211     /* listes de noeuds, faces, tetraedres ... */
212     ELEMENT *tab_ele ;
213     ELEMENT *tab_ele2 ;
214     SURFACE *tab_surf;
215     int nb_surf ;
216     FACE *front ;
217     SEGMENT *tab_segment[1000] ;
218     int nb_segment ;
219     NOEUD *tab_noeud[1000];
220     int nb_noeud ;
221     FACE *tab_face[1000] ;
222     int nb_face ;
223     TETRAEDRE *tab_tetra[1000] ;
224     int nb_tetra ;
225     OBJET *tab_objet[1000] ;
226     int nb_objet ;
227     CONNEC *tab_connec[1000] ;
228     int nb_connec ;
229     OCT *tab_oct[1000] ;
230     int nb_oct ;
231     /* surfaces enveloppes */
232     ENVEL *envel ;
233    
234     } ;
235    
236     typedef struct st_gest_mem GEST_MEM ;
237    
238     /* ************************************************************ */
239     /* ------------------------------------------------------------ */
240     /* **************** definition de macros ******************** */
241     /* operateurs sur des structures, but : ameliorer l'alloc */
242     /* ------------------------------------------------------------ */
243     /* ************************************************************ */
244     #define NB_ENT 12
245     //#ifdef __STDC__
246     #define mknom(a,b) a##b
247     //#else
248     //#define mknom(a,b) a/**/b
249     //#endif
250     #define TAB(entite) mknom(gest->tab_,entite)
251     #define NBR(entite) mknom(gest->nb_,entite)
252     /* ------------------------------------------ */
253     /* cree et retourne l'adresse de la structure */
254     /* ------------------------------------------ */
255     #define NEW_ALLOC(entite,adresse) \
256     if (TAB(entite)[NBR(entite)>>NB_ENT] == NULL ) \
257     {\
258     TAB(entite)[NBR(entite)>>NB_ENT] = (struct mknom(st_,entite) *)calloc(1<<NB_ENT,sizeof(struct mknom(st_,entite))) ; \
259     gest->size = gest->size + sizeof(struct mknom(st_,entite)) * (1<< NB_ENT) ;\
260     if (debug) printf("alloc de %d entites de taille %d :%d\n",1<<NB_ENT,sizeof(struct mknom(st_,entite)),gest->size) ;\
261     }\
262     if ( TAB(entite)[NBR(entite)>>NB_ENT] == NULL ) adresse = NULL ; \
263     else {\
264     adresse = &(TAB(entite)[NBR(entite)>>NB_ENT][NBR(entite)-(NBR(entite)>>NB_ENT) * (1<<NB_ENT)]) ;\
265     NBR(entite) ++ ;\
266     }
267     /* ---------------------------------------------------------------- */
268     /* trouve l'adresse de la structure de type entite de numero numero */
269     /* ---------------------------------------------------------------- */
270     #define ADRESSE(entite,numero,adr)\
271     if (TAB(entite)[numero>>NB_ENT] == NULL) adr = NULL ; \
272     else adr = &(TAB(entite)[numero>>NB_ENT][numero - (numero>>NB_ENT)*(1<<NB_ENT)]) ;
273    
274     #define ERREUR_ALLOC_NULL(adr) \
275     if (adr == NULL) {\
276     if (debug) printf("%s\n","Erreur alloc") ;\
277     m3d_erreur(ERR_ALLOC) ;\
278     return(NULL) ;\
279     }
280    
281     #define ERREUR_ALLOC_FAUX(adr) \
282     if (adr == NULL) {\
283     if (debug) printf("%s\n","Erreur alloc") ;\
284     m3d_erreur(ERR_ALLOC) ;\
285     return(FAUX) ;\
286     }
287     #define ERREUR_ALLOC(adr) \
288     if (adr == NULL) {\
289     m3d_erreur(ERR_ALLOC) ;\
290     return(FAUX) ;\
291     }
292    
293     /*
294     #define PREC(entite) mknom((entite),->prec)
295     #define SUIVANT(entite) mknom((entite),->suivant)
296     */
297     #define PREC(entite) (entite)->prec
298     #define SUIVANT(entite) (entite)->suivant
299    
300     /* suppression d'un objet d'une liste */
301     /* chainage suivant prec */
302     #define SUPPRIMER_LISTE(type,tete,fin,entite)\
303     {\
304     struct mknom(st_,type) *pcourant, *scourant ;\
305     pcourant = PREC(entite) ;\
306     scourant = SUIVANT(entite) ;\
307     if (pcourant== NULL) {\
308     if (scourant == NULL) { \
309     PREC(entite) = NULL ;\
310     SUIVANT(entite) = NULL ;\
311     tete = NULL ;\
312     fin = NULL ;\
313     }\
314     else {\
315     scourant->prec = NULL ;\
316     PREC(entite) = NULL ;\
317     SUIVANT(entite) = NULL ;\
318     tete = scourant ;\
319     }\
320     }\
321     else {\
322     if (scourant == NULL) \
323     {pcourant->suivant = NULL ;\
324     PREC(entite) = NULL ;\
325     SUIVANT(entite) = NULL ;\
326     fin = pcourant ;}\
327     else {pcourant->suivant = scourant ;\
328     scourant->prec = pcourant ;\
329     PREC(entite) = NULL ;\
330     SUIVANT(entite) = NULL ;}\
331     }\
332     }
333    
334     #define INSERER_FIN_LISTE(tete,fin,entite)\
335     if (tete == NULL)\
336     /* insertion dans une liste vide */\
337     {\
338     tete = entite ;\
339     fin = entite ;\
340     }\
341     /* sinon insertion en queue de liste, la tete est inchangee */\
342     else {\
343     SUIVANT(fin) = entite ;\
344     PREC(entite) = fin ;\
345     SUIVANT(entite) = NULL ;\
346     fin = entite ;\
347     }
348    
349     #define INSERER_TETE_LISTE(tete,fin,entite)\
350     PREC(entite) = NULL ;\
351     if (tete != NULL) {\
352     /* la queue est inchangee */\
353     SUIVANT(entite) = tete ;\
354     PREC(tete) = entite ;\
355     tete = entite ;\
356     }\
357     else {\
358     SUIVANT(entite) = NULL ;\
359     tete = entite ;\
360     fin = entite ;\
361     }
362    
363     #define INSERE_TETE(tete,entite)\
364     SUIVANT(entite) = tete ;\
365     tete = entite ;
366    
367     #define INIT_ALLOC(entite)\
368     {\
369     int i ;\
370     for (i=0;i<1000;i++)\
371     {\
372     TAB(entite)[i] = NULL ;\
373     NBR(entite) = 0 ;\
374     }\
375     }
376    
377     #define FREE_ALLOC(entite)\
378     {\
379     int i ;\
380     i = 0 ;\
381     while (TAB(entite)[i] !=NULL)\
382     {\
383     free(TAB(entite)[i]) ;\
384     TAB(entite)[i] = NULL ;\
385     i++ ;\
386     }\
387     }
388    
389     #define TRI_REEL(tab,nb_val,iordre)\
390     {\
391     int i ;\
392     for (i=0;i<nb_val;i++) iordre[i] = i+1 ;\
393     trirea(tab,&nb_val,iordre) ; \
394     }
395    
396     #define BOITE(n1,n2,n3,n4,vmin,vmax)\
397     vmin[0] = min(coord[x(n1->num)],coord[x(n2->num)]) ;\
398     vmin[0] = min(vmin[0],coord[x(n3->num)]) ;\
399     if (n4 != NULL) vmin[0] = min(vmin[0],coord[x(n4->num)]) ;\
400     vmin[1] = min(coord[y(n1->num)],coord[y(n2->num)]) ;\
401     vmin[1] = min(vmin[1],coord[y(n3->num)]) ;\
402     if (n4 != NULL) vmin[1] = min(vmin[1],coord[y(n4->num)]) ;\
403     vmin[2] = min(coord[z(n1->num)],coord[z(n2->num)]) ;\
404     vmin[2] = min(vmin[2],coord[z(n3->num)]) ;\
405     if (n4 != NULL) vmin[2] = min(vmin[2],coord[z(n4->num)]) ;\
406     vmax[0] = max(coord[x(n1->num)],coord[x(n2->num)]) ;\
407     vmax[0] = max(vmax[0],coord[x(n3->num)]) ;\
408     if (n4 != NULL) vmax[0] = max(vmax[0],coord[x(n4->num)]) ;\
409     vmax[1] = max(coord[y(n1->num)],coord[y(n2->num)]) ;\
410     vmax[1] = max(vmax[1],coord[y(n3->num)]) ;\
411     if (n4 != NULL) vmax[1] = max(vmax[1],coord[y(n4->num)]) ;\
412     vmax[2] = max(coord[z(n1->num)],coord[z(n2->num)]) ;\
413     vmax[2] = max(vmax[2],coord[z(n3->num)]) ;\
414     if (n4 != NULL) vmax[2] = max(vmax[2],coord[z(n4->num)]) ;
415    
416     #define NOEUD_TETRA(tetra,noeud)\
417     if (((tetra->face2)->n1 != (tetra->face1)->n1) && ((tetra->face2)->n1 != (tetra->face1)->n2) && ((tetra->face2)->n1 != (tetra->face1)->n3)) \
418     noeud = (tetra->face2)->n1 ;\
419     else if (((tetra->face2)->n2 != (tetra->face1)->n1) && ((tetra->face2)->n2 != (tetra->face1)->n2) && ((tetra->face2)->n2 != (tetra->face1)->n3)) \
420     noeud = (tetra->face2)->n2 ;\
421     else noeud = (tetra->face2)->n3 ;