PUMA
219
Portable University Model of the Atmosphere
|
00001 module mpimod 00002 use pumamod 00003 use mpi 00004 00005 integer :: mpi_itype = MPI_INTEGER4 00006 integer :: mpi_rtype = MPI_REAL4 00007 integer :: mpi_ltype = MPI_LOGICAL 00008 00009 end module mpimod 00010 ! 00011 ! interface routines to MPI: 00012 ! 00013 00014 00015 ! ================ 00016 ! SUBROUTINE MPBCI 00017 ! ================ 00018 00019 subroutine mpbci(k) ! broadcast 1 integer 00020 use mpimod 00021 00022 call mpi_bcast(k,1,mpi_itype,NROOT,myworld,mpinfo) 00023 00024 return 00025 end subroutine mpbci 00026 00027 ! ================= 00028 ! SUBROUTINE MPBCIN 00029 ! ================= 00030 00031 subroutine mpbcin(k,n) ! broadcast n integer 00032 use mpimod 00033 00034 integer k(n) 00035 00036 call mpi_bcast(k,n,mpi_itype,NROOT,myworld,mpinfo) 00037 00038 return 00039 end subroutine mpbcin 00040 00041 ! ================ 00042 ! SUBROUTINE MPBCR 00043 ! ================ 00044 00045 subroutine mpbcr(p) ! broadcast 1 real 00046 use mpimod 00047 00048 call mpi_bcast(p,1,mpi_rtype,NROOT,myworld,mpinfo) 00049 00050 return 00051 end subroutine mpbcr 00052 00053 ! ================= 00054 ! SUBROUTINE MPBCRN 00055 ! ================= 00056 00057 subroutine mpbcrn(p,n) ! broadcast n real 00058 use mpimod 00059 00060 real p(n) 00061 00062 call mpi_bcast(p,n,mpi_rtype,NROOT,myworld,mpinfo) 00063 00064 return 00065 end subroutine mpbcrn 00066 00067 ! ================ 00068 ! SUBROUTINE MPBCL 00069 ! ================ 00070 00071 subroutine mpbcl(k) ! broadcast 1 logical 00072 use mpimod 00073 logical k 00074 00075 call mpi_bcast(k,1,mpi_ltype,NROOT,myworld,mpinfo) 00076 00077 return 00078 end subroutine mpbcl 00079 00080 ! ================= 00081 ! SUBROUTINE MPSCIN 00082 ! ================= 00083 00084 subroutine mpscin(k,n) ! scatter n integer 00085 use mpimod 00086 00087 integer k(*) 00088 00089 call mpi_scatter(k,n,mpi_itype,k,n,mpi_itype & 00090 & ,NROOT,myworld,mpinfo) 00091 00092 return 00093 end subroutine mpscin 00094 00095 ! ================= 00096 ! SUBROUTINE MPSCRN 00097 ! ================= 00098 00099 subroutine mpscrn(p,n) ! scatter n real 00100 use mpimod 00101 00102 real p(*) 00103 00104 call mpi_scatter(p,n,mpi_rtype,p,n,mpi_rtype,NROOT,myworld,mpinfo) 00105 00106 return 00107 end subroutine mpscrn 00108 00109 ! ================= 00110 ! SUBROUTINE MPSCDN 00111 ! ================= 00112 00113 subroutine mpscdn(p,n) ! scatter n double precision 00114 use mpimod 00115 00116 real (kind=8) :: p(*) 00117 00118 call mpi_scatter(p,n,MPI_REAL8,p,n,MPI_REAL8,NROOT,myworld,mpinfo) 00119 00120 return 00121 end subroutine mpscdn 00122 00123 ! ================= 00124 ! SUBROUTINE MPSCGP 00125 ! ================= 00126 00127 subroutine mpscgp(pf,pp,klev) ! scatter gridpoint fields 00128 use mpimod 00129 00130 real pf(NUGP,klev) 00131 real pp(NHOR,klev) 00132 00133 do jlev = 1 , klev 00134 call mpi_scatter(pf(1,jlev),NHOR,mpi_rtype, & 00135 & pp(1,jlev),NHOR,mpi_rtype, & 00136 & NROOT,myworld,mpinfo) 00137 enddo 00138 00139 return 00140 end subroutine mpscgp 00141 00142 ! ================= 00143 ! SUBROUTINE MPGAGP 00144 ! ================= 00145 00146 subroutine mpgagp(pf,pp,klev) ! gather gridpoint fields 00147 use mpimod 00148 00149 real pf(NLON*NLAT,klev) 00150 real pp(NHOR,klev) 00151 00152 do jlev = 1 , klev 00153 call mpi_gather(pp(1,jlev),NHOR,mpi_rtype, & 00154 & pf(1,jlev),NHOR,mpi_rtype, & 00155 & NROOT,myworld,mpinfo) 00156 enddo 00157 00158 return 00159 end subroutine mpgagp 00160 00161 ! =================== 00162 ! SUBROUTINE MPGALLGP 00163 ! =================== 00164 00165 subroutine mpgallgp(pf,pp,klev) ! gather gritpoint to all 00166 use mpimod 00167 00168 real pf(NLON*NLAT,klev) 00169 real pp(NHOR,klev) 00170 00171 do jlev = 1 , klev 00172 call mpi_allgather(pp(1,jlev),NHOR,mpi_rtype, & 00173 & pf(1,jlev),NHOR,mpi_rtype, & 00174 & myworld,mpinfo) 00175 enddo 00176 00177 return 00178 end subroutine mpgallgp 00179 00180 ! ================= 00181 ! SUBROUTINE MPSCSP 00182 ! ================= 00183 00184 subroutine mpscsp(pf,pp,klev) ! scatter spectral fields 00185 use mpimod 00186 00187 real pf(NESP,klev) 00188 real pp(NSPP,klev) 00189 00190 do jlev = 1 , klev 00191 call mpi_scatter(pf(1,jlev),NSPP,mpi_rtype & 00192 & ,pp(1,jlev),NSPP,mpi_rtype & 00193 & ,NROOT,myworld,mpinfo) 00194 enddo 00195 00196 return 00197 end subroutine mpscsp 00198 00199 ! ================= 00200 ! SUBROUTINE MPGASP 00201 ! ================= 00202 00203 subroutine mpgasp(pf,pp,klev) ! gather spectral fields 00204 use mpimod 00205 00206 real pf(NESP,klev) 00207 real pp(NSPP,klev) 00208 00209 do jlev = 1 , klev 00210 call mpi_gather(pp(1,jlev),NSPP,mpi_rtype & 00211 & ,pf(1,jlev),NSPP,mpi_rtype & 00212 & ,NROOT,myworld,mpinfo) 00213 enddo 00214 00215 return 00216 end subroutine mpgasp 00217 00218 ! ================= 00219 ! SUBROUTINE MPGACS 00220 ! ================= 00221 00222 subroutine mpgacs(pcs) ! gather cross sections 00223 use mpimod 00224 00225 real pcs(NLAT,NLEV) 00226 00227 do jlev = 1 , NLEV 00228 call mpi_gather(pcs(1,jlev),NLPP,mpi_rtype & 00229 & ,pcs(1,jlev),NLPP,mpi_rtype & 00230 & ,NROOT,myworld,mpinfo) 00231 00232 enddo 00233 return 00234 end subroutine mpgacs 00235 00236 ! =================== 00237 ! SUBROUTINE MPGALLSP 00238 ! =================== 00239 00240 subroutine mpgallsp(pf,pp,klev) ! gather spectral to all 00241 use mpimod 00242 00243 real pf(NESP,klev) 00244 real pp(NSPP,klev) 00245 00246 do jlev = 1 , klev 00247 call mpi_allgather(pp(1,jlev),NSPP,mpi_rtype & 00248 & ,pf(1,jlev),NSPP,mpi_rtype & 00249 & ,myworld,mpinfo) 00250 enddo 00251 00252 return 00253 end subroutine mpgallsp 00254 00255 ! ================ 00256 ! SUBROUTINE MPSUM 00257 ! ================ 00258 00259 subroutine mpsum(psp,klev) ! sum spectral fields 00260 use mpimod 00261 00262 real psp(NESP,klev) 00263 real tmp(NESP,klev) 00264 00265 call mpi_reduce(psp(1,1),tmp(1,1),NESP*klev,mpi_rtype,MPI_SUM & 00266 & ,NROOT,myworld,mpinfo) 00267 if (mypid == NROOT) psp = tmp 00268 00269 return 00270 end subroutine mpsum 00271 00272 ! ================== 00273 ! SUBROUTINE MPSUMSC 00274 ! ================== 00275 00276 subroutine mpsumsc(psf,psp,klev) ! sum & scatter spectral 00277 use mpimod 00278 00279 real psf(NESP,klev) 00280 real psp(NSPP,klev) 00281 00282 do jlev = 1 , klev 00283 call mpi_reduce_scatter(psf(1,jlev),psp(1,jlev),nscatsp & 00284 & ,mpi_rtype,MPI_SUM,myworld,mpinfo) 00285 enddo 00286 00287 return 00288 end subroutine mpsumsc 00289 00290 ! ==================== 00291 ! SUBROUTINE MPSUMR 00292 ! ==================== 00293 00294 subroutine mpsumr(pr,kdim) ! sum kdim reals 00295 use mpimod 00296 00297 real pr(kdim) 00298 real tmp(kdim) 00299 00300 call mpi_reduce(pr,tmp,kdim,mpi_rtype,MPI_SUM,NROOT,myworld,mpinfo) 00301 if (mypid == NROOT) pr = tmp 00302 00303 return 00304 end subroutine mpsumr 00305 00306 ! ==================== 00307 ! SUBROUTINE MPSUMBCR 00308 ! ==================== 00309 00310 subroutine mpsumbcr(pr,kdim) ! sum & broadcast kdim reals 00311 use mpimod 00312 00313 real :: pr(kdim) 00314 real :: tmp(kdim) 00315 00316 call mpi_allreduce(pr,tmp,kdim,mpi_rtype,MPI_SUM,myworld,mpinfo) 00317 pr = tmp 00318 00319 return 00320 end subroutine mpsumbcr 00321 00322 ! ================== 00323 ! SUBROUTINE MPSTART 00324 ! ================== 00325 00326 subroutine mpstart ! initialization 00327 use mpimod 00328 integer :: itest = 0 00329 real :: rtest = 0.0 00330 00331 if (kind(itest) == 8) mpi_itype = MPI_INTEGER8 00332 if (kind(rtest) == 8) mpi_rtype = MPI_REAL8 00333 00334 call mpi_init(mpinfo) 00335 myworld=MPI_COMM_WORLD 00336 00337 call mpi_comm_size(myworld,npro ,mpinfo) 00338 call mpi_comm_rank(myworld,mypid,mpinfo) 00339 allocate(ympname(npro)) ; ympname(:) = ' ' 00340 call mpi_get_processor_name(ympname(1),ilen,mpinfo) 00341 00342 call mpi_gather(ympname,80,MPI_CHARACTER, & 00343 ympname,80,MPI_CHARACTER, & 00344 NROOT,myworld,mpinfo) 00345 00346 return 00347 end subroutine mpstart 00348 00349 ! ================= 00350 ! SUBROUTINE MPSTOP 00351 ! ================= 00352 00353 subroutine mpstop 00354 use mpimod 00355 00356 call mpi_barrier(myworld,mpinfo) 00357 call mpi_finalize(mpinfo) 00358 00359 return 00360 end subroutine mpstop 00361 00362 ! =================== 00363 ! SUBROUTINE MPREADGP 00364 ! =================== 00365 00366 subroutine mpreadgp(ktape,p,kdim,klev) 00367 use mpimod 00368 00369 real p(kdim,klev) 00370 real z(NLON*NLAT,klev) 00371 00372 z = 0.0 00373 00374 if (mypid == NROOT) read (ktape) z(1:NLON*NLAT,:) 00375 00376 if (kdim == NHOR) then 00377 call mpscgp(z,p,klev) 00378 else 00379 if (mypid == NROOT) p = z 00380 endif 00381 00382 return 00383 end subroutine mpreadgp 00384 00385 ! ==================== 00386 ! SUBROUTINE MPWRITEGP 00387 ! ==================== 00388 00389 subroutine mpwritegp(ktape,p,kdim,klev) 00390 use mpimod 00391 00392 real p(kdim,klev) 00393 real z(NLON*NLAT,klev) 00394 00395 if (kdim == NHOR) then 00396 call mpgagp(z,p,klev) 00397 if (mypid == NROOT) write(ktape) z(1:NLON*NLAT,:) 00398 else 00399 if (mypid == NROOT) write(ktape) p(1:NLON*NLAT,:) 00400 endif 00401 00402 return 00403 end subroutine mpwritegp 00404 00405 00406 ! ===================== 00407 ! SUBROUTINE MPWRITEGPH 00408 ! ===================== 00409 00410 subroutine mpwritegph(ktape,p,kdim,klev,ihead) 00411 use mpimod 00412 00413 real p(kdim,klev) 00414 real z(NLON*NLAT,klev) 00415 ! 00416 real(kind=4) :: zp(kdim,klev) 00417 real(kind=4) :: zz(NLON*NLAT,klev) 00418 ! 00419 00420 integer ihead(8) 00421 00422 if (kdim == NHOR) then 00423 call mpgagp(z,p,klev) 00424 if (mypid == NROOT) then 00425 write(ktape) ihead 00426 zz(:,:)=z(:,:) 00427 write(ktape) zz(1:NLON*NLAT,:) 00428 endif 00429 else 00430 if (mypid == NROOT) then 00431 write(ktape) ihead 00432 zp(:,:)=p(:,:) 00433 write(ktape) zp(1:NLON*NLAT,:) 00434 endif 00435 endif 00436 00437 return 00438 end subroutine mpwritegph 00439 00440 ! =================== 00441 ! SUBROUTINE MPREADSP 00442 ! =================== 00443 00444 subroutine mpreadsp(ktape,p,kdim,klev) 00445 use mpimod 00446 00447 real p(kdim,klev) 00448 real z(NESP,klev) 00449 00450 z = 0.0 00451 if (mypid == NROOT) read(ktape) z(1:NRSP,:) 00452 if (kdim == NSPP) then 00453 call mpscsp(z,p,klev) 00454 else 00455 if (mypid == NROOT) p = z 00456 endif 00457 00458 return 00459 end subroutine mpreadsp 00460 00461 ! ==================== 00462 ! SUBROUTINE MPWRITESP 00463 ! ==================== 00464 00465 subroutine mpwritesp(ktape,p,kdim,klev) 00466 use mpimod 00467 00468 real p(kdim,klev) 00469 real z(NESP,klev) 00470 00471 if (kdim == NSPP) then 00472 call mpgasp(z,p,klev) 00473 if (mypid == NROOT) write(ktape) z(1:NRSP,:) 00474 else 00475 if (mypid == NROOT) write(ktape) p(1:NRSP,:) 00476 endif 00477 00478 return 00479 end subroutine mpwritesp 00480 00481 00482 ! =================== 00483 ! SUBROUTINE MPI_INFO 00484 ! =================== 00485 00486 subroutine mpi_info(nprocess,npid) ! get nproc and pid 00487 use mpimod 00488 00489 myworld=MPI_COMM_WORLD 00490 00491 call mpi_comm_size(myworld,nprocess,mpinfo) 00492 call mpi_comm_rank(myworld,npid,mpinfo) 00493 00494 return 00495 end subroutine mpi_info 00496 00497 00498 ! ================== 00499 ! SUBROUTINE MPGETSP 00500 ! ================== 00501 00502 subroutine mpgetsp(yn,p,kdim,klev) 00503 use mpimod 00504 00505 character (len=*) :: yn 00506 real :: p(kdim,klev) 00507 real :: z(NESP,klev) 00508 00509 z(:,:) = 0.0 00510 if (mypid == NROOT) call get_restart_array(yn,z,NRSP,NESP,klev) 00511 call mpscsp(z,p,klev) 00512 00513 return 00514 end subroutine mpgetsp 00515 00516 00517 ! ================== 00518 ! SUBROUTINE MPGETGP 00519 ! ================== 00520 00521 subroutine mpgetgp(yn,p,kdim,klev) 00522 use mpimod 00523 00524 character (len=*) :: yn 00525 real :: p(kdim,klev) 00526 real :: z(NUGP,klev) 00527 00528 if (mypid == NROOT) call get_restart_array(yn,z,NUGP,NUGP,klev) 00529 call mpscgp(z,p,klev) 00530 00531 return 00532 end subroutine mpgetgp 00533 00534 00535 ! ================== 00536 ! SUBROUTINE MPPUTSP 00537 ! ================== 00538 00539 subroutine mpputsp(yn,p,kdim,klev) 00540 use mpimod 00541 00542 character (len=*) :: yn 00543 real :: p(kdim,klev) 00544 real :: z(NESP,klev) 00545 00546 call mpgasp(z,p,klev) 00547 if (mypid == NROOT) call put_restart_array(yn,z,NRSP,NESP,klev) 00548 00549 return 00550 end subroutine mpputsp 00551 00552 00553 ! ================== 00554 ! SUBROUTINE MPPUTGP 00555 ! ================== 00556 00557 subroutine mpputgp(yn,p,kdim,klev) 00558 use mpimod 00559 00560 character (len=*) :: yn 00561 real :: p(kdim,klev) 00562 real :: z(NUGP,klev) 00563 00564 call mpgagp(z,p,klev) 00565 if (mypid == NROOT) call put_restart_array(yn,z,NUGP,NUGP,klev) 00566 00567 return 00568 end subroutine mpputgp 00569 00570 00571 !! =================== 00572 !! SUBROUTINE MPSURFGP 00573 !! =================== 00574 ! 00575 ! subroutine mpsurfgp(yn,p,kdim,klev) 00576 ! use mpimod 00577 ! 00578 ! character (len=*) :: yn 00579 ! real :: p(kdim,klev) 00580 ! real :: z(NUGP,klev) 00581 ! 00582 ! if (mypid == NROOT) call get_surf_array(yn,z,NUGP,NUGP,klev,iread) 00583 ! call mpbci(iread) 00584 ! if (iread == 1) call mpscgp(z,p,klev) 00585 ! 00586 ! return 00587 ! end subroutine mpsurfgp 00588 ! 00589 ! 00590 !! ===================== 00591 !! SUBROUTINE MPSURFYEAR 00592 !! ===================== 00593 ! 00594 ! subroutine mpsurfyear(yn,p,kdim,kmon) 00595 ! use mpimod 00596 ! 00597 ! character (len=*) :: yn 00598 ! real :: p(kdim,kmon) 00599 ! real :: z(NUGP,kmon) 00600 ! 00601 ! if (mypid == NROOT) call get_surf_year(yn,z,NUGP,kmon,iread) 00602 ! call mpbci(iread) 00603 ! if (iread == 1) call mpscgp(z,p,kmon) 00604 ! 00605 ! return 00606 ! end subroutine mpsurfyear 00607 ! 00608 ! 00609 !! ===================== 00610 !! SUBROUTINE MP3DYEAR 00611 !! ===================== 00612 ! 00613 ! subroutine mp3dyear(yn,p,kdim,klev,kmon) 00614 ! use mpimod 00615 ! 00616 ! character (len=*) :: yn 00617 ! real :: p(kdim,klev,kmon) 00618 ! real :: z(NUGP,klev,kmon) 00619 ! 00620 ! if (mypid == NROOT) call get_3d_year(yn,z,NUGP,klev,kmon,iread) 00621 ! call mpbci(iread) 00622 ! if (iread == 1) call mpscgp(z,p,klev*kmon) 00623 ! 00624 ! return 00625 ! end subroutine mp3dyear 00626 00627 00628 ! =================== 00629 ! SUBROUTINE MPMAXVAL 00630 ! =================== 00631 00632 subroutine mpmaxval(p,kdim,klev,pmax) 00633 use mpimod 00634 00635 real :: p(kdim,klev) 00636 00637 zmax = maxval(p(:,:)) 00638 call mpi_allreduce(zmax,pmax,1,mpi_rtype,MPI_MAX,myworld,mpinfo) 00639 00640 return 00641 end subroutine mpmaxval 00642 00643 00644 ! =================== 00645 ! SUBROUTINE MPSUMVAL 00646 ! =================== 00647 00648 subroutine mpsumval(p,kdim,klev,psum) 00649 use mpimod 00650 00651 real :: p(kdim,klev) 00652 00653 zsum = sum(p(:,:)) 00654 call mpi_allreduce(zsum,psum,1,mpi_rtype,MPI_SUM,myworld,mpinfo) 00655 00656 return 00657 end subroutine mpsumval 00658 00659 ! Some dummy declarations that are use in multirun mode only 00660 00661 subroutine mrdiff(p,d,n) 00662 real :: p(n) 00663 real :: d(n) 00664 return 00665 end 00666 00667 subroutine mrsum(k) ! sum up 1 integer 00668 return 00669 end 00670 00671 subroutine mrbci(k) ! broadcast 1 integer 00672 return 00673 end 00674 00675 subroutine mrdimensions 00676 return 00677 end 00678 00679