Télécharger kres2.eso

Retour à la liste

Numérotation des lignes :

kres2
  1. C KRES2 SOURCE GOUNAND 25/03/24 21:15:04 12216
  2. SUBROUTINE KRES2()
  3. IMPLICIT REAL*8 (A-H,O-Z)
  4. IMPLICIT INTEGER (I-N)
  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, ECME
  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 HISTORIQUE : 06/04/04 : Renumérotation (double mult.)
  60. C HISTORIQUE : 06/04/04 : Scaling
  61. C HISTORIQUE : 08/04/04 : ajout ILUTP
  62. C HISTORIQUE : 09/02/06 : ajout nb prod matrice vecteur (NMATVEC)
  63. C simplification du code
  64. C HISTORIQUE : 22/02/06 : Gros nettoyage au niveau de l'entrée des
  65. C arguments (Nouvelle syntaxe)
  66. C HISTORIQUE : 08/2011 : En vue de la suppression de l'objet MATRIK
  67. C on utilise l'assemblage de RESOU.
  68. C HISTORIQUE : 04/2019 : remplacement de NOEL par NELIM
  69. C Idéalement, il faudrait reprendre ce que Pierre a fait dans
  70. C RESOU dans les fiches 10[0,1]?? et avec MREM.En vue de la
  71. C suppression de l'objet MATRIK
  72. C
  73. C***********************************************************************
  74. C Prière de PRENDRE LE TEMPS de compléter les commentaires
  75. C en cas de modification de ce sous-programme afin de faciliter
  76. C la maintenance !
  77. C***********************************************************************
  78.  
  79. -INC PPARAM
  80. -INC CCOPTIO
  81. -INC SMLREEL
  82. POINTEUR LRES.MLREEL
  83. -INC SMLENTI
  84. POINTEUR LNMV.MLENTI
  85. -INC SMCHPOI
  86. POINTEUR KCLIM.MCHPOI
  87. POINTEUR KSMBR.MCHPOI
  88. POINTEUR MCHINI.MCHPOI
  89. POINTEUR MCHSOL.MCHPOI
  90. -INC SMTABLE
  91. POINTEUR MTINV.MTABLE
  92. POINTEUR KTIME.MTABLE
  93. DIMENSION ITTIME(4)
  94. CHARACTER*16 CHARI
  95. * MATRIK est la matrice à inverser
  96. * MAPREC est la matrice dont on utilise le préconditionnement
  97. * MATASS est la matrice dont on utilise l'assemblage
  98. * pour préconditionner celui de MATRIK
  99. POINTEUR MAPREC.MATRIK
  100. POINTEUR MATASS.MATRIK
  101. *STAT -INC SMSTAT
  102. C
  103. CHARACTER*8 TYPE
  104. * Paramètre m du GMRES(m)
  105. INTEGER RESTRT
  106. INTEGER ITER
  107. REAL*8 BRTOL,RESID,RXMILU,RXILUP
  108. *
  109. REAL*8 XLFIL,XDTOL
  110. INTEGER KPREC
  111. INTEGER NMATRI
  112. INTEGER IP,KTYPI
  113. INTEGER IMPINV,IIMPR
  114. CHARACTER*4 MRENU,MMULAG
  115. LOGICAL LRIG,LTIME,LDETR,LDEPE,LASRIG,LDMULT,LOGII
  116. INTEGER IMPR,IRET
  117. C
  118. C Lecture des arguments et mise à défaut des optionnels ()
  119. C
  120. C MATRIK : La matrice lue en entrée au format MATRIK
  121. C MTINV : L'éventuelle table d'inversion (obsolète)
  122. C IMPR : Niveau d'impression solveur direct
  123. C KCLIM : Chpoint éventuel de conditions aux limites de Dirichlet
  124. C KSMBR : Chpoint second membre
  125. C KTYPI : Type de méthode de résolution
  126. C MATASS : Matrice utilisée pour préconditionner l'assemblage
  127. C MAPREC : Matrice utilisée pour préconditionner la construction du
  128. C préconditionneur
  129. C MRENU : Type de renumérotation
  130. C MMULAG : Placement des multiplicateurs de Lagrange
  131. C ISCAL : Scaling de la matrice
  132. C IOUBL : Oubli des matrices élémentaires ?
  133. C IMPINV : Niveau d'impression solveur itératif
  134. C MCHINI : Chpoint estimation de l'inconnue
  135. C ITER : Nombre maxi d'itérations à effectuer
  136. C RESID : Norme L2 maxi du résidu
  137. C BRTOL : Breakdown tolerance pour les méthodes de type Bi-CG
  138. C IRSTRT : Paramètre m de redémarrage pour GMRES
  139. C KPREC : Type du préconditionneur
  140. C RXMILU : Paramètre de relaxation pour MILU(0)
  141. C RXILUP : Paramètre de filtre pour ILU(0)-PV
  142. C XLFIL : Paramètre de remplissage pour les préconditionneurs ILUT
  143. C XDTOL : Drop tolerance pour les préconditionneurs ILUT
  144. C XSPIV : Sensibilité du pivoting pour les préconditionneurs ILUTP
  145. C LBCG : le l dans BicgStab(l)
  146. C ICALRS : façon de calculer le résidu
  147. C METASS : Méthode d'assemblage
  148. C LTIME : construit une table avec des statistiques temporelles
  149. C si vrai
  150. C LDEPE : élimine les dépendances si VRAI
  151. C et matrice d'entrée RIGIDITE
  152. C IDDOT : 0 => utilisation du produit scalaire normal dans les
  153. C méthodes itératives
  154. C 1 => utilisation du produit scalaire compensé
  155. * IMPR=4
  156.  
  157. IVALI=0
  158. XVALI=0.D0
  159. LOGII=.FALSE.
  160. IRETI=0
  161. XVALR=0.D0
  162. *inutile IOBRE=0
  163. IRETR=0
  164.  
  165. IMPR=0
  166. * WRITE(IOIMP,*) 'coucou kres2'
  167. *
  168. *STAT CALL INMSTA(MSTAT,IMPR)
  169. *
  170. CALL PRKRES(MATRIK,MTINV,IMPR,KCLIM,KSMBR,KTYPI,MATASS,MAPREC,
  171. $ MRENU,MMULAG,ISCAL,INORMU,IOUBL,IMPINV,MCHINI,ITER,RESID
  172. $ ,BRTOL,IRSTRT,KPREC,RXMILU,RXILUP,XLFIL,XDTOL,XSPIV,LBCG
  173. $ ,ICALRS,METASS,LTIME,LDEPE,MRIGID,IDDOT,IMVEC,IORINC,IRET)
  174. IF (IRET.NE.0) GOTO 9999
  175. IMPR=MAX(IMPR,IMPINV)
  176. *
  177. * Préparation de la matrice et du second membre
  178. * suivant les cas
  179. *
  180. * LASRIG=.TRUE. on utilise l'assemblage de RESOU
  181. * LDMULT=.TRUE. on dédouble les multiplicateurs de Lagrange
  182. * LDEPE=.TRUE. On fait l'élimination des relations
  183.  
  184. LASRIG=(METASS.EQ.6)
  185. * Pour l'instant, il faut toujours dédoubler les multiplicateurs
  186. * quand on assemble avec l'assemblage de RESO car le traitement des
  187. * multiplicateurs dans ldmt1 l'impose (simple multiplicateur non
  188. * prévu)
  189. LDMULT=LASRIG
  190. MRIGI0=MRIGID
  191. * Nouveau separm gère la non élimination avec NOEL
  192. IF (LDEPE) then
  193. * noel=0
  194. nelim=1
  195. ELSE
  196. * noel=1
  197. nelim=0
  198. ENDIF
  199. *dbg write(ioimp,*) 'LASRIG=',LASRIG
  200. *dbg write(ioimp,*) 'LDMULT=',LDMULT
  201. *dbg write(ioimp,*) 'LDEPE=',LDEPE,' NELIM=',NELIM
  202.  
  203. IF (MRIGID.NE.0) THEN
  204. *old IF (LDEPE) THEN
  205. IF (IMPR.GT.1) THEN
  206. WRITE(IOIMP,*)
  207. $ '*** ELIMINATION DES MULTIPLICATEURS DE LAGRANGE (LX) KRES6'
  208. ENDIF
  209. KSMBR0=KSMBR
  210. *old CALL KRES6(MRIGID,KSMBR,IDEPE,
  211. *old $ MRIGIC,KSMBRC,KSMBR0,KSMBR1)
  212. CALL KRES6(MRIGID,KSMBR,LDMULT,NELIM,
  213. $ MRIGIC,KSMBRC,KSMBR1)
  214. IF (IERR.NE.0) RETURN
  215. MRIGID=MRIGIC
  216. KSMBR=KSMBRC
  217. *old ENDIF
  218. IF (LASRIG) THEN
  219. IF (IMPR.GT.1) THEN
  220. WRITE(IOIMP,*)
  221. $ '*** ASSEMBLAGE RENUMEROTATION RESOU KRES8'
  222. ENDIF
  223. * Gestion de la normalisation
  224. NORICO=NORINC
  225. NORIDO=NORIND
  226. IF (ISCAL.EQ.0) THEN
  227. NORINC=0
  228. NORIND=0
  229. ELSE
  230. NORINC=-1
  231. NORIND=0
  232. ENDIF
  233. * Gestion de la renumérotation
  234. NUCROO=NUCROU
  235. IF (MRENU.EQ.'RIEN') THEN
  236. NUCROU=-1
  237. * La renumérotation sera en fait Reverse Cuthill-McKee dans NUMOPT
  238. ELSEIF (MRENU.EQ.'SLOA') THEN
  239. NUCROU=1
  240. * La renumérotation sera en fait Nested Dissection dans NUMOPT
  241. ELSEIF (MRENU.EQ.'GIPR'.OR.MRENU.EQ.'GIBA') THEN
  242. NUCROU=2
  243. ELSE
  244. WRITE(IOIMP,*) 'MRENU=',MRENU
  245. CALL ERREUR(5)
  246. RETURN
  247. ENDIF
  248. CALL KRES8(MRIGID,KSMBR,INORMU,
  249. $ KTYPI,ITER,RESID,ICALRS,IRSTRT,LBCG,BRTOL,IDDOT,IMVEC,
  250. $ IORINC,
  251. $ KPREC,RXMILU,RXILUP,XLFIL,XDTOL,XSPIV,
  252. $ KTIME,LTIME,
  253. $ MCHSOL,LRES,LNMV,ICVG,IMPR)
  254. IF (IERR.NE.0) RETURN
  255. * Gestion de la normalisation
  256. NORINC=NORICO
  257. NORIND=NORIDO
  258. NUCROU=NUCROO
  259. IF (LTIME) CALL ECROBJ('TABLE ',KTIME)
  260. IF (MTINV.NE.0) THEN
  261. CALL ECME(MTINV,'CVGOK',ICVG)
  262. IF (LRES.NE.0) CALL ECMO(MTINV,'CONVINV','LISTREEL',LRES)
  263. IF (LNMV.NE.0) CALL ECMO(MTINV,'NMATVEC','LISTENTI',LNMV)
  264. ENDIF
  265. * On décondense si nécessaire
  266. *
  267. * write (6,*) ' resou - mchsol '
  268. * call ecchpo(mchsol,0)
  269. * call mucpri(mchsol,mrigid,iresi)
  270. * write (6,*) ' kres - iresi '
  271. * call ecchpo(iresi,0)
  272. * WRITE(IOIMP,*) 'Avant KRES7'
  273. *old IF (MRIGI0.NE.0.AND.LDEPE) THEN
  274. IF (MRIGI0.NE.0) THEN
  275. MSOLC=MCHSOL
  276. IDTARG=1
  277. * On détruit MSOLC dans KRES7
  278. CALL KRES7(MSOLC,MRIGI0,KSMBR0,KSMBR1,IDTARG,
  279. $ MCHSOL)
  280. IF (IERR.NE.0) RETURN
  281. ENDIF
  282. CALL ECROBJ('CHPOINT ',MCHSOL)
  283. RETURN
  284. ELSE
  285. IF (IMPR.GT.1) THEN
  286. WRITE(IOIMP,*)
  287. $ '*** TRANSFORMATION RIGIDITE -> MATRIK'
  288. WRITE(IOIMP,*)
  289. $ '*** ASSEMBLAGE RENUMEROTATION KRES2'
  290. ENDIF
  291. CALL ECROBJ('RIGIDITE',MRIGID)
  292. CALL RIMA
  293. IF (IERR.NE.0) GOTO 9999
  294. CALL MACHI2(1)
  295. IF (IERR.NE.0) GOTO 9999
  296. CALL LIROBJ('MATRIK',MATRIK,1,IRET)
  297. IF(IRET.EQ.0) GOTO 9999
  298.  
  299. * Changement des noms d'inconnues du second membre
  300. IF (KSMBR.NE.0) THEN
  301. CALL ECROBJ('CHPOINT ',KSMBR)
  302. CALL MACHI2(1)
  303. CALL LIROBJ('CHPOINT ',KSMBR,1,IRET)
  304. IF (IRET.EQ.0) GOTO 9999
  305. ENDIF
  306. ENDIF
  307. * write (6,*) ' le vecteur 2'
  308. * call ecchpo(ksmbr,0)
  309. * write (6,*) ' la matrice 2'
  310. * call ecrent(5)
  311. * call ecmatk(matrik)
  312. ENDIF
  313. *
  314. SEGACT MATRIK
  315. NMATRI=IRIGEL(/2)
  316. IF(NMATRI.EQ.0)THEN
  317. C% Résolution impossible : la matrice de RIGIDITE est vide
  318. CALL ERREUR(727)
  319. RETURN
  320. ENDIF
  321. SEGDES MATRIK
  322. IF (MATASS.EQ.0) MATASS=MATRIK
  323. IF (MAPREC.EQ.0) MAPREC=MATRIK
  324. * WRITE(IOIMP,*) 'Sortie de prkres'
  325. * WRITE(IOIMP,*) 'IOUBL=',IOUBL
  326. C
  327. IF (LTIME) THEN
  328. CALL CRTABL(KTIME)
  329. call timespv(ittime,oothrd)
  330. ITI1=(ITTIME(1)+ITTIME(2))/10
  331. ELSE
  332. KTIME=0
  333. ENDIF
  334. *STAT CALL PRMSTA('Lectures',MSTAT,IMPR)
  335. *
  336. C
  337. C Assemblage proprement dit
  338. C
  339. IIMPR=0
  340. CALL KRES3(MATRIK,MATASS,MRENU,MMULAG,METASS,
  341. $ KTIME,LTIME,
  342. $ IIMPR,IRET)
  343. * Gestion du CTRL-C
  344. if (ierr.NE.0) return
  345. IF (IRET.NE.0) GOTO 9999
  346. *! WRITE(IOIMP,*) 'Aprés assemblage'
  347. *! CALL ECRENT(5)
  348. *! CALL ECROBJ('MATRIK ',MATRIK)
  349. *! CALL PRLIST
  350. IF (LTIME) THEN
  351. call timespv(ittime,oothrd)
  352. ITI2=(ITTIME(1)+ITTIME(2))/10
  353. ENDIF
  354. *STAT CALL PRMSTA('Assemblage',MSTAT,IMPR)
  355. *
  356. * "Oubli" des valeurs des matrice élémentaires
  357. * On met les tableaux de LIZAFM à 0 => à MENAGE de les supprimmer
  358. * si besoin est.
  359. *
  360. IOUBD=MOD(IOUBL,10)
  361. *! WRITE(IOIMP,*) 'IOUBD=',IOUBD
  362. IF (IOUBD.EQ.1) THEN
  363. CALL OUBIMA(MATRIK,IMPR,IRET)
  364. IF (IRET.NE.0) GOTO 9999
  365. IF (IMPR.GT.2) THEN
  366. WRITE(IOIMP,*) 'Oubli des mat. elem.'
  367. ENDIF
  368. ELSEIF (IOUBD.EQ.2) THEN
  369. call ooohor(0)
  370. SEGACT MATRIK*MOD
  371. LDETR=.FALSE.
  372. NMATE=IRIGEL(/2)
  373. DO IMATE=1,NMATE
  374. IMATRI=IRIGEL(4,IMATE)
  375. SEGACT IMATRI*MOD
  376. NSOUM =LIZAFM(/1)
  377. NTOTIN=LIZAFM(/2)
  378. DO ITOTIN=1,NTOTIN
  379. DO ISOUM=1,NSOUM
  380. IZAFM=LIZAFM(ISOUM,ITOTIN)
  381. IF (IZAFM.NE.0) THEN
  382. LDETR=.TRUE.
  383. SEGSUP IZAFM
  384. LIZAFM(ISOUM,ITOTIN)=0
  385. ENDIF
  386. ENDDO
  387. ENDDO
  388. SEGDES IMATRI
  389. ENDDO
  390. IF (IMPR.GT.2.AND.LDETR) THEN
  391. WRITE(IOIMP,*) 'Destruction des mat. elem.'
  392. ENDIF
  393. ELSEIF (IOUBD.NE.0) THEN
  394. WRITE(IOIMP,*) 'IOUBL=',IOUBL, ' non prevu'
  395. GOTO 9999
  396. ENDIF
  397. *STAT CALL PRMSTA('Oubli',MSTAT,IMPR)
  398. *! WRITE(IOIMP,*) 'Aprés oubli'
  399. C
  400. C Méthode directe
  401. C
  402. IF (KTYPI.EQ.1) THEN
  403. CALL KRES4(MATRIK,KCLIM,KSMBR,
  404. $ ISCAL,
  405. $ MCHSOL,
  406. $ IMPR,IRET)
  407. if (ierr.ne.0) return
  408. IF (IRET.NE.0) GOTO 9999
  409. *STAT CALL PRMSTA('Methode directe',MSTAT,IMPR)
  410. C
  411. C Méthodes itératives
  412. C
  413. ELSE
  414. C
  415. CALL KRES5(MATRIK,KCLIM,KSMBR,KTYPI,
  416. $ MCHINI,ITER,RESID,
  417. $ BRTOL,IRSTRT,LBCG,ICALRS,
  418. $ MAPREC,KPREC,
  419. $ RXMILU,RXILUP,XLFIL,XDTOL,XSPIV,
  420. $ ISCAL,
  421. $ KTIME,LTIME,IDDOT,IMVEC,IORINC,
  422. $ MCHSOL,LRES,LNMV,ICVG,
  423. $ IMPR,IRET)
  424. if (ierr.ne.0) return
  425. IF (IRET.NE.0) GOTO 9999
  426. *STAT CALL PRMSTA('Methode iterative',MSTAT,IMPR)
  427. IF (MTINV.NE.0) THEN
  428. CALL ECME(MTINV,'CVGOK',ICVG)
  429. CALL ECMO(MTINV,'CONVINV','LISTREEL',LRES)
  430. CALL ECMO(MTINV,'NMATVEC','LISTENTI',LNMV)
  431. ENDIF
  432. ENDIF
  433. IF (LTIME) THEN
  434. call timespv(ittime,oothrd)
  435. ITI3=(ITTIME(1)+ITTIME(2))/10
  436. CHARI='KRES ASS+RENU'
  437. CALL ECCTAB(KTIME,'MOT ',IVALI,XVALI,CHARI,LOGII,IRETI,
  438. $ 'ENTIER ',ITI2-ITI1,XVALR,CHARR,LOGIR,IRETR)
  439. CHARI='KRES PRE+RESO'
  440. CALL ECCTAB(KTIME,'MOT ',IVALI,XVALI,CHARI,LOGII,IRETI,
  441. $ 'ENTIER ',ITI3-ITI2,XVALR,CHARR,LOGIR,IRETR)
  442. CHARI='KRES TOTAL '
  443. CALL ECCTAB(KTIME,'MOT ',IVALI,XVALI,CHARI,LOGII,IRETI,
  444. $ 'ENTIER ',ITI3-ITI1,XVALR,CHARR,LOGIR,IRETR)
  445. SEGDES KTIME
  446. CALL ECROBJ('TABLE ',KTIME)
  447. ENDIF
  448. IOUBE=IOUBL/10
  449. *! WRITE(IOIMP,*) 'IOUBE=',IOUBE
  450. IF (IOUBE.GE.1) THEN
  451. call ooohor(0)
  452. SEGACT MATRIK*MOD
  453. IF (IOUBE.EQ.2) THEN
  454. PMORS=KIDMAT(4)
  455. IF (PMORS.NE.0) THEN
  456. IF (IMPR.GT.2) THEN
  457. WRITE(IOIMP,*) 'Destruction du profil morse'
  458. ENDIF
  459. SEGSUP PMORS
  460. KIDMAT(4)=0
  461. ENDIF
  462. ENDIF
  463. IZA=KIDMAT(5)
  464. IF (IZA.NE.0) THEN
  465. IF (IMPR.GT.2) THEN
  466. WRITE(IOIMP,*) 'Destruction des valeurs morses'
  467. ENDIF
  468. SEGSUP IZA
  469. KIDMAT(5)=0
  470. ENDIF
  471. PMORS=KIDMAT(6)
  472. IF (PMORS.NE.0) THEN
  473. IF (IMPR.GT.2) THEN
  474. WRITE(IOIMP,*) 'Destruction du profil du precon'
  475. ENDIF
  476. SEGSUP PMORS
  477. KIDMAT(6)=0
  478. ENDIF
  479. IZA=KIDMAT(7)
  480. IF (IZA.NE.0) THEN
  481. IF (IMPR.GT.2) THEN
  482. WRITE(IOIMP,*) 'Destruction des valeurs du precon'
  483. ENDIF
  484. SEGSUP IZA
  485. KIDMAT(7)=0
  486. ENDIF
  487. SEGDES MATRIK
  488. ELSEIF (IOUBE.NE.0) THEN
  489. WRITE(IOIMP,*) 'IOUBL=',IOUBL, ' non prevu'
  490. GOTO 9999
  491. ENDIF
  492. *
  493. * On décondense si nécessaire
  494. *
  495. * write (6,*) ' resou - mchsol '
  496. * call ecchpo(mchsol,0)
  497. * call mucpri(mchsol,mrigid,iresi)
  498. * write (6,*) ' kres - iresi '
  499. * call ecchpo(iresi,0)
  500. * WRITE(IOIMP,*) 'Avant KRES7'
  501. IF (MRIGI0.NE.0.AND.LDEPE) THEN
  502. MSOLC=MCHSOL
  503. * CALL KRES7(MRIGID,IDEPE,KSMBR0,KSMBR1,
  504. * $ MSOLC,
  505. * $ MCHSOL)
  506. CALL KRES7(MSOLC,MRIGI0,KSMBR0,KSMBR1,1,
  507. $ MCHSOL)
  508. IF (IERR.NE.0) RETURN
  509. ENDIF
  510. CALL ECROBJ('CHPOINT ',MCHSOL)
  511. *STAT CALL SUMSTA(MSTAT,IMPR)
  512. *
  513. * Normal termination
  514. *
  515. RETURN
  516. *
  517. * Format handling
  518. *
  519. *
  520. * Error handling
  521. *
  522. 9999 CONTINUE
  523. WRITE(IOIMP,*) 'An error was detected in kres2.eso'
  524. * 153 2
  525. * Opération illicite dans ce contexte
  526. CALL ERREUR(153)
  527. RETURN
  528. *
  529. * End of KRES2
  530. *
  531. END
  532.  
  533.  

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