Télécharger injection.dgibi

Retour à la liste

Numérotation des lignes :

  1. * fichier : injection.dgibi
  2. ************************************************************************
  3. ************************************************************************
  4. ************************************************************************
  5. * NOM : Injection Air/Air (Air chaud dans Air froid sur 6
  6. * secondes) dans une cavité carrée 2D plan
  7. *
  8. * DESCRIPTION : Cas test du modèle asymptotique à bas nombre de Mach
  9. * pour un écoulement laminaire, confiné, monoespèce mais
  10. * faisant intervenir la thermique
  11. *
  12. *
  13. *
  14. * LANGAGE : GIBIANE-CAST3M
  15. * AUTEUR : Etienne STUDER (CEA/DEN/DM2S/SFME/LTMF)
  16. * mél : mp1@semt2.smts.cea.fr
  17. **********************************************************************
  18. * VERSION : v1, 15/12/2005, version initiale
  19. * HISTORIQUE : v1, ??/??/2003, création
  20. * HISTORIQUE :
  21. * HISTORIQUE :
  22. ************************************************************************
  23. * Toutes les procédures sont incluses dans le fichier
  24. *
  25. *
  26. *
  27. ************************************************************************
  28. interact= FAUX ;
  29. 'OPTION' 'DIME' 2 'ELEM' QUA8 'ISOV' 'SULI' ;
  30. 'SI' ('NON' interact) ;
  31. 'OPTION' 'TRAC' 'PS' ;
  32. 'SINON' ;
  33. 'OPTION' 'TRAC' 'X' ;
  34. 'FINSI' ;
  35.  
  36. ************************************************************************
  37. * PROCEDURES LOCALES
  38. ************************************************************************
  39.  
  40. * Calcul d'un terme source VECTEUR Explicite
  41. 'DEBPROC' TSOUR ;
  42. 'ARGUMENT' RVX * 'TABLE' ;
  43. RV = RVX . 'EQEX' ;
  44. LSTINC = RVX.'LISTINCO' ;
  45. NOMV = 'EXTRAIRE' LSTINC 1 ;
  46. ROG = RVX.'ARG1' ;
  47. RGX = 'EXCO' (RV.'INCO'. ROG) 'UX' ;
  48. RGY = 'EXCO' (RV.'INCO'. ROG) 'UY' ;
  49. RGX = RGX '*' ('DOMA' (RVX.'DOMZ') 'XXDIAGSI') ;
  50. RGY = RGY '*' ('DOMA' (RVX.'DOMZ') 'XXDIAGSI') ;
  51. RGX = 'NOMC' RGX '1UN' 'NATURE' 'DISCRET' ;
  52. RGY = 'NOMC' RGY '2UN' 'NATURE' 'DISCRET' ;
  53. ROG = RGX 'ET' RGY ;
  54. AS AF = 'KOPS' 'MATRIK' ;
  55. 'FINPROC' ROG AF ;
  56.  
  57. 'DEBP' EXEC ;
  58. *
  59. ************************ K E X IMPL ************************************
  60. *
  61. * Résolution implicite d'un problème de mécanique des fluides décrit
  62. * par l'intermédiaire d'une table structurée par l'opérateur EQEX.
  63. *
  64. * E/ rv : Table décrivant le problème à résoudre
  65. * /S rv : A l'indice INCO, table contenant la solution au temps final
  66. *
  67. * L'utilisateur peut par l'intermédiaire de procédures personnelles
  68. * stocker les résultats aux temps intermédiaires, faire des tracés, ...
  69. *
  70. ************************ K E X E C ************************************
  71. *
  72. 'ARGUMENT' rv*'TABLE ' ;
  73. *
  74. *========================
  75. * TESTS et INITIALISATION
  76. *========================
  77. *
  78. *
  79. * On n'accepte que les formulations en modèle NAVIER_STOKES
  80. 'SI' ('EGA' (rv . 'NAVISTOK') 0) ;
  81. 'MESS' 'Pas de modèle NAVIER_STOKES' ;
  82. 'ERRE' 5 ;
  83. 'FINSI' ;
  84. *
  85. dimope = 'DIME' (rv . 'LISTOPER') ;
  86. *
  87. * Coefficient de relaxation mis à 1 si non initialisé
  88. 'SI' ('NON' ('EXISTE' rv 'OMEGA')) ;
  89. 'MESS' '** WARNING ** OMEGA NON défini --> OMEGA=1' ;
  90. omeg = 1.D0 ;
  91. 'SINON' ;
  92. omeg = rv . 'OMEGA' ;
  93. 'FINSI' ;
  94. *
  95. * Création de la table pour historique
  96. 'SI' ('NON' ('EXISTE' rv 'HIST')) ;
  97. rv . 'HIST' = 'TABLE' ;
  98. 'FINSI' ;
  99. *
  100. * IMPKRES : Niveau d'impression pour KRES
  101. * IMPTCRR : Fréquence d'impression pour TCRR (affichage résidu)
  102. IMPKRES = 0 ;
  103. IMPTCRR = RV . 'IMPR' ;
  104. *
  105. *========================
  106. * RESOLUTION suivant EQEX
  107. *========================
  108. *
  109. *
  110. * 1
  111. * Boucle en temps
  112. *----------------
  113. ITMA = rv . 'ITMA' ;
  114. 'SI' ('<EG' ITMA 1) ; ITMA = 1 ; 'FINSI' ;
  115. 'REPETER' bloc1 ITMA ;
  116. *
  117. * 2
  118. * Boucle de point fixe interne à un pas de temps
  119. *-----------------------------------------------
  120. 'REPETER' bloci (rv . 'NITER') ;
  121. * st mat : Second membre et matrice des opérateurs DFDT
  122. * sf mau : Second membre er matrice des opérateurs sauf DFDT
  123. st mat = 'KOPS' 'MATRIK' ;
  124. sf mau = 'KOPS' 'MATRIK' ;
  125. *
  126. * 3
  127. * Boucle sur les opérateurs de EQEX
  128. *----------------------------------
  129. 'REPETER' bloc2 dimope ;
  130. nomper = 'EXTRAIRE' &bloc2 (rv . 'LISTOPER') ;
  131. notable = 'MOT' ('TEXTE' ('CHAINE' &bloc2 nomper)) ;
  132. 'SI' ('EGA' nomper 'DFDT ') ;
  133. msi mai = ('TEXTE' nomper) (rv . notable) ;
  134. mat = mat 'ET' mai ;
  135. st = st 'ET' msi ;
  136. 'SINO' ;
  137. msi mai = ('TEXTE' nomper) (rv . notable) ;
  138. mau = mau 'ET' mai ;
  139. sf = sf 'ET' msi ;
  140. 'FINSI' ;
  141. 'FIN' bloc2 ;
  142. *--------------------------------------
  143. * Fin Boucle sur les opérateurs de EQEX
  144. * 3
  145. *
  146. s2 = sf 'ET' st ;
  147. ma1 = mau 'ET' mat ;
  148. * ES preconditionnement MATRICE
  149. 'SI' ('EXISTE' rv 'ESPR') ;
  150. ma1 = ma1 'ET' (rv.'ESPR') ;
  151. 'FINSI' ;
  152. *
  153. 'SI' ('EXISTE' rv 'CLIM') ;
  154. s1 = rv . 'CLIM' ;
  155. 'SINON' ;
  156. s1 = ' ' ;
  157. 'FINSI' ;
  158. rv . 'S2' = s2 ;
  159.  
  160. rv . 'METHINV' . 'MATASS' = ma1 ;
  161. rv . 'METHINV' . 'MAPREC' = ma1 ;
  162. res = 'KRES' ma1 'TYPI' (rv . 'METHINV')
  163. 'CLIM' s1
  164. 'SMBR' s2
  165. 'IMPR' IMPKRES ;
  166. eps = 'TCRR' res omeg (rv . 'INCO') 'IMPR' IMPTCRR ;
  167. 'MENAGE' ;
  168. 'FIN' bloci ;
  169. *---------------------------------------------------
  170. * Fin Boucle de point fixe interne à un pas de temps
  171. * 2
  172. *
  173. irt = 0 ;
  174. 'SI' ('EGA' (rv . 'ITMA') 0) ;
  175. irt = 'TCNM' rv 'NOUP';
  176. 'SINON' ;
  177. irt = 'TCNM' rv ;
  178. 'FINSI' ;
  179. 'MENAGE' ;
  180. 'SI' ('EGA' irt 1) ;
  181. 'MESSAGE' ' Temps final atteint : ' (rv . 'PASDETPS' . 'TPS') ;
  182. 'QUITTER' bloc1 ;
  183. 'FINSI' ;
  184. 'FIN' bloc1 ;
  185. *--------------------
  186. * Fin Boucle en temps
  187. * 1
  188. *
  189. ************************ K E X IMPL ************************************
  190. 'FINPROC' ;
  191.  
  192. *************************************************************************
  193. * DEFINITION DU PROBLEME
  194. *
  195. *************************************************************************
  196.  
  197. * le raffinement: nombre de points dans une dimension
  198. * en 1D la verticale et en 2D on prend le même nombre
  199. * dans les deux directions
  200. RAF = 20 ;
  201.  
  202. * Choix de la méthode d'inversion
  203. * TYPI = 1 : méthode directe
  204. * TYPI = 3 : méthode itérative
  205. TYPI = 1 ;
  206.  
  207. * Le cas test peut être traité en 1D ou 2D plan
  208. UND = VRAI ;
  209. DEUXD = FAUX ;
  210.  
  211. * Dimensions géométriques
  212. * L : largeur de l'injection au bas de la cavité (centrée)
  213. * H : hauteur de la cavité
  214. * LL : largeur de la cavité (=L en 1D)
  215. L = 0.2 ;
  216. H = 7.0 ;
  217.  
  218. 'SI' (UND) ;
  219. LL = L ;
  220. 'SINON' ;
  221. 'SI' (DEUXD) ;
  222. LL = 3.0 ;
  223. 'SINON' ;
  224. 'MESSAGE' '==> Cas imprévu !!!!!!!!' ;
  225. 'ERREUR' 5 ;
  226. 'FINSI' ;
  227. 'FINSI' ;
  228.  
  229. * Caractéristiques du gaz
  230. * GAMA : coeffcient isentropique
  231. * RAIR : constante du gaz (R/M) avec R constante des gaz et M la masse
  232. * molaire
  233. * PRANDTL : nombre de Prandtl
  234. GAMA = 1.4 ;
  235. RAIR = 287.0 ;
  236. PRANDTL = 0.71 ;
  237.  
  238. * Conditions initiales
  239. * P0 : pression (Pa)
  240. * T0 : température (K)
  241. * R0 : masse volumique (kg/m3)
  242. P0 = 1.E5 ;
  243. T0 = 300.0 ;
  244. R0 = P0 '/' (RAIR '*' T0) ;
  245.  
  246. * conditions à l'injection
  247. * MPH : débit surfacique kg/m²/s
  248. * TH : température (K)
  249. * V0 : vitesse à l'injection (m/s)
  250. MPH = 1.0 ;
  251. TH = 600.0 ;
  252. V0 = MPH '*' RAIR '*' TH '/' P0 ;
  253.  
  254. * Choix du cas de calcul par l'intermédiaire de booléens
  255. * C1 : Reynolds = 40 sans prise en compte de la gravité
  256. * C2 : idem C1 mais avec ajout de la flottabilité
  257. * C3 : idem C2 mais avec Reynolds = 800
  258. C1 = VRAI ;
  259. C2 = FAUX ;
  260. C3 = FAUX ;
  261. 'SI' (C1) ;
  262. REYN = 40.0 ;
  263. GRAV = 0.0 ;
  264. 'FINSI' ;
  265. 'SI' (C2) ;
  266. REYN = 40.0 ;
  267. GRAV = 9.81 ;
  268. 'FINSI' ;
  269. 'SI' (C3) ;
  270. REYN = 800.0 ;
  271. GRAV = 9.81 ;
  272. 'FINSI' ;
  273.  
  274. * Calcul des propriétés physiques
  275. * MU : viscosité kg/m/s
  276. * LAMBDA : conductivité thermique W/m/K
  277. * CPAIR et CVAIR : chaleur spécifique J/kg/K
  278. MU = MPH '*' L '/' REYN ;
  279. LAMBDA = MU '*' GAMA '*' RAIR '/' PRANDTL '/' (GAMA '-' 1.0) ;
  280. CPAIR = GAMA '*' RAIR '/' (GAMA '-' 1.0) ;
  281. CVAIR = CPAIR '-' RAIR ;
  282.  
  283. ********************************************************************
  284. * LE MAILLAGE
  285. ********************************************************************
  286. LS2 = L '/' 2.0 ;
  287. LLS2 = LL '/' 2.0 ;
  288. LLS4 = LL '/' 4.0 ;
  289. HS4 = H '/' 4.0 ;
  290. HS2 = H '/' 2.0 ;
  291. 3HS4 = H '*' 3.0 '/' 4.0 ;
  292.  
  293.  
  294. 'DENS' (H '/' ('FLOTTANT' RAF)) ;
  295.  
  296. 'SI' (UND) ;
  297. P1 = (LS2 '*' (-1.0)) 0.0 ;
  298. P2 = 0.0 0.0 ;
  299. P3 = (LS2) 0.0 ;
  300. P1P2 = 'DROIT' 1 P1 P2 ;
  301. P2P3 = 'DROIT' 1 P2 P3 ;
  302. P4 = 0.0 H ;
  303. P5 = (LS2 '*' (-1.0)) H ;
  304. P6 = LS2 H ;
  305. P2P4 = 'DROIT' P2 P4 ;
  306. P1P5 = 'DROIT' P1 P5 ;
  307. P3P6 = 'DROIT' P3 P6 ;
  308. P6P4 = 'DROIT' 1 P6 P4 ;
  309. P4P5 = 'DROIT' 1 P4 P5 ;
  310. BAS = P1P2 'ET' P2P3 ;
  311. HAUT = P6P4 'ET' P4P5 ;
  312. LMIL = P2P4 ;
  313. CDRO = P3P6 ;
  314. CGAU = 'INVERSE' P1P5 ;
  315. VTF = BAS P3P6 HAUT ('INVERSE' P1P5) 'DALLER' 'PLAN' ;
  316. 'FINSI' ;
  317. 'SI' (DEUXD) ;
  318. P1 = (LLS2 '*' (-1.0)) 0.0 ;
  319. P2 = 0. 0. ;
  320. P1A = (LLS4 '*' (-1.0)) 0.0 ;
  321. P1B = (LS2 '*' (-1.0)) 0.0 ;
  322. N1 = 'ENTIER' (LS2 '/' LL '*' RAF) ;
  323. N2 = 'ENTIER' ((LLS4 '-' LS2) '/' LL '*' RAF) ;
  324. N3 = ('ENTIER' (RAF '/' 2.0)) '-' N1 '-' N2 ;
  325. P1P1A = 'DROIT' N3 P1 P1A ;
  326. P1AP1B = 'DROIT' N2 P1A P1B ;
  327. P1BP2 = 'DROIT' N1 P1B P2 ;
  328. P3 = LLS2 0.0 ;
  329. P3B = LLS4 0.0 ;
  330. P3A = LS2 0.0 ;
  331. P2P3A = 'DROIT' N1 P2 P3A ;
  332. P3AP3B = 'DROIT' N2 P3A P3B ;
  333. P3BP3 = 'DROIT' N3 P3B P3 ;
  334. INJE = P1BP2 'ET' P2P3A ;
  335. BAS = P1P1A 'ET' P1AP1B 'ET' P3AP3B 'ET' P3BP3 ;
  336. NY = 'ENTIER' (RAF '/' 4.0) ;
  337. P4 = 0.0 H ;
  338. HAUT = 'INVERSE' ('PLUS' (BAS 'ET' INJE) P4) ;
  339. P5 = 'POIN' HAUT 'PROC' ('PLUS' P1 P4 ) ;
  340. P6 = 'POIN' HAUT 'PROC' ('PLUS' P3 P4) ;
  341. CDRO = 'DROIT' ('ENTIER' RAF) P3 P6 ;
  342. CGAU = 'DROIT' ('ENTIER' RAF) P5 P1 ;
  343. VTF = (BAS 'ET' INJE) CDRO HAUT CGAU 'DALLER' 'PLAN' ;
  344. 'FINSI' ;
  345.  
  346. * Choix des éléments Vitesse/Pression
  347. DISCR = 'QUAF' ;
  348. KPRES = 'CENTREP1' ;
  349. * Choix du décentrement des termes convectifs
  350. KSUPG = 'CENTREE' ;
  351. * Choix de la condensation ou non de la matrice masse
  352. *KMASS = 'EF' ;
  353. KMASS = 'EFM1' ;
  354. * Choix de l'ordre du schéma en temps
  355. * BDF1 : Euler
  356. * BDF2 : schéma à 3 pas (ordre 2)
  357. KBDF = 'BDF1' ;
  358. *KBDF = 'BDF2' ;
  359.  
  360. * positionnement des constantes pour le schéma en temps BDF1 (EULER) ou
  361. * BDF2
  362. 'SI' ('EGA' KBDF 'BDF1') ;
  363. AN = 1.0 ;
  364. ANM1 = (-1.0) ;
  365. ANM2 = 0.0 ;
  366. ADT = 1.0 ;
  367. 'SINON' ;
  368. 'SI' ('EGA' KBDF 'BDF2') ;
  369. AN = 3.0 ;
  370. ANM1 = (-4.0) ;
  371. ANM2 = 1.0 ;
  372. ADT = 2.0 ;
  373. 'SINON' ;
  374. 'MESSAGE' '==> ERREUR indice KBDF' ;
  375. 'ERREUR' 5 ;
  376. 'FINSI' ;
  377. 'FINSI' ;
  378. * Precision maillage
  379. EPS0 = 1.D-7 ;
  380.  
  381. * Definition des objets modèles
  382.  
  383. MVTF = 'CHANGER' VTF 'QUAF' ;
  384. $VTF = 'MODELISER' MVTF 'NAVIER_STOKES' DISCR ;
  385.  
  386. * Traitement de la geometrie
  387. *-----------------------------------------------------------
  388. * Creation de l'enveloppe du maillage FLUIDE
  389. GEO = TABLE ;
  390. GEO.'MVTF' = MVTF ;
  391. GEO.'$VTF' = $VTF ;
  392. GEO.'EPSI' = EPS0 ;
  393.  
  394. GEO.'MHAUT' = HAUT ;
  395. GEO.'$HAUT' = 'MODE' GEO.'MHAUT' 'NAVIER_STOKES' DISCR ;
  396. TOTO = 'DOMA' GEO.'$HAUT' 'CENTRE' ;
  397.  
  398. GEO.'MBAS' = BAS ;
  399. GEO.'$BAS' = 'MODE' GEO.'MBAS' 'NAVIER_STOKES' DISCR ;
  400. TOTO1 = 'DOMA' GEO.'$BAS' 'CENTRE' ;
  401.  
  402. 'SI' (DEUXD) ;
  403. GEO.'MINJE' = INJE ;
  404. GEO.'$INJE' = 'MODE' GEO.'MINJE' 'NAVIER_STOKES' DISCR ;
  405. TOTO2 = 'DOMA' GEO.'$INJE' 'CENTRE' ;
  406. 'FINSI' ;
  407.  
  408. GEO.'MCDRO' = CDRO ;
  409. GEO.'$CDRO' = 'MODE' GEO.'MCDRO' 'NAVIER_STOKES' DISCR ;
  410. TUTU1 = 'DOMA' GEO.'$CDRO' 'CENTRE' ;
  411.  
  412. GEO.'MCGAU' = CGAU ;
  413. GEO.'$CGAU' = 'MODE' GEO.'MCGAU' 'NAVIER_STOKES' DISCR ;
  414. TUTU2 = 'DOMA' GEO.'$CGAU' 'CENTRE' ;
  415.  
  416. 'SI' (UND) ;
  417. GEO.'MMENVF' = CGAU 'ET' BAS 'ET' CDRO 'ET' HAUT ;
  418. 'FINSI' ;
  419. 'SI' (DEUXD) ;
  420. GEO.'MMENVF' = CGAU 'ET' BAS 'ET' CDRO 'ET' HAUT 'ET' INJE ;
  421. 'FINSI' ;
  422.  
  423. GEO.'$MENVF' = 'MODE' GEO.'MMENVF' 'NAVIER_STOKES' DISCR ;
  424. TUTU3 = 'DOMA' GEO.'$MENVF' 'CENTRE' ;
  425.  
  426. DOMF = 'DOMA' GEO.'$VTF' 'FACE' ;
  427. 'SI' (UND) ;
  428. 'ELIMINATION' DOMF (TOTO 'ET' TOTO1
  429. 'ET' TUTU1 'ET' TUTU2 'ET' TUTU3) EPS0 ;
  430. 'FINSI' ;
  431. 'SI' (DEUXD) ;
  432. 'ELIMINATION' DOMF (TOTO 'ET' TOTO1 'ET' TOTO2
  433. 'ET' TUTU1 'ET' TUTU2 'ET' TUTU3) EPS0 ;
  434. 'FINSI' ;
  435. *
  436. * Récupération des maillages a partir des MODELES
  437.  
  438. GEO.'BAS' = 'DOMA' GEO.'$BAS' 'MAILLAGE' ;
  439. GEO.'HAUT' = 'DOMA' GEO.'$HAUT' 'MAILLAGE' ;
  440. GEO.'VTF' = 'DOMA' GEO.'$VTF' 'MAILLAGE' ;
  441. GEO.'CGAU' = 'DOMA' GEO.'$CGAU' 'MAILLAGE' ;
  442. GEO.'CDRO' = 'DOMA' GEO.'$CDRO' 'MAILLAGE' ;
  443. 'SI' (DEUXD) ;
  444. GEO.'INJE' = 'DOMA' GEO.'$INJE' 'MAILLAGE' ;
  445. 'FINSI' ;
  446. GEO.'MENVF' = 'DOMA' GEO.'$MENVF' 'MAILLAGE' ;
  447.  
  448. * Maillage pour la condition aux limites en pression
  449. * Il faut un point à l'intérieur
  450.  
  451. LVTF = 'CHAN' 'LINEAIRE' GEO.'VTF' ;
  452. PI0 = LVTF 'POIN' 'PROC' ( 0. H ) ;
  453. MELTI = (GEO.'MVTF') 'ELEM' 'CONTENANT' PI0 ;
  454. $ELTI = 'MODE' MELTI 'NAVIER_STOKES' DISCR ;
  455. MTPI = 'DOMA' $ELTI KPRES ;
  456. MTPI = 'ELEM' ('LECT' 1) MTPI ;
  457. 'ELIM' (MTPI 'ET' ('DOMA' GEO.'$VTF' KPRES )) EPS0 ;
  458.  
  459. * Maillage pour ouvrir en haut pour faire fuir et equilibrer
  460. * (test de la pertinence de la contrainte de divergence
  461. GEO.'HAUT1' = 'CHAN' GEO.'HAUT' 'POI1' ;
  462. 'SI' (UND) ;
  463. GEO.'HAUT1' = 'DIFF' GEO.'HAUT1' ('ELEM' GEO.'HAUT1' 2) ;
  464. 'FINSI' ;
  465. 'SI' (DEUXD) ;
  466. GEO.'HAUT1' = 'DIFF' GEO.'HAUT1' ('ELEM' GEO.'HAUT1' ('ENTIER' (RAF
  467. '/' 2.0))) ;
  468. 'FINSI' ;
  469.  
  470. * Fonction de forme en 2D pour l'injection: profil parabolique
  471. *
  472. 'SI' (DEUXD) ;
  473. XX = 'COORDONNEE' 1 GEO.'INJE' ;
  474. FX = XX '*' XX '*' (-1.0) '+' (L '**' 2.0 '/' 4.0) '*' (6.0 '/' (L
  475. '*' L)) ;
  476. 'FINSI' ;
  477.  
  478. * MTPI = 'ELEM' ('LECT' 1) ('DOMA' GEO.'$VTF' KPRES) ;
  479.  
  480. * Definition du modèle physique table RV
  481. *
  482. * - Definition de constantes intervenant dans les equations
  483. * FCPRECI = Frequence dans le methode iterative
  484. * FCPRECT= idem
  485.  
  486. FCPRECI = 1 ;
  487. FCPRECT = 1 ;
  488.  
  489. *---------------------------------------------------------------
  490. *- EQUATIONS DE LA QDM
  491. * Formulation non conservative
  492. * Decentrement KSUPG
  493. *
  494. RV = 'EQEX' 'NITER' 1 'OMEGA' 1.0 'ITMA' 0
  495. 'OPTI' 'EF' 'IMPL' KSUPG
  496. 'ZONE' (GEO.'$VTF') 'OPER' 'LAPN' 'MU' 'INCO' 'UN'
  497. 'OPTI' 'EF' 'IMPL' KSUPG
  498. 'ZONE' (GEO.'$VTF') 'OPER' 'KONV' 'RHO' 'UN' 'MU' 'INCO' 'UN'
  499. ;
  500. RV = 'EQEX' RV
  501. 'OPTI' 'EF' 'CENTREE' 'IMPL'
  502. 'ZONE' (GEO.'$VTF') 'OPER' TSOUR 'ROG' 'INCO' 'UN'
  503. ;
  504. 'SI' ('EGA' KBDF 'BDF2') ;
  505. RV = 'EQEX' RV
  506. 'OPTI' KMASS 'CENTREE' 'BDF2'
  507. 'ZONE' (GEO.'$VTF') 'OPER' 'DFDT' 'RHO' 'UNM' 'UNMM' 'DT'
  508. 'INCO' 'UN'
  509. ;
  510. 'SINON' ;
  511. RV = 'EQEX' RV
  512. 'OPTI' KMASS 'CENTREE'
  513. 'ZONE' (GEO.'$VTF') 'OPER' 'DFDT' 'RHO' 'UNM' 'DT' 'INCO' 'UN'
  514. ;
  515. 'FINSI' ;
  516. * Conditions aux limites QDM
  517. 'SI' (UND) ;
  518. RV = 'EQEX' RV 'CLIM'
  519. 'UN' 'UIMP' GEO.'CDRO' 0.
  520. 'UN' 'UIMP' GEO.'CGAU' 0.
  521. 'UN' 'UIMP' GEO.'HAUT' 0.
  522. 'UN' 'VIMP' GEO.'HAUT1' 0.
  523. 'UN' 'UIMP' GEO.'BAS' 0.
  524. 'UN' 'VIMP' GEO.'BAS' V0
  525. ;
  526. 'FINSI' ;
  527. 'SI' (DEUXD) ;
  528. RV = 'EQEX' RV 'CLIM'
  529. 'UN' 'UIMP' GEO.'CDRO' 0.
  530. 'UN' 'VIMP' GEO.'CDRO' 0.
  531. 'UN' 'UIMP' GEO.'CGAU' 0.
  532. 'UN' 'VIMP' GEO.'CGAU' 0.
  533. 'UN' 'UIMP' GEO.'BAS' 0.
  534. 'UN' 'VIMP' GEO.'BAS' 0.
  535. ;
  536. RV = 'EQEX' RV 'CLIM'
  537. 'UN' 'UIMP' GEO.'HAUT' 0.
  538. 'UN' 'VIMP' GEO.'HAUT1' 0.
  539. 'UN' 'UIMP' GEO.'INJE' 0.
  540. 'UN' 'VIMP' GEO.'INJE' (FX '*' V0)
  541. ;
  542. 'FINSI' ;
  543.  
  544.  
  545. * Inconnues
  546. RV.'INCO' = TABLE 'INCO' ;
  547. RV.'INCO'.'UN'= 'KCHT' GEO.'$VTF' 'VECT' 'SOMMET' (0.0 0.0) ;
  548. RV.'INCO'.'ROG'= 'KCHT' GEO.'$VTF' 'VECT' 'SOMMET' (0.0 0.0) ;
  549. RV.'INCO'.'MU'= MU ;
  550. RV.'INCO'.'UNM' = 'COPIER' RV.'INCO'.'UN' ;
  551. RV.'INCO'.'UNMM' = 'COPIER' RV.'INCO'.'UN' ;
  552.  
  553. *
  554. * Gestion du transitoire
  555. * TFIN : temps final en secondes
  556. * CCFL : condition CFL utilisée
  557. *
  558.  
  559. TFIN = 6.0 ;
  560.  
  561. CCFL = 2.0 ;
  562. 'SI' (DEUXD) ;
  563. DT0 = CCFL '*' R0 '*' (LL '/' ('FLOTTANT' RAF)) '/' MPH ;
  564. 'FINSI' ;
  565. 'SI' (UND) ;
  566. DT0 = CCFL '*' R0 '*' (H '/' ('FLOTTANT' RAF)) '/' MPH ;
  567. 'FINSI' ;
  568. NPAS = 'ENTIER' (TFIN '/' DT0) ;
  569.  
  570. RV.'INCO'.'DT' = ('FLOTTANT' (TFIN '/' NPAS)) ;
  571.  
  572. 'MESSAGE' 'Pas de temps choisi DT = ' RV.'INCO'.'DT' ;
  573.  
  574. RV.'METHINV'.TYPINV = TYPI ;
  575. RV.'METHINV'.TYPRENU = 'SLOANE' ;
  576. RV.'METHINV'.IMPINV = 1 ;
  577. RV.'METHINV'.NITMAX = 3000 ;
  578. RV.'METHINV'.PRECOND = 5 ;
  579. RV.'METHINV'.SCALING = 1 ;
  580. RV.'METHINV'.OUBMAT = 1 ;
  581. RV.'METHINV'.ILUTLFIL = 2.0 ;
  582. RV.'METHINV'.'FCPRECI' = FCPRECI ;
  583. RV.'METHINV'.'FCPRECT' = FCPRECT ;
  584. RV.'METHINV'.RESID = 1.e-15 ;
  585.  
  586.  
  587. RV = 'EQEX' RV
  588. 'OPTI' 'EF' KSUPG 'IMPL' KPRES
  589. 'ZONE' (GEO.'$VTF') 'OPER' 'KBBT' 'CF' 'INCO' 'UN' 'PRES'
  590. ;
  591. RV.'INCO'.'CF' = 1.0 ;
  592.  
  593. * Terme source pour la divergence: DSRC
  594. * volume fermé: on impose la pression en un point MTPI
  595. * ou on ouvre le contour en 1 point
  596. RV = 'EQEX' RV
  597. 'OPTI' 'EF' 'IMPL' 'INCOD' KPRES
  598. 'ZONE' (GEO.'$VTF') 'OPER' 'FIMP' 'DSRC' 'INCO' 'PRES'
  599. ;
  600. RV.'INCO'.'DSRC' = 0.0 ;
  601.  
  602. *RV = 'EQEX' RV
  603. * 'CLIM' 'PRES' 'TIMP' MTPI 0.
  604. *;
  605.  
  606. *-------------------------------------------------------------
  607. * Definition de l'inconnue de PRESSION
  608. *
  609. RV.'INCO'.'PRES' = 'KCHT' GEO.'$VTF' 'SCAL' KPRES 0.0 ;
  610.  
  611.  
  612. *------------------------------------------------------------
  613. * Equation en temperature: TABLE RTF
  614. *------------------------------------------------------------
  615. RT = 'EQEX' 'ITMA' 0
  616. 'OPTI' 'EF' 'IMPL' KSUPG 'NOCONS'
  617. 'ZONE' (GEO.'$VTF') 'OPER' 'LAPN' 'ALPHA' 'INCO' 'TF'
  618. 'OPTI' 'EF' 'IMPL' KSUPG 'NOCONS'
  619. 'ZONE' (GEO.'$VTF') 'OPER' 'KONV' 'RHO' 'UN' 'ALPHA' 'INCO' 'TF'
  620. 'OPTI' 'EF' 'CENTREE' 'IMPL'
  621. 'ZONE' (GEO.'$VTF') 'OPER' 'FIMP' 'STF' 'INCO' 'TF'
  622. ;
  623. 'SI' ('EGA' KBDF 'BDF2') ;
  624. RT = 'EQEX' RT
  625. 'OPTI' KMASS 'CENTREE' 'BDF2'
  626. 'ZONE' (GEO.'$VTF') 'OPER' 'DFDT' 'RHO' 'TFNM' 'TFNMM' 'DT'
  627. 'INCO' 'TF' ;
  628. 'SINON' ;
  629. RT = 'EQEX' RT
  630. 'OPTI' KMASS 'CENTREE'
  631. 'ZONE' (GEO.'$VTF') 'OPER' 'DFDT' 'RHO' 'TFNM' 'DT'
  632. 'INCO' 'TF' ;
  633. 'FINSI' ;
  634.  
  635. 'SI' (UND) ;
  636. RT = 'EQEX' RT
  637. 'CLIM' 'TF' 'TIMP' (GEO.'BAS') TH
  638. ;
  639. 'FINSI' ;
  640. 'SI' (DEUXD) ;
  641. RT = 'EQEX' RT
  642. 'CLIM' 'TF' 'TIMP' (GEO.'INJE') TH
  643. ;
  644. 'FINSI' ;
  645.  
  646.  
  647. * Table des inconnues
  648. RT.'INCO' = RV.'INCO' ;
  649. 'SI' (DEUXD) ;
  650. TIN = 'MANUEL' 'CHPO' GEO.'INJE' 'SCAL' TH ;
  651. 'FINSI' ;
  652. 'SI' (UND) ;
  653. TIN = 'MANUEL' 'CHPO' GEO.'BAS' 'SCAL' TH ;
  654. 'FINSI' ;
  655.  
  656. RV.'INCO'.'TF' = 'KCHT' (GEO.'$VTF') 'SCAL' 'SOMMET' T0 TIN ;
  657. RV.'INCO'.'TFNM' = 'COPIER' RV.'INCO'.'TF' ;
  658. RV.'INCO'.'TFNMM' = 'COPIER' RV.'INCO'.'TF' ;
  659. RV.'INCO'.'STF' = 0.0 ;
  660. RV.'INCO'.'ALPHA' = MU '/' PRANDTL ;
  661.  
  662. RT.'METHINV'.TYPINV = TYPI ;
  663. RT.'METHINV'.TYPRENU = 'SLOANE' ;
  664. RT.'METHINV'.IMPINV = 1 ;
  665. RT.'METHINV'.NITMAX = 400 ;
  666. RT.'METHINV'.PRECOND = 5 ;
  667. RT.'METHINV'.SCALING = 1 ;
  668. RT.'METHINV'.OUBMAT = 1 ;
  669. RT.'METHINV'.ILUTLFIL = 1.5 ;
  670. RT.'METHINV'.'FCPRECI' = FCPRECI ;
  671. RT.'METHINV'.'FCPRECT' = FCPRECT ;
  672. RT.'METHINV'.RESID = 1.e-15 ;
  673.  
  674. * Execution - Algorithme BAS-MACH
  675.  
  676. * Conditions initiales
  677. RV.'INCO'.'PM' = 'PROG' P0 ;
  678. RV.'INCO'.'PE' = 'PROG' P0 ;
  679. *
  680. RV.'INCO'.'RHOM' = 'PROG' R0 ;
  681. RV.'INCO'.'REM' = 'PROG' (R0 '*' CVAIR '*' T0) ;
  682. RHO = ('INVERSE' (RV.'INCO'.'TF' '*' RAIR ) ) '*' P0 ;
  683. RV.'INCO'.'RHO' = 'KCHT' (GEO.'$VTF') 'SCAL' 'SOMMET' RHO ;
  684. RV.'INCO'.'RHONM' = 'COPIER' RV.'INCO'.'RHO' ;
  685. RV.'INCO'.'RHONMM' = 'COPIER' RV.'INCO'.'RHO' ;
  686. 'MESSAGE' '==> Masse initiale (kg) : ' ('SOMT' (RV.'INCO'.'RHO'
  687. '*' ('DOMA' (GEO.'$VTF') 'XXDIAGSI'))) ;
  688.  
  689. DIAG = 'DOMA' (GEO.'$VTF') 'XXDIAGSI' ;
  690.  
  691. RV.'INCO'.'LTPS' = 'PROG' 0.0 ;
  692. RV.'INCO'.'ERRORM' = 'PROG' 0.0 ;
  693. RV.'INCO'.'ERRORE' = 'PROG' 0.0 ;
  694.  
  695. 'OPTION' 'ECHO' 1 ;
  696.  
  697. VTOTAL = 'SOMT' ('DOMA' GEO.'$VTF' 'XXDIAGSI') ;
  698. 'MESSAGE' '==> Volume Fluide (m3) : ' VTOTAL ;
  699.  
  700. 'TEMPS' 'ZERO' ;
  701.  
  702. *
  703. * Ajout d'une matrice nulle pour préconditionner le système lors de
  704. * l'utilisation de méthodes itératives
  705. *
  706. 'SI' ('EGA' TYPI 3) ;
  707. RR = 'EQEX'
  708. 'OPTI' 'EF' KSUPG 'IMPL' KPRES
  709. 'ZONE' (GEO.'$VTF') 'OPER' 'KMAB' 1.0 'INCO' 'UN' 'PRES' ;
  710. RR.'INCO' = RV.'INCO' ;
  711. A1 B1 = 'KMAB' RR.'1KMAB' ;
  712. M1 = 'KCHT' GEO.'$VTF' 'VECT' 'SOMMET' 'COMP' '1UN' '2UN' (0. 0.) ;
  713. RV.'ESPR' = 'KOPS' 'CMCT' B1 B1 M1 ;
  714. 'FINSI' ;
  715.  
  716. 'SI' (interact) ;
  717. 'TRACER' (GEO.'VTF') 'TITR' 'Maillage...' ;
  718. 'FINSI' ;
  719.  
  720.  
  721. *
  722. * Boucle en temps
  723. *
  724.  
  725. 'REPETER' BCLTPS NPAS ;
  726. NBT = 'DIME' (RV.'INCO'.'PM') ;
  727. * temps courant
  728. TCOUR = ('EXTRAIRE' (RV.'INCO'.'LTPS') NBT)
  729. '+' (RV.'INCO'.'DT') ;
  730. 'MESSAGE' '==> Instant courant : ' TCOUR ;
  731. RV.'PASDETPS'.'TPS' = TCOUR ;
  732. RV.'PASDETPS'.'NUPASDT' = (RV.'PASDETPS'.'NUPASDT') '+' 1 ;
  733. RT.'PASDETPS'.'NUPASDT' = (RT.'PASDETPS'.'NUPASDT') '+' 1 ;
  734. RV.'INCO'.'LTPS' = (RV.'INCO'.'LTPS') 'ET' ('PROG' TCOUR) ;
  735. * grandeurs 0D à l'instant précédent
  736. PMN = 'EXTRAIRE' (RV.'INCO'.'PM') NBT ;
  737. RHOMN = 'EXTRAIRE' (RV.'INCO'.'RHOM') NBT ;
  738. REMN = 'EXTRAIRE' (RV.'INCO'.'REM') NBT ;
  739. * Calcul de la masse 0D à l'instant courant
  740. RHOM = 'MAXIMUM' ('PROG' (RHOMN '+' ((RV.'INCO'.'DT') '*'L '*'
  741. MPH '/' VTOTAL)) 0.0 ) ;
  742. PEXACT = P0 '+' (GAMA '*' MPH '*' RAIR '*' TH '*' L '*' TCOUR
  743. '/' VTOTAL) ;
  744. * Energie 0D
  745. REM = REMN '+' ((RV.'INCO'.'DT') '*'
  746. (MPH '*' L '*' CPAIR '*' TH)
  747. '/' VTOTAL) ;
  748. * boucle interne
  749. ERIM = 'PROG' ;
  750. NBCLI = 100 ;
  751. 'REPETER' BCLI NBCLI ;
  752. * Calcul d'une pression garantissant la conservation de la masse pour
  753. * cela on derive l'integrale de la loi d'état par rapport au temps
  754. TFOLD = 'COPIER' (RV.'INCO'.'TF') ;
  755. IRTF = 'INVERSE' (RAIR '*' (RV.'INCO'.'TF')) ;
  756. IRTF = 'SOMT' (DIAG '*' IRTF) ;
  757. PM = RHOM '*' VTOTAL '/' IRTF ;
  758. 'MESSAGE' '==> Pression 0D (Pa) : ' PM PEXACT ;
  759. * Nouvelle condition aux limites pour la vitesse à l'injection
  760. UINJ = MPH '*' RAIR '*' TH '/' PM ;
  761. 'MESSAGE' '==> Vitesse injection : ' UINJ ;
  762. RV.'CLIM' = (RV.'CLIM') '/' V0 '*' UINJ ;
  763. V0 = UINJ ;
  764. * Correction du champ de vitesse transportant compte tenu de
  765. * cette nouvelle vitesse: RV.'INCO'.'UN'
  766. UNX = 'EXCO' (RV.'CLIM') '1UN' ;
  767. UNX = 'NOMC' 'UX' UNX 'NATURE' 'DISCRET' ;
  768. UNY = 'EXCO' (RV.'CLIM') '2UN' ;
  769. UNY = 'NOMC' 'UY' UNY 'NATURE' 'DISCRET' ;
  770. 'MESSAGE' ('MAXIMUM' ('EXCO' (RV.'INCO'.'UN') 'UY')) ;
  771. RV.'INCO'.'UN' = 'KCHT' (GEO.'$VTF') 'VECT' 'SOMMET'
  772. (RV.'INCO'.'UN') (UNX 'ET' UNY) ;
  773. 'MESSAGE' ('MAXIMUM' ('EXCO' (RV.'INCO'.'UN') 'UY')) ;
  774. * Calcul du terme source de l'équation d'énergie
  775. 'SI' ('EGA' RV.'PASDETPS'.'NUPASDT' 2) ;
  776. DPDT = (PM '-' PMN) '/' (RV.'INCO'.'DT') ;
  777. 'SINON' ;
  778. PMNN = 'EXTRAIRE' (RV.'INCO'.'PM') (NBT '-' 1) ;
  779. DPDT = ((AN '*' PM) '+' (ANM1 '*' PMN) '+' (ANM2
  780. '*' PMNN)) '/' (ADT '*' (RV.'INCO'.'DT')) ;
  781. 'FINSI' ;
  782. 'MESSAGE' '==> DPDT C.M. : ' DPDT ;
  783. RV.'INCO'.'STF' = DPDT '/' CPAIR ;
  784. * Calcul de la densité à partir de la loi d'état et du terme source
  785. * de la divergence
  786. RHO = ('INVERSE' ((RV.'INCO'.'TF') '*' RAIR ) ) '*' PM ;
  787. 'SI' ('EGA' RV.'PASDETPS'.'NUPASDT' 2) ;
  788. DRDT = RHO '-' (RV.'INCO'.'RHONM')
  789. '/' RV.'INCO'.'DT' ;
  790. 'SINON' ;
  791. DRDT = ((RHO '*' AN) '+' ((RV.'INCO'.'RHONM') '*'
  792. ANM1) '+' ((RV.'INCO'.'RHONMM') '*' ANM2))
  793. '/' (ADT '*' (RV.'INCO'.'DT')) ;
  794. 'FINSI' ;
  795. 'MESSAGE' ('MAXIMUM' DRDT) ('MINIMUM' DRDT) ;
  796. GRHO = 'KOPS' RHO (GEO.'$VTF') 'GRADS' ;
  797. DRDT = DRDT '+' ( 'PSCA' (RV.'INCO'.'UN') ('MOTS' 'UX' 'UY')
  798. GRHO ('MOTS' 'UX' 'UY')) ;
  799. 'MESSAGE' ('MAXIMUM' DRDT) ('MINIMUM' DRDT) ;
  800. DRDT = DRDT '*' ('INVERSE' RHO) ;
  801. 'MESSAGE' ('MAXIMUM' DRDT) ('MINIMUM' DRDT) ;
  802. * Compatibilité (intégrale sur le volume de Div u = 0
  803. INTD = 'SOMT' (DRDT '*' DIAG) ;
  804. * calcul de l'integrale de u.n ds sur le contour
  805. TATA = 'DBIT' (RV.'INCO'.'UN') GEO.'$MENVF' ;
  806. 'MESSAGE' '==> Compatibilité Div u : ' INTD TATA;
  807. RV.'INCO'.'DSRC' = DRDT '-' (INTD '/' VTOTAL) '+' (TATA '/'
  808. VTOTAL) ;
  809. 'MESSAGE' ('SOMT' ((RV.'INCO'.'DSRC') '*' DIAG)) ;
  810. *
  811. RV.'INCO'.'RHO' = 'KCHT' (GEO.'$VTF') 'SCAL' 'SOMMET' RHO ;
  812. * résolution de la QDM et de la temperature
  813. * Pour le cas BDF2, le premier pas doit être traité en BDF1
  814. PAS1 = 'EGA' &BCLI 1 ;
  815. 'SI' (('EGA' KBDF 'BDF2') 'ET' PAS1) ;
  816. 'SI' ('EGA' RV.'PASDETPS'.'NUPASDT' 2) ;
  817. RV.'4DFDT'.'KOPT'.'ISCHT' = 0 ;
  818. RT.'4DFDT'.'KOPT'.'ISCHT' = 0 ;
  819. 'OUBLIER' (RT.'4DFDT') 'ARG3' ;
  820. 'OUBLIER' (RT.'4DFDT') 'ARG4' ;
  821. 'OUBLIER' (RV.'4DFDT') 'ARG3' ;
  822. 'OUBLIER' (RV.'4DFDT') 'ARG4' ;
  823. RT.'4DFDT'.'ARG3' = 'DT' ;
  824. RV.'4DFDT'.'ARG3' = 'DT' ;
  825. RT.'4DFDT'.'IARG' = 3 ;
  826. RV.'4DFDT'.'IARG' = 3 ;
  827. 'FINSI' ;
  828. 'SI' ('EGA' RV.'PASDETPS'.'NUPASDT' 3) ;
  829. RV.'4DFDT'.'KOPT'.'ISCHT' = 1 ;
  830. RT.'4DFDT'.'KOPT'.'ISCHT' = 1 ;
  831. RT.'4DFDT'.'ARG4' = 'DT' ;
  832. RV.'4DFDT'.'ARG4' = 'DT' ;
  833. RT.'4DFDT'.'ARG3' = 'TFNMM' ;
  834. RV.'4DFDT'.'ARG3' = 'UNMM' ;
  835. RT.'4DFDT'.'IARG' = 4 ;
  836. RV.'4DFDT'.'IARG' = 4 ;
  837. 'FINSI' ;
  838. 'FINSI' ;
  839. RV.'METHINV'.'XINIT' = ('NOMC' '1UN' ('EXCO' (RV.'INCO'.'UN')
  840. 'UX') 'NATURE' 'DISCRET') '+'
  841. ('NOMC' '2UN' ('EXCO' (RV.'INCO'.'UN')
  842. 'UY') 'NATURE' 'DISCRET') '+'('NOMC' 'PRES' (RV.'INCO'.'PRES')
  843. 'NATURE' 'DISCRET') ;
  844. 'TEMPS' 'ZERO' ;
  845. EXEC RV ;
  846. TABTPS = TEMP 'NOEC';
  847. tcpu = TABTPS.'TEMPS_CPU'.'INITIAL';
  848. 'MESSAGE' 'Temps CPU = ' tcpu ' ms' ;
  849. RT.'METHINV'.'XINIT' = 'NOMC' 'TF' (RV.'INCO'.'TF') 'NATURE'
  850. 'DISCRET' ;
  851. 'TEMPS' 'ZERO' ;
  852. EXEC RT ;
  853. TABTPS = TEMP 'NOEC';
  854. tcpu = TABTPS.'TEMPS_CPU'.'INITIAL';
  855. 'MESSAGE' 'Temps CPU = ' tcpu ' ms' ;
  856. RES1 = 'MAXIMUM' ('ABS' ((RV.'INCO'.'TF') '-' TFOLD)) ;
  857. 'MESSAGE' '==> Residu RES : ' RES1 ;
  858. 'SI' ( interact) ;
  859. VUN1 = RV.'INCO'.'UN' ;
  860. VVN1 = 'VECT' VUN1 0.05 'UX' 'UY' 'JAUN' ;
  861. 'TRACE' RV.'INCO'.'TF' (GEO.'VTF') ('CONT' (GEO.'VTF'))
  862. VVN1 'TITR' 'Champ de vitesse' ;
  863. 'FINSI' ;
  864. * Calcul du terme source de la QDM
  865. RHOP = (-1.) '*' (R0 '-' RHO) ;
  866. * Ajout du terme de compressibilité dans Div Tau
  867. TAU = 'KCHT' (GEO.'$VTF') 'SCAL' 'SOMMET' ((RV.'INCO'.'DSRC')
  868. '*' (RV.'INCO'.'MU') '*' (-1.0) '/' 3.0 ) ;
  869. 'MESSAGE' ('MAXIMUM' TAU) ('MINIMUM' TAU) ;
  870. DIVT = 'KOPS' TAU GEO.'$VTF' 'GRADS' ;
  871. DIVTX = ('EXCO' DIVT 'UX') ;
  872. DIVTY = ('EXCO' DIVT 'UY') ;
  873. 'MESSAGE' ('MAXIMUM' DIVTX) ('MINIMUM' DIVTX) ;
  874. 'MESSAGE' ('MAXIMUM' DIVTY) ('MINIMUM' DIVTY) ;
  875. RHOP = 'KCHT' (GEO.'$VTF') 'SCAL' 'SOMMET' RHOP ;
  876. ROGX = 'KCHT' (GEO.'$VTF') 'SCAL' 'SOMMET' DIVTX ;
  877. ROGY = 'KCHT' (GEO.'$VTF') 'SCAL' 'SOMMET'
  878. (DIVTY '-' (GRAV '*' RHOP));
  879. ROGX = 'NOMC' 'UX' ROGX 'NATU' 'DISCRET' ;
  880. ROGY = 'NOMC' 'UY' ROGY 'NATU' 'DISCRET' ;
  881. RV.'INCO'.'ROG'= 'KCHT' (GEO.'$VTF') 'VECT' 'SOMMET'
  882. (ROGX 'ET' ROGY) ;
  883. * Affichage des bilans masse et énergie 0D/MultiD
  884. MAS0D = RHOM '*' VTOTAL ;
  885. ENR0D = REM '*' VTOTAL ;
  886. MMD = 'SOMT' (DIAG '*' RHO) ;
  887. EMD = 'SOMT' (DIAG '*' RHO '*' CVAIR '*' (RV.'INCO'.'TF')) ;
  888. 'MESSAGE' 'Bilan Masse 0D: ' MAS0D ' MultiD: ' MMD ;
  889. 'MESSAGE' 'Bilan Energie 0D: ' ENR0D ' MultiD: ' EMD ;
  890. ERRORM = (MMD '-' MAS0D) '/' MAS0D '*' 100.0 ;
  891. ERRORE = (EMD '-' ENR0D) '/' ENR0D '*' 100.0 ;
  892. ERIM = ERIM 'ET' ('PROG' ERRORE) ;
  893. 'SI' ('<' RES1 1.E-4) ;
  894. 'QUITTER' BCLI ;
  895. 'FINSI' ;
  896. 'FIN' BCLI ;
  897. 'LISTE' ERIM ;
  898. * sauvegarde et avancement en temps
  899. RV.'INCO'.'RHOM' = RV.'INCO'.'RHOM' 'ET' ('PROG' RHOM) ;
  900. RV.'INCO'.'REM' = RV.'INCO'.'REM' 'ET' ('PROG' REM) ;
  901. RV.'INCO'.'PM' = RV.'INCO'.'PM' 'ET' ('PROG' PM) ;
  902. RV.'INCO'.'PE' = RV.'INCO'.'PE' 'ET' ('PROG' PEXACT) ;
  903. RV.'INCO'.'ERRORM' = RV.'INCO'.'ERRORM' 'ET'
  904. ('PROG' ERRORM) ;
  905. RV.'INCO'.'ERRORE' = RV.'INCO'.'ERRORE' 'ET'
  906. ('PROG' ERRORE) ;
  907. * mise a jour des champ multiD pour DFDT
  908. RV.'INCO'.'UNMM' = 'COPIER' RV.'INCO'.'UNM' ;
  909. RV.'INCO'.'TFNMM' = 'COPIER' RV.'INCO'.'TFNM' ;
  910. RV.'INCO'.'RHONMM' = 'COPIER' RV.'INCO'.'RHONM' ;
  911. RV.'INCO'.'UNM' = 'COPIER' RV.'INCO'.'UN' ;
  912. RV.'INCO'.'RHONM' = 'COPIER' RV.'INCO'.'RHO' ;
  913. RV.'INCO'.'TFNM' = 'COPIER' RV.'INCO'.'TF' ;
  914. 'SI' ('>EG' TCOUR TFIN) ;
  915. 'QUITTER' BCLTPS ;
  916. 'FINSI' ;
  917. 'FIN' BCLTPS ;
  918. 'TEMPS' ;
  919. 'TEMPS' 'PLACE' ;
  920. *
  921. * Creation d'une droite pour profil
  922. *
  923. mt = GEO.'VTF' ;
  924. m1 = 'CHANGER' mt 'LIGNE' ;
  925. pt1 = 'POIN' mt 'PROC' P2 ;
  926. pt2 = 'POIN' mt 'PROC' P4 ;
  927. p1 = 'POIN' mt 'DROIT' pt1 pt2 0.0001 ;
  928. LMIL = 'ELEM' m1 'APPUYE' p1 ;
  929. *
  930. * Profils sur LMIL
  931. * Température (K)
  932.  
  933. *==========================================================
  934. * Test de non régression
  935. *==========================================================
  936.  
  937. EVVT = 'EVOL' 'VERT' 'CHPO' (RV.'INCO'.'TF') LMIL ;
  938. TREF = 'PROG' 600.00 606.09 612.13 618.67
  939. 625.08 632.00 638.60 645.44 651.59
  940. 657.31 661.76 664.78 665.90 664.56 660.85
  941. 653.95 644.61 631.93 617.26 599.82
  942. 581.27 561.06 540.84 520.36 500.90
  943. 482.52 465.88 451.28 438.73 428.64
  944. 420.45 414.55 410.11 407.34 405.46
  945. 404.50 403.94 403.74 403.65 403.63 403.62 ;
  946. TCOUR = 'EXTRAIRE' EVVT 'ORDO' ;
  947. DT = 'MAXIMUM' (TCOUR '-' TREF) 'ABS' ;
  948. 'MESSAGE' 'Erreur max sur le profil de température = ' DT ;
  949. DTMAX = 0.01 ;
  950. TEST = DT < DTMAX ;
  951. 'SI' TEST ;
  952. 'ERREUR' 0 ;
  953. 'SINON' ;
  954. 'ERREUR' 5 ;
  955. 'FINSI' ;
  956.  
  957. *
  958. * End of dgibi file INJECTION AIR/AIR
  959. *
  960. *
  961. 'FIN' ;
  962.  
  963.  
  964.  
  965.  
  966.  
  967.  
  968.  
  969.  

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