Télécharger yns.eso

Retour à la liste

Numérotation des lignes :

yns
  1. C YNS SOURCE FANDEUR 22/01/03 21:15:58 11136
  2. SUBROUTINE YNS
  3. PARAMETER (NTB=1)
  4. IMPLICIT INTEGER(I-N)
  5. IMPLICIT REAL*8 (A-H,O-Z)
  6. C***********************************************************************
  7. C
  8. C CET OPERATEUR DISCRETISE LES EQUATIONS DE NAVIER STOKES
  9. C EN 2D SUR LES ELEMENTS QUA4 ET TRI3 PLAN OU AXI
  10. C EN 3D SUR LES ELEMENTS CUB8 ET PRI6
  11. C LES OPERATEURS SONT "SOUS-INTEGRES"
  12. C
  13. C SYNTAXE :
  14. C ---------
  15. C 1/ Cas incompressible
  16. C
  17. C du/dt + u Grad u = nu Lapl u - 1/ro Grad p < + S >
  18. C < + g beta (T-Tref) >
  19. C
  20. C 'OPER' 'NS' nu 'INCO' UN :
  21. C 'OPER' 'NS' nu s 'INCO' UN :
  22. C 'OPER' 'NS' nu gb tn tref 'INCO' UN :
  23. C
  24. C
  25. C 2/ Cas compressible
  26. C
  27. C dG/dt + u Grad G + G Div u = mu Lapl u - Grad p < + S >
  28. C
  29. C 'OPER' 'NS' mu un 'INCO' GN :
  30. C 'OPER' 'NS' mu un s 'INCO' GN :
  31. C
  32. C
  33. C nu,mu viscosité cinématique resp. dynamique
  34. C FLOTTANT ou CHPOINT SCAL CENTRE
  35. C s source volumique de qdm
  36. C POINT ou CHPOINT VECT CENTRE
  37. C gb coéfficient de flottabilité (g*beta où g est l'accéllération
  38. C de la pesanteur et beta le coéfficient de dilatabilité)
  39. C POINT ou CHPOINT VECT CENTRE
  40. C tn Champ de température CHPOINT SCAL SOMMET
  41. C tref température de référence
  42. C FLOTTANT ou CHPOINT SCAL SOMMET
  43. C
  44. C Champ de vitesse -> VITESS
  45. C un Champ de vitesse transportant -> UTRANS
  46. C CHPOINT VECT SOMMET
  47. C gn Champ de vitesse massique (transporté) -> IZTU1 (Inconnue)
  48. C CHPOINT VECT SOMMET
  49. C
  50. C************************************************************************
  51.  
  52.  
  53. -INC PPARAM
  54. -INC CCOPTIO
  55. -INC CCGEOME
  56. -INC SIZFFB
  57. POINTEUR IZF1.IZFFM,IZH2.IZHR
  58. PARAMETER (LRV=64)
  59. SEGMENT PETROV
  60. REAL*8 WT(LRV,NP,NPG,KDIM),WS(LRV,NP,NPG,KDIM),HK(LRV,IDIM,NP,NPG)
  61. REAL*8 UIL(LRV,IDIM,NPG),DUIL(LRV,IDIM,NPG)
  62. REAL*8 PGSK(LRV,NPG),RPGK(LRV,NPG),AIRE(LRV),ANUK(LRV),COEFK(LRV)
  63. REAL*8 AJK(LRV,IDIM,IDIM,NPG)
  64. ENDSEGMENT
  65.  
  66. SEGMENT SOURCE
  67. REAL*8 SO(NBEL,NPG,KSDIM)
  68. ENDSEGMENT
  69.  
  70. -INC SMCHAML
  71. -INC SMCOORD
  72. -INC SMLENTI
  73. POINTEUR IPADI.MLENTI,IPADU.MLENTI,IPADS.MLENTI,IPADF.MLENTI
  74. POINTEUR IPADQ.MLENTI
  75. -INC SMELEME
  76. POINTEUR MELEMC.MELEME,IGEOM0.MELEME,MELEMS.MELEME
  77. POINTEUR MELEMI.MELEME,MELEP1.MELEME,MELEMK.MELEME
  78. -INC SMCHPOI
  79. POINTEUR MCHPIN.MCHPOI
  80. POINTEUR IZTU1.MPOVAL,IZGG1.MPOVAL
  81. POINTEUR VITESS.MPOVAL,UTRANS.MPOVAL
  82. POINTEUR VISCO.MPOVAL,IZTGG2.MPOVAL
  83. POINTEUR IZTGG3.MPOVAL,IZTGG4.MPOVAL
  84. POINTEUR IZVOL.MPOVAL,IZTCO.MPOVAL
  85.  
  86. POINTEUR IPM.IZAFM
  87.  
  88. SEGMENT IMATRS
  89. INTEGER LIZAFS(NBSOUS,NBME)
  90. ENDSEGMENT
  91. POINTEUR IPMS.IZAFM,IPS1.IZAFM,IPS2.IZAFM,IPS3.IZAFM
  92.  
  93. -INC SMLMOTS
  94. POINTEUR LINCO.MLMOTS
  95. LOGICAL IKMMDZ
  96. CHARACTER*8 NOMZ,NOMI,TYPE,TYPC,NOM0,NOMA,NOM,MTYP,CHAI
  97. CHARACTER*4 NOM4(3)
  98. CHARACTER*8 LTAB(NTB)
  99. DIMENSION KTAB(NTB),IXV(9),RO(1)
  100. SAVE IPAS
  101. DATA LTAB/'KIZX '/,IPAS/0/
  102. C*****************************************************************************
  103. CNS
  104. C write(6,*)' DEBUT NS '
  105. segact mcoord
  106. C Nouvelle directive EQUA de EQEX
  107. MTYP=' '
  108. ro(1)=1.D0
  109. CALL QUETYP(MTYP,0,IRET)
  110. IF(IRET.EQ.0)THEN
  111. C% On attend un des objets : %m1:8 %m9:16 %m17:24 %m25:32 %m33:40
  112. MOTERR( 1: 8) = 'CHAI '
  113. MOTERR( 9:16) = 'MMODEL '
  114. MOTERR(17:24) = 'TABLE '
  115. CALL ERREUR(472)
  116. RETURN
  117. ENDIF
  118.  
  119. IF(MTYP.EQ.'MMODEL')THEN
  120. CALL YTCLSF('TCLS ','NS ')
  121. RETURN
  122. ELSEIF(MTYP.EQ.'MOT ')THEN
  123. CALL LIRCHA(CHAI,1,IRET)
  124. CALL YTCLSF(CHAI,'NS ')
  125. RETURN
  126. ENDIF
  127. C Fin Nouvelle directive EQUA de EQEX
  128. C
  129. C- Récupération de la table xNS type KIZX (pointeur MTABX)
  130. C
  131. C write(6,*)'Récupération de la table xNS type KIZX pointeur MTABX'
  132. CALL LIROBJ('TABLE',MTABX,1,IRET)
  133. IF(IRET.EQ.0)THEN
  134. C% On ne trouve pas d'objet de type %m1:8
  135. MOTERR( 1: 8) = 'TABLE'
  136. CALL ERREUR(37)
  137. RETURN
  138. ENDIF
  139. C A tout hazard on regarde si une table ne peut en cacher une autre
  140. C write(6,*)' A tout hazard on regarde si une table ne peut en',
  141. C &' cacher une autre !!!!!'
  142.  
  143. MTYP='TABLE'
  144. CALL QUETYP(MTYP,0,IRET)
  145. IF(IRET.EQ.1)THEN
  146. CALL LIROBJ('TABLE',MTABX,1,IRET)
  147. IF(IRET.EQ.0)THEN
  148. C Indice %m1:8 : Indice %m9:16 non trouvé dans la table %m17:24
  149. MOTERR( 1: 8) = ' MTABX '
  150. MOTERR( 9:16) = ' MTABX '
  151. MOTERR(17:24) = ' KIZX '
  152. CALL ERREUR(786)
  153. RETURN
  154. ENDIF
  155. C write(6,*)' EH BEN C EST LE CAS une KIZX // !!!!!!!!!!'
  156. cc call pplist(MTABX)
  157. TYPE=' '
  158. CALL ACMO(MTABX,'DOMZ',TYPE,MMDZ)
  159. CALL LEKMOD(MMDZ,MTABZ,INEFMD)
  160. TYPE='TABLE'
  161. CALL ECMO(MTABX,'TDOMZ',TYPE,MTABZ)
  162. ELSE
  163. TYPE=' '
  164. CALL ACMO(MTABX,'DOMZ',TYPE,MMDZ)
  165. IF(TYPE.EQ.'MMODEL')THEN
  166. CALL LEKMOD(MMDZ,MTABZ,INEFMD)
  167. TYPE='TABLE'
  168. CALL ECMO(MTABX,'TDOMZ',TYPE,MTABZ)
  169. ELSE
  170. MTABZ=MMDZ
  171. TYPE='TABLE'
  172. CALL ECMO(MTABX,'TDOMZ',TYPE,MTABZ)
  173. ENDIF
  174. ENDIF
  175.  
  176.  
  177.  
  178. C.......................................................................
  179. C
  180. C- Récupération de la table RV type EQEX (pointeur MTAB1)
  181. C
  182. C write(6,*)'Récupération de la table RV type EQEX (pointeur MTAB1)'
  183. CALL LEKTAB(MTABX,'EQEX',MTAB1)
  184. IF(MTAB1.EQ.0)THEN
  185. C Indice %m1:8 : Indice %m9:16 non trouvé dans la table %m17:24
  186. MOTERR( 1: 8) = ' EQEX '
  187. MOTERR( 9:16) = ' EQEX '
  188. MOTERR(17:24) = ' KIZX '
  189. CALL ERREUR(786)
  190. RETURN
  191. ENDIF
  192. CALL ACME(MTAB1,'NAVISTOK',NASTOK)
  193. IF(NASTOK.EQ.0)THEN
  194. CALL ZNS(MTABX,MTAB1)
  195. RETURN
  196. ENDIF
  197. C
  198. C- Récupération de la table INCO (pointeur KINC)
  199. C
  200. C write(6,*)'Récupération de la table INCO (pointeur KINC)'
  201. CALL LEKTAB(MTAB1,'INCO',KINC)
  202. IF(KINC.EQ.0)THEN
  203. C Indice %m1:8 : Indice %m9:16 non trouvé dans la table %m17:24
  204. MOTERR( 1: 8) = ' INCO '
  205. MOTERR( 9:16) = ' INCO '
  206. MOTERR(17:24) = ' EQEX '
  207. CALL ERREUR(786)
  208. RETURN
  209. ENDIF
  210. C.......................................................................
  211.  
  212. CALL ACMM(MTABX,'NOMZONE',NOMZ)
  213.  
  214. MTYP='MMODEL'
  215. CALL QUETYP(MTYP,0,IRET)
  216. IF(IRET.EQ.1)THEN
  217. CALL LIROBJ('MMODEL',MMDZ,1,IRET)
  218. CALL LIROBJ('TABLE',MTABZ,1,IRET)
  219. ELSE
  220. C On récupère le MODEL partitionné ou non
  221. TYPE=' '
  222. CALL ACMO(MTABX,'DOMZ',TYPE,MMDZ)
  223. IF(TYPE.NE.'MMODEL')THEN
  224. C Indice %m1:8 : Indice %m9:16 non trouvé dans la table %m17:24
  225. MOTERR( 1: 8) = ' DOMZ '
  226. MOTERR( 9:16) = ' DOMZ '
  227. MOTERR(17:24) = ' KIZX '
  228. CALL ERREUR(786)
  229. RETURN
  230. ENDIF
  231. ENDIF
  232. C*****************************************************************************
  233. C OPTIONS
  234. C KIMPL = 0 -> EXPL 1 -> IMPL 2 -> SEMI
  235. C KFORM = 0 -> SI 1 -> EF 2 -> VF 3 -> EFMC
  236. C IDCEN = 0-> rien 1-> CENTREE 2-> SUPGDC 3-> SUPG 4-> TVISQUEU 5-> CNG
  237.  
  238. IAXI=0
  239. IF(IFOMOD.EQ.0)IAXI=2
  240. C
  241. C- Récupération de la table des options KOPT (pointeur KOPTI)
  242. C
  243. CALL LEKTAB(MTABX,'KOPT',KOPTI)
  244. IF (KOPTI.EQ.0) THEN
  245. C Indice %m1:8 : Indice %m9:16 non trouvé dans la table %m17:24
  246. MOTERR( 1: 8) = ' KOPT '
  247. MOTERR( 9:16) = ' KOPT '
  248. MOTERR(17:24) = ' KIZX '
  249. CALL ERREUR(786)
  250. RETURN
  251. ENDIF
  252.  
  253. C write(6,*)' Avant les options '
  254. CALL ACME(KOPTI,'MTRMASS ',MMPG)
  255. IPG=0
  256. IF(MMPG.EQ.3)IPG=1
  257. CALL ACME(KOPTI,'IDCEN',IDCEN)
  258. KDIM=1
  259. IF(IDCEN.EQ.2)KDIM=IDIM
  260. CALL ACME(KOPTI,'IKOMP',IKOMP)
  261. CALL ACME(KOPTI,'ALE',KKALE)
  262. CALL ACME(KOPTI,'KIMPL',KIMPL)
  263. CALL ACME(KOPTI,'KPOIN',KPRE)
  264. CALL ACME(KOPTI,'KFORM',KFORM)
  265. CALL ACME(KOPTI,'IDIV',IDIV)
  266. CALL ACMF(KOPTI,'CMD',CMD)
  267.  
  268. IF(KFORM.NE.0.AND.KFORM.NE.1)THEN
  269. C Option %m1:8 incompatible avec les données
  270. MOTERR( 1: 8) = 'EF/EFM1 '
  271. CALL ERREUR(803)
  272. RETURN
  273. ENDIF
  274. CALL ACME(KOPTI,'KMACO',KMACO)
  275. CALL ACMF(KOPTI,'AIMPL',AIMPL)
  276. AG=AIMPL
  277. AD=AIMPL-1.D0
  278. IF (IERR.NE.0) RETURN
  279.  
  280. C write(6,*)' Apres les options '
  281. C*****************************************************************************
  282. C
  283. C- Récupération de la table DOMAINE associée au domaine local
  284. C
  285.  
  286. C E/ MMODEL : Pointeur de la table contenant l'information cherchée
  287. C /S IPOINT : Pointeur sur la table DOMAINE
  288. C /S INEFMD : Type formulation INEFMD=1 LINE,=2 MACRO,=3 QUADRATIQUE
  289. C INEFMD=4 LINB
  290.  
  291. c?? CALL LEKMOD(MMDZ,MTABZ,INEFMD)
  292.  
  293. CALL LEKTAB(MTABZ,'MAILLAGE',MELEME)
  294. CALL LEKTAB(MTABZ,'SOMMET',MELEMS)
  295. CALL LEKTAB(MTABZ,'CENTRE',MELEMC)
  296. CALL LEKTAB(MTABZ,'MACRO1',MELEMK)
  297. IF (IERR.NE.0) RETURN
  298.  
  299. MELEMQ=MELEMC
  300. MELEP1=MELEMC
  301. IF(KPRE.EQ.3)THEN
  302. CALL LEKTAB(MTABZ,'CENTREP0',MELEMQ)
  303. MELEP1=MELEMQ
  304. ELSEIF(KPRE.EQ.4)THEN
  305. CALL LEKTAB(MTABZ,'CENTREP1',MELEMQ)
  306. CALL LEKTAB(MTABZ,'ELTP1NC ',MELEP1)
  307. ENDIF
  308.  
  309. C*************************************************************************
  310. C VERIFICATIONS SUR LES INCONNUES
  311. C
  312. C- Récupération du nombre d'inconnues et du nom de l'inconnue NOMI
  313. C
  314. TYPE='LISTMOTS'
  315. CALL ACMO(MTABX,'LISTINCO',TYPE,LINCO)
  316. IF (IERR.NE.0) RETURN
  317. SEGACT LINCO
  318. NBINC=LINCO.MOTS(/2)
  319. IF(NBINC.NE.1)THEN
  320. C Indice %m1:8 : contient plus de %i1 %m9:16
  321. MOTERR( 1:8) = 'LISTINCO'
  322. INTERR(1) = 1
  323. MOTERR(9:16) = ' MOTS '
  324. CALL ERREUR(799)
  325. RETURN
  326. ENDIF
  327.  
  328. NOMI=LINCO.MOTS(1)
  329. NOMA=NOMI
  330. DO 15 I=1,IDIM
  331. WRITE(NOM4(I),FMT='(I1,A3)')I,NOMI(1:3)
  332. 15 CONTINUE
  333. C
  334. C- Récupération de l'inconnue
  335. C
  336. TYPE=' '
  337. CALL ACMO(KINC,NOMI,TYPE,MCHPOI)
  338. IF(TYPE.NE.'CHPOINT ')THEN
  339. C Indice %m1:8 : ne contient pas un objet de type %m9:16
  340. MOTERR( 1: 8) = 'INC '//NOMI
  341. MOTERR( 9:16) = 'CHPOINT '
  342. CALL ERREUR(800)
  343. RETURN
  344. ELSE
  345. CALL LICHTL(MCHPOI,IZTU1,TYPC,MELEMI)
  346. MCHPIN=MCHPOI
  347. NINKO = IZTU1.VPOCHA(/2)
  348. NPTI = IZTU1.VPOCHA(/1)
  349. IF (NINKO.NE.IDIM) THEN
  350. C Indice %m1:8 : Le %m9:16 n'a pas le bon nombre de composantes
  351. MOTERR( 1: 8) = 'INC '//NOMI
  352. MOTERR( 9:16) = 'CHPOINT '
  353. CALL ERREUR(784)
  354. RETURN
  355. ENDIF
  356. C write(6,*)' MCHPOI,MELEMI=',MCHPOI,MELEMI
  357. C On fait pointer ces deux tableaux sur le champ U inconu (tjs présent) pour
  358. C eviter de les enlever lors de l'appel FORTRAN si les options sont absentes
  359. UTRANS=IZTU1
  360. IKW=0
  361. VITESS=IZTU1
  362. ENDIF
  363.  
  364. C*****************************************************************************
  365. C Le domaine de definition est donne par le SPG de la premiere inconnue
  366. C Les inconnues suivantes devront posseder ce meme pointeur
  367. C On verifie que les points de la zone sont tous inclus dans ce SPG
  368.  
  369. CALL KRIPAD(MELEMI,IPADI)
  370. IPADS=IPADI
  371. NPTS=NPTI
  372. IF(MELEMI.NE.MELEMS)THEN
  373. CALL KRIPAD(MELEMS,IPADS)
  374. NPTS=MELEMS.NUM(/2)
  375. ENDIF
  376. IPADU=IPADI
  377. NPTU=NPTI
  378. IPADF=IPADI
  379. NPTF=NPTS
  380.  
  381. IF(IPAS.EQ.0)THEN
  382. IPAS=1
  383. CALL VERPAD(IPADI,MELEME,IRET)
  384. IF(IRET.NE.0)THEN
  385. C Indice %m1:8 : L'objet %m9:16 n'a pas le bon support géométrique
  386. MOTERR(1: 8) = 'INC '//NOMI
  387. MOTERR(9:16) = 'CHPOINT '
  388. CALL ERREUR(788)
  389. RETURN
  390. ENDIF
  391. ENDIF
  392.  
  393. C*****************************************************************************
  394.  
  395.  
  396. C--Cas non conservatif et non conservatif EF (KFORM=1)
  397. IF(KFORM.EQ.1)THEN
  398. IHV=1
  399.  
  400. C WRITE(6,*)' IDCEN=',IDCEN
  401. IF(IDCEN.EQ.4.OR.IDCEN.EQ.5)THEN
  402. CALL LEKTAB(MTAB1,'PASDETPS',MTABT)
  403. C CALL LEKTAB(MTABZ,'PASDETPS',MTABT)
  404. C write(6,*)' LEKTAB rend MTABT=',MTABT,' <==============='
  405. IF(MTABT.EQ.0)THEN
  406. C write(6,*)' On cree une table pasdetps !!!!!!!!!!!!!!!!!!!!!!'
  407. CALL CRTABL(MTABT)
  408. CALL ECMM(MTABT,'SOUSTYPE','PASDETPS')
  409. C CALL ECMO(MTABZ,'PASDETPS','TABLE ',MTABT)
  410. CALL ECMO(MTAB1,'PASDETPS','TABLE ',MTABT)
  411. DT=1.E30
  412. DTP=1.E30+DT
  413. IPAT=1
  414. CALL ECME(MTABT,'NUPASDT',IPAT)
  415. CALL ECMF(MTABT,'DELTAT',DTP)
  416. ELSE
  417. CALL ACMF(MTABT,'DELTAT',DTP)
  418. ENDIF
  419.  
  420. CMT=CMD
  421. DT=0.D0
  422.  
  423. C IDCEN 1-> CENTREE 2-> SUPGDC 3-> SUPG 4-> TVISQUEU 5-> CNG
  424. TYPE=' '
  425. CALL ACMO(MTABT,'DELTAT',TYPE,IENT)
  426. IF(TYPE.NE.'ENTIER')THEN
  427. WRITE(6,*)' Opérateur NS '
  428. WRITE(6,*)' On reclame un pas de temps'
  429. C Option %m1:8 incompatible avec les données
  430. MOTERR( 1: 8) = ' EF '
  431. CALL ERREUR(803)
  432. RETURN
  433. ELSE
  434. CALL ACMF(MTABT,'DELTAT',DT)
  435. CMT=DT
  436. ENDIF
  437. ENDIF
  438.  
  439. CALL YCLS('NS ',MTABX,MTABZ,IHV,MCHPIN,NOMI,MATRIK,MCHPO1)
  440.  
  441. CALL ECROBJ('MATRIK',MATRIK)
  442. CALL ECROBJ('CHPOINT',MCHPO1)
  443.  
  444. C write(6,*)' FIN YNSSSSSSSSSSSSSSS'
  445. RETURN
  446. ENDIF
  447.  
  448.  
  449. C*************************************************************************
  450. C Lecture des coefficients
  451. C Type du coefficient :
  452. C IK1=0 CHPOINT IK1=1 scalaire IK1=2 vecteur
  453.  
  454. CALL ACME(MTABX,'IARG',IARG)
  455. IF(IKOMP.EQ.0)THEN
  456. IF(IARG.NE.1.AND.IARG.NE.2.AND.IARG.NE.3.AND.IARG.NE.4.
  457. & AND.IARG.NE.5.AND.IARG.NE.6)THEN
  458. WRITE(6,*)' Opérateur NS : option incompressible '
  459. WRITE(6,*)' Nombre d''arguments incorrect : ',IARG
  460. WRITE(6,*)' On attend 1 , 2 , 3 , 4 , 5 ou 6'
  461. C Indice %m1:8 : nombre d'arguments incorrect
  462. MOTERR(1:8) = 'IARG '
  463. CALL ERREUR(804)
  464. RETURN
  465. ENDIF
  466. ELSEIF(IKOMP.EQ.1.AND.KFORM.NE.1)THEN
  467. IF(IARG.NE.2.AND.IARG.NE.3)THEN
  468. WRITE(6,*)' Opérateur NS : option conservative'
  469. WRITE(6,*)' Nombre d''arguments incorrect : ',IARG
  470. WRITE(6,*)' On attend 2 ou 3 '
  471. C Indice %m1:8 : nombre d'arguments incorrect
  472. MOTERR(1:8) = 'IARG '
  473. CALL ERREUR(804)
  474. RETURN
  475. ENDIF
  476. ENDIF
  477.  
  478. IBOU=0
  479. IF(KKALE.NE.0)IARG=IARG-1
  480.  
  481. IF(IKOMP.EQ.0)THEN
  482. C--Cas non conservatif
  483. IXV(1)=MELEMC
  484. IXV(2)=1
  485. IXV(3)=0
  486. CALL LEKCOF('Opérateur NS :',
  487. & MTABX,KINC,1,IXV,IZTG1,VISCO,NPT1,NC1,IK1,IRET)
  488. IF(IRET.EQ.0)RETURN
  489.  
  490. IK2=-1
  491. IKG=-1
  492. IK3=-1
  493. IK4=-1
  494. IZTGG2=VISCO
  495. IZTGG3=VISCO
  496. IZTGG4=VISCO
  497.  
  498. C 2ème coef -> source S
  499. IF(IARG.GE.2)THEN
  500. C?? IXV(1)=-MELEMQ
  501. IXV(1)=-MELEMC
  502. IXV(2)=0
  503. IXV(3)=1
  504. IXV(4)=-MELEMS
  505. c? IRET=3
  506. c? IF(KPRE.EQ.5)THEN
  507. IRET = 4
  508. c? ENDIF
  509. CALL LEKCOF('Opérateur NS :',
  510. & MTABX,KINC,2,IXV,IZTG2,IZTGG2,NELG,NC2,IKG,IRET)
  511. IF(IRET.EQ.0)RETURN
  512.  
  513. IF(IKG.EQ.2)IKG=1
  514.  
  515.  
  516. IF(IARG.EQ.4)THEN
  517. C Cas incompressible Boussinesq
  518. IBOU = 1
  519. IXV(1)=MELEMS
  520. IXV(2)=0
  521. IXV(3)=0
  522. CALL LEKCOF('Opérateur NS :',
  523. & MTABX,KINC,3,IXV,IZTG3,IZTGG3,NPT3,NC3,IK3,IRET)
  524. IF(IRET.EQ.0)RETURN
  525.  
  526. IXV(1)=MELEMS
  527. IXV(2)=1
  528. IXV(3)=0
  529. CALL LEKCOF('Opérateur NS :',
  530. & MTABX,KINC,4,IXV,IZTG4,IZTGG4,NPT4,NC4,IK4,IRET)
  531. IF(IRET.EQ.0)RETURN
  532.  
  533. ENDIF
  534. ENDIF
  535.  
  536.  
  537. C--Cas compressible
  538. ELSEIF(IKOMP.EQ.1)THEN
  539. C 1er coef : mu , viscosité dynamique
  540. IXV(1)=MELEMC
  541. IXV(2)=1
  542. IXV(3)=0
  543. CALL LEKCOF('Opérateur NS :',
  544. & MTABX,KINC,1,IXV,IZTG1,VISCO,NPT1,NC1,IK1,IRET)
  545. IF(IRET.EQ.0)RETURN
  546.  
  547. IK2=-1
  548. IKG=-1
  549. IK3=-1
  550. IK4=-1
  551. IZTGG2=VISCO
  552. IZTGG3=VISCO
  553. IZTGG4=VISCO
  554.  
  555. IF(IARG.EQ.3)THEN
  556. C 3ème coef : rog
  557. IXV(1)=-MELEMC
  558. IXV(2)=0
  559. IXV(3)=1
  560. IXV(4)=-MELEMS
  561. IRET=3
  562. IF(KPRE.EQ.5)THEN
  563. IRET = 4
  564. ENDIF
  565. CALL LEKCOF('Opérateur NS :',
  566. & MTABX,KINC,3,IXV,IZTG2,IZTGG2,NELG,NC2,IKG,IRET)
  567. IF(IRET.EQ.0)RETURN
  568.  
  569. IF(IKG.EQ.2)IKG=1
  570.  
  571. ENDIF
  572.  
  573. C 2ème coef : un , champ de vitesse transportant
  574. IXV(1)=-MELEMS
  575. IXV(2)=0
  576. IXV(3)=0
  577. CALL LEKCOF('Opérateur NS :',
  578. & MTABX,KINC,2,IXV,MUTRAN,UTRANS,NPTU,NC3,IKW,IRET)
  579. IF(IRET.EQ.0)RETURN
  580. IPADU=IPADS
  581. NPTU=NPTS
  582. VITESS=UTRANS
  583. NPTF=NPTS
  584.  
  585. ENDIF
  586.  
  587. IF(KKALE.NE.0)THEN
  588.  
  589. IF(IKOMP.EQ.1)THEN
  590. WRITE(6,*)' Operateur NS : option compressible en ALE'
  591. WRITE(6,*)' Cas non prévu pour l''instant: '
  592. C Option %m1:8 incompatible avec les données
  593. MOTERR( 1: 8) = ' EF '
  594. CALL ERREUR(803)
  595. RETURN
  596. ENDIF
  597.  
  598. C n+ unième coefficient champ de vitesse transportant pour l'option ALE
  599. IXV(1)=-MELEMS
  600. IXV(2)=0
  601. IXV(3)=0
  602. CALL LEKCOF('Opérateur NS :',
  603. & MTABX,KINC,(IARG+1),IXV,MUTRAN,UTRANS,NPTW,NCW,IKW,IRET)
  604. IF(IRET.EQ.0)RETURN
  605. IPADU=IPADS
  606.  
  607. ENDIF
  608.  
  609. C write(6,*)' Operateur NS : Fin lecture Arguments '
  610. C Fin lecture Arguments ************************************************
  611.  
  612. IF(IPADS.EQ.0)CALL KRIPAD(MELEMS,IPADS)
  613.  
  614. CALL LEKTAB(MTAB1,'PASDETPS',MTABT)
  615. C CALL LEKTAB(MTABZ,'PASDETPS',MTABT)
  616. C write(6,*)' LEKTAB rend MTABT=',MTABT,' <==============='
  617. IF(MTABT.EQ.0)THEN
  618. C write(6,*)' On cree une table pasdetps !!!!!!!!!!!!!!!!!!!!!!'
  619. CALL CRTABL(MTABT)
  620. CALL ECMM(MTABT,'SOUSTYPE','PASDETPS')
  621. C CALL ECMO(MTABZ,'PASDETPS','TABLE ',MTABT)
  622. CALL ECMO(MTAB1,'PASDETPS','TABLE ',MTABT)
  623. DT=1.E30
  624. DTP=1.E30+DT
  625. IPAT=1
  626. CALL ECME(MTABT,'NUPASDT',IPAT)
  627. ELSE
  628. CALL ACMF(MTABT,'DELTAT',DTP)
  629. ENDIF
  630.  
  631. C*********** FORMULATIONS **********
  632.  
  633. IF(KFORM.EQ.0)THEN
  634. C Formulation EFM1
  635. IF(INEFMD.NE.1.AND.INEFMD.NE.2)THEN
  636. C Option %m1:8 incompatible avec les données
  637. MOTERR( 1: 8) = ' EFM1 '
  638. CALL ERREUR(803)
  639. RETURN
  640. ENDIF
  641.  
  642. IF(KIMPL.NE.0)THEN
  643. WRITE(6,*)' Operateur NS '
  644. C Option %m1:8 incompatible avec les données
  645. MOTERR( 1: 8) = ' EF '
  646. CALL ERREUR(803)
  647. RETURN
  648. ENDIF
  649.  
  650. CALL LEKTAB(MTABZ,'XXDXDY',MCHPOI)
  651. IF(MCHPOI.EQ.0)GO TO 90
  652. CALL LICHTL(MCHPOI,IZTCO,TYPC,IGEOM)
  653. NELZ=IZTCO.VPOCHA(/1)
  654.  
  655. CALL LEKTAB(MTABZ,'MATESI',MATRIK)
  656. IF(MATRIK.EQ.0)GO TO 90
  657. SEGACT MATRIK
  658.  
  659. IMATRI=IRIGEL(4,1)
  660. SEGACT IMATRI
  661.  
  662. CALL LEKTAB(MTABZ,'XXPSOML',MCHELM)
  663. IF(MCHELM.EQ.0)GO TO 90
  664. SEGACT MCHELM
  665.  
  666. CALL LEKTAB(MTABZ,'XXVOLUM',MCHVOL)
  667. IF(MCHVOL.EQ.0)GO TO 90
  668. CALL LICHTL(MCHVOL,IZVOL,TYPC,IGEOM)
  669.  
  670. NC=IZTU1.VPOCHA(/2)
  671. TYPE='SOMMET'
  672. CALL KRCHPT(TYPE,MELEMS,NC,IZG1,NOM4)
  673.  
  674. CALL LICHTM(IZG1,IZGG1,TYPC,IGEOM)
  675.  
  676. SEGACT MELEME
  677. NUTOEL=0
  678. DT=1.E30
  679.  
  680.  
  681. DO 1 L=1,MAX(1,LISOUS(/1))
  682. IPT1=MELEME
  683. IF(LISOUS(/1).NE.0)IPT1=LISOUS(L)
  684. SEGACT IPT1
  685. IZAFM=LIZAFM(L,1)
  686. IPM1=IZAFM
  687. SEGACT IZAFM
  688. IF(IAXI.NE.0)THEN
  689. IPM1=LIZAFM(L,2)
  690. SEGACT IPM1
  691. ENDIF
  692.  
  693. NOM0=NOMS(IPT1.ITYPEL)//' '
  694. CALL KALPBG(NOM0,'FONFORM0',IZFFM)
  695. SEGACT IZFFM*MOD
  696.  
  697. MCHAML=ICHAML(L)
  698. SEGACT MCHAML
  699. MELVAL=IELVAL(1)
  700. SEGACT MELVAL
  701.  
  702. IF(IMACHE(L).NE.IPT1)THEN
  703. goto 90
  704. ENDIF
  705.  
  706. NP =IPT1.NUM(/1)
  707. NBEL=IPT1.NUM(/2)
  708. IES=IDIM
  709.  
  710. CALL YCVI(AM,IPM1.AM,VELCHE,
  711. & IPT1.NUM,NBEL,NUTOEL,NPTS,IES,NP,IAXI,IKOMP,IARG,
  712. & VISCO.VPOCHA,IK1,IZTGG2.VPOCHA,IKG,NELG,IZTGG3.VPOCHA,IK3,
  713. & IZTGG4.VPOCHA,IK4,IPADS.LECT,
  714. & UTRANS.VPOCHA,IPADU.LECT,NPTU,
  715. & IZTU1.VPOCHA,IZGG1.VPOCHA,IPADI.LECT,
  716. & VITESS.VPOCHA,IPADF.LECT,NPTF,
  717. & IZVOL.VPOCHA,IZTCO.VPOCHA,NELZ,IDCEN,IPG,
  718. & DTM1,DT,DTT1,DTT2,NUEL,DIAEL,IZFFM.FN)
  719.  
  720. SEGDES IZAFM,IPT1
  721. IF(IAXI.NE.0)SEGDES IPM1
  722. NUTOEL=NUTOEL+NBEL
  723.  
  724. 1 CONTINUE
  725. SEGDES MELEME
  726. SEGDES IZTCO
  727.  
  728.  
  729. IF(DT.LT.DTP)THEN
  730. CALL ECMF(MTABT,'DELTAT',DT)
  731. CALL ECMM(MTABT,'OPER','NS')
  732. CALL ECMM(MTABT,'ZONE',NOMZ)
  733. CALL ECMF(MTABT,'DTCONV',DTT1)
  734. CALL ECMF(MTABT,'DTDIFU',DTT2)
  735. CALL ECMF(MTABT,'DIAEL',DIAEL)
  736. CALL ECME(MTABT,'NUEL',NUEL)
  737. ENDIF
  738.  
  739. SEGDES VISCO,IZTGG2,IZTGG3,IZTGG4
  740. SEGDES UTRANS,IZTU1
  741. SEGDES IZGG1
  742. SEGDES IZVOL,IZTCO
  743. SEGDES LINCO
  744. SEGDES MATRIK,IMATRI
  745.  
  746. NRIGE=7
  747. NKID =9
  748. NKMT =7
  749. NMATRI=0
  750. SEGINI MATRIK
  751. CALL ECROBJ('MATRIK',MATRIK)
  752. CALL ECROBJ('CHPOINT',IZG1)
  753.  
  754. C*************************************************************************
  755. ELSE IF(KFORM.EQ.1)THEN
  756. C CAS FORMULATION EF
  757.  
  758. IF(KIMPL.EQ.0)THEN
  759. C Option %m1:8 incompatible avec les données
  760. MOTERR( 1: 8) = ' EF '
  761. CALL ERREUR(803)
  762. RETURN
  763. ENDIF
  764.  
  765. DT=0.D0
  766. IF(IDCEN.EQ.4)THEN
  767. TYPE=' '
  768. CALL ACMO(MTABT,'DELTAT',TYPE,IENT)
  769. IF(TYPE.NE.'ENTIER')THEN
  770. WRITE(6,*)' Opérateur NS '
  771. WRITE(6,*)' On reclame un pas de temps'
  772. C Option %m1:8 incompatible avec les données
  773. MOTERR( 1: 8) = ' EF '
  774. CALL ERREUR(803)
  775. RETURN
  776. ELSE
  777. CALL ACMF(MTABT,'DELTAT',DT)
  778. ENDIF
  779. ENDIF
  780.  
  781. IHV=1
  782. NUTOEL=0
  783. NINKO=VITESS.VPOCHA(/2)
  784. SEGACT MELEME
  785.  
  786. TYPE=' '
  787. CALL ACMO(MTABX,'MATELM',TYPE,MATRIK)
  788. IF(TYPE.EQ.'MATRIK'.AND.KMACO.NE.0)THEN
  789. SEGACT MATRIK
  790. NMATRI=IRIGEL(/2)
  791. MELEME=IRIGEL(1,1)
  792. SEGACT MELEME
  793. IMATRI=IRIGEL(4,1)
  794. SEGACT IMATRI
  795. NBME=LIZAFM(/2)
  796. NINKO=NBME
  797. MELEMS=KSPGP
  798.  
  799.  
  800. ELSE
  801.  
  802. NRIGE=7
  803. NKID =9
  804. NKMT =7
  805. NMATRI=1
  806. SEGINI MATRIK
  807. IRIGEL(1,1)=MELEME
  808. IRIGEL(2,1)=MELEME
  809.  
  810. IRIGEL(7,1)=2
  811. NBME=NINKO
  812.  
  813.  
  814. C --- CALCUL DU NOMBRE DE PAQUETS DE LRV ELEMENTS ---
  815. CALL KPLRVM(MELEME,LRV,NBSOUS)
  816. SEGACT MELEME
  817.  
  818. SEGINI IMATRI,IMATRS
  819. IRIGEL(4,1)=IMATRI
  820. KSPGP=MELEMS
  821. KSPGD=MELEMS
  822. IF(NBME.EQ.1)THEN
  823. LISPRI(1)=NOMI(1:4)//' '
  824. LISDUA(1)=NOMA(1:4)//' '
  825. ELSE
  826. DO 102 I=1,NBME
  827. WRITE(NOM,FMT='(I1,A7)')I,NOMI(1:7)
  828. LISPRI(I)=NOM(1:4)//' '
  829. WRITE(NOM,FMT='(I1,A7)')I,NOMA(1:7)
  830. LISDUA(I)=NOM(1:4)//' '
  831. 102 CONTINUE
  832. ENDIF
  833.  
  834. LL=0
  835. NUTOEL=0
  836. DO 101 L=1,MAX(1,LISOUS(/1))
  837. IPT1=MELEME
  838. IF(LISOUS(/1).NE.0)IPT1=LISOUS(L)
  839. SEGACT IPT1
  840. NOM0=NOMS(IPT1.ITYPEL)//' '
  841. CALL KALPBG(NOM0,'FONFORM ',IZFFM)
  842. SEGACT IZFFM*MOD
  843. IZHR=KZHR(1)
  844. SEGACT IZHR*MOD
  845. NES=GR(/1)
  846. NPG=GR(/3)
  847.  
  848. NP = IPT1.NUM(/1)
  849. MP = NP
  850. NBELEM=IPT1.NUM(/2)
  851.  
  852. NNN=MOD(NBELEM,LRV)
  853. IF(NNN.EQ.0) NPACK=NBELEM/LRV
  854. IF(NNN.NE.0) NPACK=1+(NBELEM-NNN)/LRV
  855.  
  856. DO 701 KPACK=1,NPACK
  857.  
  858. C --- POUR CHAQUE PAQUET DE LRV ELEMENTS ou moins
  859. KDEB=1+(KPACK-1)*LRV
  860. KFIN=MIN(NBELEM,KDEB+LRV-1)
  861.  
  862. NBEL=KFIN-KDEB+1
  863. LL=LL+1
  864.  
  865. SEGINI IPM1,IPS1
  866. LIZAFM(LL,1)=IPM1
  867. LIZAFS(LL,1)=IPS1
  868. IPM2=IPM1
  869. IPM3=IPM1
  870. IPS2=IPS1
  871. IPS3=IPS1
  872. IF(NBME.GE.2)THEN
  873. IF(IDCEN.EQ.2)SEGINI IPM2,IPS2
  874. LIZAFM(LL,2)=IPM2
  875. LIZAFS(LL,2)=IPS2
  876. ENDIF
  877. IF(NBME.GE.3)THEN
  878. IF(IDCEN.EQ.2)SEGINI IPM3,IPS3
  879. LIZAFM(LL,3)=IPM3
  880. LIZAFS(LL,3)=IPS3
  881. ENDIF
  882.  
  883.  
  884. KITT=2
  885. KJTT=IK1
  886.  
  887. SEGINI PETROV
  888.  
  889. C write(6,*)' YNS appel a XCONV '
  890. CALL XCONV(FN,GR,PG,XYZ,HR,PGSQ,RPG,AJ,
  891. & NES,IDIM,NP,NPG,IAXI,AG,AD,IDIV,CMD,IKOMP,LRV,
  892. & WT,WS,HK,PGSK,RPGK,AIRE,AJK,UIL,DUIL,COEFK,ANUK,
  893. & RO,1,UTRANS.VPOCHA,NPTU,IPADU.LECT,VISCO.VPOCHA,IK1,
  894. & IPT1.NUM(1,KDEB),NBEL,NUTOEL,XCOOR,
  895. & IPM1.AM,IPM2.AM,IPM3.AM,
  896. & IPS1.AM,IPS2.AM,IPS3.AM,
  897. & NINKO,IDCEN,DT,
  898. & IZTU1.VPOCHA,NPTI,IPADI.LECT)
  899. SEGSUP PETROV
  900.  
  901. C write(6,*)' YNS retour XCONV '
  902. CALL XLAPL(FN,GR,PG,XYZ,HR,PGSQ,RPG,NES,IDIM,NP,NPG,IAXI,
  903. & VISCO.VPOCHA,VISCO.VPOCHA,VISCO.VPOCHA,KITT,KJTT,IK1,
  904. & IPT1.NUM(1,KDEB),NBEL,NUTOEL,XCOOR,AIMPL,IKOMP,
  905. & IPM1.AM,IPM2.AM,IPM3.AM,
  906. & IPS1.AM,IPS2.AM,IPS3.AM,
  907. & NINKO,IHV,IARG,VISCO.VPOCHA)
  908. C write(6,*)' YNS retour XLAPL '
  909.  
  910. NUTOEL=NUTOEL+NBEL
  911. SEGDES IPM1,IPS1
  912. IF(NBME.GE.2.AND.IDCEN.EQ.2)SEGDES IPM2,IPS2
  913. IF(NBME.GE.3.AND.IDCEN.EQ.2)SEGDES IPM3,IPS3
  914.  
  915. 701 CONTINUE
  916. SEGDES IPT1
  917. 101 CONTINUE
  918.  
  919. ENDIF
  920.  
  921.  
  922. IF(KIMPL.EQ.2.OR.KIMPL.EQ.0.OR.IARG.GT.1)THEN
  923. NAT=2
  924. NSOUPO=1
  925. SEGACT MELEMS
  926. N=MELEMS.NUM(/2)
  927. NC=NINKO
  928. SEGINI MCHPO1,MSOUP1,MPOVA1
  929. MCHPO1.IFOPOI=IFOUR
  930. MCHPO1.MOCHDE=' '
  931. MCHPO1.MTYPOI='SMBR'
  932. MCHPO1.JATTRI(1)=2
  933. MCHPO1.IPCHP(1)=MSOUP1
  934. DO 177 N=1,NINKO
  935. MSOUP1.NOCOMP(N)=LISDUA(N)
  936. 177 CONTINUE
  937. MSOUP1.IGEOC=MELEMS
  938. MSOUP1.IPOVAL=MPOVA1
  939. ENDIF
  940.  
  941. IF(IARG.EQ.2.OR.IARG.EQ.4)THEN
  942. CALL KRIPAD(MELEMQ,IPADQ)
  943. IF(IARG.EQ.2)THEN
  944. IKAS=2
  945. ELSEIF(IARG.EQ.4)THEN
  946. IKAS=3
  947. ENDIF
  948.  
  949.  
  950.  
  951. IF(INEFMD.EQ.2.AND.KPRE.EQ.2)MELEMK=MELEME
  952. IF(KPRE.EQ.5)MELEMK=MELEME
  953. SEGACT MELEMK
  954.  
  955. NUTOEL=0
  956. DO 1102 L=1,MAX(1,MELEMK.LISOUS(/1))
  957. IPT1=MELEMK
  958. IF(MELEMK.LISOUS(/1).NE.0)IPT1=MELEMK.LISOUS(L)
  959. SEGACT IPT1
  960.  
  961.  
  962. C write(6,*)' KPRE.=',KPRE,' IBOU=',IBOU
  963. IF(KPRE.EQ.5)THEN
  964.  
  965. NOM0=NOMS(IPT1.ITYPEL)//' '
  966.  
  967. CALL KALPBG(NOM0,'FONFORM ',IZFFM)
  968.  
  969.  
  970. SEGACT IZFFM*MOD
  971. IZHR=KZHR(1)
  972. IZH2=KZHR(2)
  973. SEGACT IZHR*MOD
  974. NES=GR(/1)
  975. NPG=GR(/3)
  976. IZF1=KTP(1)
  977. NP = IPT1.NUM(/1)
  978. NBEL=IPT1.NUM(/2)
  979. NPT =MPOVA1.VPOCHA(/1)
  980.  
  981. KSDIM=IDIM
  982. SEGINI SOURCE
  983.  
  984. IF(IBOU.EQ.0)THEN
  985. C 1ere CONFIG ROG seul
  986. IF(IKG.EQ.0)THEN
  987. C SOURCE CENTRE
  988. DO K=1,NBEL
  989. DO N=1,IDIM
  990. DO LG=1,NPG
  991. SO(K,LG,N)=IZTGG2.VPOCHA(K+NUTOEL,N)
  992. ENDDO
  993. ENDDO
  994. ENDDO
  995.  
  996. ELSEIF(IKG.EQ.1)THEN
  997. C SOURCE CONSTANTE
  998. DO K=1,NBEL
  999. DO N=1,IDIM
  1000. DO LG=1,NPG
  1001. SO(K,LG,N)=IZTGG2.VPOCHA(1,N)
  1002. ENDDO
  1003. ENDDO
  1004. ENDDO
  1005.  
  1006. ELSEIF(IKG.EQ.4)THEN
  1007. C SOURCE SOMMET
  1008. DO K=1,NBEL
  1009. DO N=1,IDIM
  1010. DO LG=1,NPG
  1011. U0=0.D0
  1012. DO 3547 I=1,NP
  1013. I1=IPADS.LECT(IPT1.NUM(I,K))
  1014. U0=U0+FN(I,LG)*IZTGG2.VPOCHA(I1,N)
  1015. 3547 CONTINUE
  1016. SO(K,LG,N)=U0
  1017. ENDDO
  1018. ENDDO
  1019. ENDDO
  1020.  
  1021. ENDIF
  1022.  
  1023. ELSEIF(IBOU.EQ.1)THEN
  1024. C 2eme CONFIG ROG*(T-T0) Boussinesq
  1025.  
  1026. IF(IKG.EQ.0)THEN
  1027. C SOURCE CENTRE
  1028. DO K=1,NBEL
  1029. DO N=1,IDIM
  1030. DO LG=1,NPG
  1031. U0=0.D0
  1032. DO 4548 I=1,NP
  1033. I1=IPADS.LECT(IPT1.NUM(I,K))
  1034. I0=(1-I1)*IK4 + I1
  1035. U0=U0-FN(I,LG)*(IZTGG2.VPOCHA(K+NUTOEL,N)*
  1036. & (IZTGG3.VPOCHA(I1,1)-IZTGG4.VPOCHA(I0,1)))
  1037. 4548 CONTINUE
  1038. SO(K,LG,N)=U0
  1039. ENDDO
  1040. ENDDO
  1041. ENDDO
  1042.  
  1043. ELSEIF(IKG.EQ.1)THEN
  1044. C SOURCE CONSTANTE
  1045.  
  1046. DO K=1,NBEL
  1047. DO N=1,IDIM
  1048. DO LG=1,NPG
  1049. U0=0.D0
  1050. DO 4549 I=1,NP
  1051. I1=IPADF.LECT(IPT1.NUM(I,K))
  1052. I0=(1-I1)*IK4 + I1
  1053. U0=U0-FN(I,LG)*(IZTGG2.VPOCHA(1,N)
  1054. &*(IZTGG3.VPOCHA(I1,1)-IZTGG4.VPOCHA(I0,1)))
  1055. 4549 CONTINUE
  1056. SO(K,LG,N)=U0
  1057. ENDDO
  1058. ENDDO
  1059. ENDDO
  1060.  
  1061. ELSEIF(IKG.EQ.4)THEN
  1062. C SOURCE SOMMET
  1063. DO K=1,NBEL
  1064. DO N=1,IDIM
  1065. DO LG=1,NPG
  1066. U0=0.D0
  1067. DO 4547 I=1,NP
  1068. I1=IPADS.LECT(IPT1.NUM(I,K))
  1069. I0=(1-I1)*IK4 + I1
  1070. U0=U0-FN(I,LG)*(IZTGG2.VPOCHA(I1,N)*
  1071. & (IZTGG3.VPOCHA(I1,1)-IZTGG4.VPOCHA(I0,1)))
  1072. 4547 CONTINUE
  1073. SO(K,LG,N)=U0
  1074. ENDDO
  1075. ENDDO
  1076. ENDDO
  1077.  
  1078. ENDIF
  1079. ENDIF
  1080.  
  1081. C write(6,*)' Appel a YSOUR'
  1082. CALL YSOUR(FN,GR,PG,XYZ,HR,PGSQ,RPG,
  1083. & NES,IDIM,NP,NPG,IAXI,IPT1.NUM,SO,IDIM,IPADS.LECT,
  1084. & NBEL,NUTOEL,XCOOR,MPOVA1.VPOCHA,NPT)
  1085. C write(6,*)' Retour YSOUR'
  1086.  
  1087. SEGSUP IZFFM,IZF1,IZHR,IZH2
  1088. SEGSUP SOURCE
  1089.  
  1090. ELSE
  1091.  
  1092. IF(INEFMD.EQ.3)THEN
  1093. C Quadratiques
  1094. IF(KPRE.EQ.2)NOM0=NOMS(IPT1.ITYPEL)//'PRP0'
  1095. IF(KPRE.EQ.3)NOM0=NOMS(IPT1.ITYPEL)//'PRP0'
  1096. IF(KPRE.EQ.4)NOM0=NOMS(IPT1.ITYPEL)//'PRP1'
  1097. ELSEIF(INEFMD.EQ.2)THEN
  1098. C Macro
  1099. IF(KPRE.EQ.2)NOM0=NOMS(IPT1.ITYPEL)//' '
  1100. IF(KPRE.EQ.3)NOM0=NOMS(IPT1.ITYPEL)//'MCP0'
  1101. IF(KPRE.EQ.4)NOM0=NOMS(IPT1.ITYPEL)//'MCP1'
  1102. ELSE
  1103. IF(KPRE.EQ.2)NOM0=NOMS(IPT1.ITYPEL)//' '
  1104. ENDIF
  1105.  
  1106. CALL KALPBG(NOM0,'FONFORM ',IZFFM)
  1107.  
  1108.  
  1109. SEGACT IZFFM*MOD
  1110. IZHR=KZHR(1)
  1111. IZH2=KZHR(2)
  1112. SEGACT IZHR*MOD
  1113. NES=GR(/1)
  1114. NPG=GR(/3)
  1115. IZF1=KTP(1)
  1116. SEGACT IZF1*MOD
  1117. MP1=IZF1.FN(/1)
  1118. NP = IPT1.NUM(/1)
  1119. MP = NP
  1120. NBEL=IPT1.NUM(/2)
  1121. NELG=IZTGG2.VPOCHA(/1)
  1122. NPT =MPOVA1.VPOCHA(/1)
  1123.  
  1124. SEGACT MELEP1
  1125. C write(6,*)' Appel a XSOUR '
  1126. CALL XSOUR(FN,IZF1.FN,GR,PG,XYZ,HR,PGSQ,RPG,
  1127. & NES,IDIM,NP,MP1,NPG,IAXI,IPT1.NUM,IKAS,KPRE,
  1128. & IZTGG2.VPOCHA,IKG,NELG,IPADQ.LECT,MELEP1.NUM,
  1129. & IZTGG3.VPOCHA,IK3,IZTGG4.VPOCHA,IK4,IPADS.LECT,
  1130. & NBEL,NUTOEL,XCOOR,MPOVA1.VPOCHA,NPT)
  1131. C write(6,*)' Retour XSOUR '
  1132.  
  1133. SEGSUP IZFFM,IZF1,IZHR,IZH2
  1134.  
  1135. ENDIF
  1136.  
  1137. NUTOEL=NUTOEL+NBEL
  1138. 1102 CONTINUE
  1139.  
  1140. ENDIF
  1141.  
  1142. IF(KIMPL.EQ.2.OR.KIMPL.EQ.0)THEN
  1143.  
  1144. SEGACT MELEME
  1145. LL=0
  1146.  
  1147. DO 1533 L=1,MAX(1,LISOUS(/1))
  1148. IPT1=MELEME
  1149. IF(LISOUS(/1).NE.0)IPT1=LISOUS(L)
  1150. SEGACT IPT1
  1151.  
  1152. NP=IPT1.NUM(/1)
  1153. NBELEM=IPT1.NUM(/2)
  1154.  
  1155. NNN=MOD(NBELEM,LRV)
  1156. IF(NNN.EQ.0) NPACK=NBELEM/LRV
  1157. IF(NNN.NE.0) NPACK=1+(NBELEM-NNN)/LRV
  1158.  
  1159. DO 1534 KPACK=1,NPACK
  1160.  
  1161. C --- POUR CHAQUE PAQUET DE LRV ELEMENTS ou moins
  1162. KDEB=1+(KPACK-1)*LRV
  1163. KFIN=MIN(NBELEM,KDEB+LRV-1)
  1164. NBEL=KFIN-KDEB+1
  1165. LL=LL+1
  1166.  
  1167. DO 2 N=1,NINKO
  1168. IPMS=LIZAFS(LL,N)
  1169. SEGACT IPMS
  1170. DO 12 K=1,NBEL
  1171. KI=KDEB+K-1
  1172.  
  1173. DO 13 J=1,NP
  1174. UU=0.D0
  1175. IU=IPADS.LECT(IPT1.NUM(J,KI))
  1176. DO 14 I=1,NP
  1177. IK=IPADI.LECT(IPT1.NUM(I,KI))
  1178. UU=UU+IPMS.AM(K,I,J)*IZTU1.VPOCHA(IK,N)
  1179. 14 CONTINUE
  1180. MPOVA1.VPOCHA(IU,N)=MPOVA1.VPOCHA(IU,N)+UU
  1181. 13 CONTINUE
  1182. 12 CONTINUE
  1183. SEGDES IPMS
  1184. 2 CONTINUE
  1185.  
  1186. 1534 CONTINUE
  1187. SEGDES IPT1
  1188. 1533 CONTINUE
  1189. ENDIF
  1190.  
  1191. SEGDES IPM1,IPM2,IPM3
  1192. IPS=IPS1
  1193. SEGSUP IPS1
  1194. IF(IPS2.NE.IPS)SEGSUP IPS2
  1195. IF(IPS3.NE.IPS)SEGSUP IPS3
  1196.  
  1197. CALL ECROBJ('MATRIK',MATRIK)
  1198.  
  1199. IF(KIMPL.EQ.2.OR.KIMPL.EQ.0.OR.IARG.GT.1)THEN
  1200.  
  1201. CALL ECROBJ('CHPOINT',MCHPO1)
  1202. ELSE
  1203. NAT=2
  1204. NSOUPO=0
  1205. SEGINI MCHPOI
  1206. JATTRI(1)=2
  1207. CALL ECROBJ('CHPOINT',MCHPOI)
  1208. ENDIF
  1209.  
  1210. SEGDES IMATRI
  1211. SEGDES MELEME,MATRIK
  1212. IF(IK1.EQ.0)THEN
  1213. SEGDES VISCO
  1214. ENDIF
  1215.  
  1216. C*************************************************************************
  1217. ELSE IF(KFORM.EQ.2)THEN
  1218. C CAS FORMULATION VF
  1219. WRITE(6,*)' FORMULATION VF NON DISPONIBLE '
  1220. ENDIF
  1221. C*************************************************************************
  1222.  
  1223. IPDI=IPADI
  1224. SEGSUP IPADI
  1225. IF(IPADS.NE.IPDI)SEGSUP IPADS
  1226. IPAS=1
  1227. C write(6,*)' FIN de NS'
  1228. RETURN
  1229. 90 CONTINUE
  1230. WRITE(6,*)' Interuption anormale de NS'
  1231. C Option %m1:8 incompatible avec les données
  1232. MOTERR( 1: 8) = ' EF '
  1233. CALL ERREUR(803)
  1234. RETURN
  1235. 1002 FORMAT(10(1X,1PE11.4))
  1236. 1008 FORMAT( 8(1X,1PE11.4))
  1237. 1001 FORMAT(20(1X,I5))
  1238. END
  1239.  
  1240.  
  1241.  
  1242.  
  1243.  
  1244.  
  1245.  
  1246.  
  1247.  
  1248.  
  1249.  
  1250.  
  1251.  
  1252.  
  1253.  
  1254.  
  1255.  
  1256.  
  1257.  
  1258.  
  1259.  
  1260.  
  1261.  
  1262.  
  1263.  
  1264.  
  1265.  
  1266.  
  1267.  
  1268.  
  1269.  
  1270.  
  1271.  
  1272.  
  1273.  
  1274.  
  1275.  
  1276.  
  1277.  
  1278.  

© Cast3M 2003 - Tous droits réservés.
Mentions légales