Télécharger kres11.eso

Retour à la liste

Numérotation des lignes :

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

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