Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Use threading to speed up MERRA2 aerosol interpolation and NoahMP initialization #605

Merged
merged 6 commits into from
Apr 15, 2021
88 changes: 68 additions & 20 deletions physics/GFS_phys_time_vary.fv3.F90
Original file line number Diff line number Diff line change
Expand Up @@ -305,8 +305,8 @@ subroutine GFS_phys_time_vary_init (
call setindxh2o (im, xlat_d, jindx1_h, jindx2_h, ddy_h)
endif

!> - Call setindxaer() to initialize aerosols data
!$OMP section
!> - Call setindxaer() to initialize aerosols data
if (iaerclm) then
call setindxaer (im, xlat_d, jindx1_aer, &
jindx2_aer, ddy_aer, xlon_d, &
Expand All @@ -317,8 +317,8 @@ subroutine GFS_phys_time_vary_init (
jamin=min(minval(jindx1_aer), jamin)
jamax=max(maxval(jindx2_aer), jamax)
endif
!$OMP section

!$OMP section
!> - Call setindxci() to initialize IN and CCN data
if (iccn == 1) then
call setindxci (im, xlat_d, jindx1_ci, &
Expand Down Expand Up @@ -376,10 +376,14 @@ subroutine GFS_phys_time_vary_init (
!$OMP end sections

!$OMP end parallel

if (errflg/=0) return

if (iaerclm) then
call read_aerdataf (iamin, iamax, jamin, jamax, me,master,iflip, &
idate,errmsg,errflg)
endif
call read_aerdataf (iamin, iamax, jamin, jamax, me, master, iflip, &
idate, errmsg, errflg)
if (errflg/=0) return
end if

if (lsm == lsm_noahmp) then
if (all(tvxy < zero)) then
Expand Down Expand Up @@ -431,15 +435,34 @@ subroutine GFS_phys_time_vary_init (
smoiseq(:,:) = missing_value
zsnsoxy(:,:) = missing_value

imn = idate(2)

!$OMP parallel do num_threads(nthrds) default(none) &
!$OMP shared(im,lsoil,con_t0c,landfrac,tsfcl,tvxy,tgxy,tahxy) &
!$OMP shared(snowd,canicexy,canliqxy,canopy,eahxy,cmxy,chxy) &
!$OMP shared(fwetxy,sneqvoxy,weasd,alboldxy,qsnowxy,wslakexy) &
!$OMP shared(taussxy,albdvis,albdnir,albivis,albinir,emiss) &
!$OMP shared(waxy,wtxy,zwtxy,imn,vtype,xlaixy,xsaixy,lfmassxy) &
!$OMP shared(stmassxy,rtmassxy,woodxy,stblcpxy,fastcpxy) &
!$OMP shared(isbarren_table,isice_table,isurban_table) &
!$omp shared(iswater_table,laim_table,sla_table,bexp_table) &
!$omp shared(stc,smc,slc,tg3,snowxy,tsnoxy,snicexy,snliqxy) &
!$omp shared(zsnsoxy,STYPE,SMCMAX_TABLE,SMCWLT_TABLE,zs,dzs) &
!$omp shared(DWSAT_TABLE,DKSAT_TABLE,PSISAT_TABLE,smoiseq) &
!$OMP shared(smcwtdxy,deeprechxy,rechxy,errmsg,errflg) &
!$OMP private(vegtyp,masslai,masssai,snd,dzsno,dzsnso,isnow) &
!$OMP private(soiltyp,bexp,smcmax,smcwlt,dwsat,dksat,psisat,ddz)
do ix=1,im
if (landfrac(ix) >= drythresh) then
tvxy(ix) = tsfcl(ix)
tgxy(ix) = tsfcl(ix)
tahxy(ix) = tsfcl(ix)

if (snowd(ix) > 0.01_kind_phys .and. tsfcl(ix) > con_t0c ) tvxy(ix) = con_t0c
if (snowd(ix) > 0.01_kind_phys .and. tsfcl(ix) > con_t0c ) tgxy(ix) = con_t0c
if (snowd(ix) > 0.01_kind_phys .and. tsfcl(ix) > con_t0c ) tahxy(ix) = con_t0c
if (snowd(ix) > 0.01_kind_phys .and. tsfcl(ix) > con_t0c ) then
tvxy(ix) = con_t0c
tgxy(ix) = con_t0c
tahxy(ix) = con_t0c
end if

canicexy(ix) = 0.0_kind_phys
canliqxy(ix) = canopy(ix)
Expand All @@ -463,14 +486,12 @@ subroutine GFS_phys_time_vary_init (
albinir(ix) = 0.2_kind_phys
emiss(ix) = 0.95_kind_phys


waxy(ix) = 4900.0_kind_phys
wtxy(ix) = waxy(ix)
zwtxy(ix) = (25.0_kind_phys + 2.0_kind_phys) - waxy(ix) / 1000.0_kind_phys / 0.2_kind_phys

vegtyp = vtype(ix)
if (vegtyp == 0) vegtyp = 7
imn = idate(2)

if ((vegtyp == isbarren_table) .or. (vegtyp == isice_table) .or. (vegtyp == isurban_table) .or. (vegtyp == iswater_table)) then

Expand Down Expand Up @@ -552,7 +573,6 @@ subroutine GFS_phys_time_vary_init (
else
errmsg = 'Error in GFS_phys_time_vary.fv3.F90: Problem with the logic assigning snow layers in Noah MP initialization'
errflg = 1
return
endif

! Now we have the snowxy field
Expand Down Expand Up @@ -628,6 +648,9 @@ subroutine GFS_phys_time_vary_init (
endif

enddo ! ix
!$OMP end parallel do

if (errflg/=0) return

deallocate(dzsno)
deallocate(dzsnso)
Expand Down Expand Up @@ -748,6 +771,20 @@ subroutine GFS_phys_time_vary_timestep_init (
return
end if

!$OMP parallel num_threads(nthrds) default(none) &
!$OMP shared(kdt,nsswr,lsswr,clstp,imfdeepcnv,cal_pre,random_clds) &
!$OMP shared(fhswr,fhour,seed0,cnx,cny,nrcm,wrk,rannie,rndval) &
!$OMP shared(rann,im,isc,jsc,imap,jmap,ntoz,me,idate,jindx1_o3,jindx2_o3) &
!$OMP shared(ozpl,ddy_o3,h2o_phys,jindx1_h,jindx2_h,h2opl,ddy_h,iaerclm,master) &
!$OMP shared(levs,prsl,iccn,jindx1_ci,jindx2_ci,ddy_ci,iindx1_ci,iindx2_ci) &
!$OMP shared(ddx_ci,in_nm,ccn_nm,do_ugwp_v1,jindx1_tau,jindx2_tau,ddy_j1tau) &
!$OMP shared(ddy_j2tau,tau_amf) &
!$OMP private(iseed,iskip,i,j,k)

!$OMP sections

!$OMP section

!--- switch for saving convective clouds - cnvc90.f
!--- aka Ken Campana/Yu-Tai Hou legacy
if ((mod(kdt,nsswr) == 0) .and. (lsswr)) then
Expand All @@ -764,6 +801,8 @@ subroutine GFS_phys_time_vary_timestep_init (
clstp = 0100
endif

!$OMP section

!--- random number needed for RAS and old SAS and when cal_pre=.true.
! imfdeepcnv < 0 when ras = .true.
if ( (imfdeepcnv <= 0 .or. cal_pre) .and. random_clds ) then
Expand All @@ -789,29 +828,23 @@ subroutine GFS_phys_time_vary_timestep_init (

endif ! imfdeepcnv, cal_re, random_clds

!$OMP section
!> - Call ozinterpol() to make ozone interpolation
if (ntoz > 0) then
call ozinterpol (me, im, idate, fhour, &
jindx1_o3, jindx2_o3, &
ozpl, ddy_o3)
endif

!$OMP section
!> - Call h2ointerpol() to make stratospheric water vapor data interpolation
if (h2o_phys) then
call h2ointerpol (me, im, idate, fhour, &
jindx1_h, jindx2_h, &
h2opl, ddy_h)
endif

!> - Call aerinterpol() to make aerosol interpolation
if (iaerclm) then
call aerinterpol (me, master, im, idate, fhour, &
jindx1_aer, jindx2_aer, &
ddy_aer, iindx1_aer, &
iindx2_aer, ddx_aer, &
levs, prsl, aer_nm)
endif

!$OMP section
!> - Call ciinterpol() to make IN and CCN data interpolation
if (iccn == 1) then
call ciinterpol (me, im, idate, fhour, &
Expand All @@ -821,13 +854,28 @@ subroutine GFS_phys_time_vary_timestep_init (
levs, prsl, in_nm, ccn_nm)
endif

!$OMP section
!> - Call cires_indx_ugwp to read monthly-mean GW-tau diagnosed from FV3GFS-runs that resolve GW-activ
if (do_ugwp_v1) then
call tau_amf_interp(me, master, im, idate, fhour, &
jindx1_tau, jindx2_tau, &
ddy_j1tau, ddy_j2tau, tau_amf)
endif

!$OMP end sections
!$OMP end parallel

!> - Call aerinterpol() to make aerosol interpolation
if (iaerclm) then
! aerinterpol is using threading inside, don't
! move into OpenMP parallel section above
call aerinterpol (me, master, nthrds, im, idate, &
fhour, jindx1_aer, jindx2_aer,&
ddy_aer, iindx1_aer, &
iindx2_aer, ddx_aer, &
levs, prsl, aer_nm)
endif

!> - Call gcycle() to repopulate specific time-varying surface properties for AMIP/forecast runs
if (nscyc > 0) then
if (mod(kdt,nscyc) == 1) THEN
Expand Down
2 changes: 1 addition & 1 deletion physics/aerclm_def.F
Original file line number Diff line number Diff line change
Expand Up @@ -13,7 +13,7 @@ module aerclm_def
real (kind=kind_phys), allocatable, dimension(:,:,:,:) :: aer_pres
real (kind=kind_phys), allocatable, dimension(:,:,:,:,:) :: aerin

data aer_time/15.5, 45., 74.5, 105., 135.5, 166., 196.5,
data aer_time/15.5, 45., 74.5, 105., 135.5, 166., 196.5,
& 227.5, 258., 288.5, 319., 349.5, 380.5/

data specname /'DU001','DU002','DU003','DU004','DU005',
Expand Down
63 changes: 39 additions & 24 deletions physics/aerinterp.F90
Original file line number Diff line number Diff line change
Expand Up @@ -43,15 +43,23 @@ SUBROUTINE read_aerdata (me, master, iflip, idate, errmsg, errflg)
endif
!
!! ===================================================================
!! check if all necessary files exist
!! ===================================================================
do imon = 1, timeaer
write(mn,'(i2.2)') imon
fname=trim("aeroclim.m"//mn//".nc")
inquire (file = fname, exist = file_exist)
if (.not. file_exist) then
errmsg = 'Error in read_aerdata: file ' // trim(fname) // ' not found'
errflg = 1
return
endif
enddo
!
!! ===================================================================
!! fetch dim spec and lat/lon from m01 data set
!! ===================================================================
fname=trim("aeroclim.m"//'01'//".nc")
inquire (file = fname, exist = file_exist)
if (.not. file_exist) then
errmsg = 'Error in read_aerdata: file ' // trim(fname) // ' not found'
errflg = 1
return
endif
call nf_open(fname , nf90_NOWRITE, ncid)

vname = trim(specname(1))
Expand Down Expand Up @@ -117,27 +125,16 @@ SUBROUTINE read_aerdataf (iamin, iamax, jamin, jamax, &
endif

! allocate local working arrays
if (.not. allocated(buff)) then
allocate (buff(lonsaer, latsaer, levsw))
allocate (pres_tmp(lonsaer,levsw))
endif
if (.not. allocated(buffx)) then
allocate (buffx(lonsaer, latsaer, levsw,1))
endif
allocate (buff(lonsaer, latsaer, levsw))
allocate (pres_tmp(lonsaer, levsw))
allocate (buffx(lonsaer, latsaer, levsw, 1))

!! ===================================================================
!! loop thru m01 - m12 for aer/pres array
!! ===================================================================
do imon = 1, timeaer
write(mn,'(i2.2)') imon
fname=trim("aeroclim.m"//mn//".nc")
inquire (file = fname, exist = file_exist)
if (.not. file_exist) then
errmsg = 'Error in read_aerdata: file ' // trim(fname) // ' not found'
errflg = 1
return
endif

call nf_open(fname , nf90_NOWRITE, ncid)

! ====> construct 3-d pressure array (Pa)
Expand Down Expand Up @@ -259,7 +256,7 @@ END SUBROUTINE setindxaer
!**********************************************************************
!**********************************************************************
!
SUBROUTINE aerinterpol(me,master,npts,IDATE,FHOUR,jindx1,jindx2, &
SUBROUTINE aerinterpol(me,master,nthrds,npts,IDATE,FHOUR,jindx1,jindx2, &
ddy,iindx1,iindx2,ddx,lev,prsl,aerout)
!
USE MACHINE, ONLY : kind_phys
Expand All @@ -270,7 +267,7 @@ SUBROUTINE aerinterpol(me,master,npts,IDATE,FHOUR,jindx1,jindx2, &
!

integer JINDX1(npts), JINDX2(npts),iINDX1(npts),iINDX2(npts)
integer me,idate(4), master
integer me,idate(4), master, nthrds
integer IDAT(8),JDAT(8)
!
real(kind=kind_phys) DDY(npts), ddx(npts),ttt
Expand Down Expand Up @@ -317,7 +314,16 @@ SUBROUTINE aerinterpol(me,master,npts,IDATE,FHOUR,jindx1,jindx2, &
tx2 = 1.0 - tx1
if (n2 > 12) n2 = n2 -12

!
#ifndef __GFORTRAN__
!$OMP parallel num_threads(nthrds) default(none) &
!$OMP shared(npts,ntrcaer,aerin,aer_pres,prsl) &
!$OMP shared(ddx,ddy,jindx1,jindx2,iindx1,iindx2) &
!$OMP shared(aerpm,aerpres,aerout,n1,n2,lev,nthrds) &
!$OMP private(l,j,k,ii,i1,i2,j1,j2,temj,temi) &
!$OMP copyin(tx1,tx2) firstprivate(tx1,tx2)

!$OMP do
#endif
DO L=1,levsaer
DO J=1,npts
J1 = JINDX1(J)
Expand All @@ -341,8 +347,12 @@ SUBROUTINE aerinterpol(me,master,npts,IDATE,FHOUR,jindx1,jindx2, &
+TEMI*DDY(j)*aer_pres(I1,J2,L,n2)+DDX(j)*TEMJ*aer_pres(I2,J1,L,n2))
ENDDO
ENDDO
#ifndef __GFORTRAN__
!$OMP end do

! don't flip, input is the same direction as GFS (bottom-up)
!$OMP do
#endif
DO J=1,npts
DO L=1,lev
if(prsl(j,L).ge.aerpres(j,1)) then
Expand Down Expand Up @@ -371,7 +381,12 @@ SUBROUTINE aerinterpol(me,master,npts,IDATE,FHOUR,jindx1,jindx2, &
endif
ENDDO !L-loop
ENDDO !J-loop
!
#ifndef __GFORTRAN__
!$OMP end do

!$OMP end parallel
#endif

RETURN
END SUBROUTINE aerinterpol

Expand Down