Salome HOME
Homard executable
[modules/homard.git] / src / tool / AP_Conversion / pcsono.F
1       subroutine pcsono ( numnp1, numnp2, typint, deraff, option,
2      >                    nbpara, carenf, carchf, nrfonc,
3      >                    hetnoe, ancnoe,
4      >                    nnoeho, nnoeca, nnosho,
5      >                    hetare, somare, filare,
6      >                    np2are,
7      >                    hettri, aretri, filtri,
8      >                    hetqua, arequa, filqua,
9      >                    tritet, cotrte, aretet,
10      >                    filtet, hettet,
11      >                    quahex, coquhe, arehex,
12      >                    filhex, hethex, fhpyte,
13      >                    facpen, cofape, arepen,
14      >                    filpen, hetpen, fppyte,
15      >                    facpyr, cofapy, arepyr,
16      >                    ulsort, langue, codret )
17 c ______________________________________________________________________
18 c
19 c                             H O M A R D
20 c
21 c Outil de Maillage Adaptatif par Raffinement et Deraffinement d'EDF R&D
22 c
23 c Version originale enregistree le 18 juin 1996 sous le numero 96036
24 c aupres des huissiers de justice Simart et Lavoir a Clamart
25 c Version 11.2 enregistree le 13 fevrier 2015 sous le numero 2015/014
26 c aupres des huissiers de justice
27 c Lavoir, Silinski & Cherqui-Abrahmi a Clamart
28 c
29 c    HOMARD est une marque deposee d'Electricite de France
30 c
31 c Copyright EDF 1996
32 c Copyright EDF 1998
33 c Copyright EDF 2002
34 c Copyright EDF 2020
35 c ______________________________________________________________________
36 c
37 c    aPres adaptation - Conversion de SOlution - NOeud
38 c     -                 -             --         --
39 c ______________________________________________________________________
40 c .        .     .        .                                            .
41 c .  nom   . e/s . taille .           description                      .
42 c .____________________________________________________________________.
43 c . numnp1 . e   .    1   . nombre de noeuds de la fonction si P1      .
44 c . numnp2 . e   .    1   . nombre de noeuds de la fonction si P2      .
45 c . typint . e   .   1    . type d'interpolation                       .
46 c .        .     .        .  0, si automatique                         .
47 c .        .     .        .  elements : 0 si intensif, sans orientation.
48 c .        .     .        .             1 si extensif, sans orientation.
49 c .        .     .        .             2 si intensif, avec orientation.
50 c .        .     .        .             3 si extensif, avec orientation.
51 c .        .     .        .  noeuds : 1 si degre 1                     .
52 c .        .     .        .           2 si degre 2                     .
53 c .        .     .        .           3 si iso-P2                      .
54 c . deraff . e   .    1   . vrai, s'il y a eu du deraffinement en      .
55 c .        .     .        . passant de l'iteration n a n+1 ; faux sinon.
56 c . option . e   .    1   . option du traitement                       .
57 c .        .     .        . -1 : Pas de changement dans le maillage    .
58 c .        .     .        .  0 : Adaptation complete                   .
59 c .        .     .        .  1 : Modification de degre                 .
60 c . nbpara . e   .   1    . nombre de parametres a enregistrer par fonc.
61 c . carenf .   s .nbpara* . caracteristiques entieres des fonctions :  .
62 c .        .     .  nnfopa.  1 : 0, pour une fonction ancienne isolee  .
63 c .        .     .        .      1, pour une ancienne associee a une   .
64 c .        .     .        .         autre fonction                     .
65 c .        .     .        .      -1, pour une nouvelle fonction        .
66 c .        .     .        .  2 : typcha                                .
67 c .        .     .        .  3 : typgeo                                .
68 c .        .     .        .  4 : typass                                .
69 c .        .     .        .  5 : ngauss                                .
70 c .        .     .        .  6 : nnenmx                                .
71 c .        .     .        .  7 : nnvapr                                .
72 c .        .     .        .  8 : carsup                                .
73 c .        .     .        .  9 : nbtafo                                .
74 c .        .     .        . 10 : anvale                                .
75 c .        .     .        . 11 : anvalr                                .
76 c .        .     .        . 12 : anobch                                .
77 c .        .     .        . 13 : anprpg                                .
78 c .        .     .        . 14 : anlipr                                .
79 c .        .     .        . 15 : npenmx                                .
80 c .        .     .        . 16 : npvapr                                .
81 c .        .     .        . 17 : apvale                                .
82 c .        .     .        . 18 : apvalr                                .
83 c .        .     .        . 19 : apobch                                .
84 c .        .     .        . 20 : apprpg                                .
85 c .        .     .        . 21 : apvatt                                .
86 c .        .     .        . 22 : apvane                                .
87 c .        .     .        . 23 : antyas                                .
88 c .        .     .        . 24 : aptyas                                .
89 c .        .     .        . 25 : numero de la 1ere fonction associee   .
90 c .        .     .        . 26 : numero de la 2nde fonction associee   .
91 c . carchf . es  .nbpara* . caracteristiques caracteres des fonctions :.
92 c .        .     .  nnfopa.  1 : nom de la fonction                    .
93 c .        .     .        .  2 : nom de la fonction n associee         .
94 c .        .     .        .  3 : nom de la fonction p associee         .
95 c .        .     .        .  4 : obpcan                                .
96 c .        .     .        .  5 : obpcap                                .
97 c .        .     .        .  6 : obprof                                .
98 c .        .     .        .  7 : oblopg                                .
99 c .        .     .        .  8 : si aux points de Gauss, nom de la     .
100 c .        .     .        .      fonction n ELNO correspondante        .
101 c .        .     .        .  9 : si aux points de Gauss, nom de la     .
102 c .        .     .        .      fonction p ELNO correspondante        .
103 c . nrfonc . e   .   1    . numero de la fonction principale           .
104 c . hetnoe . e   . nbnoto . historique de l'etat des noeuds            .
105 c . ancnoe . e   . nbnoto . ancien numero de noeud si deraffinement    .
106 c . nnoeho . e   . renoto . numero des noeuds en entree pour homard    .
107 c . nnoeca . e   . renoto . numero des noeuds en entree dans le calcul .
108 c . nnosho . e   . rsnoto . numero des noeuds en sortie pour homard    .
109 c . hetare . e   . nbarto . historique de l'etat des aretes            .
110 c . somare . e   .2*nbarto. numeros des extremites d'arete             .
111 c . filare . e   . nbarto . premiere fille des aretes                  .
112 c . np2are . e   . nbarto . numero des noeuds p2 milieux d'aretes      .
113 c . hettri . e   . nbtrto . historique de l'etat des triangles         .
114 c . aretri . e   .nbtrto*3. numeros des 3 aretes des triangles         .
115 c . filtri . e   . nbtrto . premier fils des triangles                 .
116 c . hetqua . e   . nbquto . historique de l'etat des quadrangles       .
117 c . arequa . e   .nbquto*4. numeros des 4 aretes des quadrangles       .
118 c . filqua . e   . nbquto . premier fils des quadrangles               .
119 c . tritet . e   .nbtecf*4. numeros des 4 triangles des tetraedres     .
120 c . cotrte . e   .nbtecf*4. code des 4 triangles des tetraedres        .
121 c . aretet . e   .nbteca*6. numeros des 6 aretes des tetraedres        .
122 c . filtet . e   . nbteto . premier fils des tetraedres                .
123 c . hettet . e   . nbteto . historique de l'etat des tetraedres        .
124 c . quahex . e   .nbhecf*6. numeros des 6 quadrangles des hexaedres    .
125 c . coquhe . e   .nbhecf*6. code des 6 quadrangles des hexaedres       .
126 c . arehex . e   .nbheca12. numeros des 12 aretes des hexaedres        .
127 c . filhex . e   . nbheto . premier fils des hexaedres                 .
128 c . hethex . e   . nbheto . historique de l'etat des hexaedres        .
129 c . fhpyte . e   .2*nbheco. fhpyte(1,j) = numero de la 1ere pyramide   .
130 c .        .     .        . fille de l'hexaedre k tel que filhex(k) =-j.
131 c .        .     .        . fhpyte(2,j) = numero du 1er tetraedre      .
132 c .        .     .        . fils de l'hexaedre k tel que filhex(k) = -j.
133 c . facpen . e   .nbpecf*5. numeros des faces des pentaedres           .
134 c . cofape . e   .nbpecf*5. codes des faces des pentaedres             .
135 c . filpen . e   . nbpeto . premier fils des pentaedres                .
136 c . hetpen . e   . nbpeto . historique de l'etat des pentaedres        .
137 c . arepen . e   .nbpeca*9. numeros des 9 aretes des pentaedres        .
138 c . fppyte . e   .2*nbpeco. fppyte(1,j) = numero de la 1ere pyramide   .
139 c .        .     .        . fille du pentaedre k tel que filpen(k) =-j .
140 c .        .     .        . fppyte(2,j) = numero du 1er tetraedre      .
141 c .        .     .        . fils du pentaedre k tel que filpen(k) = -j .
142 c . facpyr . e   .nbpycf*5. numeros des 5 faces des pyramides          .
143 c . cofapy . e   .nbpycf*5. codes des 5 faces des pyramides            .
144 c . arepyr . e   .nbpyca*8. numeros des 8 aretes des pyramides         .
145 c . ulsort . e   .   1    . numero d'unite logique de la liste standard.
146 c . langue . e   .    1   . langue des messages                        .
147 c .        .     .        . 1 : francais, 2 : anglais                  .
148 c . codret . es  .    1   . code de retour des modules                 .
149 c .        .     .        . 0 : pas de probleme                        .
150 c .        .     .        . 1 : probleme                               .
151 c ______________________________________________________________________
152 c
153 c====
154 c 0. declarations et dimensionnement
155 c====
156 c
157 c 0.1. ==> generalites
158 c
159       implicit none
160       save
161 c
162       character*6 nompro
163       parameter ( nompro = 'PCSONO' )
164 c
165 #include "nblang.h"
166 c
167 c 0.2. ==> communs
168 c
169 #include "envex1.h"
170 c
171 #include "gmreel.h"
172 #include "gmenti.h"
173 c
174 #include "envca1.h"
175 #include "nomber.h"
176 #include "nombno.h"
177 #include "nombar.h"
178 #include "nombtr.h"
179 #include "nombqu.h"
180 #include "nombte.h"
181 #include "nombhe.h"
182 #include "nombpe.h"
183 #include "nombpy.h"
184 #include "nombsr.h"
185 c
186 c 0.3. ==> arguments
187 c
188       integer numnp1, numnp2
189       integer typint
190       integer option
191 c
192       integer nbpara
193       integer carenf(nbpara,*)
194       integer nrfonc
195 c
196       integer hetnoe(nbnoto), ancnoe(nbnoto)
197       integer nnoeho(renoto), nnoeca(renoto)
198       integer nnosho(rsnoto)
199       integer hetare(nbarto), somare(2,nbarto), filare(nbarto)
200       integer np2are(nbarto)
201       integer hettri(nbtrto), aretri(nbtrto,3), filtri(nbtrto)
202       integer hetqua(nbquto), arequa(nbquto,4), filqua(nbquto)
203       integer tritet(nbtecf,4),  cotrte(nbtecf,4), aretet(nbteca,6)
204       integer filtet(nbteto), hettet(nbteto)
205       integer quahex(nbhecf,6), coquhe(nbhecf,6), arehex(nbheca,12)
206       integer hethex(nbheto), filhex(nbheto), fhpyte(2,nbheco)
207       integer facpen(nbpecf,5), cofape(nbpecf,5), arepen(nbpeca,9)
208       integer filpen(nbpeto), hetpen(nbpeto), fppyte(2,nbpeco)
209       integer facpyr(nbpycf,5), cofapy(nbpycf,5), arepyr(nbpyca,8)
210 c
211       character*8 carchf(nbpara,*)
212 c
213       logical deraff
214 c
215       integer ulsort, langue, codret
216 c
217 c 0.4. ==> variables locales
218 c
219       integer iaux, jaux
220       integer codre1, codre2
221       integer codre0
222 c
223       integer typfon, typcha, typgeo, typass
224       integer ngauss, nnenmx, nnvapr, carsup, nbtafo
225       integer n1vale, n1valr, n1obpr, n1obch, n1lipr
226       integer npenmx, npvapr
227       integer p1vale, p1valr, p1obpr, p1obch, p1vatt
228       integer p1vane, p1tyas
229       integer adpcan, adpcap
230       integer nrfon2, nrfon3
231       integer typprf, typin0
232 c
233       integer nbfop1, nbfop2
234       integer pvap1h, pvap2h
235 c
236       character*8 nofonc, obpcan, obpcap, obprof
237       character*8 oblopg
238       character*8 nvap1h, nvap2h
239 c
240       integer nbmess
241       parameter ( nbmess = 120 )
242       character*80 texte(nblang,nbmess)
243 c
244 c 0.5. ==> initialisations
245 c ______________________________________________________________________
246 c
247 c====
248 c 1. initialisations
249 c====
250 c
251 #include "impr01.h"
252 #include "impr03.h"
253 c
254 #ifdef _DEBUG_HOMARD_
255       write (ulsort,texte(langue,1)) 'Entree', nompro
256       call dmflsh (iaux)
257 #endif
258 c
259       codret = 0
260 c
261 #include "esimpr.h"
262 c
263       texte(1,4) = '(5x,''Fonctions '',a2,'' :'')'
264       texte(1,5) =
265      > '(5x,''. nombre de noeuds dans leur definition     : '',i10)'
266       texte(1,6) =
267      > '(5x,''. nombre de noeuds dans le maillage initial : '',i10)'
268       texte(1,7) =
269      > '(5x,''. nombre de valeurs du profil : '',i10)'
270       texte(1,8) = '(''... Premiere(s) valeur(s) : '',5i10)'
271       texte(1,9) = '(''... Dernieres valeurs     : '',5i10)'
272       texte(1,10) = '(''. Interpolation '',a)'
273 c
274       texte(2,4) = '(5x,''Fonctions '',a2,'' :'')'
275       texte(2,5) =
276      > '(5x,''. number of nodes in their definition : '',i10)'
277       texte(2,6) =
278      > '(5x,''. number of nodes in the initial mesh : '',i10)'
279       texte(2,7) =
280      > '(5x,''. length of profile : '',i10)'
281       texte(2,8) = '(''... First value(s) : '',5i10)'
282       texte(2,9) = '(''... Last value(s)  : '',5i10)'
283       texte(2,10) = '(''. '',a,'' interpolation '')'
284 c
285 #ifdef _DEBUG_HOMARD_
286       write (ulsort,90002) 'option', option
287 #endif
288 c
289 c====
290 c 2. grandeurs utiles
291 c====
292 c 2.1. ==> recuperation
293 c
294       if ( codret.eq.0 ) then
295 c
296       iaux = nrfonc
297 #ifdef _DEBUG_HOMARD_
298       write (ulsort,texte(langue,3)) 'PCFOR2', nompro
299 #endif
300       call pcfor2 ( nbpara, carenf, carchf,
301      >              iaux,
302      >              typfon, typcha, typgeo, typass,
303      >              ngauss, nnenmx, nnvapr, carsup, nbtafo,
304      >              n1vale, n1valr, n1obpr, n1obch, n1lipr,
305      >              npenmx, npvapr,
306      >              p1vale, p1valr, p1obpr, p1obch, p1vatt,
307      >              p1vane, p1tyas,
308      >              nrfon2, nrfon3,
309      >              nofonc,
310      >              obpcan, obpcap, obprof, adpcan, adpcap,
311      >              oblopg,
312      >              ulsort, langue, codret )
313 c
314       endif
315 cgn      write(ulsort,*) 'apres pcfor2'
316 cgn      write(ulsort,90002) 'carsup', carsup
317 cgn      write(ulsort,90002) 'nnvapr', nnvapr
318 cgn      write(ulsort,90002) 'nbtafo', nbtafo
319 c
320 c 2.2. ==> type de profil
321 #ifdef _DEBUG_HOMARD_
322       write (ulsort,90002) '2.2. type de profil ;  codret', codret
323 #endif
324 c
325       if ( codret.eq.0 ) then
326 c
327       if ( degre.eq.2 .and. nnvapr.gt.0 ) then
328 #ifdef _DEBUG_HOMARD_
329       write (ulsort,texte(langue,3)) 'PCSPRN', nompro
330 #endif
331         call pcsprn ( typprf, numnp1,
332      >                hetnoe, nnoeho,
333      >                nnvapr, imem(n1lipr) )
334 c
335       else
336 c
337         typprf = 0
338 c
339       endif
340 c
341       endif
342 c
343 c 2.3. ==> grandeurs deduites
344 c          . Si on n'a rien de special sur le profil, on est fidele
345 c            au degre
346 c          . Sinon, c'est une fonction de degre 1
347 c
348       if ( codret.eq.0 ) then
349 c
350       if ( typprf.eq.0 ) then
351 c
352         if ( degre.eq.1 ) then
353           nbfop1 = nbtafo
354           nbfop2 = 0
355         else
356           nbfop1 = 0
357           nbfop2 = nbtafo
358         endif
359 c
360       else
361 c
362         nbfop1 = nbtafo
363         nbfop2 = 0
364 c
365       endif
366 c
367 #ifdef _DEBUG_HOMARD_
368       write (ulsort,90002) 'nbfop1, nbfop2, typint',
369      >                      nbfop1, nbfop2, typint
370       write (ulsort,texte(langue,7)) nnvapr
371       if ( nnvapr.gt.0 ) then
372         write (ulsort,texte(langue,8))
373      >         (imem(iaux),iaux=n1lipr,n1lipr+min(4,nnvapr-1))
374         if ( nnvapr.gt.5 ) then
375           write (ulsort,texte(langue,9))
376      >           (imem(iaux),iaux=n1lipr+nnvapr-5,n1lipr+nnvapr-1)
377         endif
378       endif
379 #endif
380 c
381       endif
382 c
383 c 2.4. ==> verification des coherences de taille
384 #ifdef _DEBUG_HOMARD_
385       write (ulsort,90002) '2.4. verification ;  codret', codret
386 #endif
387 c
388       if ( codret.eq.0 ) then
389 c
390       if ( nnvapr.le.0 ) then
391 c
392       if ( nbfop1.ne.0 ) then
393         if ( numnp1.ne.reno1i ) then
394           write (ulsort,texte(langue,4)) 'P1'
395           write (ulsort,texte(langue,5)) numnp1
396           write (ulsort,texte(langue,6)) reno1i
397           write (ulsort,texte(langue,7)) nnvapr
398           codret = 4
399         endif
400       endif
401 c
402       if ( nbfop2.ne.0 ) then
403         if ( numnp2.ne.renoto ) then
404           write (ulsort,texte(langue,4)) 'P2'
405           write (ulsort,texte(langue,5)) numnp2
406           write (ulsort,texte(langue,6)) renoto
407           write (ulsort,texte(langue,7)) nnvapr
408           codret = 4
409         endif
410       endif
411 c
412       endif
413 c
414       endif
415 c
416 c 2.5. ==> type d'interpolation
417 #ifdef _DEBUG_HOMARD_
418       write (ulsort,90002) '2.5. type interpolation ;  codret', codret
419 #endif
420 c
421       if ( codret.eq.0 ) then
422 c
423       typin0 = -1
424 c
425       if ( option.eq.1 ) then
426         if ( nbfop1.ne.0 ) then
427           typin0 = 5
428         else
429           typin0 = 4
430         endif
431       elseif ( typint.eq.0 ) then
432         if ( nbfop1.ne.0 ) then
433           typin0 = 1
434         else
435           typin0 = 2
436         endif
437       elseif ( typint.eq.1 ) then
438         typin0 = 1
439         if ( nbfop1.eq.0 ) then
440           write (ulsort,texte(langue,100+typin0))
441           write (ulsort,texte(langue,117)) 1, nbfop1
442           codret = 251
443         endif
444       elseif ( typint.eq.2 .or. typint.eq.3 ) then
445         typin0 = typint
446         if ( nbfop2.eq.0 ) then
447           write (ulsort,texte(langue,100+typin0))
448           write (ulsort,texte(langue,117)) 2, nbfop2
449           codret = 252
450         endif
451       endif
452 c
453       endif
454 c
455 #ifdef _DEBUG_HOMARD_
456       write (ulsort,90002) 'typin0', typin0
457       if ( typin0.ge.0 .and. typin0.le.5 ) then
458         write (ulsort,texte(langue,100+typin0))
459         write (ulsort,texte(langue,117)) 1, nbfop1
460         write (ulsort,texte(langue,117)) 2, nbfop2
461       endif
462 #endif
463 c
464 c====
465 c 3. interpolation des variables aux noeuds
466 c    remarque : si les fonctions sont inexistantes dans l'une des
467 c               categories, on alloue quand meme les tableaux. les
468 c               longueurs sont nulles donc on ne perd pas de place.
469 c               la lisibilite du programme compense le peu de temps cpu
470 c               necessaire a cela.
471 c====
472 #ifdef _DEBUG_HOMARD_
473       write (ulsort,90002) '3. redistribution ;  codret', codret
474 #endif
475 c
476 c 3.1. ==> allocations des tableaux intermediaires pour les fonctions
477 c          aux noeuds
478 c
479       if ( codret.eq.0 ) then
480 c
481       if ( nbfop1.ne.0 ) then
482         iaux = nbfop1 * max(renoto,nbnoto)
483       else
484         iaux = 0
485       endif
486       call gmalot ( nvap1h, 'reel    ', iaux, pvap1h, codre1 )
487       iaux = nbfop2 * max(renoto,nbnoto)
488       call gmalot ( nvap2h, 'reel    ', iaux, pvap2h, codre2 )
489 c
490       codre0 = min ( codre1, codre2 )
491       codret = max ( abs(codre0), codret,
492      >               codre1, codre2 )
493 c
494       endif
495 c
496 c 3.2. ==> redistribution des valeurs dans la numerotation homard
497 c
498 cgn       write (*,*) 'fonctions P1'
499 cgn       write (*,92010)
500 cgn     >(rmem(iaux),iaux=n1valr,n1valr-1+nbfop1*max(nnvapr,renoto))
501 cgn       write (*,*) 'fonctions P2'
502 cgn      write (*,92010)
503 cgn     >(rmem(iaux),iaux=n1valr,n1valr-1+nbfop2*max(nnvapr,renoto))
504 c
505       if ( codret.eq.0 ) then
506 c
507 #ifdef _DEBUG_HOMARD_
508       write (ulsort,texte(langue,3)) 'PCSRHO', nompro
509 #endif
510       call pcsrho ( nbfop1, nbfop2, numnp1, numnp2,
511      >              deraff, option,
512      >              hetnoe, ancnoe,
513      >              nnoeho, nnoeca,
514      >              nnvapr, imem(n1lipr), imem(adpcan), imem(adpcap),
515      >              rmem(n1valr), rmem(n1valr),
516      >              rmem(pvap1h), rmem(pvap2h),
517      >              ulsort, langue, codret )
518 c
519       endif
520 cgn      write (*,*)'apres pcsrho'
521 cgn      call gmprsx (nompro, nvap1h )
522 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=2073,2073)
523 cgn      write(*,92010) (rmem(pvap1h-1+iaux),iaux=1,nbnoto)
524 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=1,nbnoto)
525 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
526 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,nbfop2*nbnoto)
527 c
528 c====
529 c 4. Interpolation p1 des variables aux noeuds
530 c====
531 c
532       if ( codret.eq.0 ) then
533 c
534       if ( typin0.eq.1 ) then
535 c
536 #ifdef _DEBUG_HOMARD_
537       write (ulsort,90002) '4. ==> Interpolation p1 ; codret', codret
538 #endif
539 c
540         write (ulsort,texte(langue,10)) 'P1'
541 c
542 c 4.1. ==> interpolation p1 pour les aretes decoupees
543 c
544 #ifdef _DEBUG_HOMARD_
545       write (ulsort,texte(langue,3)) 'PCS1AR', nompro
546 #endif
547         call pcs1ar ( nbfop1, imem(adpcap),
548      >                hetare, somare, filare,
549      >                rmem(pvap1h) )
550 cgn      write (*,*)'apres pcs1ar'
551 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=2073,2073)
552 cgn      write(*,92010) (rmem(pvap1h-1+iaux),iaux=2073,2073)
553 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
554 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
555 c
556 c 4.2. ==> interpolation p1 pour les quadrangles decoupes
557 c
558         if ( nbquto.ne.0 ) then
559 c
560 #ifdef _DEBUG_HOMARD_
561       write (ulsort,texte(langue,3)) 'PCS1QU', nompro
562 #endif
563           call pcs1qu ( nbfop1, imem(adpcap),
564      >                  somare,
565      >                  hetqua, arequa, filqua,
566      >                  rmem(pvap1h) )
567 cgn      write (*,*)'apres pcs1qu'
568 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
569 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
570 c
571         endif
572 c
573 c 4.3. ==> interpolation p1 pour les hexaedres decoupes
574 c
575         if ( nbheto.ne.0 ) then
576 c
577 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
578 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
579 #ifdef _DEBUG_HOMARD_
580       write (ulsort,texte(langue,3)) 'PCS1HE', nompro
581 #endif
582           call pcs1he ( nbfop1, imem(adpcap),
583      >                  somare,
584      >                  aretri, arequa,
585      >                  tritet, cotrte, aretet,
586      >                  quahex, coquhe, arehex,
587      >                  filhex, hethex, fhpyte,
588      >                  facpyr, cofapy, arepyr,
589      >                  rmem(pvap1h) )
590 c
591 cgn      write (*,*)'apres pcs1he'
592 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
593 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
594         endif
595 c
596 c 4.3. ==> interpolation p1 pour les pentaedres decoupes
597 c
598         if ( nbpeto.ne.0 ) then
599 c
600 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
601 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
602 #ifdef _DEBUG_HOMARD_
603       write (ulsort,texte(langue,3)) 'PCS1PE', nompro
604 #endif
605           call pcs1pe ( nbfop1, imem(adpcap),
606      >                  somare,
607      >                  aretri, arequa,
608      >                  tritet, cotrte,
609      >                  facpen, cofape, arepen,
610      >                  filpen, hetpen, fppyte,
611      >                  facpyr, cofapy, arepyr,
612      >                  rmem(pvap1h) )
613 c
614 cgn      write (*,*)'apres pcs1pe'
615 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
616 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
617         endif
618 c
619 c====
620 c 5. Interpolation p2 des variables aux noeuds
621 c    Attention a respecter l'ordre dans l'enchainement des appels
622 c====
623 c
624       elseif ( typin0.eq.2 ) then
625 c
626 #ifdef _DEBUG_HOMARD_
627       write (ulsort,90002) '5. Interpolation p2 ; codret', codret
628 #endif
629 c
630         write (ulsort,texte(langue,10)) 'P2'
631 c
632 c 5.1. ==> pour les aretes decoupees
633 c
634 #ifdef _DEBUG_HOMARD_
635       write (ulsort,texte(langue,3)) 'PCS2AR', nompro
636 #endif
637         call pcs2ar ( nbfop2, imem(adpcap), rmem(pvap2h),
638      >                hetare, somare, np2are, filare )
639 cgn       write (*,*)'apres pcs2ar'
640 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=1,nbnoto)
641 cgn       write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,nbfop2*nbnoto-1)
642 cgn       write (*,92010)
643 cgn     > (rmem(pvap2h-1+iaux),iaux=nbfop2*nbnoto-1,nbfop2*nbnoto-1)
644 c
645 c 5.2. ==> pour les tetraedres decoupes
646 c
647         if ( nbtema.ne.0 ) then
648 c
649 #ifdef _DEBUG_HOMARD_
650       write (ulsort,texte(langue,3)) 'PCS2TE', nompro
651 #endif
652           call pcs2te ( nbfop2, imem(adpcap), rmem(pvap2h),
653      >                  tritet, cotrte, aretet,
654      >                  hettet, filtet,
655      >                  somare, np2are,
656      >                  aretri )
657 c
658         endif
659 c
660 c 5.3. ==> quadrangles et hexaedres decoupes
661 c          Remarque : avec les hexaedres, il faut faire deux passages
662 c                     pour gerer les raffinements sur deux niveaux
663 c                     Tant pis si des calculs sont faits deux fois.
664 c
665         jaux = 1
666         if ( nbheto.ne.0 ) then
667           jaux = 2
668         endif
669 c
670         do 53 , iaux = 1 , jaux
671 c
672 c 5.3.1. ==> pour les quadrangles jaux
673 c
674           if ( nbquto.ne.0 ) then
675 c
676 #ifdef _DEBUG_HOMARD_
677       write (ulsort,texte(langue,3)) 'PCS2QU', nompro
678 #endif
679             call pcs2qu ( nbfop2, imem(adpcap), rmem(pvap2h),
680      >                    hetqua, arequa, filqua,
681      >                    somare, np2are,
682      >                    aretri )
683 cgn       write(ulsort,*)'apres pcs2qu'
684 cgn       write(ulsort,91011) (imem(adpcap-1+iaux),iaux=1,nbnoto)
685 cgn       write(ulsort,92010)(rmem(pvap2h-1+iaux),iaux=1,nbnoto)
686 c
687           endif
688 c
689 c 5.3.2. ==> pour les hexaedres decoupes
690 c
691           if ( nbheto.ne.0 ) then
692 c
693 cgn      write (*,*)'avant pcs2he'
694 cgn      call gmprsx ( 'avant pcs2he',nvap2h )
695 cgn      write(*,91011)(imem(adpcap-1+iaux),iaux=1,nbnoto)
696 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
697 #ifdef _DEBUG_HOMARD_
698       write (ulsort,texte(langue,3)) 'PCS2HE', nompro
699 #endif
700             call pcs2he ( nbfop2, imem(adpcap), rmem(pvap2h),
701      >                    hethex, quahex, coquhe, arehex,
702      >                    filhex, fhpyte,
703      >                    somare, np2are,
704      >                    aretri,
705      >                    hetqua, arequa, filqua,
706      >                    tritet, cotrte, aretet,
707      >                    facpyr, cofapy, arepyr,
708      >                    ulsort, langue, codret )
709 c
710 cgn      write (*,*)'apres pcs2he'
711 cgn      write(*,91011)(imem(adpcap-1+iaux),iaux=1,nbnoto)
712 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,nbfop2*nbnoto)
713 cgn      call gmprsx ( 'apres pcs2he',nvap2h )
714 cgn      call gmprsx ( 'apres pcs2he',obpcap)
715           endif
716 c
717    53   continue
718 c
719 c 5.4. ==> pour les pentaedres decoupes
720 c
721         if ( nbpeto.ne.0 ) then
722 c
723 cgn      call gmprsx ( 'avant pcs2pe',nvap2h )
724 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
725 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
726 #ifdef _DEBUG_HOMARD_
727       write (ulsort,texte(langue,3)) 'PCS2PE', nompro
728 #endif
729           call pcs2pe ( nbfop2, imem(adpcap), rmem(pvap2h),
730      >                  hetpen, facpen, cofape, filpen, fppyte,
731      >                  somare, np2are,
732      >                  aretri, arequa,
733      >                  tritet, cotrte,
734      >                  facpyr, cofapy,
735      >                  ulsort, langue, codret )
736 c
737 cgn      write (*,*)'apres pcs2pe'
738 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
739 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,nbfop2*nbnoto)
740 cgn      call gmprsx ( 'apres pcs2he',nvap2h )
741 cgn      call gmprsx ( 'apres pcs2he',obpcap)
742         endif
743 cgn      call gmprsx ( 'au final',nvap2h )
744 cgn      call gmprsx ( 'au final',obpcap)
745 c
746 c 5.5. ==> pour les triangles decoupes
747 c
748         if ( nbtrma.ne.0 ) then
749 c
750 #ifdef _DEBUG_HOMARD_
751       write (ulsort,texte(langue,3)) 'PCS2TR', nompro
752 #endif
753           call pcs2tr ( nbfop2, imem(adpcap), rmem(pvap2h),
754      >                  hettri, aretri, filtri,
755      >                  somare, np2are )
756 cgn       write (*,*)'apres pcs2tr'
757 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=1,nbnoto)
758 cgn       write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,nbnoto)
759 cgn       write (*,92010)
760 cgn     > (rmem(pvap2h-1+iaux),iaux=nbfop2*nbnoto,nbfop2*nbnoto)
761 c
762         endif
763 c
764 c====
765 c 6. Interpolation iso-p2 des variables aux noeuds
766 c    Attention a respecter l'ordre dans l'enchainement des appels
767 c====
768 c
769       elseif ( typin0.eq.3 ) then
770 c
771 #ifdef _DEBUG_HOMARD_
772       write (ulsort,90002) '6. Interpolation iso-p2 ; codret', codret
773 #endif
774 c
775         write (ulsort,texte(langue,10)) 'iso-P2'
776 cgn      write (*,*)'au debut'
777 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
778 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,2*nbnoto)
779 c
780 c 6.1. ==> pour les aretes decoupees
781 c
782 #ifdef _DEBUG_HOMARD_
783       write (ulsort,texte(langue,3)) 'PCSIAR', nompro
784 #endif
785         call pcsiar ( nbfop2, imem(adpcap), rmem(pvap2h),
786      >                hetare, somare, np2are, filare )
787 cgn      write (*,*)'apres pcsiar'
788 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
789 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,2*nbnoto)
790 c
791 c 6.2. ==> pour les tetraedres decoupes
792 c
793         if ( nbtema.ne.0 ) then
794 c
795 #ifdef _DEBUG_HOMARD_
796       write (ulsort,texte(langue,3)) 'PCSITE', nompro
797 #endif
798           call pcsite ( nbfop2, imem(adpcap), rmem(pvap2h),
799      >                  tritet, cotrte, aretet,
800      >                  hettet, filtet,
801      >                  somare, np2are,
802      >                  aretri )
803 c
804         endif
805 c
806 c 6.3. ==> quadrangles et hexaedres decoupes
807 c          Remarque : avec les hexaedres, il faut faire deux passages
808 c                     pour gerer les raffinements sur deux niveaux
809 c                     Tant pis si des calculs sont faits deux fois.
810 c
811         jaux = 1
812         if ( nbheto.ne.0 ) then
813           jaux = 2
814         endif
815 c
816         do 63 , iaux = 1 , jaux
817 c
818 c 6.3.1. ==> pour les quadrangles decoupes
819 c
820         if ( nbquto.ne.0 ) then
821 c
822 #ifdef _DEBUG_HOMARD_
823       write (ulsort,texte(langue,3)) 'PCSIQU', nompro
824 #endif
825           call pcsiqu ( nbfop2, imem(adpcap), rmem(pvap2h),
826      >                  hetqua, arequa, filqua,
827      >                  somare, np2are,
828      >                  aretri )
829 cgn      write (*,*)'apres pcsiqu'
830 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
831 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,2*nbnoto)
832 c
833         endif
834 c
835 c 6.3.2. ==> pour les hexaedres decoupes
836 c
837         if ( nbheto.ne.0 ) then
838 c
839 #ifdef _DEBUG_HOMARD_
840       write (ulsort,texte(langue,3)) 'PCSIHE', nompro
841 #endif
842           call pcsihe ( nbfop2, imem(adpcap), rmem(pvap2h),
843      >                  hethex, quahex, coquhe, arehex,
844      >                  filhex, fhpyte,
845      >                  somare, np2are,
846      >                  aretri,
847      >                  arequa,
848      >                  tritet, cotrte, aretet,
849      >                  facpyr, cofapy, arepyr )
850 c
851 cgn      write (*,*)'apres pcs2he'
852 cgn      write(*,91011)(imem(adpcap-1+iaux),iaux=1,nbnoto)
853 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,nbfop2*nbnoto)
854 cgn      call gmprsx ( 'apres pcs2he',nvap2h )
855 cgn      call gmprsx ( 'apres pcs2he',obpcap)
856 c
857         endif
858 c
859    63   continue
860 c
861 c 6.4. ==> pour les pentaedres decoupes
862 c
863         if ( nbpeto.ne.0 ) then
864 c
865 cgn      call gmprsx ( 'avant pcsipe',nvap2h )
866 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
867 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*nbnoto)
868 #ifdef _DEBUG_HOMARD_
869       write (ulsort,texte(langue,3)) 'PCSIPE', nompro
870 #endif
871           call pcsipe ( nbfop2, imem(adpcap), rmem(pvap2h),
872      >                  hetpen, facpen, cofape, filpen, fppyte,
873      >                  somare, np2are,
874      >                  aretri, arequa,
875      >                  tritet, cotrte,
876      >                  facpyr, cofapy,
877      >                  ulsort, langue, codret )
878 c
879         endif
880 c
881 c 6.5. ==> pour les triangles decoupes
882 c
883         if ( nbtrma.ne.0 ) then
884 c
885 #ifdef _DEBUG_HOMARD_
886       write (ulsort,texte(langue,3)) 'PCSITR', nompro
887 #endif
888           call pcsitr ( nbfop2, imem(adpcap), rmem(pvap2h),
889      >                  hettri, aretri, filtri,
890      >                  somare, np2are )
891 cgn      write (*,*)'apres pcsitr'
892 cgn      write (*,*)(imem(adpcap-1+iaux),iaux=1,nbnoto)
893 cgn      write (*,92010)(rmem(pvap2h-1+iaux),iaux=1,2*nbnoto)
894 c
895         endif
896 c
897 c====
898 c 7. Interpolation des variables aux noeuds P1 vers P2
899 c    Les nouveaux noeuds P2 sont tous des milieux d'aretes (cf mmcnp2)
900 c====
901 c
902       elseif ( typin0.eq.4 ) then
903 c
904 #ifdef _DEBUG_HOMARD_
905       write (ulsort,90002) '7. modification P1 vers P2 ; codret', codret
906 #endif
907 c
908         write (ulsort,texte(langue,10)) 'P1'
909 c
910 #ifdef _DEBUG_HOMARD_
911       write (ulsort,texte(langue,3)) 'PCSMAR', nompro
912 #endif
913         call pcsmar ( nbfop2, imem(adpcap),
914      >                somare, np2are,
915      >                rmem(pvap2h) )
916 cgn      write (*,*)'apres pcsmar'
917 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=2073,2073)
918 cgn      write(*,92010) (rmem(pvap2h-1+iaux),iaux=1,nbnoto)
919 c
920 c====
921 c 8. Interpolation des variables aux noeuds P2 vers P1
922 c    Rien n'est a faire puisque la copie des valeurs sur les noeuds P1
923 c    a eu lieu au depart
924 c====
925 c
926       elseif ( typin0.eq.5 ) then
927 c
928 #ifdef _DEBUG_HOMARD_
929       write (ulsort,90002) '8. modification P2 vers P1 ; codret', codret
930 #endif
931 c
932         write (ulsort,texte(langue,10)) 'P1'
933 c
934 c====
935 c 9. Type inconnu
936 c====
937 c
938       else
939 c
940 #ifdef _DEBUG_HOMARD_
941       write (ulsort,90002) '9. type inconnu ; codret', codret
942 #endif
943 c
944         write (ulsort,texte(langue,109))
945         codret = 90
946 c
947       endif
948 c
949       endif
950 c
951 c====
952 c 10. redistribution des valeurs dans la numerotation du calcul
953 c====
954 #ifdef _DEBUG_HOMARD_
955       write (ulsort,90002) '10. ==> redistribution ; codret', codret
956 #endif
957 cgn       write (ulsort,90002)'avant pcsrc0 : nbfop1, nbfop2',
958 cgn     >                            nbfop1, nbfop2
959 cgn      do 888 , iaux = 1 , nbnoto
960 cgn        if ( imem(adpcap-1+iaux).eq.0 ) then
961 cgn          write(ulsort,*) 'Noeud', iaux
962 cgn        endif
963 cgn  888 continue
964 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=1,nbnoto)
965 cgn      write (*,92010)(rmem(pvap1h-1+iaux),iaux=1,nbfop1*rsnoto-1)
966 c
967       if ( nbfop1.ne.0 ) then
968 c
969         if ( codret.eq.0 ) then
970 c
971 #ifdef _DEBUG_HOMARD_
972       write (ulsort,texte(langue,3)) 'PCSRC0_d1', nompro
973 #endif
974         call pcsrc0 ( nbfop1, rsnoto,
975      >                imem(adpcap), nnosho,
976      >                rmem(pvap1h), rmem(p1valr) )
977 c
978 cgn       write (*,*)'apres pcsrc0'
979 cgn      write(*,91011) (imem(adpcap-1+iaux),iaux=2073,2073)
980 cgn      write(*,92010) (rmem(p1valr-1+iaux),iaux=2073,2073)
981 cgn      write(*,92010)(rmem(iaux),iaux=p1valr,p1valr-1+nbfop1*rsnoto-1)
982         endif
983 c
984       endif
985 c
986       if ( nbfop2.ne.0 ) then
987 c
988         if ( codret.eq.0 ) then
989 c
990 #ifdef _DEBUG_HOMARD_
991       write (ulsort,texte(langue,3)) 'PCSRC0_d2', nompro
992 #endif
993         call pcsrc0 ( nbfop2, rsnoto,
994      >                imem(adpcap), nnosho,
995      >                rmem(pvap2h), rmem(p1valr) )
996 cgn       write (*,*)nbnoto, rsnoto
997 cgn       write (*,*)'apres pcsrc0'
998 cgn       write (*,92010)(rmem(iaux),iaux=p1valr,p1valr-1+nbfop2*rsnoto)
999 cgn       write (*,92010)
1000 cgn     > (rmem(p1valr-1+iaux),iaux=nbfop2*nbnoto-5,nbfop2*nbnoto)
1001 cgn      do 889 , iaux = 1 , nbnoto
1002 cgn        if ( imem(adpcap-1+iaux).eq.0 ) then
1003 cgn          write(ulsort,*) 'Noeud', iaux
1004 cgn        endif
1005 cgn  889 continue
1006 c
1007         endif
1008 c
1009       endif
1010 c
1011 c====
1012 c 11. menage
1013 c====
1014 c
1015       if ( codret.eq.0 ) then
1016 c
1017 #ifdef _DEBUG_HOMARD_
1018       write (ulsort,texte(langue,3)) 'GMLBOJ', nompro
1019 #endif
1020 c
1021       call gmlboj ( nvap1h, codre1 )
1022       call gmlboj ( nvap2h, codre2 )
1023 c
1024       codre0 = min ( codre1, codre2 )
1025       codret = max ( abs(codre0), codret,
1026      >               codre1, codre2 )
1027 c
1028       endif
1029 c
1030 c====
1031 c 10. la fin
1032 c====
1033 c
1034       if ( codret.ne.0 ) then
1035 c
1036 #include "envex2.h"
1037 c
1038       write (ulsort,texte(langue,1)) 'Sortie', nompro
1039       write (ulsort,texte(langue,2)) codret
1040 c
1041       endif
1042 cgn      write(ulsort,90002)'debut du profil',
1043 cgn     >          (imem(adpcap+iaux),iaux=0,4)
1044 cgn      write(ulsort,90004)'debut des valeurs',
1045 cgn     >          (rmem(pvap2h+iaux),iaux=0,4)
1046 cgn        do 999 , iaux = 1 , nbnoto
1047 cgn         if ( imem(adpcap-1+iaux).eq.0 ) then
1048 cgn          write(ulsort,90002) 'oubli de', iaux
1049 cgn            endif
1050 cgn  999 continue
1051 c
1052 #ifdef _DEBUG_HOMARD_
1053       write (ulsort,texte(langue,1)) 'Sortie', nompro
1054       call dmflsh (iaux)
1055 #endif
1056 c
1057       end