]> SALOME platform Git repositories - modules/homard.git/blob - src/tool/Creation_Maillage/cmraff.F
Salome HOME
Homard executable
[modules/homard.git] / src / tool / Creation_Maillage / cmraff.F
1       subroutine cmraff ( nomail,
2      >                    indnoe, indare, indtri, indqua,
3      >                    indtet, indhex, indpen,
4      >                    lgopts, taopts, lgetco, taetco,
5      >                    ulsort, langue, codret )
6 c ______________________________________________________________________
7 c
8 c                             H O M A R D
9 c
10 c Outil de Maillage Adaptatif par Raffinement et Deraffinement d'EDF R&D
11 c
12 c Version originale enregistree le 18 juin 1996 sous le numero 96036
13 c aupres des huissiers de justice Simart et Lavoir a Clamart
14 c Version 11.2 enregistree le 13 fevrier 2015 sous le numero 2015/014
15 c aupres des huissiers de justice
16 c Lavoir, Silinski & Cherqui-Abrahmi a Clamart
17 c
18 c    HOMARD est une marque deposee d'Electricite de France
19 c
20 c Copyright EDF 1996
21 c Copyright EDF 1998
22 c Copyright EDF 2002
23 c Copyright EDF 2020
24 c ______________________________________________________________________
25 c
26 c    Creation du Maillage - RAFFinement
27 c    -           -          ----
28 c ______________________________________________________________________
29 c .        .     .        .                                            .
30 c .  nom   . e/s . taille .           description                      .
31 c .____________________________________________________________________.
32 c . nomail . e   .  ch8   . nom de l'objet contenant le maillage       .
33 c . indnoe . es  .   1    . indice du dernier noeud cree               .
34 c . indare . es  .   1    . indice de la derniere arete creee          .
35 c . indtri . es  .   1    . indice du dernier triangle cree            .
36 c . indqua . es  .   1    . indice du dernier quadrangle cree          .
37 c . indtet . es  .   1    . indice du dernier tetraedre cree           .
38 c . indhex . es  .   1    . indice du dernier hexaedre cree            .
39 c . indpen . es  .   1    . indice du dernier pentaedre cree           .
40 c . lgopts . e   .   1    . longueur du tableau des options caracteres .
41 c . taopts . e   . lgopts . tableau des options caracteres             .
42 c . lgetco . e   .   1    . longueur du tableau de l'etat courant      .
43 c . taetco . e   . lgetco . tableau de l'etat courant                  .
44 c . ulsort . e   .   1    . numero d'unite logique de la liste standard.
45 c . langue . e   .    1   . langue des messages                        .
46 c .        .     .        . 1 : francais, 2 : anglais                  .
47 c . codret . es  .    1   . code de retour des modules                 .
48 c .        .     .        . 0 : pas de probleme                        .
49 c ______________________________________________________________________
50 c
51 c====
52 c 0. declarations et dimensionnement
53 c====
54 c
55 c 0.1. ==> generalites
56 c
57       implicit none
58       save
59 c
60       character*6 nompro
61       parameter ( nompro = 'CMRAFF' )
62 c
63 #include "nblang.h"
64 c
65 c 0.2. ==> communs
66 c
67 #include "envex1.h"
68 c
69 #include "gmenti.h"
70 #include "gmreel.h"
71 c
72 #include "envca1.h"
73 #include "nombno.h"
74 #include "nombar.h"
75 #include "nombtr.h"
76 #include "nombqu.h"
77 #include "nombte.h"
78 #include "nombhe.h"
79 #include "nombpe.h"
80 #include "nouvnb.h"
81 #include "impr02.h"
82 c
83 c 0.3. ==> arguments
84 c
85       character*8 nomail
86 c
87       integer indnoe, indare, indtri, indqua
88       integer indtet, indhex, indpen
89 c
90       integer lgopts
91       character*8 taopts(lgopts)
92 c
93       integer lgetco
94       integer taetco(lgetco)
95 c
96       integer ulsort, langue, codret
97 c
98 c 0.4. ==> variables locales
99 c
100       integer codava
101       integer nretap, nrsset
102       integer iaux, jaux
103 c
104       integer codre0, codre1, codre2
105       integer pdecar, pdecfa
106       integer phetno, pcoono, pareno
107       integer phetar, psomar, pfilar, pmerar, pnp2ar
108       integer phettr, paretr, pfiltr, ppertr, pnivtr, adnmtr
109       integer phetqu, parequ, pfilqu, pperqu, pnivqu, adnmqu
110       integer phette, ptrite, pcotrt, pfilte, pperte
111       integer phethe, pquahe, pcoquh, pfilhe, pperhe, adnmhe
112       integer phetpe, pfacpe, pcofap, pfilpe, pperpe
113       integer pfamno, pcfano
114       integer pfamar, pcfaar
115       integer pfamtr, pcfatr
116       integer pfamqu, pcfaqu
117       integer pfamte
118       integer pfamhe
119       integer pfampe
120 c
121       character*6 saux
122       character*8 norenu
123       character*8 nhnoeu, nhmapo, nharet, nhtria, nhquad
124       character*8 nhtetr, nhhexa, nhpyra, nhpent
125       character*8 nhelig
126       character*8 nhvois, nhsupe, nhsups
127       character*8 ntrav1, ntrav2
128 c
129       integer nbmess
130       parameter ( nbmess = 10 )
131       character*80 texte(nblang,nbmess)
132 c
133 c 0.5. ==> initialisations
134 c ______________________________________________________________________
135 c
136 c====
137 c 1. messages
138 c====
139 c
140       codava = codret
141 c
142 c=======================================================================
143       if ( codava.eq.0 ) then
144 c=======================================================================
145 c
146 c 1.3. ==> les messages
147 c
148 #include "impr01.h"
149 c
150 #ifdef _DEBUG_HOMARD_
151       write (ulsort,texte(langue,1)) 'Entree', nompro
152       call dmflsh (iaux)
153 #endif
154 c
155       texte(1,4) = '(/,a6,'' RAFFINEMENT STANDARD DU MAILLAGE'')'
156       texte(1,5) = '(39(''=''),/)'
157       texte(1,6) =
158      > '(5x,''Nombre de '',a,'' crees :'',i10,'' ; total : '',i10)'
159 c
160       texte(2,4) = '(/,a6,'' STANDARD REFINEMENT OF MESH'')'
161       texte(2,5) = '(34(''=''),/)'
162       texte(2,6) =
163      > '(5x,''Number of new '',a,'' :'',i10,'' ; total : '',i10)'
164 c
165 c 1.4. ==> le numero de sous-etape
166 c
167       nretap = taetco(1)
168       nrsset = taetco(2) + 1
169       taetco(2) = nrsset
170 c
171       call utcvne ( nretap, nrsset, saux, iaux, codret )
172 c
173 c 1.5. ==> le titre
174 c
175       write (ulsort,texte(langue,4)) saux
176       write (ulsort,texte(langue,5))
177 c
178 #include "impr03.h"
179 c
180 c====
181 c 2. recuperation des pointeurs
182 c====
183 c
184 c 2.1. ==> structure generale
185 c
186       if ( codret.eq.0 ) then
187 c
188 #ifdef _DEBUG_HOMARD_
189         write (ulsort,texte(langue,3)) 'UTNOMH', nompro
190 #endif
191       call utnomh ( nomail,
192      >                sdim,   mdim,
193      >               degre, maconf, homolo, hierar,
194      >              rafdef, nbmane, typcca, typsfr, maextr,
195      >              mailet,
196      >              norenu,
197      >              nhnoeu, nhmapo, nharet,
198      >              nhtria, nhquad,
199      >              nhtetr, nhhexa, nhpyra, nhpent,
200      >              nhelig,
201      >              nhvois, nhsupe, nhsups,
202      >              ulsort, langue, codret)
203 c
204       endif
205 cgn      call gmprot(nompro//' -0',nharet//'.Famille.EntiFamm',1,26)
206 cgn      call gmprot(nompro//' -0', nharet//'.Famille.EntiFamm',27,118)
207 c
208 c 2.2. ==> tableaux
209 c
210       if ( codret.eq.0 ) then
211 c
212 #ifdef _DEBUG_HOMARD_
213       write (ulsort,texte(langue,3)) 'UTAD01', nompro
214 #endif
215       iaux = 7770
216       call utad01 ( iaux, nhnoeu,
217      >              phetno,
218      >              pfamno, pcfano,   jaux,
219      >              pcoono, pareno,   jaux,  jaux,
220      >              ulsort, langue, codret )
221 c
222 #ifdef _DEBUG_HOMARD_
223       write (ulsort,texte(langue,3)) 'UTAD02_ar', nompro
224 #endif
225       iaux = 7770
226       if ( degre.eq.2 ) then
227         iaux = iaux*13
228       endif
229       call utad02 ( iaux, nharet,
230      >              phetar, psomar, pfilar, pmerar,
231      >              pfamar, pcfaar,   jaux,
232      >              jaux  , pnp2ar,   jaux,
233      >                jaux,   jaux,   jaux,
234      >              ulsort, langue, codret )
235 c
236       if ( nbtrto.ne.0 ) then
237 c
238 #ifdef _DEBUG_HOMARD_
239       write (ulsort,texte(langue,3)) 'UTAD02_tr', nompro
240 #endif
241         iaux = 85470
242         if ( mod(mailet,2).eq.0 ) then
243           iaux = iaux*19
244         endif
245         call utad02 ( iaux, nhtria,
246      >                phettr, paretr, pfiltr, ppertr,
247      >                pfamtr, pcfatr,   jaux,
248      >                pnivtr,   jaux,   jaux,
249      >                adnmtr,   jaux,   jaux,
250      >                ulsort, langue, codret )
251 c
252       endif
253 c
254       if ( nbquto.ne.0 ) then
255 c
256 #ifdef _DEBUG_HOMARD_
257       write (ulsort,texte(langue,3)) 'UTAD02_qu', nompro
258 #endif
259         iaux = 85470
260         if ( mod(mailet,3).eq.0 ) then
261           iaux = iaux*19
262         endif
263         call utad02 ( iaux, nhquad,
264      >                phetqu, parequ, pfilqu, pperqu,
265      >                pfamqu, pcfaqu,   jaux,
266      >                pnivqu,   jaux,   jaux,
267      >                adnmqu,   jaux,   jaux,
268      >                ulsort, langue, codret )
269 c
270       endif
271 c
272       if ( nbteto.ne.0 ) then
273 c
274 #ifdef _DEBUG_HOMARD_
275       write (ulsort,texte(langue,3)) 'UTAD02_te', nompro
276 #endif
277         iaux = 2730
278         call utad02 ( iaux, nhtetr,
279      >                phette, ptrite, pfilte, pperte,
280      >                pfamte,   jaux,   jaux,
281      >                jaux  , pcotrt,   jaux,
282      >                  jaux,   jaux,   jaux,
283      >                ulsort, langue, codret )
284 c
285       endif
286 c
287       if ( nbheto.ne.0 ) then
288 c
289 #ifdef _DEBUG_HOMARD_
290       write (ulsort,texte(langue,3)) 'UTAD02_he', nompro
291 #endif
292         iaux = 2730
293         if ( mod(mailet,5).eq.0 ) then
294           iaux = iaux*19
295         endif
296         call utad02 ( iaux, nhhexa,
297      >                phethe, pquahe, pfilhe, pperhe,
298      >                pfamhe,   jaux,   jaux,
299      >                jaux  , pcoquh,   jaux,
300      >                adnmhe,   jaux,   jaux,
301      >                ulsort, langue, codret )
302 c
303       endif
304 c
305       if ( nbpeto.ne.0 ) then
306 c
307 #ifdef _DEBUG_HOMARD_
308       write (ulsort,texte(langue,3)) 'UTAD02_pe', nompro
309 #endif
310         iaux = 2730
311         call utad02 ( iaux, nhpent,
312      >                phetpe, pfacpe, pfilpe, pperpe,
313      >                pfampe,   jaux,   jaux,
314      >                jaux  , pcofap,   jaux,
315      >                  jaux,   jaux,   jaux,
316      >                ulsort, langue, codret )
317 c
318       endif
319 c
320       endif
321 c
322       if ( codret.eq.0 ) then
323 c
324       ntrav1 = taopts(11)
325       call gmadoj ( ntrav1, pdecar, iaux, codre1 )
326       ntrav2 = taopts(12)
327       call gmadoj ( ntrav2, pdecfa, iaux, codre2 )
328 c
329       codre0 = min ( codre1, codre2 )
330       codret = max ( abs(codre0), codret,
331      >               codre1, codre2 )
332 cgn      call gmprsx (nompro,ntrav1)
333 cgn      call gmprsx (nompro,ntrav2)
334 c
335       endif
336 c
337 c====
338 c 3. decoupage des aretes en degre 1 et degre 2
339 c====
340 c
341 #ifdef _DEBUG_HOMARD_
342       write (ulsort,90002) '3. aretes ; codret', codret
343       write (ulsort,90002) 'indare', indare
344 #endif
345 c
346 #ifdef _DEBUG_HOMARD_
347       do 3011 , iaux = 1 , nbarto
348       if ( iaux.eq.-1661 ) then
349         write (ulsort,9001) 'arete', iaux,
350      >  imem(psomar+2*iaux-2), imem(psomar+2*iaux-1),
351      >  imem(pfilar+iaux-1),imem(phetar+iaux-1),imem(pdecar+iaux)
352       endif
353  3011 continue
354  9001 format(a,i8,', som :',2i8,', fille',i8,', e/d',i3,3i2)
355 #endif
356 c
357       if ( codret.eq.0 ) then
358 c
359       if ( degre.eq.1 ) then
360 c
361 #ifdef _DEBUG_HOMARD_
362         write (ulsort,texte(langue,3)) 'CMRDA1', nompro
363 #endif
364         call cmrda1
365      >         ( rmem(pcoono), imem(phetno), imem(pareno), imem(psomar),
366      >           imem(phetar), imem(pfilar), imem(pmerar), imem(pdecar),
367      >           imem(pcfaar), imem(pfamar), imem(pfamno),
368      >           indnoe, indare,
369      >           ulsort, langue, codret )
370 c
371       elseif ( degre.eq.2 ) then
372 c
373 #ifdef _DEBUG_HOMARD_
374         write (ulsort,texte(langue,3)) 'CMRDA2', nompro
375 #endif
376         call cmrda2
377      >         ( imem(phetno), imem(psomar), imem(phetar), imem(pfilar),
378      >           imem(pmerar), imem(pdecar), imem(pnp2ar),
379      >           imem(pcfaar), imem(pfamar),
380      >           indare,
381      >           ulsort, langue, codret )
382 c
383       else
384         write(ulsort,90050) degre
385         codret = 5
386       endif
387 c
388       endif
389 #ifdef _DEBUG_HOMARD_
390 cgn      do 30012 , iaux = 1 , nbarto
391 cgn      if ( iaux.eq.478 ) then
392 cgn        write (ulsort,9001) 'arete', iaux,
393 cgn     >  imem(psomar+2*iaux-2), imem(psomar+2*iaux-1),
394 cgn     >   imem(pfilar+iaux-1),imem(phetar+iaux-1),imem(pdecar+iaux)
395 cgn      endif
396 cgn30012 continue
397 #endif
398 cgn      write(ulsort,90002) 'indare', indare
399 cgn      call gmprot(nompro//' -1',nharet//'.Famille.EntiFamm',1,26)
400 cgn      call gmprot(nompro//' -1', nharet//'.Famille.EntiFamm',27,118)
401 c
402 c====
403 c 4. decoupage des triangles
404 c====
405 c
406 #ifdef _DEBUG_HOMARD_
407       write (ulsort,90002) '4. triangles ; codret', codret
408 #endif
409 c
410       if ( codret.eq.0 ) then
411 c
412       if ( nbtrto.ne.0 ) then
413 c
414 #ifdef _DEBUG_HOMARD_
415       write (ulsort,texte(langue,3)) 'CMRDTR', nompro
416 #endif
417       call cmrdtr
418      >         ( imem(psomar), imem(phetar), imem(pfilar), imem(pmerar),
419      >           imem(paretr), imem(phettr), imem(pfiltr), imem(ppertr),
420      >           imem(pnivtr), imem(pdecfa),
421      >           imem(pfamar), imem(pfamtr),
422      >           indare, indtri,
423      >           imem(pcfatr),
424      >           ulsort, langue, codret )
425 c
426       endif
427 c
428       endif
429 c
430 c====
431 c 5. decoupage des quadrangles
432 c====
433 c
434 #ifdef _DEBUG_HOMARD_
435       write (ulsort,90002) '5. quadrangles ; codret', codret
436 #endif
437 c
438       if ( codret.eq.0 ) then
439 c
440       if ( nbquto.ne.0 ) then
441 c
442 #ifdef _DEBUG_HOMARD_
443       write (ulsort,texte(langue,3)) 'CMRDQU', nompro
444 #endif
445       call cmrdqu
446      >         ( rmem(pcoono), imem(phetno), imem(pareno),
447      >           imem(psomar), imem(phetar), imem(pfilar), imem(pmerar),
448      >           imem(parequ), imem(phetqu), imem(pfilqu), imem(pperqu),
449      >           imem(pnivqu), imem(adnmqu), imem(pdecfa),
450      >           imem(pfamno), imem(pfamar), imem(pfamqu),
451      >           indnoe, indare, indqua,
452      >           imem(pcfaqu),
453      >           ulsort, langue, codret )
454 c
455       endif
456 c
457       endif
458 c
459 c====
460 c 6. decoupage des tetraedres
461 c====
462 c
463 #ifdef _DEBUG_HOMARD_
464       write (ulsort,90002) '6. tetraedres ; codret', codret
465 #endif
466 c
467       if ( codret.eq.0 ) then
468 c
469       if ( nbteto.ne.0 ) then
470 c
471 #ifdef _DEBUG_HOMARD_
472       write (ulsort,texte(langue,3)) 'CMRDTE', nompro
473 #endif
474         call cmrdte
475      >         ( rmem(pcoono),
476      >           imem(psomar), imem(phetar), imem(pfilar), imem(pmerar),
477      >           imem(paretr), imem(phettr),
478      >           imem(pfiltr), imem(ppertr), imem(pnivtr),
479      >           imem(ptrite), imem(pcotrt), imem(phette), imem(pfilte),
480      >           imem(pperte),
481      >           imem(pfamar), imem(pfamtr), imem(pfamte),
482      >           indare, indtri, indtet,
483      >           ulsort, langue, codret )
484 c
485       endif
486 c
487       endif
488 c
489 c====
490 c 7. decoupage des hexaedres
491 c====
492 c
493 #ifdef _DEBUG_HOMARD_
494       write (ulsort,90002) '7. hexaedres ; codret', codret
495 #endif
496 c
497       if ( codret.eq.0 ) then
498 c
499       if ( nbheto.ne.0 ) then
500 c
501 #ifdef _DEBUG_HOMARD_
502       write (ulsort,texte(langue,3)) 'CMRDHE', nompro
503 #endif
504         call cmrdhe
505      >         ( rmem(pcoono), imem(phetno), imem(pareno),
506      >           imem(psomar), imem(phetar), imem(pfilar), imem(pmerar),
507      >           imem(parequ), imem(phetqu),
508      >           imem(pfilqu), imem(pperqu), imem(pnivqu),
509      >           imem(pquahe), imem(pcoquh), imem(phethe),
510      >           imem(pfilhe), imem(pperhe), imem(adnmhe),
511      >           imem(pfamno), imem(pfamar), imem(pfamqu), imem(pfamhe),
512      >           indnoe, indare, indqua, indhex,
513      >           ulsort, langue, codret )
514 c
515       endif
516 c
517       endif
518 c
519 c====
520 c 8 decoupage des pentaedres
521 c====
522 c
523 #ifdef _DEBUG_HOMARD_
524       write (ulsort,90002) '8. pentaedres ; codret', codret
525 #endif
526 c
527       if ( codret.eq.0 ) then
528 c
529       if ( nbpeto.ne.0 ) then
530 c
531 #ifdef _DEBUG_HOMARD_
532       write (ulsort,texte(langue,3)) 'CMRDPE', nompro
533 #endif
534         call cmrdpe
535      >         ( imem(psomar), imem(phetar), imem(pfilar), imem(pmerar),
536      >           imem(paretr), imem(phettr),
537      >           imem(pfiltr), imem(ppertr), imem(pnivtr),
538      >           imem(parequ), imem(phetqu),
539      >           imem(pfilqu), imem(pperqu), imem(pnivqu),
540      >           imem(pfacpe), imem(pcofap), imem(phetpe),
541      >           imem(pfilpe), imem(pperpe),
542      >           imem(pfamar), imem(pfamtr), imem(pfamqu), imem(pfampe),
543      >           indare, indtri, indqua, indpen,
544      >           ulsort, langue, codret )
545 c
546       endif
547 c
548       endif
549 c
550 c====
551 c 9. verifications des nombres d'entites crees et impressions
552 c====
553 c
554 #ifdef _DEBUG_HOMARD_
555       write (ulsort,90002) '9. verifications ; codret', codret
556 #endif
557 c
558       if ( codret.eq.0 ) then
559 c
560       write(ulsort,texte(langue,6)) mess14(langue,3,-1),
561      >                              indnoe-nbnoto, indnoe
562       write(ulsort,texte(langue,6)) mess14(langue,3,1),
563      >                              indare-nbarto, indare
564       if ( nbtrto.ne.0 ) then
565         write(ulsort,texte(langue,6)) mess14(langue,3,2),
566      >                                indtri-nbtrto, indtri
567       endif
568       if ( nbquto.ne.0 ) then
569         write(ulsort,texte(langue,6)) mess14(langue,3,4),
570      >                                indqua-nbquto, indqua
571       endif
572       if ( nbteto.ne.0 ) then
573         write(ulsort,texte(langue,6)) mess14(langue,3,3),
574      >                                indtet-nbteto, indtet
575       endif
576       if ( nbheto.ne.0 ) then
577         write(ulsort,texte(langue,6)) mess14(langue,3,6),
578      >                                indhex-nbheto, indhex
579       endif
580       if ( nbpeto.ne.0 ) then
581         write(ulsort,texte(langue,6)) mess14(langue,3,7),
582      >                                indpen-nbpeto, indpen
583       endif
584 c
585       iaux = 0
586       if ( degre.eq.1 ) then
587         if ( indnoe.ne.nouvno ) then
588           write(ulsort,90100) mess14(langue,3,-1), permno-nbnoto
589           iaux = iaux + 1
590         endif
591       endif
592       if ( indare.ne.permar ) then
593         write(ulsort,90100) mess14(langue,3,1), permar-nbarto
594         iaux = iaux + 1
595       endif
596       if ( nbtrto.ne.0 .and. indtri.ne.permtr ) then
597         write(ulsort,90100) mess14(langue,3,2), permtr-nbtrto
598         iaux = iaux + 1
599       endif
600       if ( nbquto.ne.0 .and. indqua.ne.permqu ) then
601         write(ulsort,90100) mess14(langue,3,4), permqu-nbquto
602         iaux = iaux + 1
603       endif
604       if ( nbteto.ne.0 .and. indtet.ne.permte ) then
605         write(ulsort,90100) mess14(langue,3,3), permte-nbteto
606         iaux = iaux + 1
607       endif
608       if ( nbheto.ne.0 .and. indhex.ne.permhe ) then
609         write(ulsort,90100) mess14(langue,3,6), permhe-nbheto
610         iaux = iaux + 1
611       endif
612       if ( nbpeto.ne.0 .and. indpen.ne.permpe ) then
613         write(ulsort,90100) mess14(langue,3,7), permpe-nbpeto
614         iaux = iaux + 1
615       endif
616 c
617       if ( iaux.ne.0 ) then
618         codret = 4
619       endif
620 c
621       endif
622 cgn      call gmprsx(nompro,nhnoeu)
623 c
624 90050 format(/,5x,'Le degre d''interpolation est : ',i1,
625      >       /,5x,'Seules les discretisations de degres 1 et 2 sont ',
626      >            'supportees',/)
627 90100 format(/,5x,'Nombre de ',a,' crees incorrect, ',
628      >       /,5x,'On en attendait ',i10,'  nouveaux.')
629 c
630 c====
631 c 10. la fin
632 c====
633 c
634       if ( codret.ne.0 ) then
635 c
636 #include "envex2.h"
637 c
638       write (ulsort,texte(langue,1)) 'Sortie', nompro
639       write (ulsort,texte(langue,2)) codret
640 c
641       endif
642 c
643 #ifdef _DEBUG_HOMARD_
644       write (ulsort,texte(langue,1)) 'Sortie', nompro
645       call dmflsh (iaux)
646 #endif
647 c
648 c=======================================================================
649       endif
650 c=======================================================================
651 c
652       end