Télécharger bicg.eso

Retour à la liste

Numérotation des lignes :

  1. C BICG SOURCE PV 16/11/17 21:58:18 9180
  2. SUBROUTINE BICG(KMORS,KISA,B,X,
  3. $ KPREC,INVDIA,ILUM,ILUI,
  4. $ LRES,LNMV,
  5. $ R,RTLD,Z,ZTLD,P,PTLD,Q,QTLD,
  6. $ ITER,INMV,RESID,BRTOL,ICALRS,IDDOT,IMVEC,IMPR,IRET)
  7. IMPLICIT INTEGER(I-N)
  8. IMPLICIT REAL*8 (A-H,O-Z)
  9. C***********************************************************************
  10. C NOM : BICG
  11. C DESCRIPTION :
  12. C Résolution d'un système linéaire Ax=b
  13. C par la méthode du bi-gradient conjugué préconditionné.
  14. C
  15. C LANGAGE : FORTRAN 77 + chouhia ESOPE (pour les E/S)
  16. C AUTEUR : Stéphane GOUNAND (CEA/DRN/DMT/SEMT/TTMF)
  17. C mél : gounand@semt2.smts.cea.fr
  18. C REFERENCE (bibtex-like) :
  19. C @BOOK{templates,
  20. C AUTHOR={R.Barrett, M.Berry, T.F.Chan, J.Demmel, J.Donato,
  21. C J.Dongarra, V.Eijkhout, R.Pozo, C.Romine,
  22. C H. Van der Vorst},
  23. C TITLE={Templates for the Solution of Linear Systems :
  24. C Building Blocks for Iterative Methods},
  25. C PUBLISHER={SIAM}, YEAR={1994}, ADDRESS={Philadelphia,PA} }
  26. C -> URL : http://www.netlib.org/templates/Templates.html
  27. C***********************************************************************
  28. C APPELES (BLAS 1) : GCOPY, GAXPY (subroutines)
  29. C GNRM2, GDOT (functions)
  30. C APPELES (Calcul) : GMOMV : produit Matrice Morse - Vecteur
  31. C PSOLVE : Preconditionner solve
  32. C***********************************************************************
  33. C ENTREES :
  34. C ENTREES/SORTIES :
  35. C SORTIES :
  36. C CODE RETOUR (IRET) : 0 si ok
  37. C <0 si problème (non utilisé !)
  38. C >0 si l'algorithme n'a pas convergé
  39. C***********************************************************************
  40. C VERSION : v1, 09/02/06
  41. C HISTORIQUE : v1, 09/02/06, création
  42. C HISTORIQUE :
  43. C HISTORIQUE :
  44. C***********************************************************************
  45. C Prière de PRENDRE LE TEMPS de compléter les commentaires
  46. C en cas de modification de ce sous-programme afin de faciliter
  47. C la maintenance !
  48. C***********************************************************************
  49. *
  50. * .. Includes .. a commenter si mise au point ok
  51. -INC CCOPTIO
  52. -INC CCREEL
  53. -INC SMLREEL
  54. POINTEUR LRES.MLREEL
  55. -INC SMLENTI
  56. POINTEUR LNMV.MLENTI
  57. POINTEUR KMORS.PMORS
  58. POINTEUR KISA.IZA
  59. POINTEUR INVDIA.IZA
  60. POINTEUR ILUM.PMORS
  61. POINTEUR ILUI.IZA
  62. POINTEUR B.IZA
  63. POINTEUR X.IZA
  64. * .. Work Vectors for CG
  65. POINTEUR R.IZA,P.IZA,Q.IZA,Z.IZA
  66. POINTEUR RTLD.IZA,PTLD.IZA,QTLD.IZA,ZTLD.IZA
  67. * .. Scalar Arguments ..
  68. INTEGER ITER, KPREC, IMPR, IRET
  69. REAL*8 RESID
  70. *
  71. *
  72. * .. Variables locales
  73. * .. Parameters
  74. REAL*8 ZERO , ONE
  75. PARAMETER (ZERO=0.0D0, ONE=1.0D0)
  76. * ..
  77. INTEGER MAXIT
  78. REAL*8 ALPHA,BETA,RHO,RHO1,TOL,GNRM2,OMEGA
  79. REAL*8 GDOT,BNRM2
  80. * ..
  81. * .. External subroutines and functions..
  82. EXTERNAL GAXPY,GCOPY,GDOT,GNRM2
  83. * ..
  84. * .. Executable Statements ..
  85. *
  86.  
  87. IRET = 0
  88. INMV = 0
  89. IWARN = 0
  90. IGRAN = 0
  91. MGRAN = 5
  92. MAXIT = ITER
  93. TOL = RESID
  94. IF(IMPR.GT.5) WRITE(IOIMP,*) 'Entrée dans bicg.eso'
  95. *
  96. * Set parameter tolerances.
  97. *
  98. RHOTOL = BRTOL
  99. OMETOL = BRTOL
  100. RRMAX = 1.D0/(XZPREC**0.25D0)
  101. *
  102. * Calcul du résidu initial
  103. *
  104. C r(0)=b
  105. CALL GCOPY(B,R)
  106.  
  107. C r(0)=b-Ax
  108. CALL GMOMV(IMVEC,'N',-ONE,KMORS,KISA,X,ONE,R)
  109. INMV=INMV+1
  110. *
  111. ITER = 0
  112. BNRM2 = GNRM2(B)
  113. IF(IMPR.GT.5) THEN
  114. WRITE(IOIMP,*) '||B||=',BNRM2
  115. ENDIF
  116. RNRM2 = GNRM2(R)
  117. IF (ICALRS.EQ.1) BNRM2=RNRM2
  118. IF (BNRM2.LT.XPETIT) BNRM2 = ONE
  119. RESID = RNRM2 / BNRM2
  120. IF(IMPR.GT.5) THEN
  121. WRITE(IOIMP,*) 'Initial residual=',RESID
  122. ENDIF
  123. RESMAX = RESID * RRMAX
  124. LRES.PROG(ITER+1)=RESID
  125. LNMV.LECT(ITER+1)=INMV
  126. *
  127. IF (RESID.LE.TOL) THEN
  128. IRET = 0
  129. GOTO 30
  130. ENDIF
  131. CALL GCOPY(R,RTLD)
  132. *
  133. 10 CONTINUE
  134. *
  135. * Perform (Preconditioned) BiConjugate Gradient Iteration.
  136. *
  137. ITER = ITER +1
  138. IF(IMPR.GT.7) THEN
  139. WRITE(IOIMP,*) 'ITER=',ITER
  140. ENDIF
  141. *
  142. * Preconditioner Solve
  143. C Mz(i-1)=r(i-1)
  144. CALL PSOLVE('N',KPREC,KMORS,KISA,INVDIA,ILUM,ILUI
  145. *! CALL PSOLVE(KPREC,KMORS,KISA,INVDIA,ILUM,ILUI
  146. $ ,Z,R)
  147. C M'z'(i-1)=r'(i-1)
  148. *! CALL PSOLVE(KPREC,KMORS,KISA,INVDIA,ILUM,ILUI
  149. CALL PSOLVE('T',KPREC,KMORS,KISA,INVDIA,ILUM,ILUI
  150. $ ,ZTLD,RTLD)
  151. C rho(i-1)=z(i-1)'r'(i-1)
  152. RHO = GDOT(Z,RTLD,IDDOT)
  153. C WRITE(IOIMP,*) ' RHO=',RHO
  154. C XR1 = GDOT(Z,R)
  155. C WRITE(IOIMP,*) ' XR1=',XR1
  156. C XR2 = GDOT(ZTLD,RTLD)
  157. C WRITE(IOIMP,*) ' XR2=',XR2
  158. IF (ABS(RHO).LT.RHOTOL) THEN
  159. RHO=SIGN(XZPREC,RHO)
  160. IWARN=IWARN+1
  161. ENDIF
  162. *
  163. * Compute direction vector P
  164. *
  165. IF (ITER.GT.1) THEN
  166. C beta(i-1)=rho(i-1)/rho(i-2)
  167. BETA = RHO / RHO1
  168. C p(i)=z(i-1)+beta(i-1)p(i-1)
  169. C p'(i)=z'(i-1)+beta(i-1)p'(i-1)
  170. CALL GAXPY(BETA,P,Z)
  171. CALL GCOPY(Z,P)
  172. CALL GAXPY(BETA,PTLD,ZTLD)
  173. CALL GCOPY(ZTLD,PTLD)
  174. ELSE
  175. C p(1)=z(0)
  176. C p'(1)=z'(0)
  177. CALL GCOPY(Z,P)
  178. CALL GCOPY(ZTLD,PTLD)
  179. ENDIF
  180. *
  181. * Compute scalar ALPHA (save A*p to q)
  182. *
  183. C q(i)=Ap(i)
  184. CALL GMOMV(IMVEC,'N',ONE,KMORS,KISA,P,ZERO,Q)
  185. INMV=INMV+1
  186. C q'(i)=A'p'(i)
  187. CALL GMOMV(IMVEC,'T',ONE,KMORS,KISA,PTLD,ZERO,QTLD)
  188. INMV=INMV+1
  189. C alpha(i)=rho(i-1)/(p'(i)'q(i))
  190. OMEGA = GDOT(PTLD,Q,IDDOT)
  191. C WRITE(IOIMP,*) ' OMEGA=',OMEGA
  192. C XOM1 = GDOT(P,Q)
  193. C WRITE(IOIMP,*) ' XOM1 =',XOM1
  194. C XOM2 = GDOT(PTLD,QTLD)
  195. C WRITE(IOIMP,*) ' XOM2 =',XOM2
  196. IF (ABS(OMEGA).LT.OMETOL) THEN
  197. OMEGA=SIGN(XZPREC,OMEGA)
  198. IWARN=IWARN+1
  199. ENDIF
  200. ALPHA = RHO / OMEGA
  201. *
  202. * Compute current solution vector x
  203. *
  204. C x(i)=x(i-1)+alpha(i)p(i)
  205. CALL GAXPY(ALPHA,P,X)
  206. *
  207. * Compute residual vector R, find norm,
  208. * then check for tolerance
  209. *
  210. C r(i)=r(i-1)-alpha(i)q(i)
  211. CALL GAXPY(-ALPHA,Q,R)
  212. *
  213. RESID = GNRM2(R) / BNRM2
  214. IF(IMPR.GT.7) THEN
  215. WRITE(IOIMP,*) 'Résidu=',RESID
  216. ENDIF
  217. LRES.PROG(ITER+1)=RESID
  218. LNMV.LECT(ITER+1)=INMV
  219. IF (RESID.GT.RESMAX) THEN
  220. IGRAN=IGRAN+1
  221. ELSE
  222. IGRAN=0
  223. ENDIF
  224. *
  225. * Iteration successful
  226. *
  227. IF (RESID.LE.TOL) THEN
  228. IRET = 0
  229. *
  230. * Iteration fails
  231. *
  232. * ELSE IF (ITER.EQ.MAXIT.OR.IGRAN.GT.MGRAN) THEN
  233. ELSE IF (INMV.GE.MAXIT.OR.IGRAN.GT.MGRAN) THEN
  234. IRET = 1
  235. *
  236. * Do some more
  237. *
  238. ELSE
  239. C r'(i)=r'(i-1)-alpha(i)q'(i)
  240. CALL GAXPY(-ALPHA,QTLD,RTLD)
  241. RHO1 = RHO
  242. GOTO 10
  243. ENDIF
  244. *
  245. * Normal termination
  246. *
  247. 30 CONTINUE
  248. IF (IWARN.GT.0.AND.IMPR.GT.0) THEN
  249. WRITE(IOIMP,*) 'BiCG nearly failed ',IWARN
  250. $ ,' times.'
  251. ENDIF
  252. RETURN
  253. *
  254. * End of BICG
  255. *
  256. END
  257.  
  258.  
  259.  
  260.  
  261.  
  262.  
  263.  
  264.  
  265.  
  266.  
  267.  
  268.  
  269.  

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