PCMRH.f90 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563
  1. program PCMRH
  2. !***************************************************************
  3. ! Parallel CMRH program with MPI
  4. !
  5. ! Created by S. Duminil sep. 2012
  6. !***************************************************************
  7. use mpi
  8. external SSWAP,DCOPY,DSCAL,DGEMV,DSWAP
  9. intrinsic dabs,dsqrt
  10. !********************************************
  11. ! VARIABLES
  12. !********************************************
  13. integer,parameter ::Nprocs_max=16
  14. integer ::rang,Nprocs,n,nl,code,i,j,k,choix,solchoice,iprint
  15. real(kind=8) ::init(5)
  16. integer,dimension(MPI_STATUS_SIZE) ::statut
  17. real(kind=8),allocatable,dimension(:,:) ::AL
  18. real(kind=8) ::btemp,bmax,btemp2
  19. real(kind=8),allocatable,dimension(:) ::b,sol,rs,temp,x,workl,work,tempn,c,s,xl,cl,temp2,xx,yy,valeurs
  20. integer,allocatable,dimension(:) ::p
  21. real*8 ::dnorm2,D_inv,beta,dnorm,dnrm2,tol,alpha,H_k1,dgam,r
  22. integer ::indice,idamax,k1,n_k,knl,knlp,inl,inlp
  23. real*8 dtime, delta, delta0, dt1, dt2, dt3
  24. real*8 etime, eelta, eelta0, et1, et2, et3
  25. real*8 t1, t2, t3, tm1, tm2, tm3
  26. real*8 dtime1, delta1, delta01
  27. real*8 etime1, eelta1, eelta01,l,p1,p2
  28. integer nb_ligne
  29. real(kind=8) ::temps,tempsd,tempsf,temps2
  30. integer*4,dimension(3) ::timearray
  31. real ::rand
  32. real*4 da(2),ea(2),da1(2),ea1(2)
  33. !********************************************
  34. ! MPI CONFIG
  35. !********************************************
  36. call MPI_INIT(code)
  37. call MPI_COMM_RANK(MPI_COMM_WORLD,rang,code)
  38. call MPI_COMM_SIZE(MPI_COMM_WORLD,Nprocs,code)
  39. call MPI_BARRIER(MPI_COMM_WORLD,code)
  40. call itime(timearray) ! Get the current time
  41. i = rand ( timearray(1)+timearray(2)+timearray(3) )
  42. !************************************
  43. ! Matrix Initialisation
  44. !************************************
  45. if (rang==0) then
  46. open(unit=10,file='inputfile.dat',status='old')
  47. read(10,*)init
  48. close(10)
  49. iprint=int(init(1))
  50. n=int(init(2))
  51. choix=int(init(3))
  52. solchoice=int(init(4))
  53. tol=init(5)
  54. end if
  55. tempsd=MPI_WTIME()
  56. call MPI_BCAST(n,1,MPI_INTEGER,0,MPI_COMM_WORLD,code)
  57. call MPI_BCAST(choix,1,MPI_INTEGER,0,MPI_COMM_WORLD,code)
  58. call MPI_BCAST(tol,1,MPI_DOUBLE_PRECISION,0,MPI_COMM_WORLD,code)
  59. if (mod(n,Nprocs)==0)then
  60. nl=n/Nprocs
  61. else
  62. nl=(n/Nprocs)+1
  63. end if
  64. allocate(b(n),rs(n),x(n),AL(nl,n),workl(nl),sol(n),temp(n),temp2(n))
  65. if (choix<3)then
  66. call matrice(AL,nl,n,2.d0,1.d0,choix,rang)
  67. else
  68. allocate(valeurs(n))
  69. open(unit=30,file='matrixfile.dat',status='old')
  70. do i=1,n
  71. call position(i,nl,inl,inlp)
  72. read(30,*)valeurs
  73. if (rang.eq.inl)then
  74. AL(inlp,:)=valeurs(:)
  75. endif
  76. enddo
  77. close(30)
  78. endif
  79. !********************************************
  80. ! End Matrix initialisation
  81. !********************************************
  82. !********************************************
  83. ! Solution initialisation
  84. !********************************************
  85. call MPI_BCAST(solchoice,1,MPI_INTEGER,0,MPI_COMM_WORLD,code)
  86. call solution(AL,workl,sol,nl,n,solchoice)
  87. call MPI_ALLGATHER(workl,nl,MPI_DOUBLE_PRECISION,b,nl,MPI_DOUBLE_PRECISION,MPI_COMM_WORLD,code)
  88. !********************************************
  89. ! End solution initialisation
  90. !********************************************
  91. !********************************************
  92. ! Permutation Vector initialisation
  93. ! and computation of initial residual
  94. !********************************************
  95. if (rang==0) then
  96. allocate(p(n))
  97. do i=1,n
  98. p(i)=i
  99. end do
  100. end if
  101. dnorm2=dnrm2(n,b,1)
  102. if (rang==0)then
  103. if (iprint==2)then
  104. print*,''
  105. print*,'la norme 2 de r0=',dnorm2
  106. print*,''
  107. end if
  108. end if
  109. !********************************************
  110. ! Computation of max b(i)
  111. !********************************************
  112. indice = idamax(n,b,1)
  113. beta = b(indice)
  114. D_inv = 1.d0/beta
  115. !********************************************
  116. ! First permutation
  117. !********************************************
  118. if (indice .ne.1) then
  119. call dswap(1,b(indice),1,b(1),1)
  120. call dswap(nl,AL(1,indice),1,AL(1,1),1)
  121. call position(indice,nl,inl,inlp)
  122. if (inl==0) then
  123. if (rang==inl)then
  124. call dswap(n,AL(inlp,1),nl,AL(1,1),nl)
  125. end if
  126. else
  127. if (rang==inl) then
  128. call dcopy(n,AL(inlp,1),nl,temp,1)
  129. end if
  130. call MPI_BCAST(temp,n,MPI_DOUBLE_PRECISION,inl,MPI_COMM_WORLD,code)
  131. if (rang==0)then
  132. call dcopy(n,AL(1,1),nl,temp2,1)
  133. end if
  134. call MPI_BCAST(temp2,n,MPI_DOUBLE_PRECISION,0,MPI_COMM_WORLD,code)
  135. if (rang==inl)then
  136. call dcopy(n,temp2,1,AL(inlp,1),nl)
  137. end if
  138. if (rang==0)then
  139. call dcopy(n,temp,1,AL(1,1),nl)
  140. end if
  141. end if
  142. if (rang==0)then
  143. call sswap(1,p(indice),1,p(1),1)
  144. end if
  145. end if
  146. !********************************************
  147. ! End permutation
  148. !********************************************
  149. ! call dscal(n,0.d0,x,1)
  150. call dscal(n,D_inv,b,1)
  151. rs(1)=beta
  152. dnorm=dabs(beta)
  153. allocate(tempn(nl),work(n),c(n),s(n),xl(nl),cl(nl))
  154. call dscal(nl,0.d0,xl,1)
  155. call MPI_BARRIER(MPI_COMM_WORLD,code)
  156. k=0
  157. if (rang==0) then
  158. if (iprint==2)then
  159. write(*,*)'iter., pseudo residu '
  160. end if
  161. open(unit=30,file='presidual.dat',status='unknown')
  162. write(30,*)'Iter., Pseudo Residual Norm'
  163. end if
  164. !tol=1.d-7
  165. delta0 = dtime(da)
  166. eelta0 = etime(ea)
  167. tempsf=MPI_WTIME()
  168. !********************************************
  169. ! Loop
  170. !********************************************
  171. do while((dnorm.gt.tol).and.(k.lt.n))
  172. k=k+1;
  173. n_k=n-k;
  174. k1=k+1;
  175. if (rang==0) then
  176. if (iprint==2)then
  177. write(*,10)k-1,dnorm
  178. end if
  179. write(30,10)k-1,dnorm
  180. end if
  181. 10 format(i5,2x,e12.4)
  182. !********************************************
  183. ! workl = A*b
  184. !********************************************
  185. call dgemv('N',nl,n_k+1,1.d0,AL(1,k),nl,b(k),1,0.d0,workl,1)
  186. ! ! call dscal(nl,0.d0,workl,1)
  187. ! ! do j=k1,n
  188. ! ! do i=1,nl
  189. ! ! workl(i)=workl(i)+AL(i,j)*b(j)
  190. ! ! end do
  191. ! ! end do
  192. !********************************************
  193. ! A(:,k)=b
  194. !********************************************
  195. call position(k1,nl,knl,knlp)
  196. if (rang>knl)then
  197. call dcopy(nl,b(rang*nl+1),1,AL(1,k),1)
  198. else
  199. if (rang==knl) then
  200. call dcopy(nl-knlp+1,b(k1),1,AL(knlp,k),1)
  201. end if
  202. end if
  203. !********************************************
  204. ! j Loop
  205. !********************************************
  206. do j=1,k
  207. !********************************************
  208. ! A(j,k)=workl(j)
  209. ! workl(j)=0
  210. !********************************************
  211. call position(j,nl,knl,knlp)
  212. if (rang==knl) then
  213. AL(knlp,k)=workl(knlp)
  214. workl(knlp)=0.d0
  215. alpha=-AL(knlp,k)
  216. end if
  217. call MPI_BCAST(alpha,1,MPI_DOUBLE_PRECISION,knl,MPI_COMM_WORLD,code)
  218. !********************************************
  219. ! workl=workl - alpha A(:,j)
  220. !********************************************
  221. call position(j+1,nl,knl,knlp)
  222. if (rang >knl)then
  223. call daxpy(nl,alpha,AL(1,j),1,workl,1)
  224. end if
  225. if (rang==knl) then
  226. call daxpy(nl-knlp+1,alpha,AL(knlp,j),1,workl(knlp),1)
  227. end if
  228. end do
  229. !********************************************
  230. ! End j Loop
  231. !********************************************
  232. if (k.lt.n)then
  233. call MPI_ALLGATHER(workl,nl,MPI_DOUBLE_PRECISION,b,nl,MPI_DOUBLE_PRECISION,MPI_COMM_WORLD,code)
  234. !********************************************
  235. ! Search indice i_0
  236. !********************************************
  237. indice=k+idamax(n_k,b(k1),1)
  238. H_k1=b(indice)
  239. !********************************************
  240. ! Permutation
  241. !********************************************
  242. if (indice.ne.k1) then
  243. if (rang==0) then
  244. call sswap(1,p(indice),1,p(k1),1)
  245. end if
  246. call dswap(1,b(indice),1,b(k1),1)
  247. call dswap(nl,AL(1,indice),1,AL(1,k1),1)
  248. call MPI_BARRIER(MPI_COMM_WORLD,code)
  249. call position(k1,nl,knl,knlp)
  250. call position(indice,nl,inl,inlp)
  251. if (inl==knl) then
  252. if (rang==inl)then
  253. call dswap(n,AL(inlp,1),nl,AL(knlp,1),nl)
  254. end if
  255. else
  256. if (rang==inl) then
  257. call dcopy(n,AL(inlp,1),nl,temp,1)
  258. end if
  259. call MPI_BCAST(temp,n,MPI_DOUBLE_PRECISION,inl,MPI_COMM_WORLD,code)
  260. if (rang==knl)then
  261. call dcopy(n,AL(knlp,1),nl,temp2,1)
  262. end if
  263. call MPI_BCAST(temp2,n,MPI_DOUBLE_PRECISION,knl,MPI_COMM_WORLD,code)
  264. if (rang==inl)then
  265. call dcopy(n,temp2,1,AL(inlp,1),nl)
  266. end if
  267. if (rang==knl)then
  268. call dcopy(n,temp,1,AL(knlp,1),nl)
  269. end if
  270. end if
  271. end if
  272. !********************************************
  273. ! End Permutation
  274. !********************************************
  275. if (H_k1.ne.0.d0)then
  276. !********************************************
  277. ! b=b/(b)_(i_0)
  278. !********************************************
  279. D_inv=1.d0/H_k1
  280. ! call dcopy(n,work,1,b,1)
  281. call dscal(n,D_inv,b,1)
  282. end if
  283. else
  284. H_k1=0.d0
  285. call dscal(n,0.d0,b,1)
  286. end if
  287. !********************************************
  288. ! Application of Givens rotations
  289. !********************************************
  290. if (k.gt.1)then
  291. do i=2,k
  292. call position(i-1,nl,knl,knlp)
  293. call position(i,nl,inl,inlp)
  294. if (inl==knl)then
  295. if (rang==inl)then
  296. btemp=AL(knlp,k)
  297. AL(knlp,k)=c(i-1)*btemp+s(i-1)*AL(inlp,k)
  298. AL(inlp,k)=-s(i-1)*btemp+c(i-1)*AL(inlp,k)
  299. end if
  300. else
  301. if (rang==inl)then
  302. btemp=AL(inlp,k)
  303. end if
  304. if (rang==knl)then
  305. btemp2=AL(knlp,k)
  306. end if
  307. call MPI_BCAST(btemp,1,MPI_DOUBLE_PRECISION,inl,MPI_COMM_WORLD,code)
  308. call MPI_BCAST(btemp2,1,MPI_DOUBLE_PRECISION,knl,MPI_COMM_WORLD,code)
  309. net=net+2
  310. if (rang==inl)then
  311. AL(inlp,k)=-s(i-1)*btemp2+c(i-1)*btemp
  312. end if
  313. if (rang==knl)then
  314. AL(knlp,k)=c(i-1)*btemp2+s(i-1)*btemp
  315. end if
  316. end if
  317. end do
  318. end if
  319. if (k.le.n)then
  320. call position(k,nl,knl,knlp)
  321. if (rang==knl)then
  322. btemp=AL(knlp,k)
  323. end if
  324. call MPI_BCAST(btemp,1,MPI_DOUBLE_PRECISION,knl,MPI_COMM_WORLD,code)
  325. dgam=dsqrt(btemp**2+H_k1**2)
  326. c(k)=dabs(btemp)/dgam
  327. s(k)=(H_k1*dabs(btemp))/(btemp*dgam)
  328. rs(k1)=-s(k)*rs(k)
  329. rs(k)=c(k)*rs(k)
  330. dnorm=dabs(rs(k1))
  331. if (rang==knl)then
  332. AL(knlp,k)=c(k)*btemp+s(k)*H_k1
  333. end if
  334. end if
  335. !********************************************
  336. ! End Givens rotations
  337. !********************************************
  338. end do
  339. !********************************************
  340. ! Solve H_k d_k = alpha e_1
  341. !********************************************
  342. call position(k,nl,knl,knlp)
  343. do j=k,1,-1
  344. call position(j,nl,inl,inlp)
  345. if (rang==inl) then
  346. rs(j)=rs(j)/AL(inlp,j)
  347. btemp=rs(j)
  348. end if
  349. call MPI_BCAST(btemp,1,MPI_DOUBLE_PRECISION,inl,MPI_COMM_WORLD,code)
  350. rs(j)=btemp
  351. work(j)=btemp
  352. call position(j-1,nl,inl,inlp)
  353. if (rang==inl)then
  354. call daxpy(inlp,-btemp,AL(1,j),1,rs(rang*nl+1),1)
  355. end if
  356. if (rang.lt.inl)then
  357. call daxpy(nl,-btemp,AL(1,j),1,rs(rang*nl+1),1)
  358. end if
  359. end do
  360. do j=k,1,-1
  361. btemp=rs(j)
  362. call position(k,nl,knl,knlp)
  363. call position(j+1,nl,inl,inlp)
  364. if (rang==inl)then
  365. if (rang==knl)then
  366. call daxpy((knlp-inlp+1),btemp,AL(inlp,j),1,rs(rang*nl+inlp),1)
  367. end if
  368. if (rang.lt.knl)then
  369. call daxpy(nl-inlp+1,btemp,AL(inlp,j),1,rs(rang*nl+inlp),1)
  370. end if
  371. end if
  372. if (rang.gt.inl)then
  373. if (rang==knl)then
  374. call daxpy(knlp,btemp,AL(1,j),1,rs(rang*nl+1),1)
  375. end if
  376. if (rang.lt.knl)then
  377. call daxpy(nl,btemp,AL(1,j),1,rs(rang*nl+1),1)
  378. end if
  379. end if
  380. end do
  381. !********************************************
  382. ! End Loop
  383. !********************************************
  384. !********************************************
  385. ! Update x_k=x_0 + L_k d_k
  386. !********************************************
  387. call position(k1,nl,knl,knlp)
  388. if (rang.gt.knl)then
  389. call dgemv('N',nl,k,1.d0,AL(1,1),nl,work,1,0.d0,c(rang*nl+1),1)
  390. end if
  391. if (rang==knl)then
  392. call dgemv('N',nl-knlp+1,k,1.d0,AL(knlp,1),nl,work,1,0.d0,c(rang*nl+knlp),1)
  393. end if
  394. call position(k,nl,knl,knlp)
  395. if (rang.lt.knl)then
  396. call daxpy(nl,1.d0,rs(rang*nl+1),1,xl,1)
  397. end if
  398. if (rang==knl)then
  399. call daxpy(knlp,1.d0,rs(rang*nl+1),1,xl,1)
  400. call daxpy(nl-knlp,1.d0,c(rang*nl+knlp+1),1,xl(knlp+1),1)
  401. end if
  402. if (rang.gt.knl)then
  403. call daxpy(nl,1.d0,c(rang*nl+1),1,xl,1)
  404. end if
  405. !********************************************
  406. ! Reorder the components of x_k
  407. !********************************************
  408. call MPI_GATHER(xl,nl,MPI_DOUBLE_PRECISION,x,nl,MPI_DOUBLE_PRECISION,0,MPI_COMM_WORLD,code)
  409. call MPI_BARRIER(MPI_COMM_WORLD,code)
  410. if (rang==0)then
  411. do j=1,n
  412. b(p(j))=x(j)
  413. end do
  414. end if
  415. tempsf=(MPI_WTIME()-tempsf)
  416. tempsd=(MPI_WTIME()-tempsd)
  417. !write(*,*)'rang et temps',rang,tempsf,tempsd
  418. call MPI_REDUCE(tempsf,temps,1,MPI_DOUBLE_PRECISION,MPI_MAX,0,MPI_COMM_WORLD,code)
  419. call MPI_REDUCE(tempsd,temps2,1,MPI_DOUBLE_PRECISION,MPI_MAX,0,MPI_COMM_WORLD,code)
  420. if (rang==0)then
  421. if (iprint >= 2) then
  422. print*,'Results : iter, norm residual, time'
  423. write(*,40)k,dnorm,temps,temps2
  424. end if
  425. write(30,10),k,dnorm
  426. close(30)
  427. open(unit=10,file='temps.txt',status='unknown',access='append')
  428. write(unit=10,FMT=*)'Proc. Number :',Nprocs,'Iter. Number:',k,'Time:',temps,temps2
  429. close(10)
  430. end if
  431. if (rang==0)then
  432. if (iprint==2) then
  433. write(*,*)
  434. write(*,*)' la perm. et les solutions CMRH et EXACTE sont : '
  435. write(*,*)' ______________________________________________ '
  436. if (n.lt.10)then
  437. do i=1,n
  438. write(*,20)p(i), b(i), sol(i)
  439. end do
  440. else
  441. do i=1,10
  442. write(*,20)p(i), b(i), sol(i)
  443. end do
  444. end if
  445. end if
  446. end if
  447. 30 format(A37,4(f14.7,1x))
  448. 20 format(1x,i5,2(1x,f14.8))
  449. 40 format(1x,i5,4(1x,e12.4))
  450. ! deallocate(AL,b,sol,workl,temp,tempn,work,c,s,rs,xl,x,cl)
  451. if (rang==0) then
  452. deallocate(p)
  453. end if
  454. call MPI_FINALIZE(code)
  455. end program PCMRH
  456. !********************************************
  457. ! End Parallel CMRH
  458. !********************************************
  459. !********************************************
  460. ! Matrix Choice
  461. !********************************************
  462. subroutine matrice(A,nl,n,alpha,beta,choix,rang)
  463. integer ::n,i,j,choix,rang,nl,k
  464. real(kind=8),dimension(nl,n)::A
  465. real*8 ::alpha,beta
  466. if (choix==1) then
  467. do i = 1, nl
  468. k=rang*nl+i
  469. do j = 1, k
  470. A(i,j) = (2.d0*dble(j)-1.d0)/dble(n-k+j)
  471. end do
  472. do j = k+1, n
  473. A(i,j) = (2.d0*dble(k)-1.d0)/dble(n-k+j)
  474. end do
  475. end do
  476. end if
  477. if (choix==2)then
  478. do i=1,nl
  479. k=rang*nl+i
  480. do j=1,k
  481. A(i,j)=dble(k)
  482. end do
  483. do j=k+1,n
  484. A(i,j)=dble(j)
  485. end do
  486. end do
  487. end if
  488. end subroutine matrice
  489. !********************************************
  490. ! Solution Choice
  491. !********************************************
  492. subroutine solution(A,b,sol,nl,n,solchoice)
  493. integer ::n,i,nl,solchoice
  494. real(kind=8),dimension(n,n) ::A
  495. real(kind=8),dimension(n) ::sol
  496. real(kind=8),dimension(nl) ::b
  497. real ::rand
  498. do i=1,n
  499. if (solchoice==1) then
  500. sol(i)=1.d0
  501. end if
  502. if (solchoice==2) then
  503. sol(i)=dble(n-i+1)
  504. end if
  505. if (solchoice==3)then
  506. sol(i)=rand(0)
  507. end if
  508. end do
  509. call dgemv('N',nl,n,1.d0,A,nl,sol,1,0.d0,b,1)
  510. end subroutine solution
  511. !********************************************
  512. ! Where are we ?
  513. !********************************************
  514. subroutine position(k,nl,knl,knlp)
  515. integer :: k,nl,knl,knlp
  516. knl=(k-1)/nl
  517. knlp=k-knl*nl
  518. end subroutine position