Salome HOME
Homard executable
[modules/homard.git] / src / tool / Decision / deinnu.F
1       subroutine deinnu ( nomail, nohind,
2      >                    tyconf, pilraf, pilder, nivmax, nivmin,
3      >                    typseh, typseb, seuilh, seuilb, usacmp,
4      >                    filada, diammi, nbzord, iniada,
5      >                    nbsoci,
6      >                    decare, decfac,
7      >                    povoso, voisom,
8      >                    noempo,
9      >                    somare, hetare, filare, merare,
10      >                    np2are, posifa, facare,
11      >                    aretri, hettri, filtri, pertri, nivtri,
12      >                    voltri, pypetr,
13      >                    arequa, hetqua, filqua, perqua, nivqua,
14      >                    volqua,
15      >                    tritet, hettet, filtet,
16      >                    quahex, hethex, filhex,
17      >                    facpyr, hetpyr,
18      >                    facpen, hetpen, filpen,
19      >                    tabaux,
20      >                    lgopts, taopts,
21      >                    ulsort, langue, codret )
22 c ______________________________________________________________________
23 c
24 c                             H O M A R D
25 c
26 c Outil de Maillage Adaptatif par Raffinement et Deraffinement d'EDF R&D
27 c
28 c Version originale enregistree le 18 juin 1996 sous le numero 96036
29 c aupres des huissiers de justice Simart et Lavoir a Clamart
30 c Version 11.2 enregistree le 13 fevrier 2015 sous le numero 2015/014
31 c aupres des huissiers de justice
32 c Lavoir, Silinski & Cherqui-Abrahmi a Clamart
33 c
34 c    HOMARD est une marque deposee d'Electricite de France
35 c
36 c Copyright EDF 1996
37 c Copyright EDF 1998
38 c Copyright EDF 2002
39 c Copyright EDF 2020
40 c ______________________________________________________________________
41 c
42 c traitement des DEcisions - INitialisations - Non Uniforme
43 c                --          --                -   -
44 c ______________________________________________________________________
45 c .        .     .        .                                            .
46 c .  nom   . e/s . taille .           description                      .
47 c .____________________________________________________________________.
48 c . nomail . e   .  ch8   . nom de l'objet contenant le maillage       .
49 c . nohind . e   .  ch8   . nom de l'objet contenant l'indicateur      .
50 c . tyconf . e   .   1    .  0 : conforme (defaut)                     .
51 c .        .     .        .  1 : non-conforme avec au minimum 2 aretes .
52 c .        .     .        .      non decoupees en 2                    .
53 c .        .     .        .  2 : non-conforme avec 1 seul noeud        .
54 c .        .     .        .      pendant par arete                     .
55 c .        .     .        .  3 : non-conforme fidele a l'indicateur    .
56 c .        .     .        . -1 : conforme, avec des boites pour les    .
57 c .        .     .        .      quadrangles, hexaedres et pentaedres  .
58 c .        .     .        . -2 : non-conforme avec au maximum 1 arete  .
59 c .        .     .        .      decoupee en 2 (boite pour les         .
60 c .        .     .        .      quadrangles, hexaedres et pentaedres) .
61 c . pilraf . e   .   1    . pilotage du raffinement                    .
62 c .        .     .        . -1 : raffinement uniforme                  .
63 c .        .     .        .  0 : pas de raffinement                    .
64 c .        .     .        .  1 : raffinement libre                     .
65 c .        .     .        .  2 : raff. libre homogene en type d'element.
66 c . pilder . e   .   1    . pilotage du deraffinement                  .
67 c .        .     .        . 0 : pas de deraffinement                   .
68 c .        .     .        . 1 : deraffinement libre                    .
69 c .        .     .        . -1 : deraffinement uniforme                .
70 c . nivmax . e   .   1    . niveau max a ne pas depasser en raffinement.
71 c . nivmin . e   .   1    . niveau min a ne pas depasser en deraffinemt.
72 c . typseh . e   .   1    . type de seuil haut                         .
73 c .        .     .        . 1 : absolu                                 .
74 c .        .     .        . 2 : relatif                                .
75 c .        .     .        . 3 : pourcentage d'entites                  .
76 c .        .     .        . 4 : moyenne + nh*ecart-type                .
77 c .        .     .        . 5 : cible en nombre de noeuds              .
78 c . typseb . e   .   1    . type de seuil bas                          .
79 c .        .     .        . 1 : absolu                                 .
80 c .        .     .        . 2 : relatif                                .
81 c .        .     .        . 3 : pourcentage d'entites                  .
82 c .        .     .        . 4 : moyenne - nb*ecart-type                .
83 c . seuilh . e   .   1    . borne superieure de l'erreur (absolue,     .
84 c .        .     .        . relatif, pourcentage d'entites ou nh)      .
85 c . seuilb . e   .   1    . borne inferieure de l'erreur (absolue,     .
86 c .        .     .        . relatif, pourcentage d'entites ou nb)      .
87 c . usacmp . e   .   1    . usage des composantes de l'indicateur      .
88 c .        .     .        . 0 : norme L2                               .
89 c .        .     .        . 1 : norme infinie -max des valeurs absolues.
90 c .        .     .        . 2 : valeur relative si une seule composante.
91 c . filada . e   .   1    . filtrage de l'adaptation                   .
92 c .        .     .        . 0 : pas de filtrage                        .
93 c .        .     .        . >0 : filtrage                              .
94 c . diammi . e   .   1    . diametre minimal voulu                     .
95 c . nbzord . e   .    1   . nombre de zones a raffiner/deraffiner      .
96 c . iniada . e   .   1    . initialisation de l'adaptation             .
97 c .        .     .        . 0 : on garde tout (defaut)                 .
98 c .        .     .        .-1 : reactivation des mailles ou aucun      .
99 c .        .     .        .     indicateur n'est defini                .
100 c .        .     .        . 1 : raffinement des mailles ou aucun       .
101 c .        .     .        .     indicateur n'est defini                .
102 c . nbsoci . e   .   1    . cible en nombre de sommets  (-1 si non)    .
103 c . decare .  s  .0:nbarto. decisions des aretes                       .
104 c . decfac .  s  . -nbquto. decision sur les faces (quad. + tri.)      .
105 c .        .     . :nbtrto.                                            .
106 c . povoso . e   .0:nbnoto. pointeur des voisins par sommet            .
107 c . voisom . e   . nvosom . aretes voisines de chaque noeud            .
108 c . noempo . e   . nbmpto . numeros des noeuds associes aux mailles    .
109 c . somare . e   .2*nbarto. numeros des extremites d'arete             .
110 c . hetare . e   . nbarto . historique de l'etat des aretes            .
111 c . filare . e   . nbarto . premiere fille des aretes                  .
112 c . merare . e   . nbarto . mere des aretes                            .
113 c . np2are . e   . nbarto . noeud milieux des aretes                   .
114 c . posifa . e   . nbarto . pointeur sur tableau facare                .
115 c . facare . e   . nbfaar . liste des faces contenant une arete        .
116 c . aretri . e   .nbtrto*3. numeros des 3 aretes des triangles         .
117 c . hettri . e   . nbtrto . historique de l'etat des triangles         .
118 c . filtri . e   . nbtrto . premier fils des triangles                 .
119 c . pertri . e   . nbtrto . pere des triangles                         .
120 c . nivtri . e   . nbtrto . niveau des triangles                       .
121 c . voltri . e   .2*nbtrto. numeros des 2 volumes par triangle         .
122 c .        .     .        . voltri(i,k) definit le i-eme voisin de k   .
123 c .        .     .        .   0 : pas de voisin                        .
124 c .        .     .        . j>0 : tetraedre j                          .
125 c .        .     .        . j<0 : pyramide/pentaedre dans pypetr(1/2,j).
126 c . pypetr . e   .2*lgpype. pypetr(1,j) = numero de la pyramide voisine.
127 c .        .     .        . du triangle k tel que voltri(1/2,k) = -j   .
128 c .        .     .        . pypetr(2,j) = numero du pentaedre voisin   .
129 c .        .     .        . du triangle k tel que voltri(1/2,k) = -j   .
130 c . arequa . e   .nbquto*4. numeros des 4 aretes des quadrangles       .
131 c . hetqua . e   . nbquto . historique de l'etat des quadrangles       .
132 c . filqua . e   . nbquto . premier fils des quadrangles               .
133 c . perqua . e   . nbquto . pere des quadrangles                       .
134 c . nivqua . e   . nbquto . niveau des quadrangles                     .
135 c . volqua . e   .2*nbquto. numeros des 2 volumes par quadrangle       .
136 c .        .     .        . volqua(i,k) definit le i-eme voisin de k   .
137 c .        .     .        .   0 : pas de voisin                        .
138 c .        .     .        . j>0 : hexaedre j                           .
139 c .        .     .        . j<0 : pyramide/pentaedre dans pypequ(1/2,j).
140 c . tritet . e   .nbtecf*4. numeros des 4 triangles des tetraedres     .
141 c . hettet . e   . nbteto . historique de l'etat des tetraedres        .
142 c . filtet . e   . nbteto . premier fils des tetraedres                .
143 c . quahex . e   .nbhecf*6. numeros des 6 quadrangles des hexaedres    .
144 c . hethex . e   . nbheto . historique de l'etat des hexaedres         .
145 c . filhex . e   . nbheto . premier fils des hexaedres                 .
146 c . facpyr . e   .nbpycf*5. numeros des 5 faces des pyramides          .
147 c . hetpyr . e   . nbpyto . historique de l'etat des pyramides         .
148 c . facpen . e   .nbpecf*5. numeros des faces des pentaedres           .
149 c . hetpen . e   . nbpeto . historique de l'etat des pentaedres        .
150 c . filpen . e   . nbpeto . premier fils des pentaedres                .
151 c . tabaux . a   . -nbquto. tableau auxiliaire sur les faces           .
152 c .        .     . :nbtrto.  (quad. + tri.)                            .
153 c . lgopts . e   .   1    . longueur du tableau des options caracteres .
154 c . taopts . e   . lgopts . tableau des options caracteres             .
155 c . ulsort . e   .   1    . numero d'unite logique de la liste standard.
156 c . langue . e   .    1   . langue des messages                        .
157 c .        .     .        . 1 : francais, 2 : anglais                  .
158 c . codret . es  .    1   . code de retour des modules                 .
159 c .        .     .        . 0 : pas de probleme                        .
160 c ______________________________________________________________________
161 c
162 c====
163 c 0. declarations et dimensionnement
164 c====
165 c
166 c 0.1. ==> generalites
167 c
168       implicit none
169       save
170 c
171       character*6 nompro
172       parameter ( nompro = 'DEINNU' )
173 c
174 #include "nblang.h"
175 c
176 c 0.2. ==> communs
177 c
178 #include "envex1.h"
179 c
180 #include "gmenti.h"
181 #include "gmreel.h"
182 c
183 #include "nombno.h"
184 #include "nombmp.h"
185 #include "nombar.h"
186 #include "nombtr.h"
187 #include "nombqu.h"
188 #include "nombte.h"
189 #include "nombpy.h"
190 #include "nombhe.h"
191 #include "nombpe.h"
192 #include "envca1.h"
193 c
194 c 0.3. ==> arguments
195 c
196       character*8 nomail, nohind
197 c
198       integer tyconf, pilraf, pilder, nivmax, nivmin
199       integer typseh, typseb
200       integer usacmp
201       integer nbzord
202       integer filada, iniada
203       integer nbsoci
204       integer decare(0:nbarto), decfac(-nbquto:nbtrto)
205       integer povoso(0:nbnoto), voisom(*)
206       integer noempo(nbmpto)
207       integer somare(2,nbarto)
208       integer hetare(nbarto), filare(nbarto), merare(nbarto)
209       integer np2are(nbarto)
210       integer posifa(0:nbarto), facare(nbfaar)
211       integer aretri(nbtrto,3), hettri(nbtrto)
212       integer filtri(nbtrto), pertri(nbtrto), nivtri(nbtrto)
213       integer voltri(2,nbtrto), pypetr(2,*)
214       integer arequa(nbquto,4), hetqua(nbquto)
215       integer filqua(nbquto), perqua(nbquto), nivqua(nbquto)
216       integer volqua(2,nbquto)
217       integer tritet(nbtecf,4), hettet(nbteto), filtet(nbteto)
218       integer quahex(nbhecf,6), hethex(nbheto), filhex(nbheto)
219       integer facpyr(nbpycf,5), hetpyr(nbpyto)
220       integer facpen(nbpecf,5), hetpen(nbpeto), filpen(nbpeto)
221       integer tabaux(-nbquto:nbtrto)
222 c
223       double precision seuilb, seuilh
224       double precision diammi
225 c
226       integer lgopts
227       character*8 taopts(lgopts)
228 c
229       integer ulsort, langue, codret
230 c
231 c 0.4. ==> variables locales
232 c
233       integer iaux, jaux
234 c
235       integer pcoono
236       integer adnoin, adnorn, adnosu
237       integer adarin, adarrn, adarsu
238       integer adtrin, adtrrn, adtrsu
239       integer adquin, adqurn, adqusu
240       integer adtein, adtern, adtesu
241       integer adhein, adhern, adhesu
242       integer adpyin, adpyrn, adpysu
243       integer adpein, adpern, adpesu
244       integer adzord
245       integer adtra3, adtra4
246       integer nbvnoe, nbvare
247       integer nbvtri, nbvqua
248       integer nbvtet, nbvhex, nbvpyr, nbvpen
249       integer dimcst, adcocs
250 c
251       integer codre0, codre1, codre2, codre3, codre4, codre5
252       integer codre6, codre7, codre8
253 c
254       integer typind, ncmpin
255 c
256       character*8 ncazor
257       character*8 obfigr, obfidm
258       character*8 norenu
259       character*8 nhnoeu, nhmapo, nharet, nhtria, nhquad
260       character*8 nhtetr, nhhexa, nhpyra, nhpent
261       character*8 nhelig
262       character*8 nhvois, nhsupe, nhsups
263       character*8 ntrav1, ntrav2, ntrav3, ntrav4
264       character*8 ntrano, ntraar, ntratr, ntraqu
265       character*8 ntrate, ntrahe, ntrapy, ntrape
266 c
267       integer nbmess
268       parameter ( nbmess = 10 )
269       character*80 texte(nblang,nbmess)
270 c
271 c 0.5. ==> initialisations
272 c
273 #ifdef _DEBUG_HOMARD_
274       character*1 saux01(3)
275       data saux01 / 'X', 'Y', 'Z' /
276 #endif
277 c ______________________________________________________________________
278 c
279 c====
280 c 1. messages
281 c====
282 c
283 #include "impr01.h"
284 c
285 #ifdef _DEBUG_HOMARD_
286       write (ulsort,texte(langue,1)) 'Entree', nompro
287       call dmflsh (iaux)
288 #endif
289 c
290       texte(1,4) = '(''Erreur de programmation etape 3'')'
291       texte(1,5) = '(''Coordonnee '',a,'' constante :'',g13.5)'
292       texte(1,6) = '(/,5x,''Filtrage par les groupes'')'
293       texte(1,7) = '(/,5x,''Filtrage par le diametre minimal'')'
294 c
295       texte(2,4) = '(''Programming error in stage 3'')'
296       texte(2,5) = '(''Coordinate '',a,'' constant:'',g13.5)'
297       texte(2,6) = '(/,5x,''Filtering by the groups'')'
298       texte(2,7) = '(/,5x,''Filtering by the minimal diametre'')'
299 c
300 #include "impr03.h"
301 c
302 #ifdef _DEBUG_HOMARD_
303       write (ulsort,90002) 'tyconf', tyconf
304       write (ulsort,90002) 'pilraf', pilraf
305       write (ulsort,90002) 'pilder', pilder
306       write (ulsort,90002) 'usacmp', usacmp
307       write (ulsort,90002) 'nivmax', nivmax
308       write (ulsort,90002) 'nivmin', nivmin
309       write (ulsort,90002) 'typseh', typseh
310       write (ulsort,90004) 'seuilh', seuilh
311       write (ulsort,90002) 'typseb', typseb
312       write (ulsort,90004) 'seuilb', seuilb
313       write (ulsort,90002) 'filada', filada
314       write (ulsort,90002) 'nbzord', nbzord
315       write (ulsort,90004) 'diammi', diammi
316 #endif
317 c
318       obfigr = taopts(29)
319 cgn      write (ulsort,90003) 'obfigr', obfigr
320       obfidm = taopts(28)
321 cgn      write (ulsort,90003) 'obfidm', obfidm
322 c
323 c====
324 c 2. gestion des tableaux
325 c====
326 c
327 c 2.1. ==> structure generale
328 c
329       if ( codret.eq.0 ) then
330 c
331 #ifdef _DEBUG_HOMARD_
332       write (ulsort,texte(langue,3)) 'UTNOMH', nompro
333 #endif
334       call utnomh ( nomail,
335      >                sdim,   mdim,
336      >               degre, maconf, homolo, hierar,
337      >              rafdef, nbmane, typcca, typsfr, maextr,
338      >              mailet,
339      >              norenu,
340      >              nhnoeu, nhmapo, nharet,
341      >              nhtria, nhquad,
342      >              nhtetr, nhhexa, nhpyra, nhpent,
343      >              nhelig,
344      >              nhvois, nhsupe, nhsups,
345      >              ulsort, langue, codret)
346 c
347       endif
348 c
349 c 2.2. ==> tableaux
350 c
351       if ( nbzord.ne.0 ) then
352 c
353         if ( codret.eq.0 ) then
354 c
355 #ifdef _DEBUG_HOMARD_
356       write (ulsort,texte(langue,3)) 'UTAD01', nompro
357 #endif
358         iaux = 57
359         call utad01 ( iaux, nhnoeu,
360      >                  jaux,
361      >                  jaux,   jaux,   jaux,
362      >                pcoono,   jaux,   jaux, adcocs,
363      >                ulsort, langue, codret )
364 c
365         endif
366 c
367         if ( codret.eq.0 ) then
368 c
369         call gmliat ( nhnoeu, 2, dimcst, codret )
370 c
371         endif
372 c
373 #ifdef _DEBUG_HOMARD_
374         if ( dimcst.ne.0 ) then
375           write (ulsort,texte(langue,5)) saux01(dimcst), rmem(adcocs)
376         endif
377 #endif
378 c
379       endif
380 c
381 c====
382 c 3. Decompte des nombres de valeurs pour les 'faux' indicateurs
383 c====
384 #ifdef _DEBUG_HOMARD_
385       write (ulsort,90002) '3. Decompte ; codret', codret
386 #endif
387 c
388       nbvnoe = 0
389       nbvare = 0
390       nbvtri = 0
391       nbvqua = 0
392       nbvtet = 0
393       nbvpyr = 0
394       nbvhex = 0
395       nbvpen = 0
396 c
397 c 3.1. ==> Uniforme et filtre
398 c
399       if ( ( pilraf.eq.-1 .or. pilder.eq.-1 ) .and.
400      >     ( filada.ne.0 .or. diammi.gt.0.d0 ) ) then
401 c
402         if ( codret.eq.0 ) then
403 c
404         if ( filada.ne.0 ) then
405 c
406 #ifdef _DEBUG_HOMARD_
407       write (ulsort,texte(langue,3)) 'DEINI5', nompro
408 #endif
409         call deini5 ( obfigr,
410      >                nbvnoe, nbvare,
411      >                nbvtri, nbvqua,
412      >                nbvtet, nbvhex, nbvpyr, nbvpen,
413      >                ulsort, langue, codret )
414 c
415         endif
416 c
417         endif
418         typind = 1
419 c
420 c 3.2. ==> Par zone
421 c
422       elseif ( ( pilraf.gt.0 .or. pilder.gt.0 ) .and.
423      >           nbzord.ne.0 ) then
424 c
425         if ( codret.eq.0 ) then
426 c
427         nbvare = nbarto
428         typind = 0
429 c
430         endif
431 c
432 c 3.3. ==> Cas du raffinement ou deraffinement par un indicateur
433 c
434       elseif ( ( pilraf.gt.0 .or. pilder.gt.0 ) .and.
435      >           nbzord.eq.0 ) then
436 c
437         if ( codret.eq.0 ) then
438 c
439 #ifdef _DEBUG_HOMARD_
440       write (ulsort,texte(langue,3)) 'DEINI0', nompro
441 #endif
442         call deini0 ( nohind, typind, ncmpin,
443      >                nbvnoe, nbvare,
444      >                nbvtri, nbvqua,
445      >                nbvtet, nbvhex, nbvpyr, nbvpen,
446      >                adnoin, adnorn, adnosu,
447      >                adarin, adarrn, adarsu,
448      >                adtrin, adtrrn, adtrsu,
449      >                adquin, adqurn, adqusu,
450      >                adtein, adtern, adtesu,
451      >                adhein, adhern, adhesu,
452      >                adpyin, adpyrn, adpysu,
453      >                adpein, adpern, adpesu,
454      >                ulsort, langue, codret )
455 c
456         endif
457 c
458 c 3.4. ==> Autres cas impossibles
459 c
460       elseif ( pilder.ne.-1 ) then
461 c
462         codret = 2
463         write (ulsort,texte(langue,4))
464 c
465       endif
466 #ifdef _DEBUG_HOMARD_
467       write (ulsort,90002) 'typind', typind
468       write (ulsort,90002) 'ncmpin', ncmpin
469       write (ulsort,90002)
470      >' nbvnoe, nbvare, nbvtri, nbvqua, nbvtet, nbvhex, nbvpyr, nbvpen',
471      >                     nbvnoe, nbvare,
472      >                     nbvtri, nbvqua,
473      >                     nbvtet, nbvhex, nbvpyr, nbvpen
474 #endif
475 c
476 c====
477 c 4. Allocations des eventuels tableaux entiers :
478 c    . pour une adaptation selon des zones
479 c    . pour un indicateur reel
480 c====
481 #ifdef _DEBUG_HOMARD_
482       write (ulsort,90002) '4. Allocations ; codret', codret
483 #endif
484 c
485       if ( codret.eq.0 ) then
486 c
487       if ( typind.eq.0 .or. typind.eq.3 ) then
488 c
489         if ( nbvnoe.eq.0 ) then
490           iaux = 0
491         else
492           iaux = nbnoto
493         endif
494         call gmalot ( ntrano, 'entier  ', iaux, adnoin, codre1 )
495 c
496         if ( nbvare.eq.0 ) then
497           iaux = 0
498         else
499           iaux = nbarto
500         endif
501         call gmalot ( ntraar, 'entier  ', iaux, adarin, codre2 )
502 c
503         if ( nbvtri.eq.0 ) then
504           iaux = 0
505         else
506           iaux = nbtrto
507         endif
508         call gmalot ( ntratr, 'entier  ', iaux, adtrin, codre3 )
509 c
510         if ( nbvqua.eq.0 ) then
511           iaux = 0
512         else
513           iaux = nbquto
514         endif
515         call gmalot ( ntraqu, 'entier  ', iaux, adquin, codre4 )
516 c
517         if ( nbvtet.eq.0 ) then
518           iaux = 0
519         else
520           iaux = nbteto
521         endif
522         call gmalot ( ntrate, 'entier  ', iaux, adtein, codre5 )
523 c
524         if ( nbvpyr.eq.0 ) then
525           iaux = 0
526         else
527           iaux = nbpyto
528         endif
529         call gmalot ( ntrapy, 'entier  ', iaux, adpyin, codre6 )
530 c
531         if ( nbvhex.eq.0 ) then
532           iaux = 0
533         else
534           iaux = nbheto
535         endif
536         call gmalot ( ntrahe, 'entier  ', iaux, adhein, codre7 )
537 c
538         if ( nbvpen.eq.0 ) then
539           iaux = 0
540         else
541           iaux = nbpeto
542         endif
543         call gmalot ( ntrape, 'entier  ', iaux, adpein, codre8 )
544 c
545         codre0 = min ( codre1, codre2, codre3, codre4, codre5,
546      >                 codre6, codre7, codre8 )
547         codret = max ( abs(codre0), codret,
548      >                 codre1, codre2, codre3, codre4, codre5,
549      >                 codre6, codre7, codre8 )
550 c
551       endif
552 c
553       endif
554 c
555 c====
556 c 5. Remplissage des tableaux entiers
557 c====
558 c
559 #ifdef _DEBUG_HOMARD_
560       write (ulsort,90002) '5. remplissage ; codret', codret
561 #endif
562 c
563       if ( pilraf.gt.0 .or. pilder.gt.0 ) then
564 c
565 c 5.1. ==> Cas du raffinement ou deraffinement par des zones
566 c          geometriques : on convertit en un indicateur entier sur
567 c          les aretes
568 c
569         if ( nbzord.ne.0 ) then
570 c
571 c 5.1.1. ==> Recuperation de la structure
572 c
573           if ( codret.eq.0 ) then
574 c
575           ncazor = taopts(19)
576 #ifdef _DEBUG_HOMARD_
577           call gmprsx (nompro, ncazor )
578 #endif
579 c
580           call gmadoj ( ncazor, adzord, iaux, codre1 )
581           call gmalot ( ntrav1, 'entier  ', nbnoto, adnosu, codre2 )
582           call gmalot ( ntrav2, 'entier  ', nbarto, adarsu, codre3 )
583 c
584           codre0 = min ( codre1, codre2, codre3 )
585           codret = max ( abs(codre0), codret,
586      >                   codre1, codre2, codre3 )
587 c
588           endif
589 c
590 c 5.1.2. ==> Deploiement
591 c
592           if ( codret.eq.0 ) then
593 c
594 #ifdef _DEBUG_HOMARD_
595         write (ulsort,texte(langue,3)) 'DEINZR', nompro
596 #endif
597           call deinzr ( nbzord, rmem(adzord),
598      >                  rmem(pcoono), dimcst, rmem(adcocs),
599      >                  somare, hetare,
600      >                  imem(adnosu), imem(adarsu), imem(adarin),
601      >                  ulsort, langue, codret )
602 c
603           endif
604 c
605 c 5.3. ==> Cas du raffinement ou deraffinement par un indicateur reel :
606 c          on convertit en un indicateur entier
607 c
608         else
609 c
610           if ( typind.eq.3 ) then
611 c
612             if ( codret.eq.0 ) then
613 c
614 #ifdef _DEBUG_HOMARD_
615       write (ulsort,texte(langue,3)) 'DEINRI', nompro
616 #endif
617             call deinri
618      >       ( pilraf, pilder,
619      >         typseh, typseb, seuilh, seuilb, nbsoci,
620      >         usacmp,
621      >         nbvpen, nbvpyr, nbvhex, nbvtet,
622      >         nbvqua, nbvtri, nbvare, nbvnoe,
623      >         imem(adnosu), rmem(adnorn), imem(adnoin),
624      >         imem(adarsu), rmem(adarrn), imem(adarin),
625      >         imem(adtrsu), rmem(adtrrn), imem(adtrin),
626      >         imem(adqusu), rmem(adqurn), imem(adquin),
627      >         imem(adtesu), rmem(adtern), imem(adtein),
628      >         imem(adhesu), rmem(adhern), imem(adhein),
629      >         imem(adpysu), rmem(adpyrn), imem(adpyin),
630      >         imem(adpesu), rmem(adpern), imem(adpein),
631      >         ulsort, langue, codret)
632 c
633             endif
634 c
635           endif
636 c
637         endif
638 c
639       endif
640 c
641 c====
642 c 6. Elaboration des decisions sur les faces et les aretes
643 c====
644 #ifdef _DEBUG_HOMARD_
645       write (ulsort,90002) '6. Elaboration ; codret', codret
646 #endif
647 c
648 c 6.1. ==> Cas du raffinement/deraffinement uniforme
649 c
650       if ( pilraf.eq.-1 .or. pilder.eq.-1 ) then
651 c
652         if ( codret.eq.0 ) then
653 c
654 #ifdef _DEBUG_HOMARD_
655         write (ulsort,texte(langue,3)) 'DEINUN', nompro
656 #endif
657         call deinun ( pilraf, pilder, nivmax, nivmin,
658      >                decfac, decare,
659      >                hetare,
660      >                hettri,
661      >                hetqua,
662      >                ulsort, langue, codret )
663 c
664         endif
665 c
666 c 6.2. ==> Zone ou indicateur
667 c
668       else
669 c
670         if ( codret.eq.0 ) then
671 c
672 #ifdef _DEBUG_HOMARD_
673         write (ulsort,texte(langue,3)) 'DEINII', nompro
674 #endif
675         call deinii
676      >       ( pilraf, pilder, nivmax, nivmin, iniada,
677      >         decare, decfac,
678      >         somare, hetare, filare, merare, np2are,
679      >         posifa, facare,
680      >         aretri, hettri, filtri, pertri, nivtri,
681      >         arequa, hetqua, filqua, perqua, nivqua,
682      >         tritet, hettet, filtet,
683      >         quahex, hethex, filhex,
684      >         facpyr, hetpyr,
685      >         facpen, hetpen, filpen,
686      >         nbvpen, nbvpyr, nbvhex, nbvtet,
687      >         nbvqua, nbvtri, nbvare, nbvnoe,
688      >         imem(adnosu), imem(adnoin),
689      >         imem(adarsu), imem(adarin),
690      >         imem(adtrsu), imem(adtrin),
691      >         imem(adqusu), imem(adquin),
692      >         imem(adtesu), imem(adtein),
693      >         imem(adhesu), imem(adhein),
694      >         imem(adpysu), imem(adpyin),
695      >         imem(adpesu), imem(adpein),
696      >         ulsort, langue, codret)
697 c
698         endif
699 c
700       endif
701 #ifdef _DEBUG_HOMARD_
702       if ( codret.eq.0 ) then
703 c
704       iaux = 2
705       call delist ( nomail, 'DEINII', iaux,
706      >              lgopts, taopts,
707      >              ulsort, langue, codret )
708 c
709       endif
710 #endif
711 c
712 c====
713 c 7. Filtrages
714 c====
715 #ifdef _DEBUG_HOMARD_
716       write (ulsort,90002) '7. Filtrages ; codret', codret
717 #endif
718 cgn      write(*,*)'decare'
719 cgn      write(*,91030)(decare(iaux),iaux=1,nbarto)
720 cgn      write(*,*)'decfac quad'
721 cgn      write(*,91030)(decfac(iaux),iaux=-nbquto,-1)
722 cgn      write(*,*)'decfac tria'
723 cgn      write(*,91030)(decfac(iaux),iaux=1,nbtrto)
724 c
725       if ( filada.ne.0 .or. diammi.gt.0.d0 ) then
726 c
727 c 7.1. ==> Tableaux de travail
728 c
729         if ( codret.eq.0 ) then
730 c
731         call gmalot ( ntrav3, 'entier  ', nbarto, adtra3, codre1 )
732         iaux = nbquto + 1 + nbtrto
733         call gmalot ( ntrav4, 'entier  ',   iaux, adtra4, codre2 )
734 c
735         codre0 = min ( codre1, codre2 )
736         codret = max ( abs(codre0), codret,
737      >                 codre1, codre2 )
738 c
739         endif
740 cgn      write(*,91030)(decare(iaux),iaux=1,nbarto)
741 cgn      write(*,91030)(decfac(iaux),iaux=1,nbtrto)
742 c
743 c 7.2. ==> Applications du ou des filtrages
744 c 7.2.1. ==> Filtrage par les groupes
745 c
746         if ( filada.ne.0 ) then
747 c
748 c 7.2.1.1. ==> Filtrage effectif
749 c
750           if ( codret.eq.0 ) then
751 c
752           write (ulsort,texte(langue,6))
753 c
754           iaux = 1
755 c
756 #ifdef _DEBUG_HOMARD_
757           write (ulsort,texte(langue,3)) 'DEINFI-groupes', nompro
758 #endif
759           call deinfi ( iaux, obfigr,
760      >                  decare, decfac, iniada,
761      >                  imem(adtra3), imem(adtra4),
762      >                  povoso, voisom,
763      >                  noempo,
764      >                  somare,
765      >                  aretri,
766      >                  arequa,
767      >                  tritet,
768      >                  quahex,
769      >                  facpyr,
770      >                  facpen,
771      >                  ulsort, langue, codret )
772 cgn      call gmprsx ( nompro, ntrav3 )
773 cgn      call gmprsx ( nompro, ntrav4 )
774 cgn      write(*,91030)(decare(iaux),iaux=1,nbarto)
775 cgn      write(*,91030)(decfac(iaux),iaux=1,nbtrto)
776 c
777           endif
778 c
779 c 7.2.1.2. ==> Affichage final
780 c
781           if ( codret.eq.0 ) then
782 c
783 #ifdef _DEBUG_HOMARD_
784         write (ulsort,texte(langue,3)) 'DECPTE', nompro
785 #endif
786           call decpte ( pilraf, pilder,
787      >                  decare, decfac,
788      >                  hettri, hetqua,
789      >                  tritet, hettet,
790      >                  quahex, hethex,
791      >                  facpyr, hetpyr,
792      >                  facpen, hetpen,
793      >                  ulsort, langue, codret )
794 c
795           endif
796 c
797         endif
798 c
799 c 7.2.2. ==> Filtrage par le diametre minimal
800 c
801         if ( diammi.gt.0.d0 ) then
802 c
803 c 7.2.2.1. ==> Filtrage effectif
804 c
805           if ( codret.eq.0 ) then
806 c
807           write (ulsort,texte(langue,7))
808 c
809           iaux = 0
810 c
811 #ifdef _DEBUG_HOMARD_
812           write (ulsort,texte(langue,3)) 'DEINFI-diametre', nompro
813 #endif
814           call deinfi ( iaux, obfidm,
815      >                  decare, decfac, iniada,
816      >                  imem(adtra3), imem(adtra4),
817      >                  povoso, voisom,
818      >                  noempo,
819      >                  somare,
820      >                  aretri,
821      >                  arequa,
822      >                  tritet,
823      >                  quahex,
824      >                  facpyr,
825      >                  facpen,
826      >                  ulsort, langue, codret )
827 cgn      write(*,*)'decare'
828 cgn      write(*,91030)(decare(iaux),iaux=1,nbarto)
829 cgn      write(*,*)'decfac quad'
830 cgn      write(*,91030)(decfac(iaux),iaux=-nbquto,-1)
831 cgn      write(*,*)'decfac tria'
832 cgn      write(*,91030)(decfac(iaux),iaux=1,nbtrto)
833 cgn      call gmprsx ( nompro, ntrav3 )
834 cgn      call gmprsx ( nompro, ntrav4 )
835 c
836           endif
837 c
838 c 7.2.2.2. ==> Affichage final
839 c
840           if ( codret.eq.0 ) then
841 c
842 #ifdef _DEBUG_HOMARD_
843         write (ulsort,texte(langue,3)) 'DECPTE', nompro
844 #endif
845           call decpte ( pilraf, pilder,
846      >                  decare, decfac,
847      >                  hettri, hetqua,
848      >                  tritet, hettet,
849      >                  quahex, hethex,
850      >                  facpyr, hetpyr,
851      >                  facpen, hetpen,
852      >                  ulsort, langue, codret )
853 c
854           endif
855 c
856         endif
857 c
858       endif
859 cgn      do 7777 , iaux = 1 , nbarto
860 cgn        if ( decare(iaux).ne.0 ) then
861 cgn          write (ulsort,90001) '.. arete e/d', iaux,
862 cgn     >    hetare(iaux), decare(iaux), somare(1,iaux), somare(2,iaux)
863 cgn        endif
864 cgn 7777 continue
865 c
866 c====
867 c 8. Corrections selon le mode de conformite
868 c====
869 #ifdef _DEBUG_HOMARD_
870       write (ulsort,90002) '8. correction ; codret', codret
871 #endif
872 c
873       if ( codret.eq.0 ) then
874 c
875 #ifdef _DEBUG_HOMARD_
876       write (ulsort,texte(langue,3)) 'DEINI4', nompro
877 #endif
878       call deini4 ( tyconf,
879      >              decare, decfac,
880      >              hetare, filare,
881      >              aretri, hettri, filtri,
882      >              voltri, pypetr,
883      >              arequa, hetqua,
884      >              volqua,
885      >              tritet, quahex, facpen, facpyr,
886      >              tabaux,
887      >              ulsort, langue, codret)
888 cgn      write(*,*)'decare'
889 cgn      write(*,91030)(decare(iaux),iaux=1,nbarto)
890 cgn      write(*,*)'decfac quad'
891 cgn      write(*,91030)(decfac(iaux),iaux=-nbquto,-1)
892 cgn      write(*,*)'decfac tria'
893 cgn      write(*,91030)(decfac(iaux),iaux=1,nbtrto)
894 c
895 cgn      do 8888 , iaux = 1 , nbarto
896 cgn        if ( decare(iaux).ne.0 ) then
897 cgn          write (ulsort,90001) '.. arete e/d', iaux,
898 cgn     >    hetare(iaux), decare(iaux), somare(1,iaux), somare(2,iaux)
899 cgn        endif
900 cgn 8888 continue
901       endif
902 c
903 c====
904 c 9. Menage
905 c====
906 #ifdef _DEBUG_HOMARD_
907       write (ulsort,90002) '9. Menage ; codret', codret
908 #endif
909 c
910       if ( codret.eq.0 ) then
911 c
912 c 9.1. ==> Zones
913 c
914       if ( ( pilraf.gt.0 .or. pilder.gt.0 ) .and.
915      >       nbzord.ne.0 ) then
916 c
917         call gmlboj ( ntrav1, codre1 )
918         call gmlboj ( ntrav2, codre2 )
919 c
920         codre0 = min ( codre1, codre2 )
921         codret = max ( abs(codre0), codret,
922      >                 codre1, codre2 )
923 c
924       endif
925 c
926 c 9.2. ==> Filtrage
927 #ifdef _DEBUG_HOMARD_
928       write (ulsort,90002) '9.2. Filtrage ; codret', codret
929 #endif
930 c
931       if ( filada.ne.0 .or. diammi.gt.0.d0 ) then
932 c
933         call gmlboj ( ntrav3, codre1 )
934         call gmlboj ( ntrav4, codre2 )
935 c
936         codre0 = min ( codre1, codre2 )
937         codret = max ( abs(codre0), codret,
938      >                 codre1, codre2 )
939 c
940       endif
941 c
942 c 9.3. ==> Temporaires
943 #ifdef _DEBUG_HOMARD_
944       write (ulsort,90002) '9.3. Temporaire ; codret', codret
945 #endif
946 c
947       if ( typind.eq.0 .or. typind.eq.3 ) then
948 c
949         call gmlboj ( ntrano, codre1 )
950         call gmlboj ( ntraar, codre2 )
951         call gmlboj ( ntratr, codre3 )
952         call gmlboj ( ntraqu, codre4 )
953         call gmlboj ( ntrate, codre5 )
954         call gmlboj ( ntrapy, codre6 )
955         call gmlboj ( ntrahe, codre7 )
956         call gmlboj ( ntrape, codre8 )
957 c
958         codre0 = min ( codre1, codre2, codre3, codre4, codre5,
959      >                 codre6, codre7, codre8 )
960         codret = max ( abs(codre0), codret,
961      >                 codre1, codre2, codre3, codre4, codre5,
962      >                 codre6, codre7, codre8 )
963 c
964         endif
965 c
966       endif
967 c
968 c====
969 c 10. la fin
970 c====
971 c
972       if ( codret.ne.0 ) then
973 c
974 #include "envex2.h"
975 c
976       write (ulsort,texte(langue,1)) 'Sortie', nompro
977       write (ulsort,texte(langue,2)) codret
978 c
979       endif
980 c
981 #ifdef _DEBUG_HOMARD_
982       write (ulsort,texte(langue,1)) 'Sortie', nompro
983       call dmflsh (iaux)
984 #endif
985 c
986       end