Télécharger kres11.eso

Retour à la liste

Numérotation des lignes :

kres11
  1. C KRES11 SOURCE GOUNAND 22/08/25 21:15:04 11434
  2. SUBROUTINE KRES11(KMORS,KIZA,KTYP,ISMBR,
  3. $ KTYPI,ITER,RESID,ICALRS,IRSTRT,LBCG,BRTOL,IDDOT,IMVEC,
  4. $ KPREC,RXMILU,RXILUP,XLFIL,XDTOL,XSPIV,
  5. $ KTIME,LTIME,
  6. $ INCX,LRES,LNMV,ICVG,IMPR,INODET)
  7. IMPLICIT REAL*8 (A-H,O-Z)
  8. IMPLICIT INTEGER (I-N)
  9. C***********************************************************************
  10. C NOM : KRES11
  11. C DESCRIPTION : - Construction du préconditionneur
  12. C - Appel des solveurs itératifs
  13. C
  14. C
  15. C LANGAGE : ESOPE
  16. C AUTEUR : Stéphane GOUNAND (CEA/DEN/DM2S/SFME/LTMF)
  17. C mél : gounand@semt2.smts.cea.fr
  18. C***********************************************************************
  19. C VERSION : v1, 29/08/2011, version initiale
  20. C HISTORIQUE : v1, 29/08/2011, création
  21. C HISTORIQUE :
  22. C HISTORIQUE :
  23. C***********************************************************************
  24.  
  25. -INC PPARAM
  26. -INC CCOPTIO
  27. -INC SMLENTI
  28. POINTEUR LNMV.MLENTI
  29. POINTEUR ATYP.MLENTI
  30. -INC SMLREEL
  31. POINTEUR LRES.MLREEL
  32. POINTEUR MAPREC.MATRIK
  33. POINTEUR AMORS.PMORS
  34. POINTEUR AISA.IZA
  35. -INC SMVECTD
  36. POINTEUR ISMBR.MVECTD
  37. POINTEUR INCX.MVECTD
  38. C Tableau de correspondance : entier <-> mot
  39. C pour le type d'inversion
  40. INTEGER NBTYPI,LNTYPI
  41. PARAMETER (NBTYPI=16)
  42. PARAMETER (LNTYPI=18)
  43. CHARACTER*(LNTYPI) LTYPI(NBTYPI)
  44. C Tableau de correspondance : entier <-> mot
  45. C pour le type de préconditionnement (cas d'une méthode itérative)
  46. INTEGER NBPREC,LNPREC
  47. PARAMETER (NBPREC=11)
  48. PARAMETER (LNPREC=8)
  49. CHARACTER*(LNPREC) LPREC(NBPREC)
  50. -INC SMTABLE
  51. POINTEUR KTIME.MTABLE
  52. DIMENSION ITTIME(4)
  53. CHARACTER*8 CHARI
  54. CHARACTER*1 CCOMP
  55. LOGICAL LTIME,LOGII
  56. C
  57. C Initialisation des tableaux
  58. C
  59. DATA LTYPI/'Choleski ',
  60. $ 'Conjugate Gradient',
  61. $ 'BiCGSTAB G ',
  62. $ 'BiCGSTAB(l) G ',
  63. $ 'GMRES(m) ',
  64. $ 'CGS-Neumaier ',
  65. $ 'Multigrid FCG ',
  66. $ 'Multigrid GCR(m) ',
  67. $ 'Bi-CG ',
  68. $ 'CG old ',
  69. $ 'BiCGSTAB old ',
  70. $ 'GMRES(m) old ',
  71. $ 'CGS ',
  72. $ 'BiCGSTAB ',
  73. $ 'BiCGSTAB(l) ',
  74. $ 'BiCGSTAB(l)F '/
  75. DATA LPREC/'None ',
  76. $ 'Jacobi ',
  77. $ 'D-ILU ',
  78. $ 'ILU(0) ',
  79. $ 'MILU(0) ',
  80. $ 'ILUT ',
  81. $ 'ILUT2 ',
  82. $ 'ILUTP ',
  83. $ 'ILUTP+0 ',
  84. $ 'ILU0-PV ',
  85. $ 'ILU0-PVf'/
  86.  
  87. IVALI=0
  88. XVALI=0.D0
  89. LOGII=.FALSE.
  90. IRETI=0
  91. XVALR=0.D0
  92. IOBRE=0
  93. IRETR=0
  94. *
  95. * Executable statements
  96. *
  97. * WRITE(IOIMP,*) 'Entrée dans kres10.eso'
  98. ICVG=0
  99. *Debug IF (KTYPI.EQ.1) KTYPI=3
  100. IF (KTYPI.LT.2.OR.KTYPI.GT.NBTYPI) THEN
  101. WRITE(IOIMP,*) 'KTYPI incorrect (2..',NBTYPI,') :',KTYPI
  102. GOTO 9999
  103. ENDIF
  104. IF (KPREC.LT.0.OR.KPREC.GT.NBPREC-1) THEN
  105. WRITE(IOIMP,*) 'PRECOND ',
  106. $ 'incorrect (0..',NBPREC-1,') :',KPREC
  107. GOTO 9999
  108. ENDIF
  109. C
  110. C Impressions
  111. C
  112. IF (IMPR.GT.1) THEN
  113. IF (IDDOT.EQ.0) CCOMP=' '
  114. IF (IDDOT.EQ.1) CCOMP='c'
  115. IF (KTYPI.EQ.5.OR.KTYPI.EQ.8) THEN
  116. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI),' (m=',IRSTRT,')'
  117. ELSEIF (KTYPI.EQ.11.OR.KTYPI.EQ.12) THEN
  118. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI),' (l=',LBCG,')'
  119. ELSE
  120. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI)
  121. ENDIF
  122. IF (KPREC.EQ.4) THEN
  123. WRITE(IOIMP,110) ' ',LPREC(KPREC+1),' (r=',RXMILU,')'
  124. 110 FORMAT (3A,D9.2,A)
  125. ELSEIF (KPREC.GE.5.AND.KPREC.LE.6) THEN
  126. WRITE(IOIMP,111) ' ',LPREC(KPREC+1),' (lfil=',XLFIL,
  127. $ ' droptol=',XDTOL,')'
  128. 111 FORMAT (3A,D9.2,A,D9.2,A)
  129. ELSEIF (KPREC.GE.7.AND.KPREC.LE.8) THEN
  130. WRITE(IOIMP,112) ' ',LPREC(KPREC+1),' (lfil=',XLFIL,
  131. $ ' droptol=',XDTOL,' pivtol=',XSPIV,
  132. $ ')'
  133. 112 FORMAT (3A,D9.2,A,D9.2,A,D9.2,A,I4,A)
  134. ELSEIF (KPREC.EQ.10) THEN
  135. WRITE(IOIMP,110) ' ',LPREC(KPREC+1),' (r=',RXILUP,')'
  136. ELSE
  137. WRITE(IOIMP,*) LPREC(KPREC+1)
  138. ENDIF
  139. ENDIF
  140. IF (LTIME) THEN
  141. call timespv(ittime,oothrd)
  142. ITI1=(ITTIME(1)+ITTIME(2))/10
  143. ENDIF
  144. C
  145. C - Construction du préconditionneur (repris sur kres5)
  146. C
  147. SEGACT ISMBR
  148. INC=ISMBR.VECTBB(/1)
  149. NRIGE=0
  150. NMATRI=0
  151. NKID=9
  152. NKMT=7
  153. SEGINI MAPREC
  154. MAPREC.KNTTT=INC
  155. MAPREC.KSYM=2
  156. NTT=0
  157. NPT=0
  158. NBLK=0
  159. SEGINI IDMAT
  160. MAPREC.KIDMAT(1)=IDMAT
  161. *
  162. * Gestion du CTRL-C
  163. if (ierr.NE.0) return
  164. IF (KPREC.EQ.1) THEN
  165. CALL MEIDIA(KMORS,KIZA,MAPREC,IMPR,IRET)
  166. IF (IRET.NE.0) GOTO 9999
  167. ELSEIF (KPREC.EQ.2) THEN
  168. CALL PRDILU(KMORS,KIZA,MAPREC,IMPR,IRET)
  169. IF (IRET.NE.0) GOTO 9999
  170. ELSEIF (KPREC.EQ.3) THEN
  171. CALL PRILU0(KMORS,KIZA,MAPREC,IMPR,IRET)
  172. IF (IRET.NE.0) GOTO 9999
  173. ELSEIF (KPREC.EQ.4) THEN
  174. CALL PRMILU(KMORS,KIZA,MAPREC,RXMILU,IMPR,IRET)
  175. IF (IRET.NE.0) GOTO 9999
  176. ELSEIF (KPREC.GE.5.AND.KPREC.LE.6) THEN
  177. IVARI=KPREC-5
  178. CALL PRILUT(KMORS,KIZA,MAPREC,XLFIL,XDTOL,IVARI,
  179. $ IMPR,IRET)
  180. IF (IRET.NE.0) GOTO 9999
  181. * WRITE(IOIMP,*) 'PRILUT !'
  182. * ELSEIF (KPREC.GE.7.AND.KPREC.LE.9) THEN
  183. ELSEIF (KPREC.GE.7.AND.KPREC.LE.8) THEN
  184. WRITE(IOIMP,*) 'KPREC=',KPREC, 'non dispo'
  185. GOTO 9999
  186. ELSEIF (KPREC.EQ.9.OR.KPREC.EQ.10) THEN
  187. CALL PRILUP(KMORS,KIZA,MAPREC,RXILUP,KPREC-9,IMPR,IRET)
  188. IF (IRET.NE.0) GOTO 9999
  189. ELSEIF (KPREC.NE.0) THEN
  190. WRITE(IOIMP,*) 'Erreur de programmation'
  191. GOTO 9999
  192. ENDIF
  193. IF (LTIME) THEN
  194. call timespv(ittime,oothrd)
  195. ITI2=(ITTIME(1)+ITTIME(2))/10
  196. ENDIF
  197. C MATRIK=MAPREC
  198. C SEGACT MATRIK
  199. C nkid=kidmat(/1)
  200. C WRITE(IOIMP,*) 'Apres constr precon'
  201. C WRITE(IOIMP,*) ' MAPREC KIDMAT NKID=',nkid
  202. C WRITE(IOIMP,2020) (KIDMAT(II),II=1,KIDMAT(/1))
  203. C 2020 FORMAT (20(2X,I12) )
  204. C
  205. C - Appel des solveurs itératifs
  206. C Apparemment, le segment MATRIK ne sert qu'à vérifier la symétrie
  207. C et à construire INCTYP pour le multigrille
  208. C
  209. NRIGE=0
  210. NMATRI=0
  211. NKID=0
  212. NKMT=0
  213. SEGINI MATRIK
  214. IF (KTYPI.NE.2) KSYM=2
  215. SEGACT ISMBR
  216. INC=ISMBR.VECTBB(/1)
  217. SEGINI INCX
  218. KS2B=ISMBR
  219. AMORS=KMORS
  220. AISA=KIZA
  221. ATYP=KTYP
  222. *
  223. * Gestion du CTRL-C
  224. if (ierr.NE.0) return
  225. IF (KTYPI.EQ.2) THEN
  226. CALL PRCG2(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  227. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  228. $ IMPR,IRET)
  229. ELSEIF (KTYPI.EQ.3) THEN
  230. LBCG=1
  231. CALL PRBCGG(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  232. $ ITER,INMV,RESID,KPREC,BRTOL,LBCG,ICALRS,IDDOT,IMVEC,
  233. $ IMPR,IRET)
  234. ELSEIF (KTYPI.EQ.4) THEN
  235. CALL PRBCGG(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  236. $ ITER,INMV,RESID,KPREC,BRTOL,LBCG,ICALRS,IDDOT,IMVEC,
  237. $ IMPR,IRET)
  238. ELSEIF (KTYPI.EQ.5) THEN
  239. C WRITE(IOIMP,*) 'ITER=',ITER
  240. C WRITE(IOIMP,*) 'INMV=',INMV
  241. C WRITE(IOIMP,*) 'RESID=',RESID
  242. C WRITE(IOIMP,*) 'KPREC=',KPREC
  243. C WRITE(IOIMP,*) 'IRSTRT=',IRSTRT
  244. C WRITE(IOIMP,*) 'ICALRS=',ICALRS
  245. C WRITE(IOIMP,*) 'IDDOT=',IDDOT
  246. CALL PRGMR2(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  247. $ ITER,INMV,RESID,KPREC,IRSTRT,ICALRS,IDDOT,IMVEC,
  248. $ IMPR,IRET)
  249. ELSEIF (KTYPI.EQ.6) THEN
  250. CALL PRCGSN(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  251. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  252. $ IMPR,IRET)
  253. ELSEIF (KTYPI.EQ.7) THEN
  254. * WRITE(IOIMP,*) 'Appel de pragmg'
  255. CALL PRAGMG(AMORS,ATYP,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,
  256. $ INCX,
  257. $ ITER,INMV,RESID,KPREC,1,ICALRS,IDDOT,IMVEC,KTIME,LTIME,
  258. $ IMPR,IRET)
  259. SEGSUP ATYP
  260. ELSEIF (KTYPI.EQ.8) THEN
  261. CALL PRAGMG(AMORS,ATYP,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,
  262. $ INCX,
  263. $ ITER,INMV,RESID,KPREC,IRSTRT,ICALRS,IDDOT,IMVEC,KTIME,
  264. $ LTIME,
  265. $ IMPR,IRET)
  266. SEGSUP ATYP
  267. ELSEIF (KTYPI.EQ.9) THEN
  268. CALL PRBCG(AMORS,AISA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  269. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  270. $ IMPR,IRET)
  271. ELSE
  272. WRITE(IOIMP,*) 'KTYPI=',KTYPI,' invalide (1..',NBTYPI,')'
  273. GOTO 9999
  274. ENDIF
  275. * Gestion du CTRL-C
  276. if (ierr.NE.0) return
  277. C IRET<0 : 'vrai erreur' ou breakdown (dans le cas de BiCGSTAB)
  278. C IRET>0 : l'itération n'a pas convergé mais on veut
  279. C quand meme la solution calculée
  280. ICVG=IRET
  281. IF (IRET.GT.0) THEN
  282. WRITE(IOIMP,*) 'Convergence to tolerance not achieved : ',
  283. $ 'ITER=',ITER,' INMV=',INMV,' RESID=',RESID
  284. ELSEIF (IRET.EQ.0) THEN
  285. IF (IMPR.GT.0) THEN
  286. WRITE(IOIMP,*) 'ITER=',ITER,' INMV=',INMV,' RESID=',RESID
  287. ENDIF
  288. ELSEIF (IRET.LT.0) THEN
  289. WRITE(IOIMP,*) 'Error or breakdown in iterative method'
  290. GOTO 9999
  291. ELSE
  292. WRITE(IOIMP,*) 'KTYPI=',KTYPI,' invalide (1..',NBTYPI,')'
  293. GOTO 9999
  294. ENDIF
  295. SEGSUP MATRIK
  296. IF (LTIME) THEN
  297. call timespv(ittime,oothrd)
  298. ITI3=(ITTIME(1)+ITTIME(2))/10
  299. ITP=ITI2-ITI1
  300. ITR=ITI3-ITI2
  301. CHARI='PRECONDI'
  302. CALL ECCTAB(KTIME,'MOT ',IVALI,XVALI,CHARI,LOGII,IRETI,
  303. $ 'ENTIER ',ITP,XVALR,CHARR,LOGIR,IRETR)
  304. CHARI='RESOLUTI'
  305. CALL ECCTAB(KTIME,'MOT ',IVALI,XVALI,CHARI,LOGII,IRETI,
  306. $ 'ENTIER ',ITR,XVALR,CHARR,LOGIR,IRETR)
  307. ENDIF
  308. C
  309. C Destruction du préconditionneur
  310. C
  311. SEGACT MAPREC
  312. IDMAT=MAPREC.KIDMAT(1)
  313. IF (IDMAT.NE.0) THEN
  314. SEGACT IDMAT
  315. IZA=IDIAG
  316. SEGSUP IZA
  317. SEGSUP IDMAT
  318. ENDIF
  319. PMORS=MAPREC.KIDMAT(6)
  320. IF (PMORS.NE.0) SEGSUP PMORS
  321. IZA=MAPREC.KIDMAT(7)
  322. IF (IZA.NE.0) SEGSUP IZA
  323. SEGSUP MAPREC
  324. C
  325. C Destruction de la matrice Morse
  326. C
  327. IF (INODET.EQ.0) THEN
  328. SEGSUP AMORS
  329. SEGSUP AISA
  330. ENDIF
  331. *
  332. * Normal termination
  333. *
  334. RETURN
  335. *
  336. * Format handling
  337. *
  338. *
  339. * Error handling
  340. *
  341. 9999 CONTINUE
  342. MOTERR(1:8)='KRES11 '
  343. CALL ERREUR(349)
  344. RETURN
  345. *
  346. * End of subroutine KRES11
  347. *
  348. END
  349.  
  350.  

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