Télécharger equ_chaleurVF2_dirneummixte.dgibi

Retour à la liste

Numérotation des lignes :

  1. * fichier : equ_chaleurVF2_dirneummixte.dgibi
  2. ************************************************************************
  3. ************************************************************************
  4. ***********************************************************************
  5. * NOM : equ_chaleurVF2_dirneummixte.dgibi
  6. * ___
  7. *
  8. * DESCRIPTION : Solution stationnaire de l'équation de la chaleur (2D)
  9. * ___________
  10. *
  11. * GEOMETRIE : Un carré
  12. * ----------
  13. *
  14. * y
  15. * ^ y=1
  16. * |------------
  17. * | |
  18. * | |
  19. * | |
  20. * |x = 0 |x=1
  21. * | |
  22. * | |
  23. * O-----------------------------> x
  24. * y=0
  25. *
  26. *
  27. * EQUATIONS :
  28. * ----------
  29. *
  30. * - Equations :
  31. *
  32. * div K GRAD T = f
  33. *
  34. * | |
  35. * avec K = |x+2 x |
  36. * |x y+2 |
  37. * | |
  38. *
  39. * f = e-(x+y) (3x+y)
  40. *
  41. * - Conditions aux limites :
  42. *
  43. * conditions de Dirichlet, conditions de flux,
  44. * conditions mixtes calculées à partir de la
  45. * restriction de la solution exacte sur le bord
  46. * exacte
  47. *
  48. * - Solution exacte :
  49. *
  50. * T(x,y)= e-(x+y)
  51. *
  52. *
  53. * DISCRETISATION : une méthode de Volume Finis en espace.
  54. * ______________
  55. *
  56. *
  57. *
  58. *
  59. *
  60. * Opérateurs utilisés : PENT (option VF implicite)
  61. * LAPN (option VF implicite)
  62. *
  63. * RESOLUTION : - Solveur BiCGStab
  64. * __________ - Préconditionneur ILU(0)
  65. *
  66. * TESTS EFFECTUES : Vérification de l'ordre proche de 2 en espace de la méthode
  67. * _______________ (utilisation d'une norme pseudo-L2) et de la
  68. * précision absolue sur le maillage le plus fin.
  69. *
  70. *
  71. *
  72. * LANGAGE : GIBIANE-CASTEM 2000
  73. * AUTEUR : C LE POTIER
  74. * mél : clepotier@cea.fr
  75. ************************************************************************
  76. * VERSION : v2, 02/03, version initiale
  77. * HISTORIQUE :
  78. ************************************************************************
  79. * Prière de PRENDRE LE TEMPS de compléter les commentaires
  80. * en cas de modification de ce sous-programme afin de faciliter
  81. * la maintenance !
  82. ******************************************************************
  83.  
  84. interact= FAUX ;
  85. complet = FAUX ;
  86. graph = VRAI ;
  87. logxmgr = FAUX ;
  88. *
  89. 'OPTION' 'DIME' 2 'ELEM' QUA4 'MODE' 'PLAN' ;
  90. 'OPTION' 'ISOV' 'SULI' ;
  91. 'OPTION' 'ECHO' 0 ;
  92. nbisov = 15 ;
  93. 'SI' ('NON' interact) ;
  94. 'OPTION' 'TRAC' 'PS' ;
  95. 'OPTION' 'ISOV' 'LIGNE' ;
  96. 'FINSI' ;
  97. *
  98. ** Erreur Linfini entre deux Champoints.
  99. *
  100. DEBPROC CALCERR vitp1*'CHPOINT' vit*'CHPOINT' ;
  101. AUX = (vitp1 - vit) ;
  102. err = MAXI (AUX) 'ABS' ;
  103. FINPROC err ;
  104. *
  105. ** Erreur L2 entre deux Champoints.
  106. *
  107. DEBPROC CALCERR2 vitp1*'CHPOINT' vit*'CHPOINT' vol*'CHPOINT' ;
  108. er2 = vitp1 '-' vit ;
  109. compv = 'EXTRAIRE' er2 'COMP' ;
  110. er2 = 'PSCAL' er2 er2 compv compv ;
  111. suppv = 'EXTRAIRE' vol 'MAIL' ;
  112. chpun = 'MANUEL' 'CHPO' suppv 1 'SCAL' 1 ;
  113. voltot = 'XTY' chpun ('MOTS' 'SCAL') vol ('MOTS' 'SCAL') ;
  114. error = 'XTY' er2 ('MOTS' 'SCAL') vol ('MOTS' 'SCAL') ;
  115. error = error '/' voltot ;
  116. error = error '**' 0.5 ;
  117. FINPROC error ;
  118. *
  119. * Procédure paramétrée (raffinement)
  120. * renvoyant l'erreur en norme L2 sur la température.
  121. * On calcule une solution de l'équation de Laplace
  122. * (équations de la chaleur) ;
  123. *
  124. 'DEBPROC' CALCUL nraff*'ENTIER' ;
  125. *nraff=2 ;
  126. *
  127. * titre global pour les dessins
  128. *
  129. titglob = 'CHAINE' ' ; nraff=' nraff ;
  130. *
  131. * Paramètres physiques
  132. *
  133. cv= 1.0 ;
  134. mu = 0.0 ;
  135. kappa = 1.0 ;
  136. rho = 1.0 ;
  137. *
  138. * Géométrie
  139. *
  140. L = 1.;
  141. H = 1.;
  142. pA = 0. 0. ;
  143. pB = L 0. ;
  144. pC = L H ;
  145. pD = 0. H ;
  146. *
  147. * Paramètres de la discrétisation de base
  148. *
  149. 'SI' complet ;
  150. nAB = 8 ;
  151. nBC = 7 ;
  152. nCD = 5 ;
  153. nDA = 9 ;
  154. 'SINON' ;
  155. nAB = 2 ;
  156. nBC = 3 ;
  157. nCD = 4 ;
  158. nDA = 5 ;
  159.  
  160. 'FINSI' ;
  161. *
  162. * Rotation et translation aditionnelle + bruit blanc
  163. * + raffinement
  164. *
  165. vtran = ((** 2 0.5) (* -1 (** 3 0.5))) ;
  166. orig = (0.D0 0.D0) ;
  167. arot = 11.3D0 ;
  168. lesdens = 'PROG' (1. '/' nAB) (1. '/' nCD) (1. '/' nBC) (1. '/' nDA) ;
  169. ampbruit = (0.15 * ('MINIMUM' lesdens)) ;
  170. *
  171. * Géométrie discrétisée
  172. *
  173. bas = 'DROIT' nAB pA pB ;
  174. droite = 'DROIT' nBC pB pC ;
  175. haut = 'DROIT' nCD pC pD ;
  176. gauche = 'DROIT' nDA pD pA ;
  177. pourtour = bas 'ET' droite 'ET' haut 'ET' gauche ;
  178. mt = 'SURFACE' pourtour 'PLAN' ;
  179. *mt = 'DALLER' bas droite haut gauche 'PLAN' ;
  180. *
  181. * On rajoute le bruit
  182. *
  183. pmt = 'CHANGER' mt 'POI1' ;
  184. pcnt= 'CHANGER' pourtour 'POI1' ;
  185. inmt= 'DIFF' pmt pcnt;
  186. brbl1 = 'BRUIT' 'BLAN' 'GAUS' 0.D0 ampbruit inmt ;
  187. brbl2 = 'BRUIT' 'BLAN' 'GAUS' 0.D0 ampbruit inmt ;
  188. brbl = 'ET' ('NOMC' 'UX' brbl1) ('NOMC' 'UY' brbl2) ;
  189. *'DEPLACER' mt 'PLUS' brbl;
  190. * 'DEPLACER' mt 'AFFI' 10.0 PA PD;
  191. *
  192. * On raffine nraff fois
  193. *
  194. 'SI' (nraff > 0) ;
  195. 'REPETER' bcl nraff ;
  196. mt = 'CHANGER' mt 'QUADRATIQUE' ;
  197. $mt = 'DOMA' mt 'MACRO' ;
  198. mt = ('DOMA '$mt 'MAILLAGE') ;
  199. 'MENAGE' ;
  200. 'FIN' bcl ;
  201. 'FINSI' ;
  202. *
  203. * Eventuellement, on trace le résultat
  204. *
  205. 'SI' graph ;
  206. titgeo = 'CHAINE' 'mt ' 'NBPO=' ('NBNO' mt)
  207. ' NBELEM=' ('NBEL' mt) titglob ;
  208. 'TRACER' mt 'NOEUD' 'TITRE' titgeo ;
  209. 'FINSI' ;
  210. *
  211. * Modèles
  212. *
  213. Mmt = 'CHAN' mt 'QUAF' ;
  214. 'ELIM' Mmt 1.E-5 ;
  215. MDNS = 'NAVIER_STOKES' ;
  216. $mt = 'MODE' Mmt MDNS 'LINE' ;
  217. mt = 'DOMA' $mt 'MAILLAGE' ;
  218.  
  219. modmt = 'MODELISER' ('DOMA' $mt 'MAILLAGE') 'THERMIQUE' ;
  220. mtm = ('DOMA' $mt 'MAILLAGE') ;
  221. mtc = 'COULEUR' ('DOMA' $mt 'CENTRE') 'JAUN' ;
  222. mtf = 'COULEUR' ('DOMA' $mt 'FACE') 'BLEU' ;
  223. cmt = 'DOMA' $mt 'ENVELOPP' ;
  224. Mcmt = 'CHAN' cmt 'QUAF' ;
  225. 'ELIM' (Mmt 'ET' Mcmt) 1.E-5 ;
  226. $cmt = 'MODE' Mcmt MDNS 'LINE' ;
  227.  
  228. * DEFINITION DES DIFFERENTS BORD DE LA FRONTIERE
  229. gau1 = cmt 'ELEM' 'COMPRIS' PD PA ;
  230. bas1 = cmt 'ELEM' 'COMPRIS' PA PB ;
  231. Mbas = 'CHAN' bas1 'QUAF' ;
  232. Mgau = 'CHAN' gau1 'QUAF' ;
  233. Mcvn = 'CHAN' cmt 'QUAF' ;
  234. 'ELIM' (Mmt 'ET' Mbas 'ET' Mgau 'ET' Mcvn) 1.E-5 ;
  235. $bas = 'MODE' Mbas MDNS 'LINE' ;
  236. $gau = 'MODE' Mgau MDNS 'LINE';
  237. $cvn = 'MODE' Mcvn MDNS 'LINE';
  238. cvn = 'DOMA' $cvn 'MAILLAGE' ;
  239.  
  240. cdir = 'MANUEL' 'POI1' pa ;
  241. pdir1 = 'MANUEL' 'POI1' (('DOMA' $bas 'CENTRE') 'POIN' 'PROC' pA) ;
  242. pdir2 = 'MANUEL' 'POI1' (('DOMA' $gau 'CENTRE') 'POIN' 'PROC' pA) ;
  243. cvn = 'DIFF' ('DOMA' $cvn 'CENTRE') pdir1 ;
  244. cvn = 'DIFF' cvn pdir2 ;
  245.  
  246. 'SI' graph ;
  247. 'TRACER' (mt et (cvn 'COULEUR' 'ROUG')) 'TITRE' 'Von Neumann B.C.' ;
  248. 'FINSI' ;
  249.  
  250.  
  251.  
  252. *
  253. *
  254. *
  255. * Solution exacte aux centres
  256. *
  257. xxc yyc = 'COORDONNEE' ('DOMA' $mt 'CENTRE') ;
  258. solexc = (exp ((-1.D0)*(xxc + yyc))) ;
  259. solexc = 'KCHT' $mt 'SCAL' 'CENTRE' solexc ;
  260. *
  261. * Gradient exact (aux faces)
  262. *
  263. xxf yyf = 'COORDONNEE' ('DOMA' $mt 'FACE') ;
  264. *
  265. *
  266. * Second membre de l'équation
  267. *
  268. sour = (exp ((-1.D0)*(xxc + yyc)))* ((3*xxc) + yyc) ;
  269. sour = 'KCHT' $mt 'SCAL' 'CENTRE' sour ;
  270.  
  271.  
  272. xxlim yylim = 'COORDONNEE' ('DOMA' $cmt 'CENTRE') ;
  273. solim = (exp ((-1.D0)*(xxlim + yylim))) ;
  274. solim = 'KCHT' $cmt 'SCAL' 'CENTRE' solim ;
  275.  
  276. *
  277. xxf yyf = 'COORDONNEE' ('DOMA' $mt 'FACE') ;
  278. gradtx = (-1.D0)*(2. + (2*xxf)) ;
  279. gradty = (-1.D0)*(1. + (xxf+yyf)) ;
  280. gradtx = 'NOMC' 'UX' gradtx 'NATU' 'DISCRET' ;
  281. gradty = 'NOMC' 'UY' gradty 'NATU' 'DISCRET' ;
  282. *
  283.  
  284.  
  285. CHYB2 = 'DOMA' $mt 'NORMALE' ;
  286. qchal = (gradtx 'ET' gradty) ;
  287. qchal = -1.0 * (gradtx + gradty) ;
  288. * CONDITION DE FLUX
  289. qlim2 = 'REDU' qchal (pdir1 ) ;
  290. CHYB2R = 'REDU' CHYB2 pdir1 ;
  291. MOT1 = 'MOTS' 'UX' 'UY' ;
  292. MOT2 = 'MOTS' 'UX' 'UY' ;
  293. FLU1 = 'PSCA' CHYB2R qlim2 MOT1 MOT2 ;
  294. qlim2 = NOMC(FLU1) 'FLUX';
  295.  
  296.  
  297. * CONDITION DE DIRICHLET
  298. solim = 'REDU' solim (pdir2);
  299. * DEFINITION DES PARAMETRES POUR LES CONDITIONS MIXTES
  300. qaux = 'REDU' qchal cvn ;
  301. qlim = (0.0*qaux);
  302. CHYB2R = 'REDU' CHYB2 cvn ;
  303. MOT1 = 'MOTS' 'UX' 'UY' ;
  304. MOT2 = 'MOTS' 'UX' 'UY' ;
  305. FLU1 = 'PSCA' CHYB2R qlim MOT1 MOT2 ;
  306. qlim = NOMC(FLU1) 'FLUX';
  307.  
  308. AUX = 1.D0 + (0.0D0*xxf);
  309. XPAR1 = 'NOMC' 'PAR1' AUX 'NATU' 'DISCRET' ;
  310. XPAR1 = 'REDU' XPAR1 cvn ;
  311. qchal = (-1.0) * (gradtx + gradty) ;
  312. qaux = 'REDU' qchal cvn ;
  313. FLU1 = 'PSCA' CHYB2R qaux MOT1 MOT2 ;
  314. XPAR2 = 'NOMC' FLU1 'PAR2';
  315. VAL = XPAR1 + XPAR2 + QLIM;
  316.  
  317.  
  318.  
  319.  
  320.  
  321. *
  322. * Calcul du tenseur
  323. *
  324. K11 = 'KOPS' xxc '+' 2. ;
  325. K22 = 'KOPS' yyc '+' 1. ;
  326. K21 = xxc ;
  327.  
  328. K11 = NOMC(K11) 'K11' 'NATU' 'DISCRET';
  329. K22 = NOMC(K22) 'K22' 'NATU' 'DISCRET';
  330. K21 = NOMC(K21) 'K21' 'NATU' 'DISCRET';
  331. PERM = (K11 'ET' K22) ;
  332. PERM = PERM 'ET' K21;
  333.  
  334. *
  335. * Mise en place du calcul numérique
  336. *
  337. * Mise en place du calcul numérique
  338. *
  339. * équation de Laplace
  340. *
  341. * on utilise une méthode de Newton pour résoudre :
  342. * - \Delta T = 0 (\Delta opérateur laplacien)
  343. * - avec mélange de conditions aux limites (Dirichlet,Neuman,Mixte)
  344. *
  345. * T_0 : estimation initiale de la solution
  346. * On a T_1 = T_0 - {\Delta'}^{-1} (\Delta T_0)
  347. *
  348. * L'opérateur 'LAPN' 'VF' nous donne la matrice \Delta' et le
  349. * résidu \Delta T_0.
  350. * On n'inverse bien évidemment pas \Delta' mais on résoud le système:
  351. * \Delta' IncT = \Delta T_0
  352. * => IncT = {\Delta'}^{-1} (\Delta T_0)
  353. *
  354. * La méthode de Newton doit converger en un pas (on vérifie que le
  355. * résidu (\Delta T_1) est nul après le premier pas.
  356. *
  357. *
  358. tmp = 0.0 ;
  359. t0 = 'KCHT' $mt 'SCAL' 'CENTRE' 0.D0 ;
  360.  
  361. gradt0 mchamt = 'PENT' $mt 'FACE' 'MPFA' t0
  362. 'DISPDIF' PERM 'TIMP' solim 'QIMP' qlim2 'MIXT' val;
  363.  
  364. jaco chpres dt = 'LAPN' 'VF' 'CLAUDEIS' 'IMPL'
  365. $mt t0 gradt0 mchamt 'QIMP' qlim2 'MIXT' val 'TIMP' solim ;
  366.  
  367.  
  368.  
  369. mamat = 'KOPS' 'MULT' -1.000D0 jaco ;
  370. matot = mamat ;
  371.  
  372. AUX = ('EXCO' chpres 'RETN' 'RETN') - (NOMC(sour) 'RETN');
  373. rv = 'EQEX' ;
  374. rv . 'METHINV' . 'TYPINV' = 3 ;
  375. rv . 'METHINV' . 'PRECOND' = 3 ;
  376. rv . 'METHINV' . 'MATASS' = matot ;
  377. rv . 'METHINV' . 'MAPREC' = matot ;
  378. res = 'KRES' matot 'TYPI' (rv . 'METHINV')
  379. 'SMBR'AUX
  380. 'IMPR' 0 ;
  381. t1 = t0 '+' ('EXCO' res 'RETN' 'SCAL') ;
  382. gradt1 = 'PENT' $mt 'FACE' 'MPFA' t1
  383. 'DISPDIF' PERM 'TIMP' solim 'QIMP' qlim2 'MIXT' val 'GRADGEO' mchamt;
  384. jacbid chpres1 dt = 'LAPN' 'VF' 'CLAUDEIS' 'IMPL'
  385. $mt t0 gradt1 mchamt 'QIMP' qlim2 'MIXT' val 'TIMP' solim ;
  386. AUX = ('EXCO' chpres1 'RETN' 'RETN') - (NOMC(sour) 'RETN');
  387. res = 'KRES' matot 'TYPI' (rv . 'METHINV')
  388. 'SMBR'AUX
  389. 'IMPR' 0 ;
  390. mres1 = 'MAXIMUM' res 'ABS' ;
  391. 'MESSAGE' ('CHAINE' 'Maxi. mres1 =' mres1) ;
  392. 'SI' ('>' mres1 1.e-5) ;
  393. 'MESSAGE' 'La méthode de Newton na pas converge en un pas.'
  394. 'ERREUR' 5 ;
  395. 'FINSI' ;
  396. *
  397. *
  398. * Résultats
  399. *
  400. 'SI' graph ;
  401. *
  402. * solutions exactes
  403. *
  404. tn = solexc ;
  405. chm_tnex = 'KCHA' $mt 'CHAM' tn ;
  406. titt = 'CHAINE' 'Température exacte' titglob ;
  407. 'TRACER' chm_tnex modmt nbisov 'TITRE' titt ;
  408. *
  409. * graphe de convergence de la méthode itérative
  410. *
  411. conver = (rv . 'METHINV' . 'CONVINV') ;
  412. dimcon = 'DIME' conver ;
  413. labs = 'PROG' 1.D0 'PAS' 1.D0 dimcon ;
  414. lord = ('LOG' conver) '/' ('LOG' 10.D0) ;
  415. evtot = 'EVOL' 'MANU' 'ITER' labs 'RESID' lord ;
  416. titev = 'CHAINE' 'Historique de convergence' titglob ;
  417. 'DESSIN' evtot 'TITR' titev 'LEGE' ;
  418. *
  419. * solutions calculées
  420. *
  421. tn = t1 ;
  422. chm_tn = 'KCHA' $mt 'CHAM' tn ;
  423. titt = 'CHAINE' 'Température calculée' titglob ;
  424. 'TRACER' chm_tn modmt nbisov 'TITRE' titt ;
  425. *
  426. * erreur
  427. *
  428. titt = 'CHAINE' 'Abs (Température calculée - Température exacte)'
  429. titglob ;
  430. 'TRACER' ('ABS' (chm_tnex '-' chm_tn)) modmt nbisov 'TITRE' titt ;
  431. 'FINSI' ;
  432. *
  433. * Calcul des erreurs par rapport à la solution analytique
  434. *
  435. vol = 'DOMA' $mt 'VOLUME' ;
  436. echloc = (('MESURE' mt) '/' ('NBEL' mt)) ** 0.5 ;
  437. tn = t1 ;
  438. errlit = CALCERR tn solexc ;
  439. errl2t = CALCERR2 tn solexc vol ;
  440. mres1 = 'MAXIMUM' solexc 'ABS' ;
  441. 'MESSAGE' ('CHAINE' 'Maxi. solexc =' mres1) ;
  442. 'MESSAGE' '-------------------------------------------------' ;
  443. 'MESSAGE' ('CHAINE' 'Erreur en norme LINF : ' errlit) ;
  444. 'MESSAGE' ('CHAINE' 'Erreur en norme L2 : ' errl2t) ;
  445. 'MESSAGE' '-------------------------------------------------' ;
  446. *'OPTION' 'DONN' 5 ;
  447. 'FINPROC' echloc errl2t ;
  448. * Fin de la procédure de calcul *
  449. *___________________________________________________________*
  450. *-----------------------------------------------------------
  451. * Paramètres du calcul
  452. *
  453. * lraff : nb raffinement du maillage (à chaque fois, on découpe
  454. * les éléments en quatre).
  455. 'SI' complet ;
  456. lraff = 'LECT' 0 PAS 1 3 ;
  457. 'SINON' ;
  458. lraff = 'LECT' 0 PAS 1 1 ;
  459. 'FINSI' ;
  460. *
  461. *-----------------------------------------------------------*
  462. * Boucles sur les différents paramètres du calcul *
  463. ordok = VRAI ;
  464. precok = VRAI ;
  465. * ordre théorique en espace de la méthode
  466. ordtth = 2 ;
  467. * erreur L2 max pour la solution (raffinement 2, complet=FAUX)
  468. errtmax = 4.D-3 ;
  469. *
  470. * Boucle sur les raffinements
  471. *
  472. lh = 'PROG' ;
  473. lerl2t = 'PROG' ;
  474. 'REPETER' iraff ('DIME' lraff) ;
  475. raff = 'EXTRAIRE' lraff &iraff ;
  476. h errt = CALCUL raff ;
  477. lh = 'ET' lh ('PROG' h ) ;
  478. lerl2t = 'ET' lerl2t ('PROG' errt) ;
  479. 'FIN' iraff ;
  480. * 'TEMPS' 'PLACE';
  481. lh = lh '/' ('EXTRAIRE' lh 1) ;
  482. lh = ('LOG' lh) '/' ('LOG' 10.D0) ;
  483. lerl2t = ('LOG' lerl2t) '/' ('LOG' 10.D0) ;
  484. tl2 = 'EVOL' 'MANU' 'Long. carac.' lh
  485. 'Err. tempér.' lerl2t ;
  486.  
  487. *
  488. * Recherche des ordres de convergence
  489. *
  490. cpt tlipol = @POMI tl2 1 'IDEM' ;
  491. ordt = cpt . 1 ;
  492. 'MESSAGE' ('CHAINE' 'ordre température=' ordt) ;
  493. *
  494. * Tracés graphiques
  495. *
  496. 'SI' graph ;
  497. tableg = 'TABLE' ;
  498. tableg . 'TITRE' = 'TABLE' ;
  499. tableg . 1 = 'MARQ CROI NOLI' ;
  500. tableg . 'TITRE' . 1 = 'Erreur L2 calculée' ;
  501. tableg . 2 = 'TIRR' ;
  502. tableg . 'TITRE' . 2 = 'Erreur L2 interpolée' ;
  503. titdest= 'CHAINE' 'Err L2 Tempér., ordre=' ordt ;
  504. 'DESSIN' (tl2 'ET' tlipol) 'TITRE' titdest tableg ;
  505. 'FINSI' ;
  506. *
  507. * Tests des ordres de convergence
  508. * Tests des erreurs L2 sur le maillage le plus fin dans le
  509. * cas complet=faux
  510. *
  511. ordok = 'ET' ordok ('EGA' ('ENTIER' ('+' ordt 0.5)) ordtth) ;
  512. 'SI' ('EGA' complet FAUX) ;
  513. precok = 'ET' precok ('<' errt errtmax) ;
  514. 'FINSI' ;
  515. 'FINSI' ;
  516. 'SI' ('NON' ordok) ;
  517. 'MESSAGE' 'On n_obtient pas un ordre de convergence correct' ;
  518. 'ERREUR' 5 ;
  519. 'FINSI' ;
  520. 'SI' ('NON' precok) ;
  521. 'MESSAGE' 'On n_obtient pas une précision correcte' ;
  522. 'ERREUR' 5 ;
  523. 'FINSI' ;
  524. 'SI' interact ;
  525. 'OPTION' 'ECHO' 1 ;
  526. 'OPTION' 'DONN' 5 ;
  527. 'SI' logxmgr ;
  528. * Sortie pour xmgr
  529. 'OPTION' ECHO 0 ;
  530. 'OPTION' 'IMPR' 'file.txt' ;
  531. ndim = 'DIME' lh ;
  532. 'REPETER' BL1 ndim ;
  533. 'MESSAGE' ('EXTRAIRE' lh &BL1) ' ' ('EXTRAIRE' lerl2t &BL1) ;
  534. 'FIN' BL1 ;
  535. lh1 = 'EXTRAIRE' tlipol 'ABSC' ;
  536. lerr = 'EXTRAIRE' tlipol 'ORDO' ;
  537. ndim = 'DIME' lh1 ;
  538. 'OPTION' 'IMPR' 'file2.txt' ;
  539. 'REPETER' BL1 ndim ;
  540. 'MESSAGE' ('EXTRAIRE' lh1 &BL1) ' ' ('EXTRAIRE' lerr &BL1) ;
  541. 'FIN' BL1 ;
  542. 'OPTION' 'IMPR' 'poubelle.txt' ;
  543. 'FIN' ;
  544. 'FINSI' ;
  545. 'FINSI' ;
  546.  
  547. 'MESSAGE' 'Tout s_est bien passé' ;
  548. 'FIN' ;
  549.  
  550.  
  551.  
  552.  
  553.  
  554.  

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