Sca_nangi : max # of angles Sca_npf : # of phase functions Sca_nskip : number of first data record lines to be skipped. Sca_anci : the number of ancillary data (definition differs between binary and text formats)
real(R_), save, allocatable :: ang(:,:) !(Sca_nangi, Sca_npf)
real(R_), save, allocatable :: phs(:,:) !(Sca_nangi, Sca_npf)
integer, save, allocatable :: na(:) !(Sca_npf)
DSET ^les05_0045.sca *OPTIONS LITTLE_ENDIAN *OPTIONS BIG_ENDIAN OPTIONS template TITLE les05_0045.sca UNDEF 1.0e+37 XDEF 203 LINEAR 1 1 YDEF 1 LINEAR 1 1 ZDEF 60 LINEAR 1 1 TDEF 1 LINEAR 00:00Z00JAN0001 1YR VARS 2 ang 1 99 Angles (deg.) phs 60 99 Phase function ENDVARS
$ ../bin/bin_dump les05_0045.sca.ctl | less
!####
subroutine mcarUtl__binSca_read(iu, nangi, nanci, nphs, ang, phs, nskip, mbswap)
!
!* Read in phase functions from a binary file *
integer, intent(in) :: iu ! file unit index
integer, intent(in) :: nangi ! # of angles
integer, intent(in) :: nanci ! # of ancillary data
integer, intent(in) :: nphs ! # of phase functions
real(R_), intent(out) :: ang(:) ! angles (degrees)
real(R_), intent(out) :: phs(:,:) ! phase functions
integer, intent(in), optional :: nskip ! # of record lines to be skipped
integer, intent(in), optional :: mbswap ! flag for byte swapping (0=no, 1=yes)
real(R4_) :: wrk(nangi + nanci) !(AUTO)
integer :: iphs, mswap, irec
! Initialize
irec = 1
mswap = 0
if (present(nskip)) irec = nskip + 1
if (present(mbswap) .and. mbswap == 1) mswap = 4
! Read in
call bin_read_o1R4(iu, wrk, irec=irec, mswap=mswap)
ang(1:nangi) = wrk(1:nangi)
do iphs = 1, nphs
irec = irec + 1
call bin_read_o1R4(iu, wrk, irec=irec, mswap=mswap)
phs(1:nangi, iphs) = wrk(1:nangi)
end do
end subroutine mcarUtl__binSca_read
nangi --> Sca_nangi nphs --> Sca_npf mbswap --> Wld_mbswap
%mdlphs
!####
subroutine mcarUtl__txtSca_read(iu, nanci, na, ang, phs)
!
!* Get phase functions from a text-format file *
integer, intent(in) :: iu ! file unit index
integer, intent(in) :: nanci ! # of ancillary data blocks
integer, intent(out) :: na(:) ! # of angles
real(R_), intent(out) :: ang(:,:) ! angles (degrees)
real(R_), intent(out) :: phs(:,:) ! phase functions
!// Note: Number of ancillary data lines should be nanci*(np + 1).
integer :: ia, ipf, nangi, npf, np, ip, naa, ianci, ierr
! Header
nangi = size(ang, 1)
npf = size(ang, 2)
! Data
ipf = 1
loop_trial: do
call skipToSign(iu, '%mdlphs', ierr)
call err_read(ierr, iu, 'mcarUtl__txtSca_read: %mdlphs is not found.')
read (iu, *, err=1, end=1)
read (iu, *, err=1, end=1) naa
call check_iI('mcarUtl__txtSca_read: naa', naa, 1, nangi)
read (iu, *, err=1, end=1) np
do ianci = 1, nanci
do ip = 1, np + 1
read (iu, *, err=1, end=1)
end do
end do
do ip = 1, np
read (iu, *, err=1, end=1)
na(ipf) = naa
do ia = 1, na(ipf)
read (iu, *, err=1, end=1) ang(ia, ipf), phs(ia, ipf)
end do
ipf = ipf + 1
if (ipf > npf) exit loop_trial
end do
end do loop_trial
return
1 call err_read(1, iu, 'mcarUtl__txtSca_read: ang & phs at ipf='//num2str_AN(ipf))
end subroutine mcarUtl__txtSca_read
Atm_nx, Atm_ny : # of volume elements along x and y axes Atm_nz : # of all layers Atm_nz3 : # of horizontally-inhomogeneous layers Atm_nkd : # of k-distribution terms Atm_np3d : # of scattering particle polydispersions in horizontally-inhomogeneous layers
[real 3-D array, {((dat(ix, iy, iz3), ix=1, Atm_nx), iy=1, Atm_ny), iz3=1, Atm_nz3}] when Atm_mtprof=0
[real 3-D array, {((dat(ix, iy, iz3), ix=1, Atm_nx), iy=1, Atm_ny), iz3=1, Atm_nz3+1}] when Atm_mtprof=1
air_temp(ix, iy, iz) = Atm_tmp1d(iz) + Atm_tmpa3d(ix, iy, iz3), where iz = Atm_iz3l + iz3 - 1.
[real 4-D array, {(((dat(ix, iy, iz3, ikd), ix=1, Atm_nx), iy=1, Atm_ny), iz3=1, Atm_nz3), ikd=1, Atm_nkd}]
abs_coef(ix, iy, iz, ikd) = Atm_abs1d(iz, ikd) + Atm_absg3d(ix, iy, iz3, ikd)
[real 4-D array, {(((dat(ix, iy, iz3, ip3d), ix=1, Atm_nx), iy=1, Atm_ny), iz3=1, Atm_nz3), ip3d=1, Atm_np3d}]
[integer 4-D array, {(((dat(ix, iy, iz3, ip3d), ix=1, Atm_nx), iy=1, Atm_ny), iz3=1, Atm_nz3), ip3d=1, Atm_np3d}]
Atm_apf3d <= -2 : isotropic scattering phase function
-2 < Atm_apf3d <= -1 : Rayleigh scattering phase function
-1 < Atm_apf3d < 1 : Henyey-Greenstein phase function with asymmetry parameter of Atm_apf1d
1 <= Atm_apf3d : a tabulated phase function (given by the database file, see the namelist &link_anchor(mcarSca_nml_init, page=User's Guide v0.10/4b){mcarSca_nml_init})
data set 1 data set 2 ...
!####
subroutine mcarAtm__read(mbswap) !* Import data from file *
!
!// Note this procedure touches only layers iz = [iz3l:iz3u].
integer, intent(in) :: mbswap
integer :: ip, ikd, iz, iz3
real(R_), allocatable :: omgp3d(:,:,:,:)
! Initialize
allocate (omgp3d(Atm_nx, Atm_ny, Atm_nz3, Atm_np3d))
! Open & read in
if (Atm_mfmt == 0) then ! text
call mcarUtl__txtAtm_read(Atm_iud, Atm_nx, Atm_ny, Atm_nz3, Atm_nkd, Atm_np3d, &
& Atm_mtprof, Atm_tmpa3d(:,:,Atm_iz3l:), Atm_abst3d(:,:,Atm_iz3l:,:), &
& Atm_apfp3d(:,:,Atm_iz3l:,:), Atm_scap3d(:,:,Atm_iz3l:,:), omgp3d)
else ! binary
call mcarUtl__binAtm_read(Atm_iud, Atm_nx, Atm_ny, Atm_nz3, Atm_nkd, Atm_np3d, &
& Atm_mtprof, Atm_tmpa3d(:,:,Atm_iz3l:), Atm_abst3d(:,:,Atm_iz3l:,:), &
& Atm_apfp3d(:,:,Atm_iz3l:,:), Atm_scap3d(:,:,Atm_iz3l:,:), omgp3d, mbswap, Atm_idread)
end if
!// Atm_scap3d() represents tentatively the extinction coefficients, here.
!// Atm_abst3d() represents tentatively the gaseous absorption coefficients, here.
! Scale gaseous absorption coefficients
if (Atm_fabs3d /= 1.0_R_) Atm_abst3d(:,:, Atm_iz3l:Atm_iz3u, :) = &
& Atm_abst3d(:,:, Atm_iz3l:Atm_iz3u, :) * Atm_fabs3d
! Bext & Omega --> Bsca & Babs
do ip = 1, Atm_np3d
do iz = Atm_iz3l, Atm_iz3u
iz3 = iz - Atm_iz3l + 1
do ikd = 1, Atm_nkd
Atm_abst3d(:, :, iz, ikd) = Atm_abst3d(:, :, iz, ikd) + Atm_scap3d(:, :, iz, ip) &
& * (1.0_R_ - omgp3d(:, :, iz3, ip)) * Atm_fext3d(ip)
end do
Atm_scap3d(:, :, iz, ip) = Atm_scap3d(:, :, iz, ip) * omgp3d(:, :, iz3, ip) * Atm_fext3d(ip)
!// Now, Atm_abst3d() represents total absorption coefficients.
!// Now, Atm_scap3d() represents scattering coefficients.
end do
end do
deallocate (omgp3d)
end subroutine mcarAtm__read
| data array | # of streams |
| Atm_tmpa3d(:,:,:) | Atm_nz3 or Atm_nz3+1 |
| Atm_abst3d(:,:,:,:) | Atm_nz3 * Atm_nkd |
| Atm_extp3d(:,:,:,1) | Atm_nz3 |
| Atm_omgp3d(:,:,:,1) | Atm_nz3 |
| Atm_apfp3d(:,:,:,1) | Atm_nz3 |
| Atm_extp3d(:,:,:,2) | Atm_nz3 |
| Atm_omgp3d(:,:,:,2) | Atm_nz3 |
| Atm_apfp3d(:,:,:,2) | Atm_nz3 |
| ... | Atm_nz3 |
DSET ^les05_0045.atm *OPTIONS LITTLE_ENDIAN *OPTIONS BIG_ENDIAN OPTIONS template TITLE les05_0045.atm UNDEF 1.0e+37 XDEF 60 LINEAR 1 1 YDEF 60 LINEAR 1 1 ZDEF 22 LINEAR 1 1 TDEF 1 LINEAR 00:00Z00JAN0001 1YR VARS 5 tmpa3d 22 99 Temperature perturbation (K) abst3d 21 99 Absorption coefficient perturbation (/m) extp3d 21 99 Extinction coefficient (/m) omgp3d 21 99 Single scattering albedo apfp3d 21 99 Phase function specification parameter ENDVARS
$ ../bin/bin_dump les05_0045.atm.ctl | less
!####
subroutine mcarUtl__binAtm_read(iu, nx, ny, nz3, nkd, np3d, mtprof, &
& tmpa3d, absg3d, apfp3d, extp3d, omgp3d, mbswap, idloc)
!
!* Read in 3-D atmosphere model parameters from a binary file *
integer, intent(in) :: iu
integer, intent(in) :: mtprof
integer, intent(in) :: nkd
integer, intent(in) :: np3d
integer, intent(in) :: nx, ny
integer, intent(in) :: nz3
real(R_), intent(out) :: tmpa3d(:, :, :), absg3d(:, :, :, :)
real(R_), intent(out) :: extp3d(:, :, :, :), omgp3d(:, :, :, :)
real(R_), intent(out) :: apfp3d(:, :, :, :)
integer, intent(in), optional :: mbswap ! flag for byte swapping (0=no, 1=yes)
integer, intent(in), optional :: idloc ! dataset location index
real(R4_) :: wrk(nx, ny) !(AUTO)
integer :: ikd, ip3d, iz, iz3max, mswap, irec
! Initialize
if (nz3 <= 0) return
if (mtprof == 0) then
iz3max = nz3
else
iz3max = nz3 + 1
end if
mswap = 0
if (present(mbswap) .and. mbswap == 1) mswap = 4
if (present(idloc)) then
irec = (iz3max + nz3 * (nkd + 3*np3d)) * (idloc - 1) ! = (record index to be read) - 1
call fileRec_set(iu, irec)
end if
! Read in
do iz = 1, iz3max ! temperature
call bin_read_o2R4(iu, wrk, mswap=mswap)
tmpa3d(1:nx, 1:ny, iz) = wrk(1:nx, 1:ny)
end do
do ikd = 1, nkd ! gaseous absorption coefficient
do iz = 1, nz3
call bin_read_o2R4(iu, wrk, mswap=mswap)
absg3d(1:nx, 1:ny, iz, ikd) = wrk(1:nx, 1:ny)
end do
end do
do ip3d = 1, np3d
do iz = 1, nz3 ! extinction coefficients
call bin_read_o2R4(iu, wrk, mswap=mswap)
extp3d(1:nx, 1:ny, iz, ip3d) = wrk(1:nx, 1:ny)
end do
do iz = 1, nz3 ! single scattering albedo
call bin_read_o2R4(iu, wrk, mswap=mswap)
omgp3d(1:nx, 1:ny, iz, ip3d) = wrk(1:nx, 1:ny)
end do
do iz = 1, nz3 ! phase function specification parameter
call bin_read_o2R4(iu, wrk, mswap=mswap)
apfp3d(1:nx, 1:ny, iz, ip3d) = wrk(1:nx, 1:ny)
end do
end do
end subroutine mcarUtl__binAtm_read
%mdla3d
!####
subroutine mcarUtl__txtAtm_read(iu, nx, ny, nz3, nkd, np3d, mtprof, &
& tmpa3d, absg3d, apfp3d, extp3d, omgp3d)
!
!* Read in 3-D atmosphere model parameters from a text file *
integer, intent(in) :: iu
integer, intent(in) :: mtprof
integer, intent(in) :: nkd
integer, intent(in) :: np3d
integer, intent(in) :: nx, ny
integer, intent(in) :: nz3
real(R_), intent(out) :: tmpa3d(:, :, 1:), absg3d(:, :, :, :)
real(R_), intent(out) :: extp3d(:, :, :, :), omgp3d(:, :, :, :)
real(R_), intent(out) :: apfp3d(:, :, :, :)
integer :: ierr, ikd, ip3d, ix, iy, iz, iz3max
! Initialize
if (nz3 <= 0) return
if (mtprof == 0) then
iz3max = nz3
else
iz3max = nz3 + 1
end if
iz = 0
call skipToSign(iu, '%mdla3d', ierr)
call err_read(ierr, iu, 'mcarUtl__txtAtm_read: %mdla3d is not found.')
! Temperature
read (iu, *, err=1, end=9)
do iz = 1, iz3max
read (iu, *, err=1) ((tmpa3d(ix, iy, iz), ix = 1, nx), iy = 1, ny)
end do
! Gaseous absorption coefficient
do ikd = 1, nkd
read (iu, *, err=2, end=9)
do iz = 1, nz3
read (iu, *, err=2) ((absg3d(ix, iy, iz, ikd), ix = 1, nx), iy = 1, ny)
end do
end do
! Cloud & aerosol parameters
do ip3d = 1, np3d
read (iu, *, err=3, end=9)
do iz = 1, nz3 ! extinction coefficients
read (iu, *, err=3) ((extp3d(ix, iy, iz, ip3d), ix = 1, nx), iy = 1, ny)
end do
read (iu, *, err=4, end=9)
do iz = 1, nz3 ! single scattering albedo
read (iu, *, err=4) ((omgp3d(ix, iy, iz, ip3d), ix = 1, nx), iy = 1, ny)
end do
read (iu, *, err=5, end=9)
do iz = 1, nz3 ! phase function specification parameter
read (iu, *, err=5) ((apfp3d(ix, iy, iz, ip3d), ix = 1, nx), iy = 1, ny)
end do
end do
! Exit
return
1 call err_read(1, iu, 'mcarUtl__txtAtm_read: tmpa3d at iz ='//num2str_AN(iz))
2 call err_read(1, iu, 'mcarUtl__txtAtm_read: absg3d at iz ='//num2str_AN(iz))
3 call err_read(1, iu, 'mcarUtl__txtAtm_read: extp3d at iz ='//num2str_AN(iz))
4 call err_read(1, iu, 'mcarUtl__txtAtm_read: omgp3d at iz ='//num2str_AN(iz))
5 call err_read(1, iu, 'mcarUtl__txtAtm_read: apfp3d at iz ='//num2str_AN(iz))
9 call err_read(-1, iu, 'mcarUtl__txtAtm_read: Unexpected EOF. Array sizes would be invalid.')
end subroutine mcarUtl__txtAtm_read
Sfc_nxb, Sfc_nyb : # of surface elements along x and y axes Sfc_npsfc : max # of surface parameters in a single surface element
real(R_), save, allocatable :: Sfc_tmps2d(:,:) !(nxb,nyb)
temp(:,:) = Sfc_tmp + Sfc_tmps2d(:,:)
integer, save, allocatable :: Sfc_jsfc2d(:,:) !(nxb,nyb)
0: black 1: Lambertian 2: Flx_diffuse-specular mixture (DSM) model 3: Rahman-Pinty-Verstraete (RPV) model 4: Li-Sparse-Ross-Thick (LSRT) model
real(R_), save, allocatable :: Sfc_psfc2d(:,:,:) !(nxb,nyb,Sfc_NPAR)
!####
subroutine mcarSfc__read(mbswap) !* Import data from file *
integer, intent(in) :: mbswap
integer :: nmax, i
! Initialize
nmax = 0
do i = 1, 4
if (Sfc_mbrdf(i) == 1) nmax = max(nmax, Sfc_npsfc(i))
end do
! Read in
Sfc_idloc = Sfc_idloc + 1
if (Sfc_mfmt == 0) then ! text
call mcarUtl__txtSfc_read(Sfc_iud, Sfc_nxb, Sfc_nyb, nmax, Sfc_tmps2d, Sfc_psfc2d, Sfc_jsfc2d)
else
if (Sfc_idloc /= Sfc_idread) Sfc_idloc = Sfc_idread
call mcarUtl__binSfc_read(Sfc_iud, Sfc_nxb, Sfc_nyb, nmax, Sfc_tmps2d, Sfc_psfc2d, &
& Sfc_jsfc2d, mbswap, Sfc_idloc)
end if
end subroutine mcarSfc__read
DSET ^missing_example_file *OPTIONS LITTLE_ENDIAN *OPTIONS BIG_ENDIAN OPTIONS template TITLE An example UNDEF 1.0e+37 XDEF 256 LINEAR 1 1 YDEF 256 LINEAR 1 1 ZDEF 5 LINEAR 1 1 * This example shows a case with Sfc_npsfc = 5, which may depends on Sfc_mbrdf TDEF 1 LINEAR 00:00Z00JAN0001 1YR VARS 3 tmps2d 1 99 Temperature (K) jsfc2d 1 99 BRDF type index psfc2d 5 99 BRDF parameters (when Sfc_npsfc = 5) ENDVARS
!####
subroutine mcarUtl__binSfc_read(iu, nxb, nyb, npsfc, tmps2d, psfc2d, jsfc2d, mbswap, idloc)
!
!* Read in surface model parameters from a binary file *
integer, intent(in) :: iu
integer, intent(in) :: nxb
integer, intent(in) :: nyb
integer, intent(in) :: npsfc
integer, intent(out) :: jsfc2d(:, :)
real(R_), intent(out) :: tmps2d(:, :), psfc2d(:, :, :)
integer, intent(in), optional :: mbswap ! flag for byte swapping (0=no, 1=yes)
integer, intent(in), optional :: idloc ! dataset location index
real(R4_) :: wrk(nxb, nyb) !(AUTO)
integer :: ipsfc, mswap, irec
! Initialize
mswap = 0
if (present(mbswap) .and. mbswap == 1) mswap = 4
if (present(idloc)) then
irec = (2 + npsfc) * (idloc - 1) ! = (record index to be read) - 1
call fileRec_set(iu, irec)
end if
! Read in
call bin_read_o2R4(iu, wrk, mswap=mswap)
tmps2d(1:nxb, 1:nyb) = wrk(1:nxb, 1:nyb)
call bin_read_o2R4(iu, wrk, mswap=mswap)
jsfc2d(1:nxb, 1:nyb) = nint(wrk(1:nxb, 1:nyb))
do ipsfc = 1, npsfc
call bin_read_o2R4(iu, wrk, mswap=mswap)
psfc2d(1:nxb, 1:nyb, ipsfc) = wrk(1:nxb, 1:nyb)
end do
end subroutine mcarUtl__binSfc_read
%mdlsfc
!####
subroutine mcarUtl__txtSfc_read(iu, nxb, nyb, npsfc, tmps2d, psfc2d, jsfc2d)
!
!* Read in surface model parameters from a text file *
integer, intent(in) :: iu
integer, intent(in) :: nxb
integer, intent(in) :: nyb
integer, intent(in) :: npsfc
integer, intent(out) :: jsfc2d(:, :)
real(R_), intent(out) :: tmps2d(:, :), psfc2d(:, :, :)
integer :: ierr, ipsfc, ixb, iyb
! Signature
call skipToSign(iu, '%mdlsfc', ierr)
call err_read(ierr, iu, 'mcarUtl__txtSfc_read: %mdlsfc is not found.')
! Surface model
ipsfc = 0
read (iu, *, err=1, end=2)
read (iu, *, err=1, end=2) ((tmps2d(ixb, iyb), ixb = 1, nxb), iyb = 1, nyb)
read (iu, *, err=1, end=2)
read (iu, *, err=1, end=2) ((jsfc2d(ixb, iyb), ixb = 1, nxb), iyb = 1, nyb)
read (iu, *, err=1, end=2)
do ipsfc = 1, npsfc
read (iu, *, err=1, end=2) ((psfc2d(ixb, iyb, ipsfc), ixb = 1, nxb), iyb = 1, nyb)
end do
return
1 call err_read( 1, iu, 'mcarUtl__txtSfc_read: ipsfc='//num2str_AN(ipsfc))
2 call err_read(-1, iu, 'mcarUtl__txtSfc_read: ipsfc='//num2str_AN(ipsfc))
end subroutine mcarUtl__txtSfc_read
2010-08-04
2010-07-26
2010-07-25
2009-08-22
2009-06-15
2009-05-10