Télécharger kresll.eso

Retour à la liste

Numérotation des lignes :

kresll
  1. C KRESLL SOURCE CB215821 20/11/25 13:33:05 10792
  2. SUBROUTINE KRESLL()
  3. IMPLICIT INTEGER(I-N)
  4. IMPLICIT REAL*8 (A-H,O-Z)
  5. C*************************************************************************
  6. C Operateur KRES
  7. C
  8. C OBJET : Resoud l'equation AX=B par différentes méthodes
  9. C * directe (factorisation Choleski)
  10. C * itérative (Gradient conjugué, BiCGSTAB,
  11. C GMRES(m)) avec différents préconditionneurs
  12. C diagonal (Jacobi), D-ILU, ILU(0) (Choleski
  13. C incomplet), MILU(0)
  14. C SYNTAXE : CHPO3 = KRES MA1 'TYPI' TAB1
  15. C 'CLIM' CHPO1
  16. C 'SMBR' CHPO2
  17. C 'IMPR' VAL1 ;
  18. C Voir la notice pour plus de précisions.
  19. C
  20. C***********************************************************************
  21. C APPELES : KRES3, KRES4, KRES5
  22. C APPELES (E/S) : LIROBJ, ECROBJ, ERREUR, LIRMOT, LIRENT, LIRTAB,
  23. C ACME, ACMO, ACMM, ACMF, ECMO
  24. C APPELES (UTIL.) : QUETYP
  25. C APPELE PAR : KOPS
  26. C***********************************************************************
  27. C***********************************************************************
  28. C HISTORIQUE : 26/10/98 : prise en compte du cas particulier
  29. C où A est diagonale. C'est en fait assez pénible
  30. C car on utilise le CHPOINT comme structure de
  31. C données pour la matrice A et les vecteurs X,B,CLIM
  32. C HISTORIQUE : 09/02/99 : on peut utiliser le préconditionnement d'une
  33. C autre matrice que celle que l'on inverse
  34. C HISTORIQUE : 01/06/99 : on modifie la partie résolution itérative
  35. C en effet, lors de l'imposition des CL. de
  36. C Dirichlet, on changeait les valeurs de la
  37. C matrice Morse.
  38. C Ceci n'est pas bon lorsqu'on veut utiliser la
  39. C matrice assemblée pour plusieurs pas de temps.
  40. C On travaille maintenant sur une copie.
  41. C HISTORIQUE : 20/12/99 : reprogrammation de l'assemblage
  42. C Le nouvel assemblage n'est, pour l'instant effectif que
  43. C lorsqu'une méthode itérative est sélectionnée (-> fabrication
  44. C d'une matrice Morse). Le nouvel assemblage est plus performant
  45. C (temps de calcul, utilisation de la mémoire) et plus général (cas
  46. C particuliers à peu près supprimés) que le précédent.
  47. C HISTORIQUE : 05/01/00 : On ne supprime plus les 0.D0 de la matrice
  48. C assemblée (cf. clmors appelé par melim). Ceci pour ne plus avoir
  49. C perte éventuelle de symétrie de la matrice assemblée. Cela permet
  50. C aussi de ne plus dupliquer la matrice assemblée.
  51. C HISTORIQUE : 13/01/00 : Mise en conformité du solveur direct avec le
  52. C nouvel assemblage.
  53. C HISTORIQUE : 22/03/00 : Rajout des préconditionneurs ILUT
  54. C HISTORIQUE : 13/04/00 : Séparation Lecture des données
  55. C Ecriture des résultats (kres2)
  56. C Assemblage kres3
  57. C Méthode directe kres4
  58. C Méthodes itératives kres5
  59. C
  60. C***********************************************************************
  61. C Prière de PRENDRE LE TEMPS de compléter les commentaires
  62. C en cas de modification de ce sous-programme afin de faciliter
  63. C la maintenance !
  64. C***********************************************************************
  65.  
  66. -INC PPARAM
  67. -INC CCOPTIO
  68. -INC SMLREEL
  69. POINTEUR LRES.MLREEL
  70. -INC SMCHPOI
  71. POINTEUR KCLIM.MCHPOI
  72. POINTEUR KSMBR.MCHPOI
  73. POINTEUR MCHINI.MCHPOI
  74. POINTEUR MCHSOL.MCHPOI
  75. -INC SMTABLE
  76. POINTEUR MTINV.MTABLE
  77. * MATRIK est la matrice à inverser
  78. * MAPREC est la matrice dont on utilise le préconditionnement
  79. * MATASS est la matrice dont on utilise l'assemblage
  80. * pour préconditionner celui de MATRIK
  81. POINTEUR MAPREC.MATRIK
  82. POINTEUR MATASS.MATRIK
  83. C
  84. CHARACTER*8 TYPE
  85. * Paramètre m du GMRES(m)
  86. INTEGER IRSTRT
  87. INTEGER ITER
  88. REAL*8 BRTOL,RESID,RXMILU
  89. *
  90. REAL*8 XLFIL,XDTOL
  91. INTEGER KPREC
  92. INTEGER NMATRI
  93. INTEGER IP,KTYPI
  94. INTEGER IMPINV,IIMPR
  95. CHARACTER*4 MRENU,MMULAG
  96. INTEGER IMPR,IRET
  97. C
  98. C Définition des options
  99. C
  100. C 'IMPR' + un entier : niveau d'impression
  101. C 'TYPI' + une table contenant les informations
  102. C sur le type de résolution
  103. C 'CLIM' + un chpoint : conditions aux limites de Dirichlet
  104. C 'SMBR' + un chpoint : second membre (forces...)
  105. C 'PREC' + une matrice dont on va utiliser le préconditionnement
  106. C
  107. INTEGER NBM
  108. PARAMETER (NBM=4)
  109. CHARACTER*4 LMOT(NBM)
  110. C Tableau de correspondance : entier <-> mot
  111. C pour le type d'inversion
  112. INTEGER NBTYPI
  113. C INTEGER LNTYPI
  114. PARAMETER (NBTYPI=5)
  115. C PARAMETER (LNTYPI=18)
  116. C CHARACTER*(LNTYPI) LTYPI(NBTYPI)
  117. LOGICAL LBID,LVAL,LTIME,LMRENU,LDEPE,LLDEPE
  118. C
  119. C Initialisation des tableaux
  120. C
  121. DATA LMOT/'IMPR','TYPI','CLIM','SMBR'/
  122. C DATA LTYPI/'Choleski ',
  123. C $ 'Conjugate Gradient',
  124. C $ 'BiCGSTAB ',
  125. C $ 'BiCGSTAB(2) ',
  126. C $ 'GMRES(m) '/
  127. C
  128. TYPE='MATRIK'
  129. CALL LIROBJ(TYPE,MATRIK,1,IRET)
  130. IF(IRET.EQ.0) GOTO 9999
  131. c
  132. SEGACT MATRIK
  133. NMATRI=IRIGEL(/2)
  134. IF(NMATRI.EQ.0)THEN
  135. C% Résolution impossible : la matrice de RIGIDITE est vide
  136. CALL ERREUR(727)
  137. RETURN
  138. ENDIF
  139. SEGDES MATRIK
  140. c Valeurs par défaut reprise de prkres
  141. KSMBR=0
  142. MTINV=0
  143. IMPR=0
  144. KCLIM=0
  145. KTYPI=1
  146. MATASS=MATRIK
  147. MAPREC=MATRIK
  148. MRENU='SLOA'
  149. LMRENU=.FALSE.
  150. MMULAG='APR2'
  151. ISCAL=0
  152. IOUBL=0
  153. IMPINV=0
  154. MCHINI=0
  155. ITER=2000
  156. RESID=1.D-10
  157. BRTOL=1.D-40
  158. IRSTRT=50
  159. KPREC=3
  160. RXMILU=1.D0
  161. RXILUP=0.5D0
  162. XLFIL=2.D0
  163. XDTOL=-1.D0
  164. XSPIV=0.1D0
  165. LBCG=4
  166. ICALRS=0
  167. METASS=5
  168. LTIME=.FALSE.
  169. KTIME=0
  170. LDEPE=.TRUE.
  171. LLDEPE=.FALSE.
  172. IDDOT=0
  173. IMVEC=2
  174. C
  175. 1 CONTINUE
  176. CALL LIRMOT(LMOT,NBM,IP,0)
  177. IF(IP.EQ.0) GOTO 2
  178. GOTO (11,12,13,14),IP
  179. C Lecture du niveau d'impression
  180. C option 'IMPR'
  181. 11 CONTINUE
  182. c write(6,*)' option IMPR'
  183. CALL LIRENT(IMPR,1,IRET)
  184. IF (IRET.EQ.0) GOTO 9999
  185. GO TO 1
  186. C Lecture de la table pour les méthodes d'inversion
  187. C option 'TYPI'
  188. 12 CONTINUE
  189. CALL LIRTAB('METHINV',MTINV,1,IRET)
  190. IF(IRET.EQ.0)RETURN
  191. CALL LIRENT(KTYPI,1,IRET)
  192. C KTYPI : Méthode d'inversion du système (cf. LTYPI)
  193. C 1 : résolution directe (Choleski)
  194. C 2 : Gradient Conjugué
  195. C 3 : Bi-Gradient Conjugué Stabilisé (BiCGSTAB)
  196. C 4 : BiCGSTAB(2)
  197. C 5 : GMRES(m)
  198. C CALL ACME(MTINV,'TYPINV',KTYPI)
  199. IF (IERR.NE.0) GOTO 9999
  200. C write(6,*)' Ktypi=',ktypi
  201. IF (KTYPI.LT.1.OR.KTYPI.GT.NBTYPI) THEN
  202. WRITE(IOIMP,*) 'L''indice TYPINV ',
  203. $ 'is out of range (1..',NBTYPI,') :',KTYPI
  204. GOTO 9999
  205. ENDIF
  206. GO TO 1
  207. C Lecture du chpoint de conditions aux limites
  208. C option CLIM
  209. 13 CONTINUE
  210. CALL QUETYP(TYPE,1,IRET)
  211. IF(IRET.EQ.0) GOTO 9999
  212. IF(TYPE.EQ.'CHPOINT')THEN
  213. CALL LIROBJ(TYPE,KCLIM,1,IRET)
  214. ELSE
  215. WRITE(IOIMP,*) 'Erreur lecture C.lim'
  216. GOTO 9999
  217. ENDIF
  218. GO TO 1
  219. C Lecture du chpoint second membre
  220. C option SMBR
  221. 14 CONTINUE
  222. CALL QUETYP(TYPE,1,IRET)
  223. IF(IRET.EQ.0)RETURN
  224. IF(TYPE.EQ.'CHPOINT')THEN
  225. CALL LIROBJ('CHPOINT',KSMBR,1,IRET)
  226. c write(6,*)' lecture 2 membre '
  227. ELSE
  228. WRITE(IOIMP,*) 'Erreur lecture second membre'
  229. GOTO 9999
  230. ENDIF
  231. GOTO 1
  232. 2 CONTINUE
  233. C write(6,*)' On a pas finit les lectures '
  234. c NAT=2
  235. c NSOUPO=0
  236. c SEGINI MCHPO1
  237. c MCHPO1.JATTRI(1)=2
  238. c CALL ECROBJ('CHPOINT',MCHPO1)
  239. c return
  240. C
  241. C Matrice dont on utilise l'assemblage
  242. C (éventuellement identique à la matrice transmise)
  243. C
  244. TYPE='MATRIK '
  245. C CALL ACMO(MTINV,'MATASS',TYPE,MATASS)
  246. C IF (IERR.NE.0) GOTO 9999
  247. CALL LIROBJ('MATRIK',MATASS,1,IRET)
  248. IF(IRET.EQ.0) GOTO 9999
  249. C - Pour l'assemblage : type de renumérotation
  250. C CALL ACMM(MTINV,'TYRENU',MRENU)
  251. C IF (IERR.NE.0) GOTO 9999
  252. CALL LIRCHA(MRENU,1,LCHAR)
  253. IF(LCHAR.EQ.0)GO TO 9999
  254. C - Pour l'assemblage : prise en compte des mult.lag
  255. C CALL ACMM(MTINV,'PCMLAG',MMULAG)
  256. C IF (IERR.NE.0) GOTO 9999
  257. CALL LIRCHA(MMULAG,1,LCHAR)
  258. IF(LCHAR.EQ.0)GO TO 9999
  259. C - Pour l'assemblage : scaling de la matrice
  260. C CALL ACME(MTINV,'SCALING',ISCAL)
  261. C IF (IERR.NE.0) GOTO 9999
  262. ISCAL=0
  263. C Niveau d'impression pour la partie résolution itérative
  264. c CALL ACME(MTINV,'IMPINV',IMPINV)
  265. c IF (IERR.NE.0) GOTO 9999
  266. c IMPR=MAX(IMPR,IMPINV)
  267. IMPR=0
  268. C
  269. C Assemblage proprement dit
  270. C
  271. IIMPR=0
  272. * SG 2016/02/09 : non à la ligne suivante : il faut que METASS soit
  273. * égale à 5 (voir remarque dans makpr2)
  274. * METASS=4
  275. CALL KRES3(MATRIK,MATASS,MRENU,MMULAG,METASS,
  276. $ KTIME,LTIME,
  277. $ IIMPR,IRET)
  278. IF (IRET.NE.0) GOTO 9999
  279. *! WRITE(IOIMP,*) 'Aprés assemblage'
  280. *! CALL ECRENT(5)
  281. *! CALL ECROBJ('MATRIK ',MATRIK)
  282. *! CALL PRLIST
  283. *
  284. * "Oubli" des valeurs des matrice élémentaires
  285. * On met les tableaux de LIZAFM à 0 => à MENAGE de les supprimmer
  286. * si besoin est.
  287. *
  288. CALL OUBIMA(MATRIK,IMPR,IRET)
  289. IF (IRET.NE.0) GOTO 9999
  290. C
  291. C Méthode directe
  292. C
  293. IF (KTYPI.EQ.1) THEN
  294. CALL KRES4(MATRIK,KCLIM,KSMBR,
  295. $ ISCAL,
  296. $ MCHSOL,
  297. $ IMPR,IRET)
  298. IF (IRET.NE.0) GOTO 9999
  299. C
  300. C Méthodes itératives
  301. C
  302. ELSEIF(KTYPI.GE.2.AND.KTYPI.LE.NBTYPI)THEN
  303. C Paramètres pour la méthode itérative proprement dite
  304. C - Champoint d'initialisation de la méthode
  305. C (i.e. estimation de l'inconnue)
  306. TYPE='CHPOINT '
  307. C CALL ACMO(MTINV,'XINIT',TYPE,MCHINI)
  308. C IF (IERR.NE.0) GOTO 9999
  309. CALL LIROBJ('CHPOINT',MCHINI,1,IRET)
  310. IF(IRET.EQ.0) GOTO 9999
  311. C - Nombre maxi d'itérations à effectuer
  312. C CALL ACME(MTINV,'NITMAX',ITER)
  313. C IF (IERR.NE.0) GOTO 9999
  314. CALL LIRENT(ITER,1,IRET)
  315. IF(IRET.EQ.0) GOTO 9999
  316. C - Norme maxi (L2 normé par le second membre) du résidu
  317. C CALL ACMF(MTINV,'RESID',RESID)
  318. C IF (IERR.NE.0) GOTO 9999
  319. CALL LIRREE(RESID,1,IRET)
  320. IF(IRET.EQ.0) GOTO 9999
  321. C - 'Breakdown tolerance' pour les méthodes de type
  322. C BiCGSTAB. Si un certain produit scalaire de vecteurs
  323. C "direction" est inférieur à cette tolérance, la
  324. C méthode s'arrete.
  325. IF (KTYPI.EQ.3.OR.KTYPI.EQ.4) THEN
  326. C CALL ACMF(MTINV,'BCGSBTOL',BRTOL)
  327. C IF (IERR.NE.0) GOTO 9999
  328. CALL LIRREE(BRTOL,1,IRET)
  329. IF(IRET.EQ.0) GOTO 9999
  330. ENDIF
  331. C - Paramètre m de redémarrage pour GMRES(m)
  332. C i.e. nombre max. de vecteurs par rapport auxquels on
  333. C orthogonalise le résidu
  334. IF (KTYPI.EQ.5) THEN
  335. C CALL ACME(MTINV,'GMRESTRT',IRSTRT)
  336. C IF (IERR.NE.0) GOTO 9999
  337. CALL LIRENT(IRSTRT,1,IRET)
  338. IF(IRET.EQ.0) GOTO 9999
  339. ENDIF
  340. C Paramètres pour les différents préconditionneurs
  341. C Matrice dont on utilise le préconditionneur
  342. C (éventuellement identique à la matrice transmise).
  343. C TYPE='MATRIK '
  344. C CALL ACMO(MTINV,'MAPREC',TYPE,MAPREC)
  345. C IF (IERR.NE.0) GOTO 9999
  346. CALL LIROBJ('MATRIK',MAPREC,1,IRET)
  347. IF(IRET.EQ.0) GOTO 9999
  348. C - Type de préconditionnement : (cf. LPREC)
  349. C 0 : pas de préconditionnement
  350. C 1 : préconditionnement par la diagonale
  351. C 2 : préconditionnement D-ILU
  352. C 3 : préconditionnement ILU(0) (Choleski)
  353. C 4 : préconditionnement MILU(0) (Choleski modifié)
  354. C 5 : préconditionnement ILUT (dual truncation)
  355. C 6 : préconditionnement ILUT2 (une variante du
  356. C précédent qui remplit mieux la mémoire et
  357. C fonctionne mieux quelquefois)
  358. C CALL ACME(MTINV,'PRECOND',KPREC)
  359. C IF (IERR.NE.0) GOTO 9999
  360. CALL LIRENT(KPREC,1,IRET)
  361. IF(IRET.EQ.0) GOTO 9999
  362. C - Paramètre de relaxation pour le préconditionnement
  363. C MILU(0) compris entre 0.D0 et 1.D0
  364. C S'il est égal à 0, on se ramène à ILU(0)
  365. C S'il est égal à 1, MILU(0) est dit non relaxé
  366. IF (KPREC.EQ.4) THEN
  367. C CALL ACMF(MTINV,'MILURELX',RXMILU)
  368. C IF (IERR.NE.0) GOTO 9999
  369. CALL LIRREE(RXMILU,1,IRET)
  370. IF(IRET.EQ.0) GOTO 9999
  371. ENDIF
  372. C - Pour une méthode ILUT, on a les deux indices suivant :
  373. C * ILUTLFIL : encombrement maximal (approximatif) du
  374. C préconditionneur, par rapport à la matrice.
  375. C * ILUTDTOL : "drop tolerance" pour le préconditionneur.
  376. C i.e. en-dessous de cette valeur relative, les
  377. C termes de la factorisation incomplète seront
  378. C oubliés.
  379. IF (KPREC.GE.5.AND.KPREC.LE.9) THEN
  380. C CALL ACMF(MTINV,'ILUTLFIL',XLFIL)
  381. C IF (IERR.NE.0) GOTO 9999
  382. CALL LIRREE(XLFIL,1,IRET)
  383. IF(IRET.EQ.0) GOTO 9999
  384. C CALL ACMF(MTINV,'ILUTDTOL',XDTOL)
  385. C IF (IERR.NE.0) GOTO 9999
  386. CALL LIRREE(XDTOL,1,IRET)
  387. IF(IRET.EQ.0) GOTO 9999
  388. ENDIF
  389. C - Pour une méthode ILUTP, on a les deux indices suivant :
  390. C * ILUTPPIV (type REEL) : sensibilité du pivoting
  391. C IF (KPREC.GE.7.AND.KPREC.LE.9) THEN
  392. C CALL ACMF(MTINV,'ILUTPPIV',XSPIV)
  393. C ENDIF
  394. CALL KRES5(MATRIK,KCLIM,KSMBR,KTYPI,
  395. $ MCHINI,ITER,RESID,
  396. $ BRTOL,IRSTRT,LBCG,ICALRS,
  397. $ MAPREC,KPREC,
  398. $ RXMILU,RXILUP,XLFIL,XDTOL,XSPIV,
  399. $ ISCAL,
  400. $ KTIME,LTIME,IDDOT,IMVEC,
  401. $ MCHSOL,LRES,LNMV,ICVG,
  402. $ IMPR,IRET)
  403. IF (IRET.NE.0) GOTO 9999
  404. CALL ECMO(MTINV,'CONVINV','LISTREEL',LRES)
  405. ELSE
  406. WRITE(IOIMP,*) 'KTYPI=',KTYPI,'out of range (1..',NBTYPI,')'
  407. GOTO 9999
  408. ENDIF
  409. CALL ECROBJ('CHPOINT ',MCHSOL)
  410. *
  411. * Normal termination
  412. *
  413. RETURN
  414. *
  415. * Format handling
  416. *
  417. *
  418. * Error handling
  419. *
  420. 9999 CONTINUE
  421. WRITE(IOIMP,*) 'An error was detected in kresll.eso'
  422. * 153 2
  423. * Opération illicite dans ce contexte
  424. CALL ERREUR(153)
  425. RETURN
  426. *
  427. * End of KRES2
  428. *
  429. END
  430.  
  431.  
  432.  
  433.  
  434.  
  435.  
  436.  
  437.  
  438.  
  439.  
  440.  
  441.  
  442.  
  443.  
  444.  
  445.  
  446.  
  447.  
  448.  
  449.  
  450.  
  451.  
  452.  
  453.  
  454.  
  455.  
  456.  
  457.  

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