Télécharger kres5.eso

Retour à la liste

Numérotation des lignes :

kres5
  1. C KRES5 SOURCE GOUNAND 25/03/24 21:15:05 12216
  2. SUBROUTINE KRES5(MATRIK,KCLIM,KSMBR,KTYPI,
  3. $ MCHINI,ITER,RESID,
  4. $ BRTOL,IRSTRT,LBCG,ICALRS,
  5. $ MAPREC,KPREC,
  6. $ RXMILU,RXILUP,XLFIL,XDTOL,XSPIV,
  7. $ ISCAL,
  8. $ KTIME,LTIME,IDDOT,IMVEC,IORINC,
  9. $ MCHSOL,LRES,LNMV,ICVG,
  10. $ IMPR,IRET)
  11. IMPLICIT REAL*8 (A-H,O-Z)
  12. IMPLICIT INTEGER (I-N)
  13. C***********************************************************************
  14. C NOM : KRES5
  15. C DESCRIPTION : Résolution d'un système par une méthode itérative
  16. C (type gradient conjugué).
  17. C
  18. C
  19. C LANGAGE : ESOPE
  20. C AUTEUR : Stéphane GOUNAND (CEA/DRN/DMT/SEMT/LTMF)
  21. C mél : gounand@semt2.smts.cea.fr
  22. C***********************************************************************
  23. C APPELES : VERMAT, MEXINI, MESMBR, MELIM, INFMAT,
  24. C MEIDIA, PRDILU, PRILU, PRMILU, PRILUT,
  25. C PRCG, PRBCGS, PRBCS2, PRGMR,
  26. C XDISP
  27. C APPELE PAR : KRES2
  28. C***********************************************************************
  29. C ENTREES : MATRIK, KCLIM, KSMBR, KTYPI, MCHINI, ITER, RESID,
  30. C BRTOL, IRSTRT, MAPREC, KPREC, RXMILU, XLFIL,
  31. C XDTOL
  32. C ENTREES/SORTIES : -
  33. C SORTIES : MCHSOL, LRES
  34. C CODE RETOUR (IRET) : = 0 si tout s'est bien passé
  35. C***********************************************************************
  36. C VERSION : v1, 14/04/2000, version initiale
  37. C HISTORIQUE : v1, 14/04/2000, création
  38. C HISTORIQUE : 06/04/04 : Scaling
  39. C HISTORIQUE : 08/04/04 : ajout ILUTP
  40. C HISTORIQUE :
  41. C***********************************************************************
  42. C Prière de PRENDRE LE TEMPS de compléter les commentaires
  43. C en cas de modification de ce sous-programme afin de faciliter
  44. C la maintenance !
  45. C***********************************************************************
  46.  
  47. -INC PPARAM
  48. -INC CCOPTIO
  49. -INC SMLENTI
  50. POINTEUR LNMV.MLENTI
  51. POINTEUR ATYP.MLENTI
  52. POINTEUR ANOD.MLENTI
  53. -INC SMLREEL
  54. POINTEUR LRES.MLREEL
  55. -INC SMCHPOI
  56. POINTEUR KCLIM.MCHPOI
  57. POINTEUR KSMBR.MCHPOI
  58. POINTEUR MCHINI.MCHPOI
  59. POINTEUR MCHSOL.MCHPOI
  60. POINTEUR MAPREC.MATRIK
  61. POINTEUR INCX.IZA
  62. POINTEUR KS2B.IZA
  63. POINTEUR KMORS.PMORS
  64. POINTEUR KISA.IZA
  65. POINTEUR AMORS.PMORS
  66. POINTEUR AISA.IZA
  67. POINTEUR NORMP.IZA
  68. POINTEUR NORMD.IZA
  69. *
  70. * .. Parameters
  71. REAL*8 ZERO ,ONE
  72. PARAMETER (ZERO=0.0D0,ONE=1.0D0)
  73. * Paramètre m du GMRES(m)
  74. INTEGER IRSTRT
  75. INTEGER ITER,IVARI
  76. REAL*8 BRTOL,RESID,RXMILU,RXILUP
  77. *
  78. REAL*8 XLFIL,XDTOL
  79. INTEGER KPREC
  80. INTEGER KTYPI
  81. C LOGICAL LCLZER
  82. LOGICAL LRACOU
  83. LOGICAL LFIRST
  84. C
  85. INTEGER IMPR,IRET
  86. C Tableau de correspondance : entier <-> mot
  87. C pour le type d'inversion
  88. INTEGER NBTYPI,LNTYPI
  89. PARAMETER (NBTYPI=11)
  90. PARAMETER (LNTYPI=18)
  91. CHARACTER*(LNTYPI) LTYPI(NBTYPI)
  92. C Tableau de correspondance : entier <-> mot
  93. C pour le type de préconditionnement (cas d'une méthode itérative)
  94. INTEGER NBPREC,LNPREC
  95. PARAMETER (NBPREC=11)
  96. PARAMETER (LNPREC=8)
  97. CHARACTER*(LNPREC) LPREC(NBPREC)
  98. -INC SMTABLE
  99. POINTEUR KTIME.MTABLE
  100. DIMENSION ITTIME(4)
  101. CHARACTER*8 CHARI
  102. CHARACTER*1 CCOMP
  103. LOGICAL LTIME,LOGII
  104. C
  105. C Initialisation des tableaux
  106. C
  107. DATA LTYPI/'Choleski ',
  108. $ 'Conjugate Gradient',
  109. $ 'BiCGSTAB G ',
  110. $ 'BiCGSTAB(l) G ',
  111. $ 'GMRES(m) ',
  112. $ 'CGS-Neumaier ',
  113. $ 'Multigrid FCG ',
  114. $ 'Multigrid GCR(m) ',
  115. $ 'Bi-CG ',
  116. $ 'AGMG FCG Stokes ',
  117. $ 'AGMG GCR(m) NS '/
  118. DATA LPREC/'None ',
  119. $ 'Jacobi ',
  120. $ 'D-ILU ',
  121. $ 'ILU(0) ',
  122. $ 'MILU(0) ',
  123. $ 'ILUT ',
  124. $ 'ILUT2 ',
  125. $ 'ILUTP ',
  126. $ 'ILUTP+0 ',
  127. $ 'ILU0-PV ',
  128. $ 'ILU0-PVf'/
  129.  
  130. IVALI=0
  131. XVALI=0.D0
  132. LOGII=.FALSE.
  133. IRETI=0
  134. XVALR=0.D0
  135. IOBRE=0
  136. IRETR=0
  137. *
  138. * Executable statements
  139. *
  140. ICVG=0
  141. IF (IMPR.GT.5) WRITE(IOIMP,*) 'Entrée dans kres5.eso'
  142. C Batterie de tests
  143. C
  144. C KTYPI : Méthode d'inversion du système (cf. LTYPI)
  145. C 1 : résolution directe (Choleski)
  146. C 2 : Gradient Conjugué
  147. C 3 : Bi-Gradient Conjugué Stabilisé (BiCGSTAB)
  148. C 4 : BiCGSTAB(2)
  149. C 5 : GMRES(m)
  150. C 6 : CGS
  151. C 7 : Multigrille symétrique
  152. C 8 : Multigrille non symétrique
  153. C 9 : BiCG
  154. IF (KTYPI.LT.2.OR.KTYPI.GT.NBTYPI) THEN
  155. WRITE(IOIMP,*) 'KTYPI incorrect (2..',NBTYPI,') :',KTYPI
  156. GOTO 9999
  157. ENDIF
  158. * Pas de préconditionnement en cas de multigrille
  159. C - Type de préconditionnement : (cf. LPREC)
  160. C 0 : pas de préconditionnement
  161. C 1 : préconditionnement par la diagonale
  162. C 2 : préconditionnement D-ILU
  163. C 3 : préconditionnement ILU(0) (Choleski)
  164. C 4 : préconditionnement MILU(0) (Choleski modifié)
  165. C 5 : préconditionnement ILUT (dual truncation)
  166. C 6 : préconditionnement ILUT2 (une variante du
  167. C précédent qui remplit mieux la mémoire et
  168. C fonctionne mieux quelquefois)
  169. C 7 : préconditionnement ILUTP
  170. C 8 : préconditionnement ILUTP version Goon
  171. C 9 : préconditionnement ILU-PV
  172. C 10 : préconditionnement ILU-PVf
  173. IF (KPREC.LT.0.OR.KPREC.GT.NBPREC-1) THEN
  174. WRITE(IOIMP,*) 'PRECOND ',
  175. $ 'incorrect (0..',NBPREC-1,') :',KPREC
  176. GOTO 9999
  177. ENDIF
  178. C
  179. C On vérifie que la matrice est correctement assemblée
  180. C
  181. CALL VERMAT(MATRIK,IMPR,IRET)
  182. IF (IRET.NE.0) GOTO 9999
  183. C
  184. C On transforme le chpoint estimation initiale de l'inconnue en vecteur
  185. C estimation initiale de l'inconnue
  186. C In MEXINI : SEGINI INCX
  187. CALL MEXINI(MATRIK,MCHINI,
  188. $ INCX,
  189. $ IMPR,IRET)
  190. IF (IRET.NE.0) GOTO 9999
  191. C
  192. C On transforme le chpoint second membre en vecteur second membre
  193. C In MESMBR : SEGINI KS2B
  194. CALL MESMBR(MATRIK,KSMBR,
  195. $ KS2B,
  196. $ IMPR,IRET)
  197. IF (IRET.NE.0) GOTO 9999
  198. C
  199. C On applique les conditions aux limites
  200. C
  201. LRACOU=(KCLIM.EQ.0)
  202. * LRACOU=.FALSE.
  203. * WRITE(IOIMP,*) 'LRACOU=',LRACOU
  204. IF (LRACOU) THEN
  205. SEGACT MATRIK
  206. AMORS=MATRIK.KIDMAT(4)
  207. AISA=MATRIK.KIDMAT(5)
  208. SEGDES MATRIK
  209. ELSE
  210. C
  211. C In MELIM : SEGINI AMORS
  212. C SEGINI AISA
  213. CALL MELIM(MATRIK,KCLIM,
  214. $ INCX,KS2B,
  215. $ AMORS,AISA,
  216. $ IMPR,IRET)
  217. IF (IRET.NE.0) GOTO 9999
  218. ENDIF
  219. C
  220. LFIRST=.FALSE.
  221. IF (ISCAL.GT.0) THEN
  222. IF (LRACOU) THEN
  223. SEGACT MATRIK
  224. NORMP=MATRIK.KKMMT(4)
  225. NORMD=MATRIK.KKMMT(5)
  226. SEGDES MATRIK
  227. IF (NORMP.EQ.0.AND.NORMD.EQ.0) THEN
  228. LFIRST=.TRUE.
  229. ENDIF
  230. ELSE
  231. NORMP=0
  232. NORMD=0
  233. ENDIF
  234. IF (NORMP.EQ.0.OR.NORMD.EQ.0) THEN
  235. C Calcul des normes primales (colonnes) et duales (lignes)
  236. C de la matrice. Norme = norme L2, soit :
  237. C {\sum_{i ou j} a_{ij}^2}^{1/2}
  238. C
  239. CALL NORMAT(AMORS,AISA,ISCAL,NORMP,NORMD,IMPR,IRET)
  240. IF (IRET.NE.0) GOTO 9999
  241. C
  242. C On norme la matrice : attention modification...
  243. C
  244. CALL NORMAM(AMORS,AISA,NORMP,NORMD,IMPR,IRET)
  245. IF (IRET.NE.0) GOTO 9999
  246. ENDIF
  247. C
  248. C On norme le second membre : attention modification...
  249. C
  250. CALL NORMV2(KS2B,NORMD,IMPR,IRET)
  251. IF (IRET.NE.0) GOTO 9999
  252. C
  253. C On dénorme l'inconnue X : attention modification...
  254. C
  255. CALL NORMV1(INCX,NORMP,IMPR,IRET)
  256. IF (IRET.NE.0) GOTO 9999
  257. ELSE
  258. NORMP=0
  259. NORMD=0
  260. ENDIF
  261. C
  262. C On donne des infos sur la matrice
  263. C
  264. CALL INFMAT(AMORS,AISA,IMPR,IRET)
  265. IF (IRET.NE.0) GOTO 9999
  266. C
  267. C ON CHOISIT LE TYPE D'INVERSION
  268. C
  269. C Méthodes de gradient conjugué
  270. C
  271. C KTYPI = 2 : Gradient conjugué (CG)
  272. C KTYPI = 3 : Bi-Gradient conjugué stabilisé (BiCGSTAB)
  273. C KTYPI = 4 : BiCGSTAB(2)
  274. C KTYPI = 5 : GMRES(m)
  275. C KTYPI = 6 : CGS
  276. C KTYPI = 7 : Algebraic Multigrid FCG
  277. C KTYPI = 8 : Algebraic Multigrid GCR(m)
  278. C KTYPI = 9 : BiCG
  279. C KTYPI = 10 : Algebraic Multigrid FCG Stokes
  280. C KTYPI = 11 : Algebraic Multigrid GCR(m) Navier-Stokes
  281. C
  282. C KPREC = 0 : pas de préconditionnement
  283. C KPREC = 1 : préconditionnement par la diagonale
  284. C KPREC = 2 : préconditionnement par factorisation incomplete
  285. C D-ILU
  286. C KPREC = 3 : préconditionnement par factorisation incomplete
  287. C ILU(0) (i.e. Crout)
  288. C KPREC = 4 : préconditionnement par factorisation incomplete
  289. C modifiée MILU(0) (i.e. Crout)
  290. C KPREC = 5 : préconditionnement par factorisation incomplete
  291. C ILUT (dual truncation strategy)
  292. C KPREC = 6 : préconditionnement ILUT2 (une variante du
  293. C précédent qui remplit mieux la mémoire et
  294. C fonctionne mieux quelquefois)
  295. C KPREC = 7 : préconditionnement ILUTP
  296. C KPREC = 8 : préconditionnement ILUTPGOO
  297. C KPREC = 9 : préconditionnement ILU-PV
  298. C KPREC = 10 : préconditionnement ILU-PV filtré
  299. IF(KTYPI.GE.2.AND.KTYPI.LE.NBTYPI)THEN
  300. IF (IMPR.GT.1) THEN
  301. IF (IDDOT.EQ.0) CCOMP=' '
  302. IF (IDDOT.EQ.1) CCOMP='c'
  303. IF (KTYPI.EQ.5.OR.KTYPI.EQ.8.OR.KTYPI.EQ.11) THEN
  304. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI),' (m=',IRSTRT,')'
  305. ELSEIF (KTYPI.EQ.4) THEN
  306. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI),' (l=',LBCG,')'
  307. ELSE
  308. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI)
  309. ENDIF
  310. IF (KPREC.EQ.4) THEN
  311. WRITE(IOIMP,110) ' ',LPREC(KPREC+1),' (r=',RXMILU,')'
  312. 110 FORMAT (3A,D9.2,A)
  313. ELSEIF (KPREC.GE.5.AND.KPREC.LE.6) THEN
  314. WRITE(IOIMP,111) ' ',LPREC(KPREC+1),' (lfil=',XLFIL,
  315. $ ' droptol=',XDTOL,')'
  316. 111 FORMAT (3A,D9.2,A,D9.2,A)
  317. ELSEIF (KPREC.GE.7.AND.KPREC.LE.8) THEN
  318. WRITE(IOIMP,112) ' ',LPREC(KPREC+1),' (lfil=',XLFIL,
  319. $ ' droptol=',XDTOL,' pivtol=',XSPIV,
  320. $ ')'
  321. 112 FORMAT (3A,D9.2,A,D9.2,A,D9.2,A,I4,A)
  322. ELSEIF (KPREC.EQ.10) THEN
  323. WRITE(IOIMP,110) ' ',LPREC(KPREC+1),' (r=',RXILUP,')'
  324. ELSE
  325. WRITE(IOIMP,*) LPREC(KPREC+1)
  326. ENDIF
  327. ENDIF
  328. IF (LTIME) THEN
  329. call timespv(ittime,oothrd)
  330. ITI1=(ITTIME(1)+ITTIME(2))/10
  331. ENDIF
  332. C
  333. C Construction (éventuelle) du préconditionneur
  334. C
  335. * Gestion du CTRL-C
  336. if (ierr.NE.0) return
  337. IF (KPREC.EQ.1) THEN
  338. CALL MEIDIA(AMORS,AISA,MAPREC,IMPR,IRET)
  339. IF (IRET.NE.0) GOTO 9999
  340. ELSEIF (KPREC.EQ.2) THEN
  341. CALL PRDILU(AMORS,AISA,MAPREC,IMPR,IRET)
  342. IF (IRET.NE.0) GOTO 9999
  343. ELSEIF (KPREC.EQ.3) THEN
  344. CALL PRILU0(AMORS,AISA,MAPREC,IMPR,IRET)
  345. IF (IRET.NE.0) GOTO 9999
  346. ELSEIF (KPREC.EQ.4) THEN
  347. CALL PRMILU(AMORS,AISA,MAPREC,RXMILU,IMPR,IRET)
  348. IF (IRET.NE.0) GOTO 9999
  349. ELSEIF (KPREC.GE.5.AND.KPREC.LE.6) THEN
  350. IVARI=KPREC-5
  351. CALL PRILUT(AMORS,AISA,MAPREC,XLFIL,XDTOL,IVARI,
  352. $ IMPR,IRET)
  353. IF (IRET.NE.0) GOTO 9999
  354. * WRITE(IOIMP,*) 'PRILUT !'
  355. * ELSEIF (KPREC.GE.7.AND.KPREC.LE.9) THEN
  356. ELSEIF (KPREC.GE.7.AND.KPREC.LE.8) THEN
  357. * Ici, on pivote les colonnes lors de la construction du
  358. * préconditionneur...
  359. * On modifiera donc IDMATP et INCX et NORMP
  360. IF (KPREC.EQ.7) THEN
  361. IVARI=0
  362. ELSEIF (KPREC.EQ.8) THEN
  363. IVARI=2
  364. ENDIF
  365. CALL PRILTP(AMORS,AISA,MAPREC,XLFIL,XDTOL,XSPIV,
  366. $ IVARI,
  367. $ INCX,NORMP,LRACOU,
  368. $ IMPR,IRET)
  369. IF (IRET.NE.0) GOTO 9999
  370. ELSEIF (KPREC.EQ.9.OR.KPREC.EQ.10) THEN
  371. CALL PRILUP(AMORS,AISA,MAPREC,RXILUP,KPREC-9,IMPR,IRET)
  372. IF (IRET.NE.0) GOTO 9999
  373. ELSEIF (KPREC.NE.0) THEN
  374. WRITE(IOIMP,*) 'Erreur de programmation'
  375. GOTO 9999
  376. ENDIF
  377. IF (LTIME) THEN
  378. call timespv(ittime,oothrd)
  379. ITI2=(ITTIME(1)+ITTIME(2))/10
  380. ENDIF
  381. C
  382. C Nettoyage des zéros stockés (utilité ?)
  383. C
  384. IF (LFIRST) THEN
  385. * IIMPR=IMPR
  386. * IMPR=6
  387. CALL CLMORS(AMORS,AISA,IMPR,IRET)
  388. IF (IRET.NE.0) GOTO 9999
  389. * IMPR=IIMPR
  390. ENDIF
  391. * SEGPRT, AMORS
  392. * SEGPRT, AISA
  393. C
  394. C Résolution itérative proprement dite
  395. C
  396. * Gestion du CTRL-C
  397. if (ierr.NE.0) return
  398. IF (LTIME) THEN
  399. call timespv(ittime,oothrd)
  400. ITI3=(ITTIME(1)+ITTIME(2))/10
  401. ENDIF
  402. INMV=0
  403. IF (KTYPI.EQ.2) THEN
  404. CALL PRCG2(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  405. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  406. $ IMPR,IRET)
  407. ELSEIF (KTYPI.EQ.3) THEN
  408. LBCG=1
  409. CALL PRBCGG(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  410. $ ITER,INMV,RESID,KPREC,BRTOL,LBCG,ICALRS,IDDOT,IMVEC,
  411. $ IMPR,IRET)
  412. ELSEIF (KTYPI.EQ.4) THEN
  413. CALL PRBCGG(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  414. $ ITER,INMV,RESID,KPREC,BRTOL,LBCG,ICALRS,IDDOT,IMVEC,
  415. $ IMPR,IRET)
  416. ELSEIF (KTYPI.EQ.5) THEN
  417. CALL PRGMR2(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  418. $ ITER,INMV,RESID,KPREC,IRSTRT,ICALRS,IDDOT,IMVEC,
  419. $ IMPR,IRET)
  420. ELSEIF (KTYPI.EQ.6) THEN
  421. CALL PRCGSN(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  422. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  423. $ IMPR,IRET)
  424. ELSEIF (KTYPI.EQ.7.OR.KTYPI.EQ.8.OR.KTYPI.EQ.10.OR.KTYPI.EQ.11)
  425. $ THEN
  426. CALL INCTYP(MATRIK,IORINC,
  427. $ ATYP,ANOD,
  428. $ IMPR,IRET)
  429. IF (IRET.NE.0) GOTO 9999
  430. IF (KTYPI.EQ.7) THEN
  431. NREST=1
  432. ELSEIF (KTYPI.EQ.8) THEN
  433. NREST=IRSTRT
  434. ELSEIF (KTYPI.EQ.10) THEN
  435. NREST=-1
  436. ELSEIF (KTYPI.EQ.11) THEN
  437. NREST=-IRSTRT
  438. ELSE
  439. CALL ERREUR(5)
  440. GOTO 9999
  441. ENDIF
  442. * WRITE(IOIMP,*) 'Appel de pragmg'
  443. CALL PRAGMG(AMORS,ATYP,ANOD,AISA,KS2B,MATRIK,MAPREC,LRES
  444. $ ,LNMV,INCX,ITER,INMV,RESID,KPREC,NREST,ICALRS,IDDOT
  445. $ ,IMVEC,KTIME,LTIME,IMPR,IRET)
  446. SEGSUP ATYP
  447. SEGSUP ANOD
  448. ELSEIF (KTYPI.EQ.9) THEN
  449. CALL PRBCG(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  450. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  451. $ IMPR,IRET)
  452. ELSE
  453. WRITE(IOIMP,*) 'KTYPI=',KTYPI,' incorrect'
  454. GOTO 9999
  455. ENDIF
  456. * Gestion du CTRL-C
  457. if (ierr.NE.0) return
  458. C IRET<0 : 'vrai erreur' ou breakdown (dans le cas de BiCGSTAB)
  459. C IRET>0 : l'itération n'a pas convergé mais on veut
  460. C quand meme la solution calculée
  461. ICVG=IRET
  462. IF (IRET.GT.0) THEN
  463. WRITE(IOIMP,*) 'Convergence to tolerance not achieved : ',
  464. $ 'ITER=',ITER,' INMV=',INMV,' RESID=',RESID
  465. ELSEIF (IRET.EQ.0) THEN
  466. IF (IMPR.GT.0) THEN
  467. WRITE(IOIMP,*) 'ITER=',ITER,' INMV=',INMV,' RESID=',RESID
  468. ENDIF
  469. ELSEIF (IRET.LT.0) THEN
  470. WRITE(IOIMP,*) 'Error or breakdown in iterative method'
  471. GOTO 9999
  472. ENDIF
  473. IF (LTIME) THEN
  474. call timespv(ittime,oothrd)
  475. ITI4=(ITTIME(1)+ITTIME(2))/10
  476. ITP=ITI2-ITI1
  477. ITR=ITI4-ITI3
  478. CHARI='PRECONDI'
  479. CALL ECCTAB(KTIME,'MOT ',IVALI,XVALI,CHARI,LOGII,IRETI,
  480. $ 'ENTIER ',ITP,XVALR,CHARR,LOGIR,IRETR)
  481. CHARI='RESOLUTI'
  482. CALL ECCTAB(KTIME,'MOT ',IVALI,XVALI,CHARI,LOGII,IRETI,
  483. $ 'ENTIER ',ITR,XVALR,CHARR,LOGIR,IRETR)
  484. ENDIF
  485. ELSE
  486. WRITE(IOIMP,*) 'KTYPI=',KTYPI,' invalide (1..',NBTYPI,')'
  487. GOTO 9999
  488. ENDIF
  489. IF (ISCAL.GT.0) THEN
  490. C
  491. C On dénorme le vecteur solution : attention modification...
  492. C
  493. CALL NORMV2(INCX,NORMP,IMPR,IRET)
  494. IF (IRET.NE.0) GOTO 9999
  495. IF (LRACOU) THEN
  496. SEGACT MATRIK*MOD
  497. MATRIK.KKMMT(4)=NORMP
  498. MATRIK.KKMMT(5)=NORMD
  499. SEGDES MATRIK
  500. ENDIF
  501. ENDIF
  502. C
  503. C Transformation du vecteur-solution en chpoint
  504. C
  505. CALL XDISP(MATRIK,INCX,MCHSOL,IMPR,IRET)
  506. IF (IRET.NE.0) GOTO 9999
  507. C
  508. C Suppression des objets temporaires
  509. C
  510. IF (LRACOU) THEN
  511. SEGACT MATRIK*MOD
  512. MATRIK.KIDMAT(4)=AMORS
  513. MATRIK.KIDMAT(5)=AISA
  514. SEGDES MATRIK
  515. ELSE
  516. SEGSUP,AMORS
  517. SEGSUP,AISA
  518. ENDIF
  519. SEGSUP INCX
  520. SEGSUP KS2B
  521. *
  522. * Normal termination
  523. *
  524. IRET=0
  525. RETURN
  526. *
  527. * Format handling
  528. *
  529. *
  530. * Error handling
  531. *
  532. 9999 CONTINUE
  533. IRET=1
  534. WRITE(IOIMP,*) 'An error was detected in kres5.eso'
  535. RETURN
  536. *
  537. * End of KRES5
  538. *
  539. END
  540.  
  541.  

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