Télécharger kres5.eso

Retour à la liste

Numérotation des lignes :

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

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