SCALE-RM
scale_cpl_phy_sfc_fixed_temp.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
29  public :: cpl_phy_sfc_fixed_temp
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35  !-----------------------------------------------------------------------------
36  !
37  !++ Private procedure
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private parameters & variables
42  !
43  logical :: initialized = .false.
44 
45  !-----------------------------------------------------------------------------
46 contains
47  !-----------------------------------------------------------------------------
50  implicit none
51  !---------------------------------------------------------------------------
52 
53  if ( initialized ) return
54 
55  log_newline
56  log_info("CPL_PHY_SFC_FIXED_TEMP_setup",*) 'Setup'
57 
58  initialized = .true.
59 
60  return
61  end subroutine cpl_phy_sfc_fixed_temp_setup
62 
63  !-----------------------------------------------------------------------------
64  subroutine cpl_phy_sfc_fixed_temp( &
65  IA, IS, IE, &
66  JA, JS, JE, &
67  TMPA, PRSA, &
68  WA, UA, VA, &
69  RHOA, QVA, LH, &
70  Z1, PBL, &
71  RHOS, PRSS, &
72  RFLXD, &
73  TMPS, WSTR, QVEF, &
74  ALBEDO, &
75  Rb, Z0M, Z0H, Z0E, &
76  calc_flag, dt, &
77  ZMFLX, XMFLX, YMFLX, &
78  SHFLX, LHFLX, QVFLX, &
79  GFLX, &
80  Ustar, Tstar, Qstar, &
81  Wstar, &
82  RLmo, &
83  U10, V10, T2, Q2 )
84  use scale_const, only: &
85  eps => const_eps, &
86  undef => const_undef, &
87  epsvap => const_epsvap, &
88  pre00 => const_pre00, &
89  rdry => const_rdry, &
90  cpdry => const_cpdry, &
91  rvap => const_rvap, &
92  stb => const_stb
93  use scale_atmos_saturation, only: &
94 ! qsat => ATMOS_SATURATION_pres2qsat_all
95  qsat => atmos_saturation_dens2qsat_all
96  use scale_bulkflux, only: &
97  bulkflux, &
98  bulkflux_diagnose_surface
99  implicit none
100 
101  integer, intent(in) :: ia, is, ie
102  integer, intent(in) :: ja, js, je
103  real(rp), intent(in) :: tmpa (ia,ja) ! temperature at the lowest atmospheric layer [K]
104  real(rp), intent(in) :: prsa (ia,ja) ! pressure at the lowest atmospheric layer [Pa]
105  real(rp), intent(in) :: wa (ia,ja) ! velocity w at the lowest atmospheric layer [m/s]
106  real(rp), intent(in) :: ua (ia,ja) ! velocity u at the lowest atmospheric layer [m/s]
107  real(rp), intent(in) :: va (ia,ja) ! velocity v at the lowest atmospheric layer [m/s]
108  real(rp), intent(in) :: rhoa (ia,ja) ! density at the lowest atmospheric layer [kg/m3]
109  real(rp), intent(in) :: qva (ia,ja) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
110  real(rp), intent(in) :: lh (ia,ja) ! latent heat at the lowest atmospheric layer [J/kg]
111  real(rp), intent(in) :: z1 (ia,ja) ! cell center height at the lowest atmospheric layer [m]
112  real(rp), intent(in) :: pbl (ia,ja) ! the top of atmospheric mixing layer [m]
113  real(rp), intent(in) :: rhos (ia,ja) ! density at the surface [kg/m3]
114  real(rp), intent(in) :: prss (ia,ja) ! pressure at the surface [Pa]
115  real(rp), intent(in) :: rflxd (ia,ja,n_rad_dir,n_rad_rgn) ! downward radiation flux at the surface (direct/diffuse,IR/near-IR/VIS) [J/m2/s]
116  real(rp), intent(in) :: tmps (ia,ja) ! surface temperature [K]
117  real(rp), intent(in) :: wstr (ia,ja) ! amount of the water storage [kg/m2]
118  real(rp), intent(in) :: qvef (ia,ja) ! efficiency of evaporation (0-1)
119  real(rp), intent(in) :: albedo (ia,ja,n_rad_dir,n_rad_rgn) ! surface albedo (direct/diffuse,IR/near-IR/VIS) (0-1)
120  real(rp), intent(in) :: rb (ia,ja) ! stomata resistance [1/s]
121  real(rp), intent(in) :: z0m (ia,ja) ! roughness length for momemtum [m]
122  real(rp), intent(in) :: z0h (ia,ja) ! roughness length for heat [m]
123  real(rp), intent(in) :: z0e (ia,ja) ! roughness length for vapor [m]
124  logical, intent(in) :: calc_flag(ia,ja) ! to decide calculate or not
125  real(dp), intent(in) :: dt ! delta time
126 
127  real(rp), intent(out) :: zmflx (ia,ja) ! z-momentum flux at the surface [kg/m/s2]
128  real(rp), intent(out) :: xmflx (ia,ja) ! x-momentum flux at the surface [kg/m/s2]
129  real(rp), intent(out) :: ymflx (ia,ja) ! y-momentum flux at the surface [kg/m/s2]
130  real(rp), intent(out) :: shflx (ia,ja) ! sensible heat flux at the surface [J/m2/s]
131  real(rp), intent(out) :: lhflx (ia,ja) ! latent heat flux at the surface [J/m2/s]
132  real(rp), intent(out) :: qvflx (ia,ja) ! water vapor flux at the surface [kg/m2/s]
133  real(rp), intent(out) :: gflx (ia,ja) ! subsurface heat flux at the surface [J/m2/s]
134  real(rp), intent(out) :: ustar (ia,ja) ! friction velocity [m/s]
135  real(rp), intent(out) :: tstar (ia,ja) ! temperature scale [K]
136  real(rp), intent(out) :: qstar (ia,ja) ! moisture scale [kg/kg]
137  real(rp), intent(out) :: wstar (ia,ja) ! convective velocity scale [m/s]
138  real(rp), intent(out) :: rlmo (ia,ja) ! inversed Obukhov length [1/m]
139  real(rp), intent(out) :: u10 (ia,ja) ! velocity u at 10m [m/s]
140  real(rp), intent(out) :: v10 (ia,ja) ! velocity v at 10m [m/s]
141  real(rp), intent(out) :: t2 (ia,ja) ! temperature at 2m [K]
142  real(rp), intent(out) :: q2 (ia,ja) ! water vapor at 2m [kg/kg]
143 
144  real(rp) :: emis ! surface longwave emission [J/m2/s]
145  real(rp) :: lwd ! surface downward longwave radiation flux [J/m2/s]
146  real(rp) :: lwu ! surface upward longwave radiation flux [J/m2/s]
147  real(rp) :: swd ! surface downward shortwave radiation flux [J/m2/s]
148  real(rp) :: swu ! surface upward shortwave radiation flux [J/m2/s]
149  real(rp) :: res ! residual
150 
151  real(rp) :: uabs ! absolute velocity [m/s]
152  real(rp) :: ra ! Aerodynamic resistance (=1/Ce) [1/s]
153 
154  real(rp) :: qvsat ! saturation water vapor mixing ratio at surface [kg/kg]
155  real(rp) :: qvs(ia,ja) ! water vapor mixing ratio at surface [kg/kg]
156  real(rp) :: rtot ! total gas constant
157  real(rp) :: qdry ! dry air mass ratio [kg/kg]
158 
159  real(rp) :: fracu10(ia,ja) ! calculation parameter for U10 [1]
160  real(rp) :: fract2 (ia,ja) ! calculation parameter for T2 [1]
161  real(rp) :: fracq2 (ia,ja) ! calculation parameter for Q2 [1]
162 
163  real(rp) :: mflux
164 
165  integer :: i, j
166  !---------------------------------------------------------------------------
167 
168  log_progress(*) 'coupler / physics / surface / FIXED-TEMP'
169 
170  ! calculate surface flux
171  !$omp parallel do &
172 #ifndef __GFORTRAN__
173  !$omp default(none) &
174  !$omp shared(IS,IE,JS,JE,EPS,UNDEF,Rdry,CPdry,bulkflux,dt, &
175  !$omp calc_flag,TMPA,QVA,QVS,LH,WA,UA,VA,Z1,PBL,PRSA,TMPS,WSTR,PRSS,RHOS,QVEF,Z0M,Z0H,Z0E,ALBEDO,RFLXD,Rb, &
176  !$omp FracU10,FracT2,FracQ2, &
177  !$omp SHFLX,LHFLX,QVFLX,GFLX,ZMFLX,XMFLX,YMFLX,Ustar,Tstar,Qstar,Wstar,RLmo,U10,V10,T2,Q2) &
178 #else
179  !$omp default(shared) &
180 #endif
181  !$omp private(qdry,Rtot,QVsat,Uabs,Ra,res,emis,LWD,LWU,SWD,SWU,MFLUX)
182  do j = js, je
183  do i = is, ie
184  if ( calc_flag(i,j) ) then
185 
186 ! qdry = 1.0_RP - QVA(i,j)
187 ! Rtot = qdry * Rdry + QVA(i,j) * Rvap
188 ! call qsat( TMPS(i,j), PRSS(i,j), qdry, QVsat )
189  call qsat( tmps(i,j), rhos(i,j), qvsat )
190 
191  qvs(i,j) = ( 1.0_rp-qvef(i,j) ) * qva(i,j) &
192  + ( qvef(i,j) ) * qvsat
193 
194  uabs = sqrt( wa(i,j)**2 + ua(i,j)**2 + va(i,j)**2 )
195 
196  call bulkflux( tmpa(i,j), tmps(i,j), & ! [IN]
197  prsa(i,j), prss(i,j), & ! [IN]
198  qva(i,j), qvs(i,j), & ! [IN]
199  uabs, z1(i,j), pbl(i,j), & ! [IN]
200  z0m(i,j), z0h(i,j), z0e(i,j), & ! [IN]
201  ustar(i,j), tstar(i,j), qstar(i,j), & ! [OUT]
202  wstar(i,j), rlmo(i,j), ra, & ! [OUT]
203  fracu10(i,j), fract2(i,j), fracq2(i,j) ) ! [OUT]
204 
205  if ( uabs < eps ) then
206  zmflx(i,j) = 0.0_rp
207  xmflx(i,j) = 0.0_rp
208  ymflx(i,j) = 0.0_rp
209  else
210  mflux = - rhos(i,j) * ustar(i,j)**2
211  zmflx(i,j) = mflux * wa(i,j) / uabs
212  xmflx(i,j) = mflux * ua(i,j) / uabs
213  ymflx(i,j) = mflux * va(i,j) / uabs
214  end if
215  shflx(i,j) = -rhos(i,j) * ustar(i,j) * tstar(i,j) * cpdry
216  qvflx(i,j) = -rhos(i,j) * ustar(i,j) * qstar(i,j) * ra / ( ra+rb(i,j) )
217  qvflx(i,j) = min( qvflx(i,j), wstr(i,j) / real(dt,rp) )
218  lhflx(i,j) = qvflx(i,j) * lh(i,j)
219 
220  emis = ( 1.0_rp-albedo(i,j,i_r_diffuse,i_r_ir) ) * stb * tmps(i,j)**4
221 
222  lwd = rflxd(i,j,i_r_diffuse,i_r_ir)
223  lwu = rflxd(i,j,i_r_diffuse,i_r_ir) * albedo(i,j,i_r_diffuse,i_r_ir) + emis
224  swd = rflxd(i,j,i_r_direct ,i_r_nir) &
225  + rflxd(i,j,i_r_diffuse,i_r_nir) &
226  + rflxd(i,j,i_r_direct ,i_r_vis) &
227  + rflxd(i,j,i_r_diffuse,i_r_vis)
228  swu = rflxd(i,j,i_r_direct ,i_r_nir) * albedo(i,j,i_r_direct ,i_r_nir) &
229  + rflxd(i,j,i_r_diffuse,i_r_nir) * albedo(i,j,i_r_diffuse,i_r_nir) &
230  + rflxd(i,j,i_r_direct ,i_r_vis) * albedo(i,j,i_r_direct ,i_r_vis) &
231  + rflxd(i,j,i_r_diffuse,i_r_vis) * albedo(i,j,i_r_diffuse,i_r_vis)
232 
233  ! calculation for residual
234  res = swd - swu + lwd - lwu - shflx(i,j) - qvflx(i,j) * lh(i,j)
235 
236  ! put residual in ground heat flux (positive for downward)
237  gflx(i,j) = res
238 
239  else ! not calculate surface flux
240  zmflx(i,j) = undef
241  xmflx(i,j) = undef
242  ymflx(i,j) = undef
243  shflx(i,j) = undef
244  lhflx(i,j) = undef
245  qvflx(i,j) = undef
246  gflx(i,j) = undef
247  ustar(i,j) = undef
248  tstar(i,j) = undef
249  qstar(i,j) = undef
250  wstar(i,j) = undef
251  rlmo(i,j) = undef
252  u10(i,j) = undef
253  v10(i,j) = undef
254  t2(i,j) = undef
255  q2(i,j) = undef
256  endif
257  enddo
258  enddo
259 
260  call bulkflux_diagnose_surface( ia, is, ie, ja, js, je, &
261  ua(:,:), va(:,:), & ! (in)
262  tmpa(:,:), qva(:,:), & ! (in)
263  tmps(:,:), qvs(:,:), & ! (in)
264  z1(:,:), z0m(:,:), z0h(:,:), z0e(:,:), & ! (in)
265  u10(:,:), v10(:,:), t2(:,:), q2(:,:), & ! (out)
266  mask = calc_flag(:,:), & ! (in)
267  fracu10 = fracu10(:,:), & ! (in)
268  fract2 = fract2(:,:), & ! (in)
269  fracq2 = fracq2(:,:) ) ! (in)
270 
271  return
272  end subroutine cpl_phy_sfc_fixed_temp
273 
scale_cpl_sfc_index::n_rad_dir
integer, parameter, public n_rad_dir
Definition: scale_cpl_sfc_index.F90:36
scale_cpl_sfc_index::i_r_direct
integer, parameter, public i_r_direct
Definition: scale_cpl_sfc_index.F90:37
scale_const::const_epsvap
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:69
scale_cpl_sfc_index::i_r_diffuse
integer, parameter, public i_r_diffuse
Definition: scale_cpl_sfc_index.F90:38
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_const::const_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:63
scale_cpl_sfc_index::i_r_ir
integer, parameter, public i_r_ir
Definition: scale_cpl_sfc_index.F90:29
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_bulkflux
module Surface bulk flux
Definition: scale_bulkflux.F90:12
scale_cpl_phy_sfc_fixed_temp::cpl_phy_sfc_fixed_temp_setup
subroutine, public cpl_phy_sfc_fixed_temp_setup
Setup.
Definition: scale_cpl_phy_sfc_fixed_temp.F90:50
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_cpl_sfc_index::i_r_nir
integer, parameter, public i_r_nir
Definition: scale_cpl_sfc_index.F90:30
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_const::const_stb
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:49
scale_cpl_sfc_index
module coupler / surface-atmospehre
Definition: scale_cpl_sfc_index.F90:11
scale_cpl_sfc_index::i_r_vis
integer, parameter, public i_r_vis
Definition: scale_cpl_sfc_index.F90:31
scale_cpl_phy_sfc_fixed_temp::cpl_phy_sfc_fixed_temp
subroutine, public cpl_phy_sfc_fixed_temp(IA, IS, IE, JA, JS, JE, TMPA, PRSA, WA, UA, VA, RHOA, QVA, LH, Z1, PBL, RHOS, PRSS, RFLXD, TMPS, WSTR, QVEF, ALBEDO, Rb, Z0M, Z0H, Z0E, calc_flag, dt, ZMFLX, XMFLX, YMFLX, SHFLX, LHFLX, QVFLX, GFLX, Ustar, Tstar, Qstar, Wstar, RLmo, U10, V10, T2, Q2)
Definition: scale_cpl_phy_sfc_fixed_temp.F90:84
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
scale_cpl_sfc_index::n_rad_rgn
integer, parameter, public n_rad_rgn
Definition: scale_cpl_sfc_index.F90:28
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
scale_cpl_phy_sfc_fixed_temp
module coupler / surface fixed temp model
Definition: scale_cpl_phy_sfc_fixed_temp.F90:12
scale_bulkflux::bulkflux
procedure(bc), pointer, public bulkflux
Definition: scale_bulkflux.F90:78