Télécharger cg2.eso

Retour à la liste

Numérotation des lignes :

cg2
  1. C CG2 SOURCE GOUNAND 25/07/15 21:15:02 12322
  2. SUBROUTINE CG2(KMORS,KISA,B,X,
  3. $ KPREC,INVDIA,ILUM,ILUI,
  4. $ LRES,LNMV,
  5. $ R,P,Q,Z,
  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 : CG2
  11. C DESCRIPTION :
  12. C Résolution d'un système linéaire Ax=b
  13. C par la méthode du 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.  
  52. -INC PPARAM
  53. -INC CCOPTIO
  54. -INC CCREEL
  55. -INC SMLREEL
  56. POINTEUR LRES.MLREEL
  57. -INC SMLENTI
  58. POINTEUR LNMV.MLENTI
  59. POINTEUR KMORS.PMORS
  60. POINTEUR KISA.IZA
  61. POINTEUR INVDIA.IZA
  62. POINTEUR ILUM.PMORS
  63. POINTEUR ILUI.IZA
  64. POINTEUR B.IZA
  65. POINTEUR X.IZA
  66. * .. Work Vectors for CG
  67. POINTEUR R.IZA,P.IZA,Q.IZA,Z.IZA
  68. * .. Scalar Arguments ..
  69. INTEGER ITER, KPREC, IMPR, IRET
  70. REAL*8 RESID
  71. *
  72. *
  73. * .. Variables locales
  74. * .. Parameters
  75. REAL*8 ZERO , ONE
  76. PARAMETER (ZERO=0.0D0, ONE=1.0D0)
  77. * ..
  78. INTEGER MAXIT
  79. REAL*8 ALPHA,BETA,RHO,RHO1,TOL,GNRM2,OMEGA
  80. REAL*8 GDOT,BNRM2
  81. * ..
  82. * .. External subroutines and functions..
  83. EXTERNAL GAXPY,GCOPY,GDOT,GNRM2
  84. * ..
  85. * .. Executable Statements ..
  86. *
  87.  
  88. IRET = 0
  89. INMV = 0
  90. IWARN = 0
  91. MWARN = 20
  92. IGRAN = 0
  93. MGRAN = 5
  94. MAXIT = ITER
  95. TOL = RESID
  96. IF(IMPR.GT.5) WRITE(IOIMP,*) 'Entrée dans cg2.eso'
  97. *
  98. * Set parameter tolerances.
  99. *
  100. RHOTOL = XPETIT**0.8D0
  101. OMETOL = RHOTOL
  102. * RRMAX = 10.D0
  103. RRMAX = 1.D0/(XZPREC**0.25D0)
  104. * write(ioimp,*) 'RHOTOL,OMETOL,RRMAX=',RHOTOL,OMETOL,RRMAX
  105. *
  106. * Calcul du résidu initial
  107. *
  108. C r(0)=b
  109. CALL GCOPY(B,R)
  110.  
  111. C r(0)=b-Ax
  112. CALL GMOMV(IMVEC,'N',-ONE,KMORS,KISA,X,ONE,R)
  113. INMV=INMV+1
  114. *
  115. ITER = 0
  116. BNRM2 = GNRM2(B)
  117. IF(IMPR.GT.5) THEN
  118. WRITE(IOIMP,*) '||B||=',BNRM2
  119. ENDIF
  120. RNRM2 = GNRM2(R)
  121. IF (ICALRS.EQ.1) BNRM2=RNRM2
  122. IF (BNRM2.LT.XPETIT) BNRM2 = ONE
  123. RESID = RNRM2 / BNRM2
  124. IF(IMPR.GT.5) THEN
  125. WRITE(IOIMP,*) 'Initial residual=',RESID
  126. ENDIF
  127. RESMAX = RESID * RRMAX
  128. LRES.PROG(ITER+1)=RESID
  129. LNMV.LECT(ITER+1)=INMV
  130. *
  131. IF (RESID.LE.TOL) THEN
  132. IRET = 0
  133. GOTO 30
  134. ENDIF
  135.  
  136. *
  137. 10 CONTINUE
  138. * Gestion du CTRL-C
  139. if (ierr.NE.0) return
  140. *
  141. * Perform (Preconditioned) Conjugate Gradient Iteration.
  142. *
  143. ITER = ITER +1
  144. IF(IMPR.GT.7) THEN
  145. WRITE(IOIMP,*) 'ITER=',ITER
  146. ENDIF
  147. *
  148. * Preconditioner Solve
  149. C Mz(i-1)=r(i-1)
  150. CALL PSOLVE('N',KPREC,KMORS,KISA,INVDIA,ILUM,ILUI
  151. $ ,Z,R)
  152. C rho(i-1)=r(i-1)'z(i-1)
  153. RHO = GDOT(R,Z,IDDOT)
  154. * WRITE(IOIMP,*) ' RHO=',RHO
  155. IF (ABS(RHO).LT.RHOTOL) THEN
  156. RHO=SIGN(RHOTOL,RHO)
  157. IWARN=IWARN+1
  158. * WRITE(IOIMP,*) ' RHOCORR,IWARN=',RHO,IWARN
  159. ENDIF
  160. *
  161. * Compute direction vector P
  162. *
  163. IF (ITER.GT.1) THEN
  164. C beta(i-1)=rho(i-1)/rho(i-2)
  165. BETA = RHO / RHO1
  166. C p(i)=z(i-1)+beta(i-1)p(i-1)
  167. CALL GAXPY(BETA,P,Z)
  168. CALL GCOPY(Z,P)
  169. ELSE
  170. C p(1)=z(0)
  171. CALL GCOPY(Z,P)
  172. ENDIF
  173. *
  174. * Compute scalar ALPHA (save A*p to q)
  175. *
  176. C q(i)=Ap(i)
  177. CALL GMOMV(IMVEC,'N',ONE,KMORS,KISA,P,ZERO,Q)
  178. INMV=INMV+1
  179. C alpha(i)=rho(i-1)/(p(i)'q(i))
  180. OMEGA = GDOT(P,Q,IDDOT)
  181. * WRITE(IOIMP,*) ' OMEGA=',OMEGA
  182. IF (ABS(OMEGA).LT.OMETOL) THEN
  183. OMEGA=SIGN(OMETOL,OMEGA)
  184. IWARN=IWARN+1
  185. * WRITE(IOIMP,*) ' OMEGACORR,IWARN=',OMEGA,IWARN
  186. ENDIF
  187. ALPHA = RHO / OMEGA
  188. *
  189. * Compute current solution vector x
  190. *
  191. C x(i)=x(i-1)+alpha(i)p(i)
  192. CALL GAXPY(ALPHA,P,X)
  193. *
  194. * Compute residual vector R, find norm,
  195. * then check for tolerance
  196. *
  197. C r(i)=r(i-1)-alpha(i)r(i)
  198. CALL GAXPY(-ALPHA,Q,R)
  199. *
  200. RESID = GNRM2(R) / BNRM2
  201. IF(IMPR.GT.7) THEN
  202. WRITE(IOIMP,*) 'Résidu=',RESID
  203. ENDIF
  204. LRES.PROG(ITER+1)=RESID
  205. LNMV.LECT(ITER+1)=INMV
  206. IF (RESID.GT.RESMAX) THEN
  207. IGRAN=IGRAN+1
  208. * WRITE(IOIMP,*) 'Résidu grand,IGRAN=',RESID,IGRAN
  209. ELSE
  210. IGRAN=0
  211. ENDIF
  212. *
  213. * Iteration successful
  214. *
  215. IF (RESID.LE.TOL) THEN
  216. IRET = 0
  217. *
  218. * Iteration fails
  219. *
  220. * ELSE IF (ITER.EQ.MAXIT.OR.IGRAN.GT.MGRAN) THEN
  221. ELSE IF (INMV.GE.MAXIT.OR.IGRAN.GT.MGRAN.OR.IWARN.GT.MWARN) THEN
  222. IRET = 1
  223. *
  224. * Do some more
  225. *
  226. ELSE
  227. RHO1 = RHO
  228. GOTO 10
  229. ENDIF
  230. *
  231. * Normal termination
  232. *
  233. 30 CONTINUE
  234. IF (IWARN.GT.0.AND.IMPR.GT.0) THEN
  235. WRITE(IOIMP,*) 'CG nearly failed ',IWARN
  236. $ ,' times.'
  237. ENDIF
  238. RETURN
  239. *
  240. * End of CG2
  241. *
  242. END
  243.  
  244.  

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