Salome HOME
Homard executable
[modules/homard.git] / src / tool / AV_Conversion / vcmnco.F
1       subroutine vcmnco ( nohman,
2      >                    nhnoeu, nharet, nhtria, nhquad, nhvois,
3      >                    noempo,
4      >                    coonoe, hetnoe, arenoe,
5      >                    coexno, nnosho, nnosca,
6      >                    somare, hetare, np2are,
7      >                    merare, filare, insoar,
8      >                    coexar, narsho, narsca,
9      >                    hettri, aretri, filtri, pertri,
10      >                    hetqua, arequa, filqua, perqua,
11      >                    coexqu, nqusho, nqusca,
12      >                    quahex, coquhe,
13      >                    ulsort, langue, codret )
14 c ______________________________________________________________________
15 c
16 c                             H O M A R D
17 c
18 c Outil de Maillage Adaptatif par Raffinement et Deraffinement d'EDF R&D
19 c
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
25 c
26 c    HOMARD est une marque deposee d'Electricite de France
27 c
28 c Copyright EDF 1996
29 c Copyright EDF 1998
30 c Copyright EDF 2002
31 c Copyright EDF 2020
32 c ______________________________________________________________________
33 c
34 c    aVant adaptation - Conversion de Maillage - Non COnformite
35 c     -                 -             -          -   --
36 c    ATTENTION : cela suppose que le rapport de non conformite est 1/2
37 c ______________________________________________________________________
38 c .        .     .        .                                            .
39 c .  nom   . e/s . taille .           description                      .
40 c .____________________________________________________________________.
41 c . nohman . e   . char*8 . nom de l'objet maillage homard iteration n .
42 c . nhnoeu . e   . char8  . nom de l'objet decrivant les noeuds        .
43 c . nharet . e   . char8  . nom de l'objet decrivant les aretes        .
44 c . nhtria . e   . char8  . nom de l'objet decrivant les triangles     .
45 c . nhquad . e   . char8  . nom de l'objet decrivant les quadrangles   .
46 c . nhvois . e   . char8  . nom de la branche Voisins                  .
47 c . noempo . es  . nbmpto . numeros des noeuds associes aux mailles    .
48 c . coonoe . e   . nbnoto . coordonnees des noeuds                     .
49 c .        .     . * sdim .                                            .
50 c . hetnoe . es  . nbnoto . historique de l'etat des noeuds            .
51 c . arenoe . es  . nbnoto . 0 pour un sommet, le numero de l'arete pour.
52 c .        .     .        . un noeud milieu                            .
53 c . coexno . es  . nbnoto*. codes de conditions aux limites portants   .
54 c .        .     . nctfno . sur les noeuds                             .
55 c . nnosho . es  . rsnoac . numero des noeuds dans HOMARD              .
56 c . nnosca . es  . rsnoto . numero des noeuds dans le calcul           .
57 c . somare . es  .2*nbarto. numeros des extremites d'arete             .
58 c . hetare . es  . nbarto . historique de l'etat des aretes            .
59 c . np2are . es  . nbarto . noeud milieux des aretes                   .
60 c . merare . es  . nbarto . mere des aretes                            .
61 c . filare . es  . nbarto . premiere fille des aretes                  .
62 c . insoar . es  . nbarto . information sur les sommets des aretes     .
63 c . coexar . es  . nbarto*. codes de conditions aux limites portants   .
64 c .        .     . nctfar . sur les aretes                             .
65 c . narsho . es  . rsarac . numero des aretes dans HOMARD              .
66 c . narsca . es  . rsarto . numero des aretes du calcul                .
67 c . hettri . e   . nbtrto . historique de l'etat des triangles         .
68 c . aretri . e   .nbtrto*3. numeros des 3 aretes des triangles         .
69 c . filtri . e   . nbtrto . premier fils des triangles                 .
70 c . pertri . e   . nbtrto . pere des triangles                         .
71 c . hetqua . e   . nbquto . historique de l'etat des quadrangles       .
72 c . arequa . e   .nbquto*4. numeros des 4 aretes des quadrangles       .
73 c . filqua . e   . nbquto . premier fils des quadrangles               .
74 c . perqua . e   . nbquto . pere des quadrangles                       .
75 c . coexqu . es  . nbquto*. codes de conditions aux limites portants   .
76 c .        .     . nctfqu . sur les quadrangles                        .
77 c . nqusho . es  . rsquac . numero des quadrangles dans HOMARD         .
78 c . nqusca . es  . rsquto . numero des quadrangles du calcul           .
79 c . quahex . es  .nbhecf*6. numeros des 6 quadrangles des hexaedres    .
80 c . coquhe . es  .nbhecf*6. codes des 6 quadrangles des hexaedres      .
81 c . ulsort . e   .   1    . numero d'unite logique de la liste standard.
82 c . langue . e   .    1   . langue des messages                        .
83 c .        .     .        . 1 : francais, 2 : anglais                  .
84 c . codret . es  .    1   . code de retour des modules                 .
85 c .        .     .        . 0 : pas de probleme                        .
86 c .        .     .        . 3 : probleme                               .
87 c ______________________________________________________________________
88 c
89 c====
90 c 0. declarations et dimensionnement
91 c====
92 c
93 c 0.1. ==> generalites
94 c
95       implicit none
96       save
97 c
98       character*6 nompro
99       parameter ( nompro = 'VCMNCO' )
100 c
101 #include "nblang.h"
102 c
103 c 0.2. ==> communs
104 c
105 #include "envex1.h"
106 #include "gmenti.h"
107 #include "gmreel.h"
108 c
109 #include "impr02.h"
110 #include "nombno.h"
111 #include "nombmp.h"
112 #include "nombar.h"
113 #include "nombtr.h"
114 #include "nombqu.h"
115 #include "nombhe.h"
116 #include "nombte.h"
117 #include "envca1.h"
118 #include "dicfen.h"
119 #include "nombsr.h"
120 c
121 c 0.3. ==> arguments
122 c
123       character*8 nohman, nhvois
124       character*8 nhnoeu, nharet, nhtria, nhquad
125 c
126       integer noempo(nbmpto)
127       integer hetnoe(nbnoto), arenoe(nbnoto)
128       integer coexno(nbnoto,nctfno)
129       integer nnosho(rsnoac), nnosca(rsnoto)
130       integer somare(2,nbarto), hetare(nbarto), np2are(nbarto)
131       integer filare(nbarto), merare(nbarto), insoar(nbarto)
132       integer coexar(nbarto,nctfar)
133       integer narsho(rsarac), narsca(rsarto)
134       integer hettri(nbtrto), aretri(nbtrto,3)
135       integer filtri(nbtrto), pertri(nbtrto)
136       integer hetqua(nbquto), arequa(nbquto,4)
137       integer filqua(nbquto), perqua(nbquto)
138       integer coexqu(nbquto,nctfqu)
139       integer nqusho(rsquac), nqusca(rsquto)
140       integer quahex(nbhecf,6), coquhe(nbhecf,6)
141 c
142       double precision coonoe(nbnoto,sdim)
143 c
144       integer ulsort, langue, codret
145 c
146 c 0.4. ==> variables locales
147 c
148       integer un
149       parameter ( un = 1 )
150       integer iaux
151       integer codre1, codre2, codre3, codre4
152       integer codre0
153       integer ppovos, pvoiso
154       integer pposif, pfacar
155       integer voarno, vofaar, vovoar, vovofa
156       integer nbanci, nbnoct, nbnocq
157       integer numead
158       integer nbgemx
159       integer conoct, conocq
160 c
161       integer adnoer, adarra, adarrb
162       integer adtrra, adtrrb
163       integer adqura, adqurb
164       integer ptrav1, ptrav2, ptrav3, ptrav4, ptrav5
165       integer ptrav6, ptrav7, ptrav8
166 c
167       character*8 nharrc, nhtrrc, nhqurc
168       character*8 ntrav1, ntrav2, ntrav3, ntrav4, ntrav5
169       character*8 ntrav6, ntrav7, ntrav8
170 c
171       integer nbmess
172       parameter ( nbmess = 10 )
173       character*80 texte(nblang,nbmess)
174 c
175 c 0.5. ==> initialisations
176 c ______________________________________________________________________
177 c
178 c====
179 c 1. preliminaires
180 c====
181 c
182 c 1.1. ==> messages
183 c
184 #include "impr01.h"
185 c
186 #ifdef _DEBUG_HOMARD_
187       write (ulsort,texte(langue,1)) 'Entree', nompro
188       call dmflsh (iaux)
189 #endif
190 c
191       texte(1,4) =
192      > '(''Nombre de paires de '',a,'' non-conformes :'',i10))'
193 c
194       texte(2,4) =
195      > '(''Number of pairs of non-conformal '',a,'' :'',i10))'
196 c
197       codret = 0
198 c
199 c====
200 c 2. Les structures
201 c====
202 c 2.1. ==> Les structures des voisins
203 c
204       if ( codret.eq.0 ) then
205 c
206       voarno = 2
207       vofaar = 2
208       vovofa = 0
209 c
210 #ifdef _DEBUG_HOMARD_
211       write (ulsort,texte(langue,3)) 'UTVOIS', nompro
212 #endif
213       call utvois ( nohman, nhvois,
214      >              voarno, vofaar, vovoar, vovofa,
215      >              ppovos, pvoiso,
216      >              nbfaar, pposif, pfacar,
217      >              ulsort, langue, codret )
218 c
219       endif
220 c
221 c====
222 c 3. Decompte et reperage des non-conformites
223 c====
224 c 3.1. ==> Tableau de travail
225 c
226       if ( codret.eq.0 ) then
227 c
228       call gmalot ( ntrav1, 'entier  ', nbarto, ptrav1, codre0 )
229       codret = max ( abs(codre0), codret )
230 c
231       endif
232 c
233 c 3.2. ==> Decompte et reperage des aretes en vis-a-vis
234 c
235       if ( codret.eq.0 ) then
236 c
237 #ifdef _DEBUG_HOMARD_
238       write (ulsort,texte(langue,3)) 'UTNC01', nompro
239 #endif
240 c
241       call utnc01 ( nbanci, nbgemx,
242      >              coonoe,
243      >              somare, merare,
244      >              aretri,
245      >              imem(ppovos), imem(pvoiso),
246      >              imem(ptrav1),
247      >              ulsort, langue, codret )
248 c
249 #ifdef _DEBUG_HOMARD_
250       write (ulsort,texte(langue,4)) mess14(langue,3,1), nbanci
251 #endif
252 c
253       endif
254 c
255 #ifdef _DEBUG_HOMARD_
256       call gmprot (nompro,ntrav1,1,nbarto)
257 #endif
258 c
259 c 3.3. ==> Enregistrement
260 c
261       if ( codret.eq.0 ) then
262 c
263       if ( nbanci.gt.0 ) then
264 c
265         maconf = 10
266         call gmecat ( nohman, 4, maconf, codret )
267 c
268       endif
269 c
270       endif
271 c
272 c====
273 c 4. Gestion des tableaux
274 c====
275 c
276 c 4.1. ==> Tableaux de travail
277 c
278       if ( nbanci.gt.0 ) then
279 c
280         if ( codret.eq.0 ) then
281 c
282         call gmalot ( ntrav2, 'entier  ', nbnoto, ptrav2, codre1 )
283         iaux = max( nbarto+1,nbnoto+1)
284         call gmalot ( ntrav4, 'entier  ', iaux, ptrav4, codre2 )
285         iaux = max(3*nbanci,nbnoto,2*nbarto,nbarto*nctfar)
286         call gmalot ( ntrav5, 'entier  ', iaux, ptrav5, codre3 )
287         iaux = nbnoto*sdim
288         call gmalot ( ntrav6, 'reel    ', iaux, ptrav6, codre4 )
289 c
290         codre0 = min ( codre1, codre2, codre3, codre4 )
291         codret = max ( abs(codre0), codret,
292      >                 codre1, codre2, codre3, codre4 )
293 c
294         endif
295 c
296 c 4.2. ==> Allocation de la memorisation des noeuds de non-conformite
297 c
298         if ( codret.eq.0 ) then
299 c
300         call gmecat ( nhnoeu, 4, nbanci, codre1 )
301         call gmaloj ( nhnoeu//'.Recollem', ' ', nbanci, adnoer, codre2 )
302 c
303         codre0 = min ( codre1, codre2 )
304         codret = max ( abs(codre0), codret,
305      >                 codre1, codre2 )
306 c
307         endif
308 c
309 c 4.3. ==> Memorisation des non-conformites en aretes
310 c
311         if ( codret.eq.0 ) then
312 c
313         iaux = 1
314 #ifdef _DEBUG_HOMARD_
315         write (ulsort,texte(langue,3)) 'VCMNC4_ar', nompro
316 #endif
317         call vcmnc4 ( iaux, nharet,
318      >                nharrc,
319      >                ulsort, langue, codret )
320 c
321         endif
322 c
323         if ( codret.eq.0 ) then
324 c
325         call gmecat ( nharrc, 1, nbanci, codre1 )
326         iaux = 2
327         call gmecat ( nharrc, 2, iaux, codre2 )
328         iaux = 2*nbanci
329         call gmaloj ( nharrc//'.ListeA', ' ', iaux, adarra, codre3 )
330         call gmaloj ( nharrc//'.ListeB', ' ', iaux, adarrb, codre4 )
331 c
332         codre0 = min ( codre1, codre2, codre3, codre4 )
333         codret = max ( abs(codre0), codret,
334      >                 codre1, codre2, codre3, codre4 )
335 c
336         endif
337 c
338 c 4.4. ==> Memorisation des non-conformites en triangles
339 c
340         if ( codret.eq.0 ) then
341 c
342         if ( nbtrto.ne.0 ) then
343 c
344           iaux = 2
345 #ifdef _DEBUG_HOMARD_
346         write (ulsort,texte(langue,3)) 'VCMNC4_tr', nompro
347 #endif
348           call vcmnc4 ( iaux, nhtria,
349      >                  nhtrrc,
350      >                  ulsort, langue, codret )
351 c
352         endif
353 c
354         endif
355 c
356 c 4.5. ==> Memorisation des non-conformites en quadrangles
357 c
358         if ( codret.eq.0 ) then
359 c
360         if ( nbquto.ne.0 ) then
361 c
362           iaux = 4
363 #ifdef _DEBUG_HOMARD_
364         write (ulsort,texte(langue,3)) 'VCMNC4_qu', nompro
365 #endif
366           call vcmnc4 ( iaux, nhquad,
367      >                  nhqurc,
368      >                  ulsort, langue, codret )
369 c
370         endif
371 c
372         endif
373 c
374       endif
375 c
376 c====
377 c 5. Renumerotation des aretes
378 c====
379 c
380       if ( codret.eq.0 ) then
381 c
382       if ( nbanci.gt.0 ) then
383 c
384 #ifdef _DEBUG_HOMARD_
385         write (ulsort,texte(langue,3)) 'VCMNC1', nompro
386 #endif
387         call vcmnc1 ( nbanci, nbgemx,
388      >                imem(adarra), imem(adarrb),
389      >                nohman, nhvois,
390      >                arenoe,
391      >                somare, hetare, np2are,
392      >                merare, filare, insoar,
393      >                coexar, narsho, narsca,
394      >                aretri, arequa,
395      >                ppovos, pvoiso,
396      >                pposif, pfacar,
397      >                imem(ptrav1), imem(ptrav4), imem(ptrav5),
398      >                ulsort, langue, codret )
399 c
400 #ifdef _DEBUG_HOMARD_
401         call gmprsx (nompro, nharrc//'.ListeA' )
402         call gmprsx (nompro, nharrc//'.ListeB' )
403 cgn        call gmprot (nompro,ntrav4,1,nbarto+1)
404 cgn        call gmprot (nompro,ntrav5,1,3*nbanci)
405 #endif
406 c
407 cgn        return
408       endif
409 c
410       endif
411 c
412 c====
413 c 6. Renumerotation des noeuds
414 c====
415 c
416       if ( codret.eq.0 ) then
417 c
418       if ( nbanci.gt.0 ) then
419 c
420 #ifdef _DEBUG_HOMARD_
421         write (ulsort,texte(langue,3)) 'VCMNC2', nompro
422 #endif
423         call vcmnc2 ( nbanci, nbgemx,
424      >                imem(adarra), imem(adarrb), imem(adnoer),
425      >                nohman, nhvois,
426      >                coonoe, hetnoe, arenoe,
427      >                coexno, nnosho, nnosca,
428      >                noempo,
429      >                somare, hetare, np2are,
430      >                merare, filare, insoar,
431      >                coexar, narsho, narsca,
432      >                aretri, arequa,
433      >                ppovos, pvoiso,
434      >                pposif, pfacar,
435      >                imem(ptrav1), imem(ptrav2),
436      >                imem(ptrav4), imem(ptrav5),
437      >                rmem(ptrav6),
438      >                ulsort, langue, codret )
439 c
440 #ifdef _DEBUG_HOMARD_
441         call gmprot (nompro,ntrav4,1,nbnoto+1)
442 #endif
443 c
444       endif
445 c
446       endif
447 c
448 c====
449 c 7. On repere chaque face du macro maillage qui est bordee par une
450 c    arete de non conformite initiale. On declare que cette face a une
451 c    mere, dont le numero est un numero fictif, ne correspondant a
452 c    aucune face possible.
453 c====
454 c
455       if ( codret.eq.0 ) then
456 c
457       if ( nbanci.gt.0 ) then
458 c
459 cgn      do iaux=1,nbquto
460 cgn      write(ulsort,*) iaux,arequa(iaux,1),arequa(iaux,2),
461 cgn     >                 arequa(iaux,3),arequa(iaux,4)
462 cgn      enddo
463 c
464         if ( codret.eq.0 ) then
465 c
466 #ifdef _DEBUG_HOMARD_
467         write (ulsort,texte(langue,3)) 'UTNC09', nompro
468 #endif
469 c
470         iaux = 0
471         call utnc09 ( nbanci, imem(adarrb), iaux,
472      >                pertri, perqua,
473      >                numead,
474      >                imem(pposif), imem(pfacar),
475      >                ulsort, langue, codret )
476 c
477         endif
478 c
479         if ( codret.eq.0 ) then
480 c
481         if ( nbtrto.gt.0 ) then
482 c
483         call gmecat ( nhtrrc, 3, numead, codre0 )
484         codret = max ( abs(codre0), codret )
485 c
486         endif
487 c
488         if ( nbquto.gt.0 ) then
489 c
490         call gmecat ( nhqurc, 3, numead, codre0 )
491         codret = max ( abs(codre0), codret )
492 c
493         endif
494 c
495         endif
496 cgn        return
497 c
498       endif
499 c
500       endif
501 c
502 c====
503 c 8. Gestion des non-conformites sur les volumes
504 c====
505 c
506       if ( codret.eq.0 ) then
507 c
508       if ( nbanci.gt.0 .and.
509      >     ( nbheto.gt.0 .or. nbteto.gt.0 ) ) then
510 c
511 c 8.1. ==> Reperage des faces dont des aretes sont non-conformes :
512 c          . les 4 aretes pour des quadrangles
513 c          . les 3 aretes pour les triangles
514 c          Les nombres nbnoct et nbnocq ne sont que des maxima : une
515 c          face peut tres bien avoir ses aretes coupees sans etre
516 c          decoupee elle-meme.
517 c
518         if ( codret.eq.0 ) then
519 c
520 #ifdef _DEBUG_HOMARD_
521         write (ulsort,texte(langue,3)) 'UTNC11', nompro
522 #endif
523 c
524         call utnc11 ( nbanci, imem(adarra),
525      >                aretri, filtri,
526      >                arequa, filqua,
527      >                filare, imem(pposif), imem(pfacar),
528      >                nbnoct, nbnocq,
529      >                ulsort, langue, codret )
530 c
531         endif
532 c
533 c 8.2. ==> Les tableaux
534 c 8.2.1. ==> Allocation de la memorisation des faces en vis-a-vis
535 c            On est dans un rapport de 4 pour 1 toujours.
536 c
537         if ( codret.eq.0 ) then
538 c
539         if ( nbquto.ne.0 ) then
540 c
541           call gmecat ( nhqurc, 1, nbnocq, codre1 )
542           iaux = 4
543           call gmecat ( nhqurc, 2, iaux, codre2 )
544           nbquri = 4*nbnocq
545           call gmaloj ( nhqurc//'.ListeA', ' ', nbquri, adqura, codre3 )
546           call gmaloj ( nhqurc//'.ListeB', ' ', nbquri, adqurb, codre4 )
547 c
548           codre0 = min ( codre1, codre2, codre3, codre4 )
549           codret = max ( abs(codre0), codret,
550      >                   codre1, codre2, codre3, codre4 )
551 c
552         endif
553 c
554         if ( nbtrto.ne.0 ) then
555 c
556           call gmecat ( nhtrrc, 1, nbnoct, codre1 )
557           iaux = 4
558           call gmecat ( nhtrrc, 2, iaux, codre2 )
559           nbtrri = 4*nbnoct
560           call gmaloj ( nhtrrc//'.ListeA', ' ', nbtrri, adtrra, codre3 )
561           call gmaloj ( nhtrrc//'.ListeB', ' ', nbtrri, adtrrb, codre4 )
562 c
563           codre0 = min ( codre1, codre2, codre3, codre4 )
564           codret = max ( abs(codre0), codret,
565      >                   codre1, codre2, codre3, codre4 )
566 c
567         endif
568 c
569         endif
570 c
571 c 8.2.2. ==> Tableaux de travail
572 c
573         if ( codret.eq.0 ) then
574 c
575         iaux = nbquto + nbtrto + 1
576         call gmalot ( ntrav3, 'entier  ', iaux, ptrav3, codre1 )
577         iaux = max( nbtrto+1,nbquto+1)
578         call gmalot ( ntrav7, 'entier  ', iaux, ptrav7, codre2 )
579         iaux =
580      >  max(5*nbnocq,5*nbnoct,4*nbquto,nbquto*nctfqu,rsquac,rsquto,
581      >      6*nbheto)
582         call gmalot ( ntrav8, 'entier  ', iaux, ptrav8, 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 8.3. ==> Repérage des faces en vis-a-vis
591 c
592         if ( codret.eq.0 ) then
593 c
594 #ifdef _DEBUG_HOMARD_
595       write (ulsort,texte(langue,3)) 'UTNC12', nompro
596 #endif
597 c
598         call utnc12 ( hettri, aretri, filtri, pertri,
599      >                hetqua, arequa, filqua, perqua,
600      >                filare, imem(pposif), imem(pfacar),
601      >                nbnocq, imem(adqura), imem(adqurb), conocq,
602      >                nbnoct, imem(adtrra), imem(adtrrb), conoct,
603      >                ulsort, langue, codret )
604 c
605 #ifdef _DEBUG_HOMARD_
606         call gmprsx (nompro, nhqurc )
607         call gmprot (nompro, nhqurc//'.ListeA', 1, 10 )
608         call gmprot (nompro, nhqurc//'.ListeA', 4*conocq, nbquri )
609         call gmprot (nompro, nhqurc//'.ListeB', 1, 10 )
610         call gmprot (nompro, nhqurc//'.ListeB', 4*conocq, nbquri )
611 #endif
612 c
613         endif
614 c
615 c 8.4 ==> Redimensionnement des tableaux
616 #ifdef _DEBUG_HOMARD_
617         write (ulsort,*) 'debut de 8.4 avec codret = ', codret
618 #endif
619 c
620        if ( codret.eq.0 ) then
621 c
622         if ( nbquto.ne.0 ) then
623 c
624           nbquri = 4*nbnocq
625           iaux = 4*conocq
626 c
627           call gmmod ( nhqurc//'.ListeA',
628      >                 adqura, nbquri, iaux, un, un, codre1 )
629           call gmmod ( nhqurc//'.ListeB',
630      >                 adqurb, nbquri, iaux, un, un, codre2 )
631           nbquri = iaux
632           nbnocq = conocq
633           call gmecat ( nhqurc, 1, nbnocq, codre3 )
634 c
635           codre0 = min ( codre1, codre2, codre3 )
636           codret = max ( abs(codre0), codret,
637      >                   codre1, codre2, codre3 )
638 c
639         endif
640 c
641         if ( nbtrto.ne.0 ) then
642 c
643           nbtrri = 4*nbnoct
644           iaux = 4*conoct
645           call gmmod ( nhtrrc//'.ListeA',
646      >                 adtrra, nbtrri, iaux, un, un, codre1 )
647           call gmmod ( nhtrrc//'.ListeB',
648      >                 adtrrb, nbtrri, iaux, un, un, codre2 )
649           nbtrri = iaux
650           nbnoct = conoct
651           call gmecat ( nhtrrc, 1, nbnoct, codre3 )
652 c
653           codre0 = min ( codre1, codre2, codre3 )
654           codret = max ( abs(codre0), codret,
655      >                   codre1, codre2, codre3 )
656 c
657         endif
658 c
659 #ifdef _DEBUG_HOMARD_
660         call gmprsx (nompro, nhqurc )
661         call gmprsx (nompro, nhqurc//'.ListeA' )
662         call gmprsx (nompro, nhqurc//'.ListeB' )
663 #endif
664 c
665         endif
666 c
667 c 8.5 ==> Renumérotation
668 c
669        if ( codret.eq.0 ) then
670 c
671 #ifdef _DEBUG_HOMARD_
672         write (ulsort,texte(langue,3)) 'VCMNC3', nompro
673 #endif
674         call vcmnc3 ( nbanci,
675      >                nbnocq, imem(adqura), imem(adqurb),
676      >                nbnoct, imem(adtrra), imem(adtrrb),
677      >                nohman, nhvois,
678      >                coonoe, hetnoe, arenoe,
679      >                coexno, nnosho, nnosca,
680      >                noempo,
681      >                somare, filare,
682      >                hettri, aretri, filtri,
683      >                hetqua, arequa, filqua, perqua,
684      >                coexqu, nqusho, nqusca,
685      >                quahex, coquhe,
686      >                ppovos, pvoiso,
687      >                pposif, pfacar,
688      >                imem(ptrav7), imem(ptrav8), rmem(ptrav6),
689      >                ulsort, langue, codret )
690 c
691         endif
692 c
693       endif
694 c
695       endif
696 c
697 c====
698 c 9. Menage
699 c====
700 #ifdef _DEBUG_HOMARD_
701         write (ulsort,*) 'debut de 9 avec codret = ', codret
702 #endif
703 c
704       if ( codret.eq.0 ) then
705 c
706       call gmlboj ( ntrav1, codre0 )
707       codret = max ( abs(codre0), codret )
708 c
709       endif
710 c
711       if ( codret.eq.0 ) then
712 c
713       if ( nbanci.gt.0 ) then
714 c
715         if ( codret.eq.0 ) then
716 c
717         call gmlboj ( ntrav2, codre1 )
718         call gmlboj ( ntrav4, codre2 )
719         call gmlboj ( ntrav5, codre3 )
720         call gmlboj ( ntrav6, codre4 )
721 c
722         codre0 = min ( codre1, codre2, codre3, codre4  )
723         codret = max ( abs(codre0), codret,
724      >                 codre1, codre2, codre3, codre4 )
725 c
726         endif
727 c
728         if ( nbheto.gt.0 .or. nbteto.gt.0 ) then
729 c
730           if ( codret.eq.0 ) then
731 c
732           call gmlboj ( ntrav7 , codre1 )
733           call gmlboj ( ntrav8 , codre2 )
734 c
735           codre0 = min ( codre1, codre2 )
736           codret = max ( abs(codre0), codret,
737      >                   codre1, codre2 )
738 c
739           endif
740 c
741         endif
742 c
743       endif
744 c
745       endif
746 c
747 c====
748 c 10. la fin
749 c====
750 c
751       if ( codret.ne.0 ) then
752 c
753 #include "envex2.h"
754 c
755       write (ulsort,texte(langue,1)) 'Sortie', nompro
756       write (ulsort,texte(langue,2)) codret
757 c
758       endif
759 c
760 #ifdef _DEBUG_HOMARD_
761       write (ulsort,texte(langue,1)) 'Sortie', nompro
762       call dmflsh (iaux)
763 #endif
764 c
765       end