Skip to content

Commit

Permalink
additional comments and clean-up for new mwRTM (GEOS_LandAssimGridCom…
Browse files Browse the repository at this point in the history
…p.F90; mwRTM_routines.F90)
  • Loading branch information
gmao-rreichle committed May 11, 2023
1 parent 3f3bb99 commit 1333041
Show file tree
Hide file tree
Showing 2 changed files with 32 additions and 34 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -2237,7 +2237,9 @@ end subroutine UPDATE_ASSIM
! ******************************************************************************

! subroutine to calculate Tb for HISTORY output

!
! IMPORTANT: hardwired mwRTM configuration for SMAP L-band Tb w/o Pellarin atm correction (RTM_ID=4)

subroutine CALC_LAND_TB(gc, import, export, clock, rc)

type(ESMF_GridComp), intent(inout) :: gc ! Gridded component
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -408,7 +408,7 @@ subroutine mwRTM_get_Tb( N_tile, freq, inc_angle, mwp, elev, &

case(3,4)

CALL mironov(freq,soilmoist(n),mwp(n)%clay,c_er)
CALL MIRONOV( freq, soilmoist(n), mwp(n)%clay, c_er )

case default

Expand Down Expand Up @@ -924,11 +924,11 @@ END SUBROUTINE DIEL_WAT

! **********************************************************************

SUBROUTINE mironov(freq,mv,clayfrac,er_r)
SUBROUTINE MIRONOV( freq, mv, clayfrac, er_r )

! Soil dielectric mixing model by Mironov et al IEEE TGRS 2009, doi:10.1109/TGRS.2008.2011631
!
! 8 May 2023: Implementation taken from SMAP L2_SM_P retrieval model; clean-up by reichle
! 8 May 2023: Implementation taken from SMAP L2_SM_P retrieval model by qliu; clean-up by reichle

IMPLICIT NONE

Expand All @@ -952,13 +952,12 @@ SUBROUTINE mironov(freq,mv,clayfrac,er_r)
REAL :: tmptaub, tmptaub2plus1, tmpdiffepsb
REAL :: tmptauu, tmptauu2plus1, tmpdiffepsu

REAL :: tmpreal, tmpepsbreal2, tmpepsbimag2
REAL :: tmpepsureal2, tmpepsuimag2
REAL :: tmpreal, tmprealb, tmprealu

! --------------------------------------------------------------------------------------

f = freq ! Section IV
C = clayfrac*100 ! Section VI
f = freq ! Section IV
C = clayfrac*100 ! Section VI

!! Mironov's regression expressions based on Curtis, Dobson, and Hallikainen datasets
!!
Expand All @@ -975,17 +974,17 @@ SUBROUTINE mironov(freq,mv,clayfrac,er_r)
!! b: bound soil water (BSW)
!! u: unbound (free) soil water (FSW)

nd = 1.634 - 0.539e-2 * C + 0.2748e-4 * C**2 ! Eqn 17
kd = 0.03952 - 0.04038e-2 * C ! Eqn 18
mvt = 0.02863 + 0.30673e-2 * C ! Eqn 19
eps0b = 79.8 - 85.4e-2 * C + 32.7e-4 * C**2 ! Eqn 20
taub = 1.062e-11 + 3.450e-12 * 1e-2 * C ! Eqn 21
sigb = 0.3112 + 0.467e-2 * C ! Eqn 22
sigu = 0.3631 + 1.217e-2 * C ! Eqn 23
eps0u = 100. ! Eqn 24
tauu = 8.5e-12 ! Eqn 25
nd = 1.634 - 0.539e-2 * C + 0.2748e-4 * C**2 ! Eqn 17
kd = 0.03952 - 0.04038e-2 * C ! Eqn 18
mvt = 0.02863 + 0.30673e-2 * C ! Eqn 19
eps0b = 79.8 - 85.4e-2 * C + 32.7e-4 * C**2 ! Eqn 20
taub = 1.062e-11 + 3.450e-12 * 1e-2 * C ! Eqn 21
sigb = 0.3112 + 0.467e-2 * C ! Eqn 22
sigu = 0.3631 + 1.217e-2 * C ! Eqn 23
eps0u = 100. ! Eqn 24
tauu = 8.5e-12 ! Eqn 25

!! Debye relaxation equations for water as a function of frequency ! Eqn 16
!! Debye relaxation equations for water as a function of frequency ! Eqn 16

tmp2PIf = 2.*MAPL_PI*f

Expand All @@ -1007,35 +1006,32 @@ SUBROUTINE mironov(freq,mv,clayfrac,er_r)

!! Refractive indices and normalized attenuation coefficients

tmpreal = 1/sqrt(2.0)
tmpreal = 1/sqrt(2.0)

tmpepsbreal2 = epsb_real**2
tmpepsbimag2 = epsb_imag**2
tmprealb = sqrt( epsb_real**2 + epsb_imag**2 )
tmprealu = sqrt( epsu_real**2 + epsu_imag**2 )

tmpepsureal2 = epsu_real**2
tmpepsuimag2 = epsu_imag**2

nb = tmpreal * sqrt( sqrt(tmpepsbreal2 + tmpepsbimag2) + epsb_real ) ! Eqn 14
kb = tmpreal * sqrt( sqrt(tmpepsbreal2 + tmpepsbimag2) - epsb_real ) ! Eqn 15
nu = tmpreal * sqrt( sqrt(tmpepsureal2 + tmpepsuimag2) + epsu_real ) ! Eqn 14
ku = tmpreal * sqrt( sqrt(tmpepsureal2 + tmpepsuimag2) - epsu_real ) ! Eqn 15
nb = tmpreal * sqrt( tmprealb + epsb_real ) ! Eqn 14
kb = tmpreal * sqrt( tmprealb - epsb_real ) ! Eqn 15
nu = tmpreal * sqrt( tmprealu + epsu_real ) ! Eqn 14
ku = tmpreal * sqrt( tmprealu - epsu_real ) ! Eqn 15

IF (mv <= mvt) THEN

nm = nd + (nb - 1.) * mv ! Eqn 12
km = kd + kb * mv ! Eqn 13
nm = nd + (nb - 1.) * mv ! Eqn 12
km = kd + kb * mv ! Eqn 13

ELSE

nm = nd + (nb - 1.) * mvt + (nu - 1.) * (mv - mvt) ! Eqn 12
km = kd + kb * mvt + ku * (mv - mvt) ! Eqn 13
nm = nd + (nb - 1.) * mvt + (nu - 1.) * (mv - mvt) ! Eqn 12
km = kd + kb * mvt + ku * (mv - mvt) ! Eqn 13

ENDIF

! complex dielectric constant of moist soil

er_r_real = nm**2 - km**2 ! Eqn 11
er_r_imag = 2. * nm * km ! Eqn 11
er_r_real = nm**2 - km**2 ! Eqn 11
er_r_imag = 2. * nm * km ! Eqn 11

er_r = CMPLX(er_r_real,er_r_imag)

Expand Down

1 comment on commit 1333041

@gmao-qliu
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The modification is fine. I tested against earlier version and got zero-diff result.

Please sign in to comment.