Skip to content
Draft
68 changes: 39 additions & 29 deletions physics/GWD/cires_tauamf_data.F90
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,7 @@
module cires_tauamf_data

use machine, only: kind_phys
use mpi_f08
!...........................................................................................
! tabulated GW-sources: GRACILE/Ern et al., 2018 and/or Resolved GWs from C384-Annual run
!...........................................................................................
Expand All @@ -20,27 +21,28 @@ module cires_tauamf_data
contains

!>
subroutine read_tau_amf(me, master, errmsg, errflg)
subroutine read_tau_amf(mpicomm, mpirank, mpiroot, errmsg, errflg)

use netcdf
integer, intent(in) :: me, master
use netcdf
use mpiutil, only: ccpp_bcast
type(MPI_Comm), intent(in) :: mpicomm
integer, intent(in) :: mpirank, mpiroot
integer :: ncid, iernc, vid, dimid, status
integer :: k

character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
!
iernc=NF90_OPEN(trim(ugwp_taufile), nf90_nowrite, ncid)
read_and_broadcast_1: if (mpirank==mpiroot) then
iernc=NF90_OPEN(trim(ugwp_taufile), nf90_nowrite, ncid)

if(iernc.ne.0) then
write(errmsg,'(*(a))') "read_tau_amf: cannot open file_limb_tab data-file ", &
trim(ugwp_taufile)
print *, 'cannot open ugwp-v1 tau-file=',trim(ugwp_taufile)
errflg = 1
return
else

endif

status = nf90_inq_dimid(ncid, "lat", DimID)
! if (status /= nf90_noerr) call handle_err(status)
Expand All @@ -50,28 +52,36 @@ subroutine read_tau_amf(me, master, errmsg, errflg)
status = nf90_inq_dimid(ncid, "days", DimID)
status = nf90_inquire_dimension(ncid, DimID, len =ntau_d2t )

if (me == master) print *, ntau_d1y, ntau_d2t, ' dimd of tau_ngw ugwp-v1 '
if (ntau_d2t .le. 0 .or. ntau_d1y .le. 0) then
print *, 'ugwp-v1 tau-file=', trim(ugwp_taufile)
print *, ' ugwp-v1: ', 'ntau_d2t=',ntau_d2t, 'ntau_d2t=',ntau_d1y
stop
endif

if (.not.allocated(ugwp_taulat)) allocate (ugwp_taulat(ntau_d1y ))
if (.not.allocated(days_limb)) allocate (days_limb(ntau_d2t))
if (.not.allocated(tau_limb)) allocate (tau_limb(ntau_d1y, ntau_d2t ))

iernc=nf90_inq_varid( ncid, 'DAYS', vid )
iernc= nf90_get_var( ncid, vid, days_limb)
iernc=nf90_inq_varid( ncid, 'LATS', vid )
iernc= nf90_get_var( ncid, vid, ugwp_taulat)
iernc=nf90_inq_varid( ncid, 'ABSMF', vid )
iernc= nf90_get_var( ncid, vid, tau_limb)

iernc=nf90_close(ncid)

endif

print *, ntau_d1y, ntau_d2t, ' dimd of tau_ngw ugwp-v1 '
if (ntau_d2t .le. 0 .or. ntau_d1y .le. 0) then
print *, 'ugwp-v1 tau-file=', trim(ugwp_taufile)
print *, ' ugwp-v1: ', 'ntau_d2t=',ntau_d2t, 'ntau_d2t=',ntau_d1y
errflg = 1
return
endif
endif read_and_broadcast_1

call ccpp_bcast(ntau_d1y, mpiroot, mpicomm, errflg)
call ccpp_bcast(ntau_d2t, mpiroot, mpicomm, errflg)

if (.not.allocated(ugwp_taulat)) allocate (ugwp_taulat(ntau_d1y ))
if (.not.allocated(days_limb)) allocate (days_limb(ntau_d2t))
if (.not.allocated(tau_limb)) allocate (tau_limb(ntau_d1y, ntau_d2t ))

read_and_broadcast_2: if (mpirank==mpiroot) then
iernc=nf90_inq_varid( ncid, 'DAYS', vid )
iernc= nf90_get_var( ncid, vid, days_limb)
iernc=nf90_inq_varid( ncid, 'LATS', vid )
iernc= nf90_get_var( ncid, vid, ugwp_taulat)
iernc=nf90_inq_varid( ncid, 'ABSMF', vid )
iernc= nf90_get_var( ncid, vid, tau_limb)
iernc=nf90_close(ncid)
endif read_and_broadcast_2

call ccpp_bcast(days_limb, mpiroot, mpicomm, errflg)
call ccpp_bcast(tau_limb, mpiroot, mpicomm, errflg)
call ccpp_bcast(ugwp_taulat, mpiroot, mpicomm, errflg)

end subroutine read_tau_amf

!>
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -107,7 +107,7 @@
standard_name = ozone_forcing
long_name = ozone forcing data
units = mixed
dimensions = (horizontal_loop_extent,vertical_dimension_of_ozone_forcing_data,number_of_coefficients_in_ozone_data)
dimensions = (horizontal_loop_extent,vertical_dimension_of_ozone_forcing_data,number_of_coefficients_in_ozone_forcing_data)
type = real
kind = kind_phys
intent = in
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -1310,7 +1310,7 @@
standard_name = ozone_forcing
long_name = ozone forcing data
units = mixed
dimensions = (horizontal_dimension,vertical_dimension_of_ozone_forcing_data,number_of_coefficients_in_ozone_data)
dimensions = (horizontal_dimension,vertical_dimension_of_ozone_forcing_data,number_of_coefficients_in_ozone_forcing_data)
type = real
kind = kind_phys
intent = inout
Expand Down
Loading