]> SALOME platform Git repositories - modules/homard.git/blob - src/tool/ES_MED/eslnum.F
Salome HOME
Homard executable
[modules/homard.git] / src / tool / ES_MED / eslnum.F
1       subroutine eslnum ( idfmed, nomamd, degre,
2      >                    nbnoto, nbelem,
3      >                    nbmapo, nbsegm, nbtria, nbtetr,
4      >                    nbquad, nbhexa, nbpent, nbpyra,
5      >                    nunoex, nuelex,
6      >                    numano, numael,
7      >                    ulsort, langue, codret )
8 c ______________________________________________________________________
9 c
10 c                             H O M A R D
11 c
12 c Outil de Maillage Adaptatif par Raffinement et Deraffinement d'EDF R&D
13 c
14 c Version originale enregistree le 18 juin 1996 sous le numero 96036
15 c aupres des huissiers de justice Simart et Lavoir a Clamart
16 c Version 11.2 enregistree le 13 fevrier 2015 sous le numero 2015/014
17 c aupres des huissiers de justice
18 c Lavoir, Silinski & Cherqui-Abrahmi a Clamart
19 c
20 c    HOMARD est une marque deposee d'Electricite de France
21 c
22 c Copyright EDF 1996
23 c Copyright EDF 1998
24 c Copyright EDF 2002
25 c Copyright EDF 2020
26 c ______________________________________________________________________
27 c
28 c  Entree-Sortie - Lecture des NUMerotations
29 c  -      -        -           ---
30 c  Par defaut, on part du principe que les elements externes sont
31 c  numerotes dans cet ordre :
32 c             tetraedres, triangles, segments, mailles-points,
33 c             quadrangles, hexaedres, pyramides, pentaedres
34 c  Voir eslmm2 pour confirmation.
35 c
36 c  Si la table de renumerotation est fournie, on ecrase la
37 c  correspondance.
38 c ______________________________________________________________________
39 c .        .     .        .                                            .
40 c .  nom   . e/s . taille .           description                      .
41 c .____________________________________________________________________.
42 c . idfmed . e   .   1    . unite logique du maillage d'entree         .
43 c . nomamd . e   . char64 . nom du maillage MED                        .
44 c . degre  . e   .   1    . degre du maillage                          .
45 c . nbnoto . e   .   1    . nombre de noeuds dans le maillage          .
46 c . nbelem . e   .   1    . nombre d'elements dans le maillage         .
47 c . nbmapo . e   .   1    . nombre de mailles-points                   .
48 c . nbsegm . e   .   1    . nombre de segments                         .
49 c . nbtria . e   .   1    . nombre de triangles                        .
50 c . nbtetr . e   .   1    . nombre de tetraedres                       .
51 c . nbhexa . e   .   1    . nombre d'hexaedres                         .
52 c . nbpyra . e   .   1    . nombre de pyramides                        .
53 c . nbpent . e   .   1    . nombre de pentaedres                       .
54 c . nuelex .  s  . nbelem . numerotation des elements en exterieur     .
55 c . nunoex .  s  . nbnoto . numerotation des noeuds en exterieur       .
56 c . numano .  s  .   1    . numero maximum des noeuds                  .
57 c . numael .  s  .   1    . numero maximum des elements                .
58 c . ulsort . e   .   1    . numero d'unite logique de la liste standard.
59 c . langue . e   .    1   . langue des messages                        .
60 c .        .     .        . 1 : francais, 2 : anglais                  .
61 c . codret . es  .    1   . code de retour des modules                 .
62 c .        .     .        . 0 : pas de probleme                        .
63 c .        .     .        . 1 : probleme                               .
64 c ______________________________________________________________________
65 c
66 c====
67 c 0. declarations et dimensionnement
68 c====
69 c
70 c 0.1. ==> generalites
71 c
72       implicit none
73       save
74 c
75       character*6 nompro
76       parameter ( nompro = 'ESLNUM' )
77 c
78 #include "nblang.h"
79 #include "consts.h"
80 c
81 c 0.2. ==> communs
82 c
83 #include "envex1.h"
84 c
85 c 0.3. ==> arguments
86 c
87       integer*8 idfmed
88       integer degre
89       integer nbnoto, nbelem
90       integer nbmapo, nbsegm, nbtria, nbtetr
91       integer nbquad, nbhexa, nbpent, nbpyra
92       integer nunoex(nbnoto), nuelex(nbelem)
93       integer numano, numael
94 c
95       character*64 nomamd
96 c
97       integer ulsort, langue, codret
98 c
99 c 0.4. ==> variables locales
100 c
101 #include "meddc0.h"
102 c
103       integer iaux, jaux
104       integer typnoe, typseg, typtri, typqua
105       integer typtet, typhex, typpyr, typpen
106       integer ibtetr, ibtria, ibsegm, ibmapo
107       integer ibquad , ibhexa, ibpyra, ibpent
108       integer codre1, codre2, codre3, codre4, codre5
109       integer codre6, codre7, codre8, codre9
110       integer codre0
111       integer ntabno, ntabpo, ntabse, ntabtr, ntabqu
112       integer ntabte, ntabhe, ntabpy, ntabpe
113       integer numdt, numit
114       integer datype, chgt, tsf
115 c
116       integer nbmess
117       parameter ( nbmess = 10 )
118       character*80 texte(nblang,nbmess)
119 c ______________________________________________________________________
120 c
121 c====
122 c 1. initialisations
123 c====
124 c
125 c 1.1. ==> les messages
126 c
127 #include "impr01.h"
128 c
129 #ifdef _DEBUG_HOMARD_
130       write (ulsort,texte(langue,1)) 'Entree', nompro
131       call dmflsh (iaux)
132 #endif
133 c
134 c====
135 c 2. grandeurs de base
136 c====
137 c
138       typnoe = 0
139       if ( degre.eq.1 ) then
140         typseg = edseg2
141         typtri = edtri3
142         typqua = edqua4
143         typtet = edtet4
144         typhex = edhex8
145         typpyr = edpyr5
146         typpen = edpen6
147       else
148         typseg = edseg3
149         typtri = edtri6
150         typqua = edqua8
151         typtet = edte10
152         typhex = edhe20
153         typpyr = edpy13
154         typpen = edpe15
155       endif
156 c
157       ibtetr = 1
158       ibtria = nbtetr + 1
159       ibsegm = nbtetr + nbtria + 1
160       ibmapo = nbtetr + nbtria + nbsegm + 1
161       ibquad = nbtetr + nbtria + nbsegm + nbmapo + 1
162       ibhexa = nbtetr + nbtria + nbsegm + nbmapo + nbquad + 1
163       ibpyra = nbtetr + nbtria + nbsegm + nbmapo + nbquad + nbhexa + 1
164       ibpent = nbtetr + nbtria + nbsegm + nbmapo + nbquad + nbhexa
165      >       + nbpyra + 1
166 c
167       numdt = ednodt
168       numit = ednoit
169       datype = edda03
170 c
171 c====
172 c 3. les renumerotations
173 c====
174 #ifdef _DEBUG_HOMARD_
175       write (ulsort,*) '3. les renumerotations ; codret = ', codret
176 #endif
177 c
178 c 3.1. ==> initialisation a la non renumerotation
179 c
180       if ( codret.eq.0 ) then
181 c
182       do 311 , iaux = 1, nbnoto
183         nunoex(iaux) = iaux
184   311 continue
185       do 312 , iaux = 1, nbelem
186         nuelex(iaux) = iaux
187   312 continue
188       numano = nbnoto
189       numael = nbelem
190 c
191       endif
192 c
193 c 3.2. ==> le nombre de noeuds et de mailles a renumeroter
194 c
195       if ( codret.eq.0 ) then
196 c
197 #ifdef _DEBUG_HOMARD_
198       write (ulsort,texte(langue,3)) 'MMHNME_NO', nompro
199 #endif
200       call mmhnme ( idfmed, nomamd, numdt, numit,
201      >              ednoeu, typnoe, datype, ednoda, chgt, tsf,
202      >              ntabno, codre1 )
203       call mmhnme ( idfmed, nomamd, numdt, numit,
204      >              edmail, edpoi1, datype, ednoda, chgt, tsf,
205      >              ntabpo, codre2 )
206       call mmhnme ( idfmed, nomamd, numdt, numit,
207      >              edmail, typseg, datype, ednoda, chgt, tsf,
208      >              ntabse, codre3 )
209       call mmhnme ( idfmed, nomamd, numdt, numit,
210      >              edmail, typtri, datype, ednoda, chgt, tsf,
211      >              ntabtr, codre4 )
212       call mmhnme ( idfmed, nomamd, numdt, numit,
213      >              edmail, typqua, datype, ednoda, chgt, tsf,
214      >              ntabqu, codre5 )
215       call mmhnme ( idfmed, nomamd, numdt, numit,
216      >              edmail, typtet, datype, ednoda, chgt, tsf,
217      >              ntabte, codre6 )
218       call mmhnme ( idfmed, nomamd, numdt, numit,
219      >              edmail, typhex, datype, ednoda, chgt, tsf,
220      >              ntabhe, codre7 )
221       call mmhnme ( idfmed, nomamd, numdt, numit,
222      >              edmail, typpyr, datype, ednoda, chgt, tsf,
223      >              ntabpy, codre8 )
224       call mmhnme ( idfmed, nomamd, numdt, numit,
225      >              edmail, typpen, datype, ednoda, chgt, tsf,
226      >              ntabpe, codre9 )
227 c
228       codre0 = min ( codre1, codre2, codre3, codre4, codre5,
229      >               codre6, codre7, codre8, codre9 )
230       codret = max ( abs(codre0), codret,
231      >               codre1, codre2, codre3, codre4, codre5,
232      >               codre6, codre7, codre8, codre9 )
233       endif
234 c
235 #ifdef _DEBUG_HOMARD_
236       if ( codret.eq.0 ) then
237  1000 format(a,' = ',10i13)
238       write (ulsort,1000) 'ntabno', ntabno
239       write (ulsort,1000) 'ntabpo', ntabpo
240       write (ulsort,1000) 'ntabse', ntabse
241       write (ulsort,1000) 'ntabtr', ntabtr
242       write (ulsort,1000) 'ntabqu', ntabqu
243       write (ulsort,1000) 'ntabte', ntabte
244       write (ulsort,1000) 'ntabhe', ntabhe
245       write (ulsort,1000) 'ntabpy', ntabpy
246       write (ulsort,1000) 'ntabpe', ntabpe
247       else
248       write (ulsort,1000) 'codrei',
249      >               codre1, codre2, codre3, codre4, codre5,
250      >               codre6, codre7, codre8, codre9
251       endif
252 #endif
253 c
254 c 3.3. ==> les tables de renumerotation
255 c
256 c 3.3.1. ==> les noeuds
257 c
258       if ( codret.eq.0 ) then
259       if ( nbnoto.gt.0 .and. ntabno.eq.nbnoto ) then
260 #ifdef _DEBUG_HOMARD_
261       write (ulsort,texte(langue,3)) 'MMHENR_NO', nompro
262 #endif
263         call mmhenr ( idfmed, nomamd, numdt, numit,
264      >                ednoeu, typnoe, nunoex,
265      >                codret )
266         if ( codret.eq.0 ) then
267           do 331 , iaux = 1, nbnoto
268             numano = max(numano,nunoex(iaux))
269   331     continue
270         endif
271       endif
272       endif
273 c
274 c 3.3.2. ==> les mailles-points
275 c
276       if ( codret.eq.0 ) then
277       if ( nbmapo.gt.0 .and. ntabpo.eq.nbmapo ) then
278 #ifdef _DEBUG_HOMARD_
279       write (ulsort,texte(langue,3)) 'MMHENR_MP', nompro
280 #endif
281         call mmhenr ( idfmed, nomamd, numdt, numit,
282      >                edmail, edpoi1, nuelex(ibmapo),
283      >                codret )
284         if ( codret.eq.0 ) then
285           jaux = ibmapo + nbmapo - 1
286           do 332 , iaux = ibmapo , jaux
287             numael = max(numael,nuelex(iaux))
288   332     continue
289         endif
290       endif
291       endif
292 c
293 c 3.3.3. ==> les segments
294 c
295       if ( codret.eq.0 ) then
296       if ( nbsegm.gt.0 .and. ntabse.eq.nbsegm ) then
297 #ifdef _DEBUG_HOMARD_
298       write (ulsort,texte(langue,3)) 'MMHENR_AR', nompro
299 #endif
300         call mmhenr ( idfmed, nomamd, numdt, numit,
301      >                edmail, typseg, nuelex(ibsegm),
302      >                codret )
303         if ( codret.eq.0 ) then
304           jaux = ibsegm + nbsegm - 1
305           do 333 , iaux = ibsegm , jaux
306             numael = max(numael,nuelex(iaux))
307   333     continue
308         endif
309       endif
310       endif
311 c
312 c 3.3.4. ==> les triangles
313 c
314       if ( codret.eq.0 ) then
315       if ( nbtria.gt.0 .and. ntabtr.eq.nbtria ) then
316 #ifdef _DEBUG_HOMARD_
317       write (ulsort,texte(langue,3)) 'MMHENR_TR', nompro
318 #endif
319         call mmhenr ( idfmed, nomamd, numdt, numit,
320      >                edmail, typtri, nuelex(ibtria),
321      >                codret )
322         if ( codret.eq.0 ) then
323           jaux = ibtria + nbtria - 1
324           do 334 , iaux = ibtria, jaux
325             numael = max(numael,nuelex(iaux))
326   334     continue
327         endif
328       endif
329       endif
330 c
331 c 3.3.5. ==> les tetraedres
332 c
333       if ( codret.eq.0 ) then
334       if ( nbtetr.gt.0 .and. ntabte.eq.nbtetr ) then
335 #ifdef _DEBUG_HOMARD_
336       write (ulsort,texte(langue,3)) 'MMHENR_TE', nompro
337 #endif
338         call mmhenr ( idfmed, nomamd, numdt, numit,
339      >                edmail, typtet, nuelex(ibtetr),
340      >                codret )
341         if ( codret.eq.0 ) then
342           jaux = ibtetr + nbtetr - 1
343           do 335 , iaux = ibtetr, jaux
344             numael = max(numael,nuelex(iaux))
345   335     continue
346         endif
347       endif
348       endif
349 c
350 c 3.3.6. ==> les quadrangles
351 c
352       if ( codret.eq.0 ) then
353       if ( nbquad.gt.0 .and. ntabqu.eq.nbquad ) then
354 #ifdef _DEBUG_HOMARD_
355       write (ulsort,texte(langue,3)) 'MMHENR_QU', nompro
356 #endif
357         call mmhenr ( idfmed, nomamd, numdt, numit,
358      >                edmail, typqua, nuelex(ibquad),
359      >                codret )
360         if ( codret.eq.0 ) then
361           jaux = ibquad + nbquad - 1
362           do 336 , iaux = ibquad, jaux
363             numael = max(numael,nuelex(iaux))
364   336     continue
365         endif
366       endif
367       endif
368 c
369 c 3.3.7. ==> les hexaedres
370 c
371       if ( codret.eq.0 ) then
372       if ( nbhexa.gt.0 .and. ntabhe.eq.nbhexa ) then
373 #ifdef _DEBUG_HOMARD_
374       write (ulsort,texte(langue,3)) 'MMHENR_HE', nompro
375 #endif
376         call mmhenr ( idfmed, nomamd, numdt, numit,
377      >                edmail, typhex, nuelex(ibhexa),
378      >                codret )
379         if ( codret.eq.0 ) then
380           jaux = ibhexa + nbhexa - 1
381           do 337 , iaux = ibhexa, jaux
382             numael = max(numael,nuelex(iaux))
383   337     continue
384         endif
385       endif
386       endif
387 c
388 c 3.3.8. ==> les pyramides
389 c
390       if ( codret.eq.0 ) then
391       if ( nbpyra.gt.0 .and. ntabpy.eq.nbpyra ) then
392 #ifdef _DEBUG_HOMARD_
393       write (ulsort,texte(langue,3)) 'MMHENR_PY', nompro
394 #endif
395         call mmhenr ( idfmed, nomamd, numdt, numit,
396      >                edmail, typpyr, nuelex(ibpyra),
397      >                codret )
398         if ( codret.eq.0 ) then
399           jaux = ibpyra + nbpyra - 1
400           do 338 , iaux = ibpyra, jaux
401             numael = max(numael,nuelex(iaux))
402   338     continue
403         endif
404       endif
405       endif
406 c
407 c 3.3.9. ==> les pentaedres
408 c
409       if ( codret.eq.0 ) then
410       if ( nbpent.gt.0 .and. ntabpe.eq.nbpent ) then
411 #ifdef _DEBUG_HOMARD_
412       write (ulsort,texte(langue,3)) 'MMHENR_PE', nompro
413 #endif
414         call mmhenr ( idfmed, nomamd, numdt, numit,
415      >                edmail, typpen, nuelex(ibpent),
416      >                codret )
417         if ( codret.eq.0 ) then
418           jaux = ibpent + nbpent - 1
419           do 339 , iaux = ibpent, jaux
420             numael = max(numael,nuelex(iaux))
421   339     continue
422         endif
423       endif
424       endif
425 c
426 c====
427 c 4. la fin
428 c====
429 c
430       if ( codret.ne.0 ) then
431 c
432 #include "envex2.h"
433 c
434       write (ulsort,texte(langue,1)) 'Sortie', nompro
435       write (ulsort,texte(langue,2)) codret
436 c
437       endif
438 c
439 #ifdef _DEBUG_HOMARD_
440       write (ulsort,texte(langue,1)) 'Sortie', nompro
441       call dmflsh (iaux)
442 #endif
443 c
444       end