1 subroutine infqen ( choix, nbenti,
7 > tritet, cotrte, aretet, hettet,
8 > quahex, coquhe, arehex, hethex,
9 > facpyr, cofapy, arepyr, hetpyr,
10 > facpen, cofape, arepen, hetpen,
11 > ntreca, nqueca, nteeca,
13 > ulsort, langue, codret )
14 c ______________________________________________________________________
18 c Outil de Maillage Adaptatif par Raffinement et Deraffinement d'EDF R&D
20 c Version originale enregistree le 18 juin 1996 sous le numero 96036
21 c aupres des huissiers de justice Simart et Lavoir a Clamart
22 c Version 11.2 enregistree le 13 fevrier 2015 sous le numero 2015/014
23 c aupres des huissiers de justice
24 c Lavoir, Silinski & Cherqui-Abrahmi a Clamart
26 c HOMARD est une marque deposee d'Electricite de France
32 c ______________________________________________________________________
34 c INformation : Qualite des ENtites
36 c ______________________________________________________________________
38 c . nom . e/s . taille . description .
39 c .____________________________________________________________________.
40 c . nbenti . e . 1 . nombre d'entites a imprimer .
41 c . coonoe . e . nbnoto . coordonnees des noeuds .
43 c . somare . e .2*nbarto. numeros des extremites d'arete .
44 c . hettri . e . nbtrto . historique de l'etat des triangles .
45 c . aretri . e .nbtrto*3. numeros des 3 aretes des triangles .
46 c . famtri . e . nbtrto . famille des triangles .
47 c . cfatri . e . nctftr*. codes des familles des triangles .
48 c . . . nbftri . 1 : famille MED .
49 c . . . . 2 : type de triangle .
50 c . . . . 3 : numero de surface de frontiere .
51 c . . . . 4 : famille des aretes internes apres raf.
52 c . . . . + l : appartenance a l'equivalence l .
53 c . hetqua . e . nbquto . historique de l'etat des quadrangles .
54 c . arequa . e .nbquto*4. numeros des 4 aretes des quadrangles .
55 c . famqua . e . nbquto . famille des quadrangles .
56 c . cfaqua . e . nctfqu*. codes des familles des quadrangles .
57 c . . . nbfqua . 1 : famille MED .
58 c . . . . 2 : type de quadrangle .
59 c . . . . 3 : numero de surface de frontiere .
60 c . . . . 4 : famille des aretes internes apres raf.
61 c . . . . 5 : famille des triangles de conformite .
62 c . . . . 6 : famille de sf active/inactive .
63 c . . . . + l : appartenance a l'equivalence l .
64 c . hettet . e . nbteto . historique de l'etat des tetraedres .
65 c . tritet . e .nbtecf*4. numeros des 4 triangles des tetraedres .
66 c . cotrte . e .nbtecf*4. code des 4 triangles des tetraedres .
67 c . aretet . e .nbteca*6. numeros des 6 aretes des tetraedres .
68 c . hethex . e . nbheto . historique de l'etat des hexaedres .
69 c . quahex . e .nbhecf*6. numeros des 6 quadrangles des hexaedres .
70 c . coquhe . e .nbhecf*6. codes des 6 quadrangles des hexaedres .
71 c . arehex . e .nbheca12. numeros des 12 aretes des hexaedres .
72 c . hetpyr . e . nbpyto . historique de l'etat des pyramides .
73 c . facpyr . e .nbpycf*5. numeros des 5 faces des pyramides .
74 c . cofapy . e .nbpycf*5. codes des faces des pyramides .
75 c . arepyr . e .nbpyca*8. numeros des 8 aretes des pyramides .
76 c . hetpen . e . nbpeto . historique de l'etat des pentaedres .
77 c . facpen . e .nbpecf*5. numeros des 5 faces des pentaedres .
78 c . cofape . e .nbpecf*5. code des 5 faces des pentaedres .
79 c . arepen . e .nbpeca*9. numeros des 9 aretes des pentaedres .
80 c . ntreca . e . * . nro des triangles dans le calcul en entree .
81 c . nqueca . e . * . nro des quadrangles dans le calcul en ent. .
82 c . nteeca . e . reteto . numero des tetraedres du code de calcul .
83 c . ulsost . e . 1 . unite logique de la sortie standard .
84 c . ulsort . e . 1 . unite logique de la sortie generale .
85 c . langue . e . 1 . langue des messages .
86 c . . . . 1 : francais, 2 : anglais .
87 c . codret . s . 1 . code de retour des modules .
88 c . . . . 0 : pas de probleme .
89 c . . . . 1 : probleme .
90 c ______________________________________________________________________
93 c 0. declarations et dimensionnement
96 c 0.1. ==> generalites
102 parameter ( nompro = 'INFQEN' )
130 double precision coonoe(nbnoto,sdim)
133 integer somare(2,nbarto)
134 integer hettri(nbtrto), aretri(nbtrto,3)
135 integer cfatri(nctftr,nbftri), famtri(nbtrto)
136 integer hetqua(nbquto), arequa(nbquto,4)
137 integer cfaqua(nctfqu,nbfqua), famqua(nbquto)
138 integer tritet(nbtecf,4), cotrte(nbtecf,4), aretet(nbteca,6)
139 integer hettet(nbteto)
140 integer quahex(nbhecf,6), coquhe(nbhecf,6), arehex(nbheca,12)
141 integer hethex(nbheto)
142 integer facpyr(nbpycf,5), cofapy(nbpycf,5), arepyr(nbpyca,8)
143 integer hetpyr(nbpyto)
144 integer facpen(nbpecf,5), cofape(nbpecf,5), arepen(nbpeca,9)
145 integer hetpen(nbpeto)
146 integer ntreca(*), nqueca(*), nteeca(*)
149 integer ulsort, langue, codret
151 c 0.4. ==> variables locales
154 integer codre1, codre2, codre3, codre4
155 integer ptrav1, ptrav2, ptrav3, ptrav4
158 integer ideb, ifin, ipas
160 integer typenh, nbenac
161 integer uldeb, ulfin, ulpas, ulecr
163 double precision vmin, vmax
165 character*8 ntrav1, ntrav2, ntrav3, ntrav4
168 parameter ( nbmess = 10 )
169 character*80 texte(nblang,nbmess)
175 c 1.1. ==> les messages
179 #ifdef _DEBUG_HOMARD_
180 write (ulsort,texte(langue,1)) 'Entree', nompro
184 texte(1,4) = '(''Quel choix : '',a,'' ?'')'
185 texte(1,5) = '(/,''Les pires '',a14,'' :'',/,25(''=''),/)'
186 texte(1,6) = '(/,''Les meilleurs '',a14,'' :'',/,30(''=''),/)'
187 texte(1,8) = '(''* Numeros | Qualite *'')'
188 texte(1,9) = '(''* HOMARD | calcul | *'')'
190 > '(''Aucune face non liee a un tetraedre dans ce maillage.'')'
192 texte(2,4) = '(''What choice : '',a,'' ?'')'
193 texte(2,5) = '(/,''Worst '',a14,'' :'',/,21(''=''),/)'
194 texte(2,6) = '(/,''Best '',a14,'' :'',/,20(''=''),/)'
195 texte(2,8) = '(''* Numbers | Quality *'')'
196 texte(2,9) = '(''* HOMARD |calculation| *'')'
198 > '(''In this mesh, all the faces are connected to tetraedra.'')'
202 10000 format (40('*'))
203 10001 format ('*',i10 ,' |',i10 ,' | ',g12.5,' *')
207 c 1.2. ==> type d'entites
209 if ( choix.eq.'tr' .or.
210 > choix.eq.'TR' ) then
214 elseif ( choix.eq.'qu' .or.
215 > choix.eq.'QU' ) then
219 elseif ( choix.eq.'te' .or.
220 > choix.eq.'TE' ) then
225 write (ulsort,texte(langue,4)) choix
229 c 1.2. ==> tableaux de travail
231 if ( codret.eq.0 ) then
234 call gmalot ( ntrav1, 'entier ', 2*ltrav1, ptrav1, codre1 )
235 call gmalot ( ntrav2, 'reel ', nbenac, ptrav2, codre2 )
236 call gmalot ( ntrav3, 'entier ', nbenac, ptrav3, codre3 )
237 call gmalot ( ntrav4, 'reel ', nbenac, ptrav4, codre4 )
239 codre0 = min ( codre1, codre2, codre3, codre4 )
240 codret = max ( abs(codre0), codret,
241 > codre1, codre2, codre3, codre4 )
245 c 1.3. ==> Preparation de l'affichage
247 uldeb = min(ulsost,ulsort)
248 ulfin = max(ulsost,ulsort)
249 ulpas = max ( 1 , ulfin-uldeb )
252 c 2. recherche des qualites globales
255 if ( codret.eq.0 ) then
257 cgn write (ulsost,90002) 'typenh',typenh
259 #ifdef _DEBUG_HOMARD_
260 write (ulsort,texte(langue,3)) 'UTB05A', nompro
262 call utb05a ( typenh,
268 > tritet, cotrte, aretet, hettet,
269 > quahex, coquhe, arehex, hethex,
270 > facpyr, cofapy, arepyr, hetpyr,
271 > facpen, cofape, arepen, hetpen,
274 > imem(ptrav1), imem(ptrav1+ltrav1),
275 > rmem(ptrav2), rmem(ptrav4),
277 > ulsort, langue, codret )
279 cgn write (ulsost,90002) 'nbeexa',nbeexa
280 cgn do 31 , iaux = 1 , nbeexa
281 cgn write (ulsost,*) imem(ptrav1+iaux-1), rmem(ptrav2+iaux-1)
290 if ( codret.eq.0 ) then
292 if ( nbeexa.ne.0 ) then
294 #ifdef _DEBUG_HOMARD_
295 write (ulsort,texte(langue,3)) 'UTTRIR', nompro
297 call uttrir ( imem(ptrav3), vmin, vmax,
298 > nbeexa, rmem(ptrav2),
299 > ulsort, langue, codret )
309 if ( codret.eq.0 ) then
311 do 40 , ulecr = uldeb , ulfin, ulpas
313 if ( nbeexa.ne.0 ) then
315 if ( nbenti.lt.0 ) then
317 ifin = max(1,nbeexa+nbenti+1)
319 write (ulecr,texte(langue,5)) mess14(langue,3,typenh)
322 ifin = min(nbeexa,nbenti)
324 write (ulecr,texte(langue,6)) mess14(langue,3,typenh)
328 write (ulecr,texte(langue,8))
329 write (ulecr,texte(langue,9))
331 if ( typenh.eq.2 ) then
332 do 41 , iaux = ideb, ifin, ipas
333 jaux = imem(ptrav3+iaux-1)
334 write (ulecr,10001) imem(ptrav1+jaux-1),
335 > ntreca(imem(ptrav1+jaux-1)), rmem(ptrav2+jaux-1)
337 elseif ( typenh.eq.3 ) then
338 do 42 , iaux = ideb, ifin, ipas
339 jaux = imem(ptrav3+iaux-1)
340 write (ulecr,10001) imem(ptrav1+jaux-1),
341 > nteeca(imem(ptrav1+jaux-1)), rmem(ptrav2+jaux-1)
343 elseif ( typenh.eq.4 ) then
344 do 43 , iaux = ideb, ifin, ipas
345 jaux = imem(ptrav3+iaux-1)
346 write (ulecr,10001) imem(ptrav1+jaux-1),
347 > nqueca(imem(ptrav1+jaux-1)), rmem(ptrav2+jaux-1)
354 write (ulecr,texte(langue,10))
366 if ( codret.eq.0 ) then
368 call gmlboj ( ntrav1, codre1 )
369 call gmlboj ( ntrav2, codre2 )
370 call gmlboj ( ntrav3, codre3 )
371 call gmlboj ( ntrav4, codre4 )
373 codre0 = min ( codre1, codre2, codre3, codre4 )
374 codret = max ( abs(codre0), codret,
375 > codre1, codre2, codre3, codre4 )
383 if ( codret.ne.0 ) then
387 write (ulsort,texte(langue,1)) 'Sortie', nompro
388 write (ulsort,texte(langue,2)) codret
392 #ifdef _DEBUG_HOMARD_
393 write (ulsort,texte(langue,1)) 'Sortie', nompro