Télécharger kres11.eso

Retour à la liste

Numérotation des lignes :

kres11
  1. C KRES11 SOURCE GOUNAND 25/04/30 21:15:10 12258
  2. SUBROUTINE KRES11(KMORS,KIZA,ISMBR,IPBLOC,
  3. $ KTYPI,ITER,RESID,ICALRS,IRSTRT,LBCG,BRTOL,IDDOT,IMVEC,
  4. $ KPREC,RXMILU,RXILUP,XLFIL,XDTOL,XSPIV,
  5. $ KTIME,LTIME,LDUMP,ISMOOT,
  6. $ INCX,LRES,LNMV,ICVG,IMPR)
  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 IPBLOC.MLENTI
  30. -INC SMLREEL
  31. POINTEUR LRES.MLREEL
  32. POINTEUR MAPREC.MATRIK
  33. POINTEUR AMORS.PMORS
  34. POINTEUR AIZA.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=11)
  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,LDUMP
  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. $ 'AGMG GCR(m) Stokes',
  69. $ 'AGMG GCR(m) NS '/
  70. c$$$ $ 'CG old ',
  71. c$$$ $ 'BiCGSTAB old ',
  72. c$$$ $ 'GMRES(m) old ',
  73. c$$$ $ 'CGS ',
  74. c$$$ $ 'BiCGSTAB ',
  75. c$$$ $ 'BiCGSTAB(l) ',
  76. c$$$ $ 'BiCGSTAB(l)F '/
  77. DATA LPREC/'None ',
  78. $ 'Jacobi ',
  79. $ 'D-ILU ',
  80. $ 'ILU(0) ',
  81. $ 'MILU(0) ',
  82. $ 'ILUT ',
  83. $ 'ILUT2 ',
  84. $ 'ILUTP ',
  85. $ 'ILUTP+0 ',
  86. $ 'ILU0-PV ',
  87. $ 'ILU0-PVf'/
  88.  
  89. IVALI=0
  90. XVALI=0.D0
  91. LOGII=.FALSE.
  92. IRETI=0
  93. XVALR=0.D0
  94. IRETR=0
  95. *
  96. * Executable statements
  97. *
  98. * WRITE(IOIMP,*) 'Entrée dans kres10.eso'
  99. ICVG=0
  100. *Debug IF (KTYPI.EQ.1) KTYPI=3
  101. IF (KTYPI.LT.2.OR.KTYPI.GT.NBTYPI) THEN
  102. WRITE(IOIMP,*) 'KTYPI incorrect (2..',NBTYPI,') :',KTYPI
  103. GOTO 9999
  104. ENDIF
  105. IF (KPREC.LT.0.OR.KPREC.GT.NBPREC-1) THEN
  106. WRITE(IOIMP,*) 'PRECOND ',
  107. $ 'incorrect (0..',NBPREC-1,') :',KPREC
  108. GOTO 9999
  109. ENDIF
  110. C
  111. C Impressions
  112. C
  113. IF (IMPR.GT.1) THEN
  114. IF (IDDOT.EQ.0) CCOMP=' '
  115. IF (IDDOT.EQ.1) CCOMP='c'
  116. IF (KTYPI.EQ.5.OR.KTYPI.EQ.8.OR.KTYPI.EQ.10.OR.KTYPI.EQ.11)
  117. $ THEN
  118. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI),' (m=',IRSTRT,')'
  119. ELSEIF (KTYPI.EQ.4) THEN
  120. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI),' (l=',LBCG,')'
  121. ELSE
  122. WRITE(IOIMP,*) CCOMP,LTYPI(KTYPI)
  123. ENDIF
  124. IF (KPREC.EQ.4) THEN
  125. WRITE(IOIMP,110) ' ',LPREC(KPREC+1),' (r=',RXMILU,')'
  126. 110 FORMAT (3A,D9.2,A)
  127. ELSEIF (KPREC.GE.5.AND.KPREC.LE.6) THEN
  128. WRITE(IOIMP,111) ' ',LPREC(KPREC+1),' (lfil=',XLFIL,
  129. $ ' droptol=',XDTOL,')'
  130. 111 FORMAT (3A,D9.2,A,D9.2,A)
  131. ELSEIF (KPREC.GE.7.AND.KPREC.LE.8) THEN
  132. WRITE(IOIMP,112) ' ',LPREC(KPREC+1),' (lfil=',XLFIL,
  133. $ ' droptol=',XDTOL,' pivtol=',XSPIV,
  134. $ ')'
  135. 112 FORMAT (3A,D9.2,A,D9.2,A,D9.2,A,I4,A)
  136. ELSEIF (KPREC.EQ.10) THEN
  137. WRITE(IOIMP,110) ' ',LPREC(KPREC+1),' (r=',RXILUP,')'
  138. ELSE
  139. WRITE(IOIMP,*) LPREC(KPREC+1)
  140. ENDIF
  141. ENDIF
  142. IF (LTIME) THEN
  143. call timespv(ittime,oothrd)
  144. ITI1=(ITTIME(1)+ITTIME(2))/10
  145. ENDIF
  146. C
  147. C - Construction du préconditionneur (repris sur kres5)
  148. C
  149. SEGACT ISMBR
  150. INC=ISMBR.VECTBB(/1)
  151. NRIGE=0
  152. NMATRI=0
  153. NKID=9
  154. NKMT=7
  155. SEGINI MAPREC
  156. MAPREC.KNTTT=INC
  157. MAPREC.KSYM=2
  158. NTT=0
  159. NPT=0
  160. NBLK=0
  161. SEGINI IDMAT
  162. MAPREC.KIDMAT(1)=IDMAT
  163. *
  164. * Gestion du CTRL-C
  165. if (ierr.NE.0) return
  166. IF (KPREC.EQ.1) THEN
  167. CALL MEIDIA(KMORS,KIZA,MAPREC,IMPR,IRET)
  168. IF (IRET.NE.0) GOTO 9999
  169. ELSEIF (KPREC.EQ.2) THEN
  170. CALL PRDILU(KMORS,KIZA,MAPREC,IMPR,IRET)
  171. IF (IRET.NE.0) GOTO 9999
  172. ELSEIF (KPREC.EQ.3) THEN
  173. CALL PRILU0(KMORS,KIZA,MAPREC,IMPR,IRET)
  174. IF (IRET.NE.0) GOTO 9999
  175. ELSEIF (KPREC.EQ.4) THEN
  176. CALL PRMILU(KMORS,KIZA,MAPREC,RXMILU,IMPR,IRET)
  177. IF (IRET.NE.0) GOTO 9999
  178. ELSEIF (KPREC.GE.5.AND.KPREC.LE.6) THEN
  179. IVARI=KPREC-5
  180. CALL PRILUT(KMORS,KIZA,MAPREC,XLFIL,XDTOL,IVARI,
  181. $ IMPR,IRET)
  182. IF (IRET.NE.0) GOTO 9999
  183. * WRITE(IOIMP,*) 'PRILUT !'
  184. * ELSEIF (KPREC.GE.7.AND.KPREC.LE.9) THEN
  185. ELSEIF (KPREC.GE.7.AND.KPREC.LE.8) THEN
  186. WRITE(IOIMP,*) 'KPREC=',KPREC, 'non dispo'
  187. GOTO 9999
  188. ELSEIF (KPREC.EQ.9.OR.KPREC.EQ.10) THEN
  189. CALL PRILUP(KMORS,KIZA,MAPREC,RXILUP,KPREC-9,IMPR,IRET)
  190. IF (IRET.NE.0) GOTO 9999
  191. ELSEIF (KPREC.NE.0) THEN
  192. WRITE(IOIMP,*) 'Erreur de programmation'
  193. GOTO 9999
  194. ENDIF
  195. IF (LTIME) THEN
  196. call timespv(ittime,oothrd)
  197. ITI2=(ITTIME(1)+ITTIME(2))/10
  198. ENDIF
  199. C MATRIK=MAPREC
  200. C SEGACT MATRIK
  201. C nkid=kidmat(/1)
  202. C WRITE(IOIMP,*) 'Apres constr precon'
  203. C WRITE(IOIMP,*) ' MAPREC KIDMAT NKID=',nkid
  204. C WRITE(IOIMP,2020) (KIDMAT(II),II=1,KIDMAT(/1))
  205. C 2020 FORMAT (20(2X,I12) )
  206. C
  207. C - Appel des solveurs itératifs
  208. C Apparemment, le segment MATRIK ne sert qu'à vérifier la symétrie
  209. C et à construire INCTYP pour le multigrille
  210. C
  211. NRIGE=0
  212. NMATRI=0
  213. NKID=0
  214. NKMT=0
  215. SEGINI MATRIK
  216. IF (KTYPI.NE.2) KSYM=2
  217. SEGACT ISMBR
  218. INC=ISMBR.VECTBB(/1)
  219. SEGINI INCX
  220. KS2B=ISMBR
  221. AMORS=KMORS
  222. AIZA=KIZA
  223. *
  224. * Gestion du CTRL-C
  225. if (ierr.NE.0) return
  226. IF (KTYPI.EQ.2) THEN
  227. CALL PRCG2(AMORS,AIZA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  228. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  229. $ IMPR,IRET)
  230. ELSEIF (KTYPI.EQ.3) THEN
  231. LBCG=1
  232. CALL PRBCGG(AMORS,AIZA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  233. $ ITER,INMV,RESID,KPREC,BRTOL,LBCG,ICALRS,IDDOT,IMVEC,
  234. $ IMPR,IRET)
  235. ELSEIF (KTYPI.EQ.4) THEN
  236. CALL PRBCGG(AMORS,AIZA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  237. $ ITER,INMV,RESID,KPREC,BRTOL,LBCG,ICALRS,IDDOT,IMVEC,
  238. $ IMPR,IRET)
  239. ELSEIF (KTYPI.EQ.5) THEN
  240. C WRITE(IOIMP,*) 'ITER=',ITER
  241. C WRITE(IOIMP,*) 'INMV=',INMV
  242. C WRITE(IOIMP,*) 'RESID=',RESID
  243. C WRITE(IOIMP,*) 'KPREC=',KPREC
  244. C WRITE(IOIMP,*) 'IRSTRT=',IRSTRT
  245. C WRITE(IOIMP,*) 'ICALRS=',ICALRS
  246. C WRITE(IOIMP,*) 'IDDOT=',IDDOT
  247. CALL PRGMR2(AMORS,AIZA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  248. $ ITER,INMV,RESID,KPREC,IRSTRT,ICALRS,IDDOT,IMVEC,
  249. $ IMPR,IRET)
  250. ELSEIF (KTYPI.EQ.6) THEN
  251. CALL PRCGSN(AMORS,AIZA,KS2B,MATRIK,MAPREC,LRES,LNMV,INCX,
  252. $ ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,IMVEC,
  253. $ IMPR,IRET)
  254. ELSEIF (KTYPI.EQ.7.OR.KTYPI.EQ.8.OR.KTYPI.EQ.10.OR.KTYPI.EQ.11)
  255. $ THEN
  256. IF (KTYPI.EQ.7) THEN
  257. NRESTS=1
  258. ELSE
  259. NRESTS=IRSTRT
  260. ENDIF
  261. * WRITE(IOIMP,*) 'Appel de pragmg'
  262. CALL PRAGMG(AMORS,AIZA,KS2B,IPBLOC,KTYPI,LRES,
  263. $ LNMV,INCX,ITER,INMV,
  264. $ RESID,KPREC,
  265. $ NRESTS,ICALRS,IDDOT,IMVEC,KTIME,LTIME,LDUMP,ISMOOT,
  266. $ IMPR,IRET)
  267. ELSEIF (KTYPI.EQ.9) THEN
  268. CALL PRBCG(AMORS,AIZA,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. *
  325. * Normal termination
  326. *
  327. RETURN
  328. *
  329. * Format handling
  330. *
  331. *
  332. * Error handling
  333. *
  334. 9999 CONTINUE
  335. MOTERR(1:8)='KRES11 '
  336. CALL ERREUR(349)
  337. RETURN
  338. *
  339. * End of subroutine KRES11
  340. *
  341. END
  342.  
  343.  

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