Salome HOME
Homard executable
[modules/homard.git] / src / tool / AP_Conversion / pcma21.F
1       subroutine pcma21 ( choixd, deltac,
2      >                    nbnoto, nbelem, famnoe, coonoe,
3      >                    famnzz, nbno3d,
4      >                    typele, noeele,
5      >                    nu3dno, famn3d, coon3d,
6      >                    ulsort, langue, codret )
7 c ______________________________________________________________________
8 c
9 c                             H O M A R D
10 c
11 c Outil de Maillage Adaptatif par Raffinement et Deraffinement d'EDF R&D
12 c
13 c Version originale enregistree le 18 juin 1996 sous le numero 96036
14 c aupres des huissiers de justice Simart et Lavoir a Clamart
15 c Version 11.2 enregistree le 13 fevrier 2015 sous le numero 2015/014
16 c aupres des huissiers de justice
17 c Lavoir, Silinski & Cherqui-Abrahmi a Clamart
18 c
19 c    HOMARD est une marque deposee d'Electricite de France
20 c
21 c Copyright EDF 1996
22 c Copyright EDF 1998
23 c Copyright EDF 2002
24 c Copyright EDF 2020
25 c ______________________________________________________________________
26 c
27 c    aPres adaptation - Conversion de MAillage - 2D/3D - phase 1
28 c     -                 -             --         -             -
29 c    Creation des noeuds supplementaires
30 c ______________________________________________________________________
31 c .        .     .        .                                            .
32 c .  nom   . e/s . taille .           description                      .
33 c .____________________________________________________________________.
34 c . choixd . e   .   1    . choix sur le calcul de delta coordonnes :  .
35 c .        .     .        . 1 : coordonnees initiales (defaut)         .
36 c .        .     .        . 2 : valeur imposee                         .
37 c .        .     .        . 3 : moyenne arithmetique des mini/maxi en  .
38 c .        .     .        .     (x,y) des mailles                      .
39 c .        .     .        . 4 : moyenne geometrique des mini/maxi en   .
40 c .        .     .        .     (x,y) des mailles                      .
41 c .        .     .        . 5 : ecart initial, divise par 2**nivsup    .
42 c . deltac . e   .   1    . valeur de delta si impose (choixd=1)       .
43 c . nbnoto . e   .   1    . nombre de noeuds du maillage externe       .
44 c . nbelem . e   .   1    . nombre d'elements du maillage externe      .
45 c . famnoe . e   . nbnoto . famille des noeuds                         .
46 c . coonoe . e   . nbnoto . coordonnees des noeuds                     .
47 c .        .     . * sdim .                                            .
48 c . famnzz . e   .   1    . famille du noeud memorisant cooinf et coosup   .
49 c . nbno3d . e   .   1    . nombre de noeuds du maillage 3d            .
50 c . typele . e   . nbelem . type des elements pour le code de calcul   .
51 c . noeele . e   . nbelem . noeuds des elements                        .
52 c .        .     .*nbmane .                                            .
53 c . nu3dno .  s  . nbnoto . numero du calcul des noeuds                .
54 c . famn3d .  s  . nbno3d . famille des noeuds  du maillage 3d         .
55 c . coon3d .  s  .nbno3d*3. coordonnees des noeuds du maillage 3d      .
56 c . ulsort . e   .   1    . numero d'unite logique de la liste standard.
57 c . langue . e   .    1   . langue des messages                        .
58 c .        .     .        . 1 : francais, 2 : anglais                  .
59 c . codret . es  .    1   . code de retour des modules                 .
60 c .        .     .        . 0 : pas de probleme                        .
61 c .        .     .        . 1 : probleme                               .
62 c ______________________________________________________________________
63 c
64 c====
65 c 0. declarations et dimensionnement
66 c====
67 c
68 c 0.1. ==> generalites
69 c
70       implicit none
71       save
72 c
73       character*6 nompro
74       parameter ( nompro = 'PCMA21' )
75 c
76 #include "nblang.h"
77 #include "consts.h"
78 c
79 c 0.2. ==> communs
80 c
81 #include "envex1.h"
82 c
83 #include "meddc0.h"
84 #include "envca1.h"
85 #include "envada.h"
86 c
87 c 0.3. ==> arguments
88 c
89       integer choixd
90       integer nbnoto, nbelem, famnzz, nbno3d
91       integer nu3dno(nbnoto)
92       integer famnoe(nbnoto), famn3d(nbno3d)
93       integer typele(nbelem), noeele(nbelem,*)
94 c
95       double precision deltac
96       double precision coon3d(nbno3d,3)
97       double precision coonoe(nbnoto,sdim)
98 c
99       integer ulsort, langue, codret
100 c
101 c 0.4. ==> variables locales
102 c
103       integer iaux, jaux, kaux
104       integer iaux1, iaux2
105 c
106       double precision cooinf, coosup
107       double precision daux, daux1, daux2
108       double precision minx, miny, maxx, maxy
109 c
110       character*1 saux01
111 c
112       integer nbmess
113       parameter ( nbmess = 20 )
114       character*80 texte(nblang,nbmess)
115 c
116 c 0.5. ==> initialisations
117 c ______________________________________________________________________
118 c
119 c====
120 c 1. messages
121 c====
122 c
123 #include "impr01.h"
124 c
125 #ifdef _DEBUG_HOMARD_
126       write (ulsort,texte(langue,1)) 'Entree', nompro
127       call dmflsh (iaux)
128 #endif
129 c
130       texte(1,4) =
131      > '(''Direction '',a1,'' : mini = '',g12.5,'' maxi = '',g12.5)'
132       texte(1,5) =
133      >'(''Nombre de noeuds attendus pour le maillage 2D :'',i10)'
134       texte(1,6) =
135      >'(''Nombre de noeuds trouves pour le maillage 2D  :'',i10)'
136       texte(1,7) = '(''Recherche du noeud de la famille '',i8)'
137       texte(1,8) = '(''Aucun noeud n''''est de la famille '',i8)'
138       texte(1,9) = '(''Impossible de retrouver cooinf et coosup.'')'
139       texte(1,10) =
140      > '(''Maille en '',a1,'' : mini = '',g12.5,'' maxi = '',g12.5)'
141       texte(1,11) = '(''Choix du calcul de delta '',a1,'' :'',i8)'
142       texte(1,12) = '(''Ce choix est inconnu.'')'
143       texte(1,13) = '(''Delta '',a1,'' initial.'')'
144       texte(1,14) = '(''Delta '',a1,'' impose.'')'
145       texte(1,15) =
146      >'(''D'',a1,'' = moyenne arithmetique des mini/maxi des mailles'')'
147       texte(1,16) =
148      > '(''D'',a1,'' = moyenne geometrique des mini/maxi des mailles'')'
149       texte(1,17) = '(''D'',a1,'' = D initial / 2**nivsup'')'
150 c
151       texte(2,4) =
152      > '(a1,''direction '','' : mini = '',g12.5,'' maxi = '',g12.5)'
153       texte(2,5) =
154      > '(''Expected number of nodes for the 2D mesh :'',i10)'
155       texte(2,6) =
156      > '(''Found number of nodes for the 2D mesh    :'',i10)'
157       texte(2,7) = '(''Searching for node with family # '',i8)'
158       texte(2,8) = '(''No node belongs to family # '',i8)'
159       texte(2,9) = '(''cooinf and coosup cannot be found.'')'
160       texte(2,10) =
161      > '(''Mesh along '',a1,'' : mini = '',g12.5,'' maxi = '',g12.5)'
162       texte(2,11) = '(''Choice for delta '',a1,'' calculation :'',i8)'
163       texte(2,12) = '(''This choice is unknown :'',i8)'
164       texte(2,13) = '(''Initial Delta '',a1,''.'')'
165       texte(2,14) = '(''Imposed Delta '',a1,''.'')'
166       texte(2,15) =
167      > '(''D'',a1,'' = arithmetic mean of mini/maxi of meshes'')'
168       texte(2,16) =
169      > '(''D'',a1,'' = geometric mean of mini/maxi of meshes'')'
170       texte(2,17) = '(''D'',a1,'' = initial / 2**nivsup'')'
171 c
172 #include "impr03.h"
173 c
174 #ifdef _DEBUG_HOMARD_
175       write (ulsort,90002) 'maextr', maextr
176 #endif
177       if ( maextr.eq.1 ) then
178         saux01 = 'X'
179         iaux1 = 2
180         iaux2 = 3
181       elseif ( maextr.eq.2 ) then
182         saux01 = 'Y'
183         iaux1 = 1
184         iaux2 = 3
185       elseif ( maextr.eq.3 ) then
186         saux01 = 'Z'
187         iaux1 = 1
188         iaux2 = 2
189       else
190         codret = 1
191       endif
192 c
193 c====
194 c 2. Quelle epaisseur ?
195 c====
196 c 2.1. ==> recuperation des cotes initiales des faces inferieures et
197 c          superieures, en examinant le noeud supplementaire :
198 c               ( x = cooinf , y = coosup )
199 c
200       if ( codret.eq.0 ) then
201 c
202 #ifdef _DEBUG_HOMARD_
203       write (ulsort,texte(langue,7)) famnzz
204 #endif
205 c
206       do 21 , iaux = 1 , nbnoto
207 c
208         if ( famnoe(iaux).eq.famnzz ) then
209 c
210           cooinf = coonoe(iaux,1)
211           coosup = coonoe(iaux,2)
212 c
213 #ifdef _DEBUG_HOMARD_
214           write (ulsort,texte(langue,4)) saux01, cooinf, coosup
215 #endif
216           goto 210
217 c
218         endif
219 c
220    21 continue
221 c
222       write (ulsort,texte(langue,8)) famnzz
223       write (ulsort,texte(langue,9))
224       codret = 12
225 c
226   210 continue
227 c
228       endif
229 c
230 c 2.2. ==> recherche eventuelle des tailles mini/maxi des mailles, selon
231 c          les axes perpendicalaires a la direction d'extrusion
232 #ifdef _DEBUG_HOMARD_
233       write (ulsort,90002) '2.2. recherche ; codret', codret
234 #endif
235 #ifdef _DEBUG_HOMARD_
236       write (ulsort,texte(langue,11)) saux01, choixd
237 #endif
238 c
239       if ( choixd.eq.3 .or. choixd.eq.4 ) then
240 c
241         if ( codret.eq.0 ) then
242 c
243 c 2.2.1. ==> initialisation des extrema pour la premiere maille trouvee
244 c
245         do 221 , iaux = 1 , nbelem
246 c
247           if ( typele(iaux).eq.edqua4 ) then
248 c
249             daux1 =
250      >    abs(coonoe(noeele(iaux,2),iaux1)-coonoe(noeele(iaux,1),iaux1))
251             daux2 =
252      >    abs(coonoe(noeele(iaux,3),iaux1)-coonoe(noeele(iaux,2),iaux1))
253           maxx = max ( daux1, daux2 )
254             minx = maxx
255 c
256             daux1 =
257      >    abs(coonoe(noeele(iaux,2),iaux2)-coonoe(noeele(iaux,1),iaux2))
258             daux2 =
259      >    abs(coonoe(noeele(iaux,3),iaux2)-coonoe(noeele(iaux,2),iaux2))
260             maxy = max ( daux1, daux2 )
261             miny = maxy
262             goto 222
263 c
264           endif
265 c
266   221   continue
267 c
268   222   continue
269 c
270 c 2.2.2. ==> parcours de toutes les mailles
271 c            on teste la non nullite au millionieme de l'ecart initial
272 c            entre le dessus et le dessous, divise par
273 c            10 puissance le niveau superieur atteint
274 c
275         daux = 1.d-6*(coosup-cooinf)/10.d0**nivsup
276 c
277         do 223 , iaux = 1 , nbelem
278 c
279           if ( typele(iaux).eq.edqua4 ) then
280 c
281             daux1 =
282      >    abs(coonoe(noeele(iaux,2),iaux1)-coonoe(noeele(iaux,1),iaux1))
283             daux2 =
284      >    abs(coonoe(noeele(iaux,3),iaux1)-coonoe(noeele(iaux,2),iaux1))
285 c
286             maxx = max ( maxx, daux1, daux2 )
287             if ( daux1.gt.daux ) then
288               minx = min ( minx, daux1 )
289             else
290               minx = min ( minx, daux2 )
291             endif
292 c
293             daux1 =
294      >    abs(coonoe(noeele(iaux,2),iaux2)-coonoe(noeele(iaux,1),iaux2))
295             daux2 =
296      >    abs(coonoe(noeele(iaux,3),iaux2)-coonoe(noeele(iaux,2),iaux2))
297             maxy = max ( maxy, daux1, daux2 )
298             if ( daux1.gt.daux ) then
299               miny = min ( miny, daux1 )
300             else
301               miny = min ( miny, daux2 )
302             endif
303 c
304           endif
305 c
306   223   continue
307 c
308 #ifdef _DEBUG_HOMARD_
309         write (ulsort,texte(langue,10)) '1', minx, maxx
310         write (ulsort,texte(langue,10)) '2', miny, maxy
311 #endif
312 c
313         endif
314 c
315       endif
316 c
317 c====
318 c 3. les noeuds de depart sont dans le plan cooinf
319 c    on cree leur correspondant dans le plan coosup
320 c====
321 #ifdef _DEBUG_HOMARD_
322       write (ulsort,90002) '3. Noeuds de depart ; codret', codret
323 #endif
324 c
325 c 3.0. ==> Calcul de l'ecart cooinf<-->coosup
326 c
327       if ( codret.eq.0 ) then
328 c
329       if ( choixd.ge.1 .and. choixd.le.5 ) then
330         write (ulsort,texte(langue,12+choixd)) saux01
331       else
332         write (ulsort,texte(langue,11)) saux01, choixd
333         write (ulsort,texte(langue,12))
334         codret = 1
335       endif
336 c
337       endif
338 c
339       if ( codret.eq.0 ) then
340 c
341       if ( choixd.eq.1 ) then
342         daux = coosup - cooinf
343       elseif ( choixd.eq.2 ) then
344         daux = deltac
345       elseif ( choixd.eq.3 ) then
346         daux = (minx+miny+maxx+maxy) * 0.25d0
347       elseif ( choixd.eq.4 ) then
348         daux = (minx*miny*maxx*maxy) ** 0.25d0
349       elseif ( choixd.eq.5 ) then
350         write (ulsort,90002) 'nivinf, nivsup',nivinf, nivsup
351         daux = ( coosup - cooinf ) / 2.d0**nivsup
352       endif
353 c
354       coosup = cooinf + daux
355       write (ulsort,texte(langue,4)) saux01, cooinf, coosup
356 c
357       jaux = 0
358 c
359       do 31 , iaux = 1 , nbnoto
360 c
361         if ( famnoe(iaux).eq.famnzz ) then
362 c
363           nu3dno(iaux) = 0
364 c
365         else
366 c
367           jaux = jaux + 1
368           coon3d(jaux,iaux1) = coonoe(iaux,1)
369           coon3d(jaux,iaux2) = coonoe(iaux,2)
370           coon3d(jaux,maextr) = cooinf
371           famn3d(jaux) = famnoe(iaux)
372           nu3dno(iaux) = jaux
373 c
374           kaux = jaux + nbnoto - 1
375           coon3d(kaux,iaux1) = coonoe(iaux,1)
376           coon3d(kaux,iaux2) = coonoe(iaux,2)
377           coon3d(kaux,maextr) = coosup
378           famn3d(kaux) = famnoe(iaux)
379 c
380         endif
381 c
382    31 continue
383 c
384       if ( kaux.ne.nbno3d ) then
385         write (ulsort,texte(langue,5)) nbno3d
386         write (ulsort,texte(langue,6)) jaux
387         codret = 3
388       endif
389 c
390       endif
391 c
392 c====
393 c 4. la fin
394 c====
395 c
396 #ifdef _DEBUG_HOMARD_
397       write (ulsort,*) '3. fin ; codret = ', codret
398 #endif
399 c
400       if ( codret.ne.0 ) then
401 c
402 #include "envex2.h"
403       write (ulsort,texte(langue,1)) 'Sortie', nompro
404       write (ulsort,texte(langue,2)) codret
405       endif
406 c
407 #ifdef _DEBUG_HOMARD_
408       write (ulsort,texte(langue,1)) 'Sortie', nompro
409       call dmflsh (iaux)
410 #endif
411 c
412       end