Télécharger test_kres_lapn.dgibi

Retour à la liste

Numérotation des lignes :

  1. * fichier : test_kres_lapn.dgibi
  2. ************************************************************************
  3. * NOM : TEST_KRES_LAPN
  4. * DESCRIPTION : Test d'options de KRES sur un laplacien
  5. * - Scaling ou pas
  6. * - solveurs itératifs : CG; BiCGStab, GMRES
  7. * - préconditionneurs ilu0, ilut, ilut + pivoting
  8. *
  9. * LANGAGE : GIBIANE-CAST3M
  10. * AUTEUR : Stéphane GOUNAND (CEA/DEN/DM2S/SFME/LTMF)
  11. * mél : gounand@semt2.smts.cea.fr
  12. **********************************************************************
  13. * VERSION : v1, 14/02/2005, version initiale
  14. * HISTORIQUE : v1, 14/02/2005, création
  15. * HISTORIQUE :
  16. * HISTORIQUE :
  17. ************************************************************************
  18. * Prière de PRENDRE LE TEMPS de compléter les commentaires
  19. * en cas de modification de ce sous-programme afin de faciliter
  20. * la maintenance !
  21. ************************************************************************
  22. *
  23. interact= FAUX ;
  24. complet = FAUX ;
  25. graph = FAUX ;
  26. ************************************************************************
  27. * NOM : MODULO
  28. * DESCRIPTION : Calcule un entier modulo un autre...
  29. *
  30. *
  31. 'DEBPROC' MODULO ;
  32. 'ARGUMENT' i*'ENTIER' j*'ENTIER' ;
  33. 'SI' ('EGA' j 0) ;
  34. 'MESSAGE' 'Impossible de faire modulo 0' ;
  35. 'ERREUR' 5 ;
  36. 'SINON' ;
  37. k=i '/' j ;
  38. mod=i '-' ( k '*'j ) ;
  39. 'RESPRO' mod ;
  40. 'FINSI' ;
  41. *
  42. * End of procedure file MODULO
  43. *
  44. 'FINPROC' ;
  45.  
  46. ************************************************************************
  47. * NOM : LOG10
  48. * DESCRIPTION : Log_10
  49. *
  50. 'DEBPROC' LOG10 ;
  51. 'ARGUMENT' fl/'FLOTTANT' ;
  52. 'ARGUMENT' lr/'LISTREEL' ;
  53. 'ARGUMENT' cp/'CHPOINT ' ;
  54. 'ARGUMENT' cm/'MCHAML ' ;
  55. 'SI' ('EXISTE' fl) ;
  56. 'RESPRO' ('/' ('LOG' fl) ('LOG' 10.D0)) ;
  57. 'FINSI' ;
  58. 'SI' ('EXISTE' lr) ;
  59. 'RESPRO' ('/' ('LOG' lr) ('LOG' 10.D0)) ;
  60. 'FINSI' ;
  61. 'SI' ('EXISTE' cp) ;
  62. 'RESPRO' ('/' ('LOG' cp) ('LOG' 10.D0)) ;
  63. 'FINSI' ;
  64. 'SI' ('EXISTE' cm) ;
  65. 'RESPRO' ('/' ('LOG' cm) ('LOG' 10.D0)) ;
  66. 'FINSI' ;
  67. *
  68. * End of procedure file LOG10
  69. *
  70. 'FINPROC' ;
  71.  
  72. ************************************************************************
  73. * NOM : FORMAR
  74. * DESCRIPTION : formate un réel de facon courte
  75. *
  76. 'DEBPROC' FORMAR ;
  77. 'ARGUMENT' fl*'FLOTTANT' ;
  78. 'ARGUMENT' vir/'ENTIER ' ;
  79. 'SI' ('NON' ('EXISTE' vir)) ;
  80. vir = 1 ;
  81. 'SINON' ;
  82. 'SI' ('<' vir 0) ;
  83. 'ERREUR' 'fournir un entier positif' ;
  84. 'FINSI' ;
  85. 'FINSI' ;
  86. 'SI' ('<' ('ABS' fl) 10.D-100) ;
  87. chfl = 'CHAINE' '0' ;
  88. 'SINON' ;
  89. *! sans le 1.D-10, ca ne fonctionne pas
  90. *! qd on entre pile poil une puissance de 10
  91. lfl = LOG10 ('ABS' fl) ;
  92. * lfl = '+' (LOG10 ('ABS' fl)) 1.D-10 ;
  93. slfl = 'SIGNE' ('ENTIER' lfl) ;
  94. 'SI' ('EGA' slfl 1) ;
  95. elfl = 'ENTIER' lfl ;
  96. 'SINON' ;
  97. elfl = '-' ('ENTIER' lfl) 1 ;
  98. 'FINSI' ;
  99. man = '/' fl ('**' 10.D0 elfl) ;
  100. *
  101. * Une verrue pour des histoires de précision...
  102. *
  103. 'SI' ('EGA' man 10.D0 ('**' 10.D0 ('*' vir -1.D0))) ;
  104. man = '/' man 10.D0 ;
  105. elfl = '+' elfl 1 ;
  106. 'FINSI' ;
  107. *
  108. sman = 'SIGNE' man ;
  109. 'SI' ('EGA' sman 1) ;
  110. fman = 'CHAINE' '(F' ('+' vir 2) '.0' vir ')' ;
  111. 'SINON' ;
  112. fman = 'CHAINE' '(F' ('+' vir 3) '.0' vir ')' ;
  113. 'FINSI' ;
  114. 'SI' ('NEG' vir 0) ;
  115. 'SI' ('NEG' elfl 0) ;
  116. chfl = 'CHAINE' 'FORMAT' fman man 'E' elfl ;
  117. 'SINON' ;
  118. chfl = 'CHAINE' 'FORMAT' fman man ;
  119. 'FINSI' ;
  120. 'SINON' ;
  121. man2 = 'ENTIER' ('+' man ('*' 0.5D0 sman)) ;
  122. 'SI' ('NEG' elfl 0) ;
  123. chfl = 'CHAINE' man2 'E' elfl ;
  124. 'SINON' ;
  125. chfl = 'CHAINE' man2 ;
  126. 'FINSI' ;
  127. 'FINSI' ;
  128. 'FINSI' ;
  129. 'RESPRO' chfl ;
  130. *
  131. * End of procedure file FORMAR
  132. *
  133. 'FINPROC' ;
  134.  
  135. *
  136. * Procédure CMAIL
  137. *
  138. 'DEBPROC' CMAIL ;
  139. 'ARGUMENT' imail*'ENTIER' ;
  140. 'ARGUMENT' iquad*'ENTIER' ;
  141. 'SI' ('EGA' imail 1) ;
  142. 'OPTION' 'DIME' 2 'ELEM' 'QUA4' ;
  143. 'SI' complet ;
  144. nmail = 75 ;
  145. 'SINON' ;
  146. nmail = 10 ;
  147. 'FINSI' ;
  148. 'SI' ('EGA' iquad 1) ;
  149. nmail = '*' nmail 2 ;
  150. 'FINSI' ;
  151. pA = 0. 0. ; pB = 1. 0. ;
  152. pC = 1. 1. ; pD = 0. 1. ;
  153. l1 = 'DROIT' nmail pA pB ;
  154. l2 = 'DROIT' nmail pB pC ;
  155. l3 = 'DROIT' nmail pC pD ;
  156. l4 = 'DROIT' nmail pD pA ;
  157. mt = 'DALLER' l1 l2 l3 l4 ;
  158. cmt = l1 'ET' l2 'ET' l3 'ET' l4 ;
  159. _mt = 'CHANGER' mt 'QUAF' ;
  160. _cmt = 'CHANGER' cmt 'QUAF' ;
  161. 'FINSI' ;
  162. 'SI' ('EGA' imail 2) ;
  163. 'OPTION' 'DIME' 2 'ELEM' 'TRI3' ;
  164. dini = 0.2D0 ;
  165. dfin = 0.004D0 ;
  166. 'SI' complet ;
  167. dini = '/' dini 4.D0 ;
  168. dfin = '/' dfin 4.D0 ;
  169. 'FINSI' ;
  170. 'SI' ('EGA' iquad 1) ;
  171. dini = '/' dini 2.D0 ;
  172. dfin = '/' dfin 2.D0 ;
  173. 'FINSI' ;
  174. pA = 0. 0. ; pB = 1. 0. ;
  175. pC = 1. 1. ; pD = 0. 1. ;
  176. l1 = 'DROIT' pA pB 'DINI' dini 'DFIN' dini ;
  177. l2 = 'DROIT' pB pC 'DINI' dini 'DFIN' dfin ;
  178. l3 = 'DROIT' pC pD 'DINI' dfin 'DFIN' dini ;
  179. l4 = 'DROIT' pD pA 'DINI' dini 'DFIN' dini ;
  180. cmt = l1 'ET' l2 'ET' l3 'ET' l4 ;
  181. mt = 'SURFACE' cmt ;
  182. _mt = 'CHANGER' mt 'QUAF' ;
  183. _cmt = 'CHANGER' cmt 'QUAF' ;
  184. 'FINSI' ;
  185. 'SI' ('EGA' imail 3) ;
  186. 'OPTION' 'DIME' 3 'ELEM' 'CUB8' ;
  187. 'SI' complet ;
  188. nmail = 15 ;
  189. 'SINON' ;
  190. nmail = 5 ;
  191. 'FINSI' ;
  192. 'SI' ('EGA' iquad 1) ;
  193. nmail = '*' nmail 2 ;
  194. 'FINSI' ;
  195. pA = 0. 0. 0. ; pB = 1. 0. 0. ;
  196. pC = 1. 1. 0. ; pD = 0. 1. 0. ;
  197. l1 = 'DROIT' nmail pA pB ;
  198. l2 = 'DROIT' nmail pB pC ;
  199. l3 = 'DROIT' nmail pC pD ;
  200. l4 = 'DROIT' nmail pD pA ;
  201. smt = 'DALLER' l1 l2 l3 l4 'PLAN' ;
  202. mt = 'VOLUME' nmail smt 'TRAN' (0. 0. 1.) ;
  203. cmt = 'ENVELOPPE' mt ;
  204. _mt = 'CHANGER' mt 'QUAF' ;
  205. _cmt = 'CHANGER' cmt 'QUAF' ;
  206. 'FINSI' ;
  207. 'ELIMINATION' ('ET' _mt _cmt) 1.D-5 ;
  208. 'SI' graph ;
  209. mtot = 'ET' _mt ('COULEUR' _cmt 'ROUG') ;
  210. titch = 'CHAINE' 'Maillage et contour nbel =' ('NBEL' mtot) ;
  211. 'TRACER' 'CACH' mtot 'TITR' titch ;
  212. 'FINSI' ;
  213. tabmod = 'TABLE' ;
  214. tabmod . 'MT' = 'TABLE' ;
  215. tabmod . 'MT' . 'QUAF' = _mt ;
  216. tabmod . 'MT' . 'LINE' = mt ;
  217. tabmod . 'CMT' = 'TABLE' ;
  218. tabmod . 'CMT' . 'QUAF' = _cmt ;
  219. tabmod . 'CMT' . 'LINE' = cmt ;
  220. 'RESPRO' tabmod ;
  221. 'FINPROC' ;
  222. *
  223. 'DEBPROC' RESOU ;
  224. 'ARGUMENT' mat/'RIGIDITE' ;
  225. 'SI' ('NON' ('EXISTE' mat) );
  226. 'ARGUMENT' mat*'MATRIK' ;
  227. 'SINON' ;
  228. lispri = 'MOTS' 'T' ;
  229. lisdua = 'MOTS' 'Q' ;
  230. mat = 'KOPS' 'CHANINCO' lispri lispri lisdua lispri
  231. ('KOPS' 'RIMA' mat) ;
  232. 'FINSI' ;
  233. 'ARGUMENT' clim*'CHPOINT' ;
  234. 'ARGUMENT' pcmlag*'ENTIER' ;
  235. 'ARGUMENT' scaling*'ENTIER' ;
  236. 'ARGUMENT' typinv*'ENTIER' ;
  237. 'ARGUMENT' typrec*'ENTIER' ;
  238. *
  239. rv = 'EQEX' ;
  240. rvc = 'TABLE' ;
  241. rvm = rv . 'METHINV' ;
  242. rvm . 'TYPINV' = typinv '+' 1 ;
  243. rvm . 'IMPINV' = 0 ;
  244. rvm . 'PCMLAG' = 'CHAINE' 'APR' ('+' pcmlag 1) ;
  245. rvm . 'SCALING' = '-' scaling 1 ;
  246. 'SI' ('EGA' typrec 1) ;
  247. rvm . 'PRECOND' = 3 ;
  248. 'SINON' ;
  249. 'SI' ('ET' ('>EG' typrec 2) ('&lt;EG' typrec 8)) ;
  250. rvm . 'ILUTLFIL' = 2. ;
  251. 'SINON' ;
  252. rvm . 'ILUTLFIL' = 3. ;
  253. 'FINSI' ;
  254. typrec = '-' typrec 1 ;
  255. typrec = MODULO ('-' typrec 1) 7 ;
  256. 'SI' ('EGA' typrec 0) ;
  257. rvm . 'PRECOND' = 5 ;
  258. 'SINON' ;
  259. typrec = '-' typrec 1 ;
  260. typrec2 = modulo typrec 2 ;
  261. rvm . 'PRECOND' = '+' typrec2 7 ;
  262. typrec3 = '/' typrec 3 ;
  263. 'SI' ('EGA' typrec3 0) ;
  264. rvm . 'ILUTPPIV' = 0. ;
  265. 'FINSI' ;
  266. 'SI' ('EGA' typrec3 1) ;
  267. rvm . 'ILUTPPIV' = 1.D-3 ;
  268. 'FINSI' ;
  269. 'SI' ('EGA' typrec3 2) ;
  270. rvm . 'ILUTPPIV' = 1.D-1 ;
  271. 'FINSI' ;
  272. 'FINSI' ;
  273. 'FINSI' ;
  274. *
  275. 'TEMPS' 'ZERO' ;
  276. solu = 'KRES' mat 'TYPI' rvm 'CLIM' clim 'IMPR' 0 ;
  277. TABTPS = TEMP 'NOEC';
  278. tcpu = TABTPS.'TEMPS_CPU'.'INITIAL';
  279. nit = 'DIME' (rvm . 'CONVINV') ;
  280. nit = 'EXTRAIRE' (rvm . 'NMATVEC') nit ;
  281. *
  282. 'RESPRO' solu nit tcpu ;
  283. 'FINPROC' ;
  284. *
  285. * Procédure CALCUL
  286. *
  287. 'DEBPROC' CALCUL ;
  288. 'ARGUMENT' imail*'ENTIER' ;
  289. 'ARGUMENT' iquad*'ENTIER' ;
  290. 'ARGUMENT' pcmlag*'ENTIER' ;
  291. 'ARGUMENT' scaling*'ENTIER' ;
  292. 'ARGUMENT' typinv*'ENTIER' ;
  293. 'ARGUMENT' typrec*'ENTIER' ;
  294. tabmod = CMAIL imail iquad ;
  295. 'SI' ('EGA' iquad 2) ;
  296. disc = 'QUAF' ;
  297. 'SINON' ;
  298. disc = 'LINE' ;
  299. 'FINSI' ;
  300. *
  301. * Calcul
  302. *
  303. * Modèle
  304. * Solution exacte
  305. mt = tabmod . 'MT' . disc ;
  306. cmt = tabmod . 'CMT' . disc ;
  307. dim3 = 'EGA' ('VALEUR' 'DIME') 3 ;
  308. xmt = 'COORDONNEE' 1 mt ;
  309. ymt = 'COORDONNEE' 2 mt ;
  310. 'SI' dim3 ;
  311. zmt = 'COORDONNEE' 3 mt ;
  312. 'FINSI' ;
  313. solex = '+' ('*' xmt ('**' 2.D0 0.5D0))
  314. ('*' ymt PI) ;
  315. 'SI' dim3 ;
  316. solex = '+' solex ('*' zmt ('**' PI 0.5D0)) ;
  317. 'FINSI' ;
  318. solex = 'NOMC' 'T' solex 'NATURE' 'DISCRET' ;
  319. * Conditions aux limites
  320. clim = 'REDU' solex cmt ;
  321. _mt = tabmod . 'MT' . 'QUAF' ;
  322. * Laplacien
  323. mlap = GLAPN _mt 'LINE' 'T' disc 'Q' disc 1.D0 ;
  324. * Masse
  325. mmas = GMASS _mt 'LINE' 'T' disc 'Q' disc 1.D0 ;
  326. * Resolution
  327. solu nit tcpu = RESOU mlap clim pcmlag scaling typinv typrec ;
  328. resit = '**' (XTMX ('-' solex solu) mmas) 0.5D0 ;
  329. 'RESPRO' nit tcpu resit ;
  330. 'FINPROC' ;
  331. *
  332. * Procédure CALCUL2 (opérateurs JPM)
  333. *
  334. 'DEBPROC' CALCUL2 ;
  335. 'ARGUMENT' imail*'ENTIER' ;
  336. 'ARGUMENT' iquad*'ENTIER' ;
  337. 'ARGUMENT' pcmlag*'ENTIER' ;
  338. 'ARGUMENT' scaling*'ENTIER' ;
  339. 'ARGUMENT' typinv*'ENTIER' ;
  340. 'ARGUMENT' typrec*'ENTIER' ;
  341. tabmod = CMAIL imail iquad ;
  342. 'SI' ('EGA' iquad 2) ;
  343. disc = 'QUAF' ;
  344. 'SINON' ;
  345. disc = 'LINE' ;
  346. 'FINSI' ;
  347. *
  348. * Calcul
  349. *
  350. * Modèle
  351. * Solution exacte
  352. mt = tabmod . 'MT' . disc ;
  353. cmt = tabmod . 'CMT' . disc ;
  354. dim3 = 'EGA' ('VALEUR' 'DIME') 3 ;
  355. xmt = 'COORDONNEE' 1 mt ;
  356. ymt = 'COORDONNEE' 2 mt ;
  357. 'SI' dim3 ;
  358. zmt = 'COORDONNEE' 3 mt ;
  359. 'FINSI' ;
  360. solex = '+' ('*' xmt ('**' 2.D0 0.5D0))
  361. ('*' ymt PI) ;
  362. 'SI' dim3 ;
  363. solex = '+' solex ('*' zmt ('**' PI 0.5D0)) ;
  364. 'FINSI' ;
  365. solex = 'NOMC' 'T' solex 'NATURE' 'DISCRET' ;
  366. * Conditions aux limites
  367. clim = 'REDU' solex cmt ;
  368. _mt = tabmod . 'MT' . 'QUAF' ;
  369. $mt = 'MODELISER' _mt 'NAVIER_STOKES' disc ;
  370. rv = 'EQEX'
  371. 'OPTI' 'EF' 'IMPL'
  372. 'ZONE' $mt 'OPER' 'LAPN' 1.D0 'INCO' 'T'
  373. 'OPTI' 'EF' 'IMPL'
  374. 'ZONE' $mt 'OPER' 'MDIA' 1.D0 'INCO' 'T'
  375. ;
  376. rv . 'INCO' = 'TABLE' 'INCO' ;
  377. rv . 'INCO' . 'T' = 'KCHT' $mt 'SCAL' 'SOMMET' 1.D0 ;
  378. toto mlap = 'LAPN' (rv . '1LAPN') ;
  379. toto mmas = 'MDIA' (rv . '2MDIA') ;
  380. * Resolution
  381. solu nit tcpu = RESOU mlap clim pcmlag scaling typinv typrec ;
  382. *resit = '**' (XTMX ('-' solex solu) mmas) 0.5D0 ;
  383. dif = '-' solex solu ;
  384. linc = 'MOTS' 'T' ;
  385. resit = '**' ('XTY' ('KOPS' mmas '*' dif) dif linc linc) 0.5D0 ;
  386. 'RESPRO' nit tcpu resit ;
  387. 'FINPROC' ;
  388. *
  389. * Test scaling
  390. *
  391. 'DEBPROC' TESTSCAL ;
  392. tabres = 'TABLE' ;
  393. pcmlag = 1 ;
  394. typinv = 2 ;
  395. typrec = 1 ;
  396. 'REPETER' iscal 3 ;
  397. tabres . &iscal = 'TABLE' ;
  398. scaling = &iscal ;
  399. 'REPETER' imail 3 ;
  400. tabres . &iscal . &imail = 'TABLE' ;
  401. 'REPETER' iquad 2 ;
  402. tabres . &iscal . &imail . &iquad = 'TABLE' ;
  403. nit tcpu resit = CALCUL2 &imail &iquad pcmlag
  404. scaling typinv typrec ;
  405. tabres . &iscal . &imail . &iquad . 'NIT' = nit ;
  406. tabres . &iscal . &imail . &iquad . 'TCPU' = tcpu ;
  407. tabres . &iscal . &imail . &iquad . 'RESIT' = resit ;
  408. 'FIN' iquad ;
  409. 'FIN' imail ;
  410. 'FIN' iscal ;
  411. *
  412. optecho = 'VALEUR' 'ECHO' ;
  413. 'OPTION' 'ECHO' 0 ;
  414. cc1 = 'CHAINE' 'Paramètres fixés : '
  415. (rcvt . 'TYPINV' . typinv) ' '
  416. (rcvt . 'TYPREC' . typrec) ;
  417. 'MESSAGE' cc1 ;
  418. 'REPETER' imail 3 ;
  419. 'REPETER' iquad 2 ;
  420. 'SAUTER' 'LIGNE' ;
  421. cc2 = 'CHAINE' (rcvt . 'MAILLAGE' . &imail) ' '
  422. (rcvt . 'ELEM' . &iquad) ;
  423. 'MESSAGE' cc2 ;
  424. cc3 = 'CHAINE' 'Nb iter.'/20 'Temps CPU'/40 'Residu reel'/60 ;
  425. 'MESSAGE' cc3 ;
  426. 'REPETER' iscal 3 ;
  427. nit = tabres . &iscal . &imail . &iquad . 'NIT' ;
  428. tcpu = tabres . &iscal . &imail . &iquad . 'TCPU' ;
  429. resit = tabres . &iscal . &imail . &iquad . 'RESIT' ;
  430. cc4 = 'CHAINE' 'SCALING : '
  431. (rcvt . 'SCALING' . &iscal)
  432. nit*30 (formar tcpu 2)*50
  433. (formar resit 4)*70 ;
  434. 'MESSAGE' cc4 ;
  435. 'FIN' iscal ;
  436. 'FIN' iquad ;
  437. 'FIN' imail ;
  438. 'OPTION' 'ECHO' optecho ;
  439. 'FINPROC' ;
  440. *
  441. * Test solveur
  442. *
  443. 'DEBPROC' TESTSOLV ;
  444. tabres = 'TABLE' ;
  445. pcmlag = 1 ;
  446. scaling = 1 ;
  447. typrec = 1 ;
  448. 'REPETER' itypinv 5 ;
  449. tabres . &itypinv = 'TABLE' ;
  450. typinv = &itypinv ;
  451. 'REPETER' imail 3 ;
  452. tabres . &itypinv . &imail = 'TABLE' ;
  453. 'REPETER' iquad 2 ;
  454. tabres . &itypinv . &imail . &iquad = 'TABLE' ;
  455. nit tcpu resit = CALCUL2 &imail &iquad pcmlag
  456. scaling typinv typrec ;
  457. tabres . &itypinv . &imail . &iquad . 'NIT' = nit ;
  458. tabres . &itypinv . &imail . &iquad . 'TCPU' = tcpu ;
  459. tabres . &itypinv . &imail . &iquad . 'RESIT' = resit ;
  460. 'FIN' iquad ;
  461. 'FIN' imail ;
  462. 'FIN' itypinv ;
  463. *
  464. optecho = 'VALEUR' 'ECHO' ;
  465. 'OPTION' 'ECHO' 0 ;
  466. cc1 = 'CHAINE' 'Paramètres fixés : '
  467. (rcvt . 'SCALING' . scaling) ' '
  468. (rcvt . 'TYPREC' . typrec) ;
  469. 'MESSAGE' cc1 ;
  470. 'REPETER' imail 3 ;
  471. 'REPETER' iquad 2 ;
  472. 'SAUTER' 'LIGNE' ;
  473. cc2 = 'CHAINE' (rcvt . 'MAILLAGE' . &imail) ' '
  474. (rcvt . 'ELEM' . &iquad) ;
  475. 'MESSAGE' cc2 ;
  476. cc3 = 'CHAINE' 'Nb iter.'/20 'Temps CPU'/40 'Residu reel'/60 ;
  477. 'MESSAGE' cc3 ;
  478. 'REPETER' itypinv 5 ;
  479. nit = tabres . &itypinv . &imail . &iquad . 'NIT' ;
  480. tcpu = tabres . &itypinv . &imail . &iquad . 'TCPU' ;
  481. resit = tabres . &itypinv . &imail . &iquad . 'RESIT' ;
  482. cc4 = 'CHAINE' 'Solveur : '
  483. (rcvt . 'TYPINV' . &itypinv)
  484. nit*30 (formar tcpu 2)*50
  485. (formar resit 4)*70 ;
  486. 'MESSAGE' cc4 ;
  487. 'FIN' itypinv ;
  488. 'FIN' iquad ;
  489. 'FIN' imail ;
  490. 'OPTION' 'ECHO' optecho ;
  491. 'FINPROC' ;
  492. *
  493. * Test preconditionneur
  494. *
  495. 'DEBPROC' TESTPREC ;
  496. tabres = 'TABLE' ;
  497. pcmlag = 1 ;
  498. scaling = 1 ;
  499. typinv = 2 ;
  500. 'REPETER' ityprec 15 ;
  501. tabres . &ityprec = 'TABLE' ;
  502. typrec = &ityprec ;
  503. 'REPETER' imail 3 ;
  504. tabres . &ityprec . &imail = 'TABLE' ;
  505. 'REPETER' iquad 2 ;
  506. * 'MESSAGE' ('CHAINE' 'typrec=' typrec) ;
  507. tabres . &ityprec . &imail . &iquad = 'TABLE' ;
  508. nit tcpu resit = CALCUL2 &imail &iquad pcmlag
  509. scaling typinv typrec ;
  510. *
  511. * Correction pour le gradient conjugué
  512. *
  513. tabres . &ityprec . &imail . &iquad . 'NIT' = nit ;
  514. tabres . &ityprec . &imail . &iquad . 'TCPU' = tcpu ;
  515. tabres . &ityprec . &imail . &iquad . 'RESIT' = resit ;
  516. 'FIN' iquad ;
  517. 'FIN' imail ;
  518. 'FIN' ityprec ;
  519. *
  520. optecho = 'VALEUR' 'ECHO' ;
  521. 'OPTION' 'ECHO' 0 ;
  522. cc1 = 'CHAINE' 'Paramètres fixés : '
  523. (rcvt . 'SCALING' . scaling) ' '
  524. (rcvt . 'TYPINV' . typinv) ;
  525. 'MESSAGE' cc1 ;
  526. 'REPETER' imail 3 ;
  527. 'REPETER' iquad 2 ;
  528. 'SAUTER' 'LIGNE' ;
  529. cc2 = 'CHAINE' (rcvt . 'MAILLAGE' . &imail) ' '
  530. (rcvt . 'ELEM' . &iquad) ;
  531. 'MESSAGE' cc2 ;
  532. cc3 = 'CHAINE' 'Nb iter.'/20 'Temps CPU'/40 'Residu reel'/60 ;
  533. 'MESSAGE' cc3 ;
  534. 'REPETER' ityprec 15 ;
  535. nit = tabres . &ityprec . &imail . &iquad . 'NIT' ;
  536. tcpu = tabres . &ityprec . &imail . &iquad . 'TCPU' ;
  537. resit = tabres . &ityprec . &imail . &iquad . 'RESIT' ;
  538. cc4 = 'CHAINE' 'Précond. : '
  539. (rcvt . 'TYPREC' . &ityprec)
  540. nit*30 (formar tcpu 2)*50
  541. (formar resit 4)*70 ;
  542. 'MESSAGE' cc4 ;
  543. 'FIN' ityprec ;
  544. 'FIN' iquad ;
  545. 'FIN' imail ;
  546. 'OPTION' 'ECHO' optecho ;
  547. 'FINPROC' ;
  548. *
  549. *
  550. * Fin des procédures
  551. * Début du programme
  552. *
  553. *
  554. *
  555. * Table de conversion texte
  556. *
  557. rcvt = 'TABLE' ;
  558. rcvt . 'MAILLAGE' = 'TABLE' ;
  559. rcvt . 'MAILLAGE' . 1 = 'CHAINE' 'Carre reg ' ;
  560. rcvt . 'MAILLAGE' . 2 = 'CHAINE' 'Tri. irreg' ;
  561. rcvt . 'MAILLAGE' . 3 = 'CHAINE' 'Cube reg' ;
  562. rcvt . 'ELEM' = 'TABLE' ;
  563. rcvt . 'ELEM' . 1 = 'CHAINE' 'Line.' ;
  564. rcvt . 'ELEM' . 2 = 'CHAINE' 'Quad.' ;
  565. rcvt . 'SCALING' = 'TABLE' ;
  566. rcvt . 'SCALING' . 1 = 'CHAINE' 'No scal.' ;
  567. rcvt . 'SCALING' . 2 = 'CHAINE' 'L2 Scaling ' ;
  568. rcvt . 'SCALING' . 3 = 'CHAINE' 'L1 Scaling ' ;
  569. rcvt . 'TYPINV' = 'TABLE' ;
  570. rcvt . 'TYPINV' . 1 = 'CHAINE' 'Conj. Grad.' ;
  571. rcvt . 'TYPINV' . 2 = 'CHAINE' 'BiCGStab ' ;
  572. rcvt . 'TYPINV' . 3 = 'CHAINE' 'BiCGStab(4)' ;
  573. rcvt . 'TYPINV' . 4 = 'CHAINE' 'GMRES(50) ' ;
  574. rcvt . 'TYPINV' . 5 = 'CHAINE' 'CGS ' ;
  575. rcvt . 'TYPREC' = 'TABLE' ;
  576. rcvt . 'TYPREC' . 1 = 'CHAINE' 'ILU(0) ' ;
  577. rcvt . 'TYPREC' . 2 = 'CHAINE' 'ILUT(2) ' ;
  578. rcvt . 'TYPREC' . 3 = 'CHAINE' 'ILUTP(2,0) ' ;
  579. rcvt . 'TYPREC' . 4 = 'CHAINE' 'ILUTP+0(2,0) ' ;
  580. rcvt . 'TYPREC' . 5 = 'CHAINE' 'ILUTP(2,-3) ' ;
  581. rcvt . 'TYPREC' . 6 = 'CHAINE' 'ILUTP+0(2,-3)' ;
  582. rcvt . 'TYPREC' . 7 = 'CHAINE' 'ILUTP(2,-1) ' ;
  583. rcvt . 'TYPREC' . 8 = 'CHAINE' 'ILUTP+0(2,-1)' ;
  584. rcvt . 'TYPREC' . 9 = 'CHAINE' 'ILUT(3) ' ;
  585. rcvt . 'TYPREC' . 10 = 'CHAINE' 'ILUTP(3,0) ' ;
  586. rcvt . 'TYPREC' . 11 = 'CHAINE' 'ILUTP+0(3,0) ' ;
  587. rcvt . 'TYPREC' . 12 = 'CHAINE' 'ILUTP(3,-3) ' ;
  588. rcvt . 'TYPREC' . 13 = 'CHAINE' 'ILUTP+0(3,-3)' ;
  589. rcvt . 'TYPREC' . 14 = 'CHAINE' 'ILUTP(3,-1) ' ;
  590. rcvt . 'TYPREC' . 15 = 'CHAINE' 'ILUTP+0(3,-1)' ;
  591. *
  592. TESTSCAL ;
  593. 'SAUTER' 10 'LIGNE' ;
  594. TESTSOLV ;
  595. 'SAUTER' 10 'LIGNE' ;
  596. TESTPREC ;
  597. 'SI' interact ;
  598. 'OPTION' 'DONN' 5 ;
  599. 'FINSI' ;
  600. *
  601. * End of dgibi file TEST_KRES
  602. *
  603. 'FIN' ;
  604.  
  605.  
  606.  
  607.  
  608.  
  609.  
  610.  
  611.  

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