Télécharger prbcg.eso

Retour à la liste

Numérotation des lignes :

prbcg
  1. C PRBCG SOURCE GOUNAND 22/08/25 21:15:08 11434
  2. SUBROUTINE PRBCG(KMORS,KISA,KS2B,MATRIK,MAPREC,
  3. $ LRES,LNMV,INCX,ITER,INMV,RESID,KPREC,BRTOL,ICALRS,IDDOT,
  4. $ IMVEC,
  5. $ IMPR,IRET)
  6. IMPLICIT INTEGER(I-N)
  7. IMPLICIT REAL*8 (A-H,O-Z)
  8. C***********************************************************************
  9. C NOM : PRBCG
  10. C DESCRIPTION :
  11. C Préparation de la résolution d'un système linéaire Ax=b
  12. C par la méthode du bigradient conjugué (préconditionnée ou non)
  13. C
  14. C
  15. C Les 4 préconditionnements disponibles pour cette méthode
  16. C sont :
  17. C * Jacobi (diagonal)
  18. C * D-ILU (Diagonal Incomplete LU factorization)
  19. C * ILU(0) (Incomplete LU factorization
  20. C of level zero ie Choleski incomplet)
  21. C * MILU(0) (Modified ILU(0)) relaxé
  22. C 5 : préconditionnement ILUT (dual truncation)
  23. C 6 : préconditionnement ILUT2 (une variante du
  24. C précédent qui remplit mieux la mémoire et
  25. C fonctionne mieux quelquefois)
  26. C
  27. C Ces préconditionneurs sont supposés déjà construits (ie
  28. C factorisés), voir les subroutines MEIDIA, PRDILU, PRILU0,
  29. C PRMILU et PRILUT.
  30. C
  31. C Ce sous-programme est en fait une interface aux sous-programmes :
  32. C cg, cgd, cgdilu, cgilu0
  33. C qui sont en Fortran presque pur (pour raison de rapidité)
  34. C et effectuent la résolution du système linéaire
  35. C proprement dite.
  36. C
  37. C LANGAGE : ESOPE
  38. C AUTEUR : Stéphane GOUNAND (CEA/DRN/DMT/SEMT/TTMF)
  39. C mél : gounand@semt2.smts.cea.fr
  40. C REFERENCE (bibtex-like) :
  41. C @BOOK{templates,
  42. C AUTHOR={R.Barrett, M.Berry, T.F.Chan, J.Demmel, J.Donato,
  43. C J.Dongarra, V.Eijkhout, R.Pozo, C.Romine,
  44. C H. Van der Vorst},
  45. C TITLE={Templates for the Solution of Linear Systems :
  46. C Building Blocks for Iterative Methods},
  47. C PUBLISHER={SIAM}, YEAR={1994}, ADDRESS={Philadelphia,PA} }
  48. C -> URL : http://www.netlib.org/templates/Templates.html
  49. C***********************************************************************
  50. C APPELES (Calcul) : CG, CGD, CGDILU, CGILU0
  51. C***********************************************************************
  52. C ENTREES : MATRIK, MAPREC, KPREC, IMPR
  53. C ENTREES/SORTIES : INCX, ITER, RESID
  54. C SORTIES : LRES, IRET
  55. C CODE RETOUR (IRET) : 0 si ok
  56. C <0 si problème
  57. C >0 si l'algorithme appelé n'a pas convergé
  58. C MATRIK : pointeur sur segment MATRIK de l'include SMMATRIK
  59. C on pioche dedans les informations nécessaires
  60. C (différents pointeurs, nb. de ddl...)
  61. C MAPREC : pointeur sur segment MATRIK de l'include SMMATRIK
  62. C on va l'utiliser comme préconditionneur
  63. C KPREC : type de préconditionnement à effectuer
  64. C 0 : pas de préconditionnement
  65. C 1 : préconditionnement Jacobi (diagonal)
  66. C 2 : préconditionnement D-ILU
  67. C 3 : préconditionnement ILU(0) (Choleski incomplet)
  68. C 4 : préconditionnement MILU(0)
  69. C (Choleski incomplet modifié)
  70. C 5 : préconditionnement ILUT (dual truncation)
  71. C 6 : préconditionnement ILUT2 (une variante du
  72. C précédent qui remplit mieux la mémoire et
  73. C fonctionne mieux quelquefois)
  74. C IMPR : niveau d'impression (0..4)
  75. C INCX : pointeur sur segment IZA de l'include SMMATRIK
  76. C C'est le vecteur des inconnues.
  77. C En entrée, contient l'estimation x(0).
  78. C En sortie, la solution trouvée (si la méthode
  79. C a convergé).
  80. C ITER : type INTEGER.
  81. C En entrée, contient le nombre maximum
  82. C d'itérations à effectuer.
  83. C En sortie, contient le nombre d'itérations
  84. C réellement effectuées.
  85. C RESID : type REAL*8.
  86. C En entrée, la valeur maximale du résidu à
  87. C convergence de l'algorithme mesurée comme suit :
  88. C norme[L2](b - A*x) / norme[L2]( b )
  89. C En sortie, la valeur effective de cette mesure.
  90. C LRES : pointeur sur segment MLREEL (include SMLREEL)
  91. C contient l'historique de la valeur de RESID en
  92. C fonction du nombre d'itérations effectuées.
  93. C
  94. C***********************************************************************
  95. C VERSION : v1, 02/04/98, version initiale
  96. C HISTORIQUE : v1, 02/04/98, création
  97. C HISTORIQUE : 09/02/99 : ajout préconditionneur <> MATRIK
  98. C HISTORIQUE : 20/12/99 : interfaçage avec le nouveau cgilu0
  99. C (factorisations incomplètes stockées au format MSR (Modified Sparse
  100. C Row) (voir la doc de Sparskit version 2+)
  101. C http://www.cs.umn.edu/Research/arpa/SPARSKIT/sparskit.html
  102. C HISTORIQUE : 22/03/00 : ajout des préconditionneurs ILUT
  103. C HISTORIQUE : 10/02/06 : Nettoyage
  104. C HISTORIQUE :
  105. C***********************************************************************
  106. C Prière de PRENDRE LE TEMPS de compléter les commentaires
  107. C en cas de modification de ce sous-programme afin de faciliter
  108. C la maintenance !
  109. C***********************************************************************
  110. *
  111. * .. Includes et pointeurs associés ..
  112.  
  113. -INC PPARAM
  114. -INC CCOPTIO
  115. -INC SMLREEL
  116. INTEGER JG
  117. POINTEUR LRES.MLREEL
  118. -INC SMLENTI
  119. POINTEUR LNMV.MLENTI
  120. POINTEUR MAPREC.MATRIK
  121. POINTEUR KMORS.PMORS
  122. POINTEUR KISA.IZA
  123. POINTEUR KS2B.IZA
  124. POINTEUR INCX.IZA
  125. POINTEUR INVDIA.IZA
  126. POINTEUR ILUM.PMORS
  127. POINTEUR ILUI.IZA
  128. * .. Work Vectors for CG
  129. POINTEUR IR.IZA,IP.IZA,IQ.IZA,IZ.IZA
  130. POINTEUR IRTLD.IZA,IPTLD.IZA,IQTLD.IZA,IZTLD.IZA
  131. * .. Scalar Arguments ..
  132. INTEGER ITER, KPREC, IMPR, IRET
  133. REAL*8 RESID
  134. INTEGER NBVA,NJA,NTT,NTTPRE
  135. * .. Executable Statements ..
  136. *
  137. IRET = 0
  138. *
  139. * On récupère les paramètres
  140. *
  141. SEGACT MATRIK
  142. SEGACT MAPREC
  143. IF (KSYM.NE.0) THEN
  144. IF (IMPR.GT.0) THEN
  145. WRITE(IOIMP,*) 'MATRIK',MATRIK,'non symétrique : ',
  146. $ 'use Bi-CGSTAB instead !'
  147. * ENDIF
  148. * IRET=-1
  149. * GOTO 9999
  150. ENDIF
  151. ENDIF
  152. * Pour le préconditionneur
  153. ILUM =MAPREC.KIDMAT(6)
  154. ILUI =MAPREC.KIDMAT(7)
  155. IDMAT=MAPREC.KIDMAT(1)
  156. SEGACT IDMAT
  157. INVDIA=IDIAG
  158. SEGDES IDMAT
  159. SEGACT KMORS
  160. NTT =KMORS.IA(/1)-1
  161. * NJA =KMORS.JA(/1)
  162. SEGACT KISA
  163. SEGACT KS2B
  164. SEGACT INCX*MOD
  165. C Paramètres des préconditionnements diagonaux et D-ILU
  166. IF (KPREC.NE.0) THEN
  167. IF (KPREC.EQ.1.OR.KPREC.EQ.2) THEN
  168. C Est-il compatible avec la matrice ?
  169. SEGACT INVDIA
  170. NTTPRE=INVDIA.A(/1)
  171. IF (NTTPRE.NE.NTT) THEN
  172. WRITE(IOIMP,*) 'La matrice et son préconditionnement'
  173. WRITE(IOIMP,*) 'ne sont pas compatibles...'
  174. WRITE(IOIMP,*) 'NTT=',NTT
  175. WRITE(IOIMP,*) 'NTTPRE=',NTTPRE
  176. IRET=-2
  177. GOTO 9999
  178. ENDIF
  179. C Paramètres des préconditionnements ILU(0), MILU(0), ILUT et ILUT2
  180. ELSEIF (KPREC.GE.3.AND.KPREC.LE.10) THEN
  181. SEGACT ILUM
  182. SEGACT ILUI
  183. NTTPRE=ILUM.IA(/1)
  184. IF (NTTPRE.NE.NTT) THEN
  185. WRITE(IOIMP,*) 'La matrice et son préconditionnement',
  186. $ 'ne sont pas compatibles...'
  187. WRITE(IOIMP,*) 'NTT=',NTT,' NTTPRE=',NTTPRE
  188. IRET=-2
  189. GOTO 9999
  190. ENDIF
  191. ELSE
  192. WRITE(IOIMP,*) 'Erreur de programmation'
  193. GOTO 9999
  194. ENDIF
  195. ENDIF
  196. C
  197. C Initialisation de l'historique de convergence
  198. C
  199. JG=ITER+1
  200. SEGINI LNMV
  201. SEGINI LRES
  202. C
  203. C
  204. C Initialisation des vecteurs de travail pour le bigradient conjugué
  205. C
  206. NBVA=NTT
  207. SEGINI IR,IRTLD,IZ,IZTLD,IP,IPTLD
  208. IQ=IZ
  209. IQTLD=IZTLD
  210. C
  211. C Appel du bigradient conjugué
  212. C
  213. CALL BICG(KMORS,KISA,KS2B,INCX,
  214. $ KPREC,INVDIA,ILUM,ILUI,
  215. $ LRES,LNMV,
  216. $ IR,IRTLD,IZ,IZTLD,IP,IPTLD,IQ,IQTLD,
  217. $ ITER,INMV,RESID,BRTOL,ICALRS,IDDOT,IMVEC,IMPR,IRET)
  218. * Gestion du CTRL-C
  219. if (ierr.NE.0) return
  220. C
  221. C Désactivation
  222. C
  223. SEGSUP IR,IRTLD,IZ,IZTLD,IP,IPTLD
  224. JG=ITER+1
  225. SEGADJ LRES
  226. SEGDES LRES
  227. SEGADJ LNMV
  228. SEGDES LNMV
  229. IF (KPREC.EQ.1.OR.KPREC.EQ.2) THEN
  230. SEGDES INVDIA
  231. ELSEIF (KPREC.GE.3.AND.KPREC.LE.9) THEN
  232. SEGDES ILUI
  233. SEGDES ILUM
  234. ENDIF
  235. SEGDES INCX
  236. SEGDES KS2B
  237. SEGDES KISA
  238. SEGDES KMORS
  239. SEGDES MAPREC
  240. SEGDES MATRIK
  241. *
  242. * Normal termination
  243. *
  244. RETURN
  245. *
  246. * Format handling
  247. *
  248. *
  249. * Error handling
  250. *
  251. 9999 CONTINUE
  252. WRITE(IOIMP,*) 'An error was detected in prbcg.eso'
  253. RETURN
  254. *
  255. * End of PRBCG.
  256. *
  257. END
  258.  
  259.  

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