mpimod.f90 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729
  1. module mpimod
  2. use pumamod
  3. ! include 'mpif.h'
  4. use mpi
  5. integer :: mpi_itype = MPI_INTEGER4
  6. integer :: mpi_rtype = MPI_REAL4
  7. integer :: mpi_ltype = MPI_LOGICAL
  8. character(len=80) ynode(NPRO) ! node names
  9. end module mpimod
  10. !
  11. ! interface routines to MPI:
  12. !
  13. ! ================
  14. ! SUBROUTINE MPBCI
  15. ! ================
  16. subroutine mpbci(k) ! broadcast 1 integer
  17. use mpimod
  18. integer :: k(*)
  19. call mpi_bcast(k,1,mpi_itype,NROOT,myworld,mpinfo)
  20. return
  21. end subroutine mpbci
  22. ! =================
  23. ! SUBROUTINE MPBCIN
  24. ! =================
  25. subroutine mpbcin(k,n) ! broadcast n integer
  26. use mpimod
  27. integer :: k(n)
  28. call mpi_bcast(k,n,mpi_itype,NROOT,myworld,mpinfo)
  29. return
  30. end subroutine mpbcin
  31. ! ================
  32. ! SUBROUTINE MPBCR
  33. ! ================
  34. subroutine mpbcr(p) ! broadcast 1 real
  35. use mpimod
  36. REAL :: p(*)
  37. call mpi_bcast(p,1,mpi_rtype,NROOT,myworld,mpinfo)
  38. return
  39. end subroutine mpbcr
  40. ! =================
  41. ! SUBROUTINE MPBCRN
  42. ! =================
  43. subroutine mpbcrn(p,n) ! broadcast n real
  44. use mpimod
  45. real :: p(n)
  46. call mpi_bcast(p,n,mpi_rtype,NROOT,myworld,mpinfo)
  47. return
  48. end subroutine mpbcrn
  49. ! ================
  50. ! SUBROUTINE MPBCL
  51. ! ================
  52. subroutine mpbcl(l) ! broadcast 1 logical
  53. use mpimod
  54. logical :: l(*)
  55. call mpi_bcast(l,1,mpi_ltype,NROOT,myworld,mpinfo)
  56. return
  57. end subroutine mpbcl
  58. ! =================
  59. ! SUBROUTINE MPSCIN
  60. ! =================
  61. subroutine mpscin(k,n) ! scatter n integer
  62. use mpimod
  63. integer :: k(n),l(n)
  64. call mpi_scatter(k,n,mpi_itype,l,n,mpi_itype,NROOT,myworld,mpinfo)
  65. k(:)=l(:)
  66. return
  67. end subroutine mpscin
  68. ! =================
  69. ! SUBROUTINE MPSCRN
  70. ! =================
  71. subroutine mpscrn(p,n) ! scatter n real
  72. use mpimod
  73. real :: p(*),l(n)
  74. call mpi_scatter(p,n,mpi_rtype,l,n,mpi_rtype,NROOT,myworld,mpinfo)
  75. p(1:n)=l(:)
  76. return
  77. end subroutine mpscrn
  78. ! =================
  79. ! SUBROUTINE MPSCDN
  80. ! =================
  81. subroutine mpscdn(p,n) ! scatter n double precision
  82. use mpimod
  83. real (kind=8) :: p(*),l(n)
  84. call mpi_scatter(p,n,MPI_REAL8,l,n,MPI_REAL8,NROOT,myworld,mpinfo)
  85. p(1:n)=l(:)
  86. return
  87. end subroutine mpscdn
  88. ! =================
  89. ! SUBROUTINE MPSCGP
  90. ! =================
  91. subroutine mpscgp(pf,pp,klev) ! scatter gridpoint fields
  92. use mpimod
  93. real :: pf(NUGP,klev)
  94. real :: pp(NHOR,klev)
  95. do jlev = 1 , klev
  96. call mpi_scatter(pf(:,jlev),NHOR,mpi_rtype, &
  97. & pp(:,jlev),NHOR,mpi_rtype, &
  98. & NROOT,myworld,mpinfo)
  99. enddo
  100. return
  101. end subroutine mpscgp
  102. ! =================
  103. ! SUBROUTINE MPGAGP
  104. ! =================
  105. subroutine mpgagp(pf,pp,klev) ! gather gridpoint fields
  106. use mpimod
  107. real :: pf(NLON*NLAT,klev)
  108. real :: pp(NHOR,klev)
  109. do jlev = 1 , klev
  110. call mpi_gather(pp(:,jlev),NHOR,mpi_rtype, &
  111. & pf(:,jlev),NHOR,mpi_rtype, &
  112. & NROOT,myworld,mpinfo)
  113. enddo
  114. return
  115. end subroutine mpgagp
  116. ! ===================
  117. ! SUBROUTINE MPGALLGP
  118. ! ===================
  119. subroutine mpgallgp(pf,pp,klev) ! gather gritpoint to all
  120. use mpimod
  121. real :: pf(NLON*NLAT,klev)
  122. real :: pp(NHOR,klev)
  123. do jlev = 1 , klev
  124. call mpi_allgather(pp(:,jlev),NHOR,mpi_rtype, &
  125. & pf(:,jlev),NHOR,mpi_rtype, &
  126. & myworld,mpinfo)
  127. enddo
  128. return
  129. end subroutine mpgallgp
  130. ! =================
  131. ! SUBROUTINE MPSCSP
  132. ! =================
  133. subroutine mpscsp(pf,pp,klev) ! scatter spectral fields
  134. use mpimod
  135. real :: pf(NESP,klev)
  136. real :: pp(NSPP,klev)
  137. do jlev = 1 , klev
  138. call mpi_scatter(pf(:,jlev),NSPP,mpi_rtype &
  139. & ,pp(:,jlev),NSPP,mpi_rtype &
  140. & ,NROOT,myworld,mpinfo)
  141. enddo
  142. return
  143. end subroutine mpscsp
  144. ! =================
  145. ! SUBROUTINE MPGASP
  146. ! =================
  147. subroutine mpgasp(pf,pp,klev) ! gather spectral fields
  148. use mpimod
  149. real :: pf(NESP,klev)
  150. real :: pp(NSPP,klev)
  151. do jlev = 1 , klev
  152. call mpi_gather(pp(:,jlev),NSPP,mpi_rtype &
  153. & ,pf(:,jlev),NSPP,mpi_rtype &
  154. & ,NROOT,myworld,mpinfo)
  155. enddo
  156. return
  157. end subroutine mpgasp
  158. ! =================
  159. ! SUBROUTINE MPGACS
  160. ! =================
  161. subroutine mpgacs(pcs) ! gather cross sections
  162. use mpimod
  163. real :: pcs(NLAT,NLEV),l(NLPP)
  164. do jlev = 1 , NLEV
  165. l(:)=pcs(1:NLPP,jlev)
  166. call mpi_gather(l,NLPP,mpi_rtype &
  167. & ,pcs(:,jlev),NLPP,mpi_rtype &
  168. & ,NROOT,myworld,mpinfo)
  169. enddo
  170. return
  171. end subroutine mpgacs
  172. ! ===================
  173. ! SUBROUTINE MPGALLSP
  174. ! ===================
  175. subroutine mpgallsp(pf,pp,klev) ! gather spectral to all
  176. use mpimod
  177. real :: pf(NESP,klev)
  178. real :: pp(NSPP,klev)
  179. do jlev = 1 , klev
  180. call mpi_allgather(pp(:,jlev),NSPP,mpi_rtype &
  181. & ,pf(:,jlev),NSPP,mpi_rtype &
  182. & ,myworld,mpinfo)
  183. enddo
  184. return
  185. end subroutine mpgallsp
  186. ! ================
  187. ! SUBROUTINE MPSUM
  188. ! ================
  189. subroutine mpsum(psp,klev) ! sum spectral fields
  190. use mpimod
  191. real :: psp(NESP*klev)
  192. real :: tmp(NESP*klev)
  193. call mpi_reduce(psp,tmp,NESP*klev,mpi_rtype,MPI_SUM &
  194. & ,NROOT,myworld,mpinfo)
  195. if (mypid == NROOT) psp = tmp
  196. return
  197. end subroutine mpsum
  198. ! ==================
  199. ! SUBROUTINE MPSUMSC
  200. ! ==================
  201. subroutine mpsumsc(psf,psp,klev) ! sum & scatter spectral
  202. use mpimod
  203. real :: psf(NESP,klev)
  204. real :: psp(NSPP,klev)
  205. do jlev = 1 , klev
  206. call mpi_reduce_scatter(psf(:,jlev),psp(:,jlev),nscatsp &
  207. & ,mpi_rtype,MPI_SUM,myworld,mpinfo)
  208. enddo
  209. return
  210. end subroutine mpsumsc
  211. ! ====================
  212. ! SUBROUTINE MPSUMR
  213. ! ====================
  214. subroutine mpsumr(pr,kdim) ! sum kdim reals
  215. use mpimod
  216. real pr(kdim)
  217. real tmp(kdim)
  218. call mpi_reduce(pr,tmp,kdim,mpi_rtype,MPI_SUM,NROOT,myworld,mpinfo)
  219. if (mypid == NROOT) pr = tmp
  220. return
  221. end subroutine mpsumr
  222. ! ====================
  223. ! SUBROUTINE MPSUMBCR
  224. ! ====================
  225. subroutine mpsumbcr(pr,kdim) ! sum & broadcast kdim reals
  226. use mpimod
  227. real pr(kdim)
  228. real tmp(kdim)
  229. call mpi_allreduce(pr,tmp,kdim,mpi_rtype,MPI_SUM,myworld,mpinfo)
  230. pr = tmp
  231. return
  232. end subroutine mpsumbcr
  233. ! ==================
  234. ! SUBROUTINE MPABORT
  235. ! ==================
  236. subroutine mpabort(ym)
  237. use mpimod
  238. character (len=* ) :: ym
  239. character (len=64) :: ystar = ' '
  240. character (len=64) :: yline = ' '
  241. character (len=64) :: yabor = 'Program aborted'
  242. character (len=64) :: ymess = ' '
  243. character (len=64) :: yhead = ' '
  244. ilmess = len_trim(ym)
  245. ilabor = len_trim(yabor)
  246. ilen = 60
  247. do j = 1 , ilen+4
  248. ystar(j:j) = '*'
  249. yline(j:j) = '-'
  250. enddo
  251. ioff = 2
  252. if (ilmess < ilen-1) ioff = ioff + (ilen - ilmess) / 2
  253. ymess(1+ioff:ilmess+ioff) = trim(ym)
  254. ioff = 2
  255. if (ilabor < ilen-1) ioff = ioff + (ilen - ilabor) / 2
  256. yhead(1+ioff:ilabor+ioff) = trim(yabor)
  257. yline(1:1) = '*'
  258. ymess(1:1) = '*'
  259. yhead(1:1) = '*'
  260. yline(2:2) = ' '
  261. ymess(2:2) = ' '
  262. yhead(2:2) = ' '
  263. j = ilen + 4
  264. yline(j:j) = '*'
  265. ymess(j:j) = '*'
  266. yhead(j:j) = '*'
  267. j = ilen + 3
  268. yline(j:j) = ' '
  269. ymess(j:j) = ' '
  270. yhead(j:j) = ' '
  271. if (mypid == NROOT) then
  272. open (44,file='Abort_Message')
  273. write(44,'(A)') trim(ystar)
  274. write(44,'(A)') trim(yhead)
  275. write(44,'(A)') trim(yline)
  276. write(44,'(A)') trim(ymess)
  277. write(44,'(A)') trim(ystar)
  278. close(44)
  279. write(nud,'(/,A)') trim(ystar)
  280. write(nud,'(A)') trim(yhead)
  281. write(nud,'(A)') trim(yline)
  282. write(nud,'(A)') trim(ymess)
  283. write(nud,'(A,/)') trim(ystar)
  284. call mpi_abort(myworld,mpinfo,mpinfo)
  285. endif
  286. stop
  287. end
  288. ! ==================
  289. ! SUBROUTINE MPSTART
  290. ! ==================
  291. subroutine mpstart(kworld) ! initialization
  292. use mpimod
  293. integer :: itest = 0
  294. real :: rtest = 0.0
  295. logical :: ltest = .true.
  296. character (80) :: myympname
  297. if (kind(itest) == 8) mpi_itype = MPI_INTEGER8
  298. if (kind(rtest) == 8) mpi_rtype = MPI_REAL8
  299. if (kworld < 0) then
  300. call mpi_init(mpinfo)
  301. myworld=MPI_COMM_WORLD
  302. else
  303. myworld = kworld
  304. endif
  305. call mpi_comm_size(myworld,nproc,mpinfo)
  306. call mpi_comm_rank(myworld,mypid,mpinfo)
  307. if (nproc .ne. NPRO .and. mypid == NROOT) then
  308. write(nud,*)'Compiled for ',NPRO,' nodes'
  309. write(nud,*)'Running on ',nproc,' nodes'
  310. call mpi_abort(myworld,mpinfo,mpinfo)
  311. endif
  312. allocate(ympname(npro)) ; ympname(:) = ' '
  313. call mpi_get_processor_name(myympname,ilen,mpinfo)
  314. call mpi_gather(myympname,80,MPI_CHARACTER, &
  315. ympname,80,MPI_CHARACTER, &
  316. NROOT,myworld,mpinfo)
  317. return
  318. end subroutine mpstart
  319. ! =================
  320. ! SUBROUTINE MPSTOP
  321. ! =================
  322. subroutine mpstop
  323. use mpimod
  324. call mpi_barrier(myworld,mpinfo)
  325. call mpi_finalize(mpinfo)
  326. return
  327. end subroutine mpstop
  328. ! ===================
  329. ! SUBROUTINE MPREADGP
  330. ! ===================
  331. subroutine mpreadgp(ktape,p,kdim,klev)
  332. use mpimod
  333. real p(kdim,klev)
  334. real z(NLON*NLAT,klev)
  335. z = 0.0
  336. if (mypid == NROOT) read (ktape) z(:,:)
  337. if (kdim == NHOR) then
  338. call mpscgp(z,p,klev)
  339. else
  340. if (mypid == NROOT) p = z
  341. endif
  342. return
  343. end subroutine mpreadgp
  344. ! ====================
  345. ! SUBROUTINE MPWRITEGP
  346. ! ====================
  347. subroutine mpwritegp(ktape,p,kdim,klev)
  348. use mpimod
  349. real p(kdim,klev)
  350. real z(NLON*NLAT,klev)
  351. if (kdim == NHOR) then
  352. call mpgagp(z,p,klev)
  353. if (mypid == NROOT) write(ktape) z(1:NLON*NLAT,:)
  354. else
  355. if (mypid == NROOT) write(ktape) p(1:NLON*NLAT,:)
  356. endif
  357. return
  358. end subroutine mpwritegp
  359. ! =====================
  360. ! SUBROUTINE MPWRITEGPH
  361. ! =====================
  362. subroutine mpwritegph(ktape,p,kdim,klev,ihead)
  363. use mpimod
  364. real p(kdim,klev)
  365. real z(NLON*NLAT,klev)
  366. !
  367. real(kind=4) :: zp(kdim,klev)
  368. real(kind=4) :: zz(NLON*NLAT,klev)
  369. !
  370. integer ihead(8)
  371. if (kdim == NHOR) then
  372. call mpgagp(z,p,klev)
  373. if (mypid == NROOT) then
  374. write(ktape) ihead
  375. zz(:,:)=z(:,:)
  376. write(ktape) zz(1:NLON*NLAT,:)
  377. endif
  378. else
  379. if (mypid == NROOT) then
  380. write(ktape) ihead
  381. zp(:,:)=p(:,:)
  382. write(ktape) zp(1:NLON*NLAT,:)
  383. endif
  384. endif
  385. return
  386. end subroutine mpwritegph
  387. ! ===================
  388. ! SUBROUTINE MPREADSP
  389. ! ===================
  390. subroutine mpreadsp(ktape,p,kdim,klev)
  391. use mpimod
  392. real p(kdim,klev)
  393. real z(NESP,klev)
  394. z = 0.0
  395. if (mypid == NROOT) read(ktape) ((z(i,j),i=1,NRSP),j=1,klev)
  396. if (kdim == NSPP) then
  397. call mpscsp(z,p,klev)
  398. else
  399. if (mypid == NROOT) p = z
  400. endif
  401. return
  402. end subroutine mpreadsp
  403. ! ====================
  404. ! SUBROUTINE MPWRITESP
  405. ! ====================
  406. subroutine mpwritesp(ktape,p,kdim,klev)
  407. use mpimod
  408. real p(kdim,klev)
  409. real z(NESP,klev)
  410. if (kdim == NSPP) then
  411. call mpgasp(z,p,klev)
  412. if (mypid == NROOT) write(ktape) ((z(i,j),i=1,NRSP),j=1,klev)
  413. else
  414. if (mypid == NROOT) write(ktape) ((z(i,j),i=1,NRSP),j=1,klev)
  415. endif
  416. return
  417. end subroutine mpwritesp
  418. ! ===================
  419. ! SUBROUTINE GET_MPI_INFO
  420. ! ===================
  421. subroutine get_mpi_info(nprocess,npid) ! get nproc and pid
  422. use mpimod
  423. myworld=MPI_COMM_WORLD
  424. call mpi_comm_size(myworld,nprocess,mpinfo)
  425. call mpi_comm_rank(myworld,npid,mpinfo)
  426. return
  427. end subroutine get_mpi_info
  428. ! ==================
  429. ! SUBROUTINE MPGETSP
  430. ! ==================
  431. subroutine mpgetsp(yn,p,kdim,klev)
  432. use mpimod
  433. character (len=*) :: yn
  434. real :: p(kdim,klev)
  435. real :: z(NESP,klev)
  436. z(:,:) = 0.0
  437. if (mypid == NROOT) call get_restart_array(yn,z,NRSP,NESP,klev)
  438. call mpscsp(z,p,klev)
  439. return
  440. end subroutine mpgetsp
  441. ! ==================
  442. ! SUBROUTINE MPGETGP
  443. ! ==================
  444. subroutine mpgetgp(yn,p,kdim,klev)
  445. use mpimod
  446. character (len=*) :: yn
  447. real :: p(kdim,klev)
  448. real :: z(NUGP,klev)
  449. if (mypid == NROOT) call get_restart_array(yn,z,NUGP,NUGP,klev)
  450. call mpscgp(z,p,klev)
  451. return
  452. end subroutine mpgetgp
  453. ! ==================
  454. ! SUBROUTINE MPPUTSP
  455. ! ==================
  456. subroutine mpputsp(yn,p,kdim,klev)
  457. use mpimod
  458. character (len=*) :: yn
  459. real :: p(kdim,klev)
  460. real :: z(NESP,klev)
  461. call mpgasp(z,p,klev)
  462. if (mypid == NROOT) call put_restart_array(yn,z,NRSP,NESP,klev)
  463. return
  464. end subroutine mpputsp
  465. ! ==================
  466. ! SUBROUTINE MPPUTGP
  467. ! ==================
  468. subroutine mpputgp(yn,p,kdim,klev)
  469. use mpimod
  470. character (len=*) :: yn
  471. real :: p(kdim,klev)
  472. real :: z(NUGP,klev)
  473. call mpgagp(z,p,klev)
  474. if (mypid == NROOT) call put_restart_array(yn,z,NUGP,NUGP,klev)
  475. return
  476. end subroutine mpputgp
  477. ! ===================
  478. ! SUBROUTINE MPSURFGP
  479. ! ===================
  480. subroutine mpsurfgp(yn,p,kdim,klev)
  481. use mpimod
  482. character (len=*) :: yn
  483. real :: p(kdim,klev)
  484. real :: z(NUGP,klev)
  485. integer :: iread(1)
  486. if (mypid == NROOT) call get_surf_array(yn,z,NUGP,klev,iread)
  487. call mpbci(iread)
  488. if (iread(1) == 1) call mpscgp(z,p,klev)
  489. return
  490. end subroutine mpsurfgp
  491. ! ===================
  492. ! SUBROUTINE MPMAXVAL
  493. ! ===================
  494. subroutine mpmaxval(p,kdim,klev,pmax)
  495. use mpimod
  496. real :: p(kdim,klev)
  497. real :: pmax(1)
  498. real zmax(1)
  499. zmax(1) = maxval(p(:,:))
  500. call mpi_allreduce(zmax,pmax,1,mpi_rtype,MPI_MAX,myworld,mpinfo)
  501. return
  502. end subroutine mpmaxval
  503. ! ===================
  504. ! SUBROUTINE MPSUMVAL
  505. ! ===================
  506. subroutine mpsumval(p,kdim,klev,psum)
  507. use mpimod
  508. real :: p(kdim,klev)
  509. real :: psum(1)
  510. real :: zsum(1)
  511. zsum = sum(p(:,:))
  512. call mpi_allreduce(zsum,psum,1,mpi_rtype,MPI_SUM,myworld,mpinfo)
  513. return
  514. end subroutine mpsumval
  515. subroutine mrsum(k) ! sum up 1 integer
  516. return
  517. end
  518. subroutine mrbci(k) ! broadcast 1 integer
  519. return
  520. end
  521. subroutine mrdiff(p,d,n)
  522. real :: p(n)
  523. real :: d(n)
  524. return
  525. end
  526. subroutine mrdimensions ! used in mpimod_multi.f90
  527. return
  528. end