23 #define _GET_VAR1 get_var1_real 25 #define _GET_VAR1 get_var1_double 105 use diag_manager_mod,
only: register_diag_field, send_data
106 use constants_mod,
only: pi=>pi_8, grav, rdgas, cp_air, kappa, cnst_radius =>
radius 107 use fms_mod,
only: write_version_number, open_namelist_file, &
108 check_nml_error, file_exist, close_file
110 use mpp_mod,
only: mpp_error, fatal, stdlog, get_unit, mpp_pe
111 use mpp_domains_mod,
only: mpp_update_domains, domain2d
112 use time_manager_mod,
only: time_type, get_time, get_date
119 use fv_mp_mod,
only: mp_reduce_sum, mp_reduce_min, mp_reduce_max, is_master
132 real(kind=R_GRID),
parameter ::
radius = cnst_radius
136 #include<file_version.h> 166 real,
allocatable::
s2c(:,:,:)
170 real(KIND=4),
allocatable,
dimension(:,:,:)::
gz3 171 real,
allocatable::
gz0(:,:)
292 k_breed,
k_trop,
p_trop,
dps_min,
kord_data,
tc_mask,
nudge_debug,
nf_ps,
nf_t,
nf_ht, &
303 subroutine fv_nwp_nudge ( Time, dt, npx, npy, npz, ps_dt, u_dt, v_dt, t_dt, q_dt, zvir, ptop, &
304 ak, bk, ts, ps, delp, ua, va, pt, nwat, q, phis, gridstruct, &
307 type(time_type),
intent(in):: time
308 integer,
intent(IN):: npx, npy
309 integer,
intent(in):: npz
310 integer,
intent(in):: nwat
311 real,
intent(in):: dt
312 real,
intent(in):: zvir, ptop
313 type(domain2d),
intent(INOUT),
target :: domain
315 real,
intent(in ),
dimension(npz+1):: ak, bk
316 real,
intent(in ),
dimension(isd:ied,jsd:jed ):: phis
317 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: pt, ua, va, delp
324 real,
intent(inout),
dimension(isd:ied,jsd:jed):: ps
326 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: u_dt, v_dt
327 real,
intent(out),
dimension(is:ie,js:je,npz):: t_dt, q_dt
328 real,
intent(out),
dimension(is:ie,js:je):: ps_dt, ts
338 real:: peln(
is:
ie,npz+1)
339 real:: pe2(
is:
ie, npz+1)
341 real,
allocatable :: ps_obs(:,:)
342 real,
allocatable,
dimension(:,:,:):: u_obs, v_obs, t_obs, q_obs
343 real,
allocatable,
dimension(:,:,:):: du_obs, dv_obs
345 integer :: seconds, days
346 integer :: i,j,k, iq, kht, n
347 real :: factor, rms, bias, co
348 real :: rdt, press(npz), profile(npz), prof_t(npz), prof_q(npz), du, dv
352 real(kind=R_GRID),
pointer,
dimension(:,:,:) :: agrid
353 real,
pointer,
dimension(:,:) :: rarea, area
355 real,
pointer,
dimension(:,:) :: sina_u, sina_v
356 real,
pointer,
dimension(:,:,:) :: sin_sg
357 real(kind=R_GRID),
pointer,
dimension(:,:,:) :: vlon, vlat
359 real,
pointer,
dimension(:,:) :: dx, dy, rdxc, rdyc
361 real(kind=R_GRID),
pointer :: da_min
363 logical,
pointer :: nested, sw_corner, se_corner, nw_corner, ne_corner
366 call mpp_error(fatal,
'==> Error from fv_nwp_nudge: module not initialized')
368 agrid => gridstruct%agrid_64
369 area => gridstruct%area
370 rarea => gridstruct%rarea
372 vlon => gridstruct%vlon
373 vlat => gridstruct%vlat
374 sina_u => gridstruct%sina_u
375 sina_v => gridstruct%sina_v
376 sin_sg => gridstruct%sin_sg
380 rdxc => gridstruct%rdxc
381 rdyc => gridstruct%rdyc
383 da_min => gridstruct%da_min
385 nested => gridstruct%nested
386 sw_corner => gridstruct%sw_corner
387 se_corner => gridstruct%se_corner
388 nw_corner => gridstruct%nw_corner
389 ne_corner => gridstruct%ne_corner
393 forecast_mode = .true.
398 call get_time (time, seconds, days)
414 press(k) = 0.5*(ak(k) + ak(k+1)) + 0.5*(bk(k)+bk(k+1))*1.e5
416 profile(k) = max(0.01, press(k)/
p_relax)
420 if( press(k) <
p_norelax ) profile(k) = 0.0
428 if ( press(k) < 10.e2 )
then 429 prof_t(k) = max(0.01, press(k)/10.e2)
438 if ( press(k) < 300.e2 )
then 439 prof_q(k) = max(0., press(k)/300.e2)
448 ptmp = ak(k+1) + bk(k+1)*1.e5
464 factor = max(1.e-5, factor)
481 call get_obs(time, dt, zvir, ak, bk, ps, ts, ps_obs, delp, pt, nwat, q, u_obs, v_obs, t_obs, q_obs, &
482 phis, ua, va, u_dt, v_dt, npx, npy, npz, factor, mask, ptop, bd, gridstruct, domain)
494 forecast_mode = .true.
501 if ( abs(ps(i,j)-ps_obs(i,j)) > 2.e2 )
then 502 ps_fac(i,j) = 2.e2 / abs(ps(i,j)-ps_obs(i,j))
518 tv(i,j) = pt(i,j,npz)*(1.+zvir*q(i,j,npz,1))
526 if(
master)
write(*,*)
'kht=', kht
528 ptmp = 0.01*
g0_sum(ps_obs,
is,
ie,
js,
je, 1, .false.,
isd,
ied,
jsd,
jed, area)
529 if(
master)
write(*,*)
'Mean PS_o=', ptmp
537 if ( abs(phis(i,j)-
gz0(i,j))/grav > 2. )
then 538 m_err(i,j) = mask(i,j)*(slp_m(i,j) - slp_n(i,j))*2.*grav/abs(phis(i,j)-
gz0(i,j))
540 m_err(i,j) = mask(i,j)*(slp_m(i,j) - slp_n(i,j))
546 call corr(slp_m, slp_n, co, area)
549 if(
master)
write(*,*)
'SLP (Pa): RMS=', rms,
' Bias=', bias
550 if(
master)
write(*,*)
'SLP correlation=',co
568 du_obs(i,j,k) = profile(k)*(u_obs(i,j,k)-ua(i,j,k))*rdt
569 dv_obs(i,j,k) = profile(k)*(v_obs(i,j,k)-va(i,j,k))*rdt
582 du_obs(i,j,k) = du_obs(i,j,k) * mask(i,j) * ps_fac(i,j)
583 dv_obs(i,j,k) = dv_obs(i,j,k) * mask(i,j) * ps_fac(i,j)
585 u_dt(i,j,k) = u_dt(i,j,k) + du_obs(i,j,k)
586 v_dt(i,j,k) = v_dt(i,j,k) + dv_obs(i,j,k)
587 ua(i,j,k) = ua(i,j,k) + du_obs(i,j,k)*dt
588 va(i,j,k) = va(i,j,k) + dv_obs(i,j,k)*dt
595 du_obs(i,j,k) = du_obs(i,j,k) * mask(i,j)
596 dv_obs(i,j,k) = dv_obs(i,j,k) * mask(i,j)
598 u_dt(i,j,k) = u_dt(i,j,k) + du_obs(i,j,k)
599 v_dt(i,j,k) = v_dt(i,j,k) + dv_obs(i,j,k)
600 ua(i,j,k) = ua(i,j,k) + du_obs(i,j,k)*dt
601 va(i,j,k) = va(i,j,k) + dv_obs(i,j,k)*dt
629 t_dt(i,j,k) = prof_t(k)*(t_obs(i,j,k)/
virq(q(i,j,k,1:
num_gas))-pt(i,j,k))*rdt*ps_fac(i,j)
631 t_dt(i,j,k) = prof_t(k)*(t_obs(i,j,k)/(1.+zvir*q(i,j,k,1))-pt(i,j,k))*rdt*ps_fac(i,j)
639 t_dt(i,j,k) = prof_t(k)*(t_obs(i,j,k)/
virq(q(i,j,k,1:
num_gas))-pt(i,j,k))*rdt
641 t_dt(i,j,k) = prof_t(k)*(t_obs(i,j,k)/(1.+zvir*q(i,j,k,1))-pt(i,j,k))*rdt
664 pe2(i,k+1) = pe2(i,k) + delp(i,j,k)
669 peln(i,k) = log(pe2(i,k))
677 h2(i,j) = h2(i,j) + (t_obs(i,j,k)-pt(i,j,k)*
virq(q(i,j,k,1:
num_gas)))*(peln(i,k+1)-peln(i,k))
679 h2(i,j) = h2(i,j) + (t_obs(i,j,k)-pt(i,j,k)*(1.+zvir*q(i,j,k,1)))*(peln(i,k+1)-peln(i,k))
684 h2(i,j) = h2(i,j) / (peln(i,npz+1)-peln(i,kht+1))
685 h2(i,j) = rdt*ps_fac(i,j)*h2(i,j)*mask(i,j)
693 t_dt(i,j,k) = h2(i,j) / (1.+zvir*q(i,j,k,1))
707 pt(i,j,k) = pt(i,j,k) + t_dt(i,j,k)*dt*mask(i,j)
715 rdt = 1./(
tau_q/factor + dt)
719 if ( press(k) >
p_wvp )
then 723 q(i,j,k,iq) = q(i,j,k,iq)*delp(i,j,k)
730 delp(i,j,k) = delp(i,j,k)*(1.-q(i,j,k,1))
731 q_dt(i,j,k) = prof_q(k)*(max(
q_min,q_obs(i,j,k))-q(i,j,k,1))*rdt*mask(i,j)
732 q(i,j,k,1) = q(i,j,k,1) + q_dt(i,j,k)*dt
733 delp(i,j,k) = delp(i,j,k)/(1.-q(i,j,k,1))
739 q(i,j,k,iq) = q(i,j,k,iq)/delp(i,j,k)
751 deallocate ( ps_obs )
754 call breed_srf_winds(time, dt, npz, u_obs, v_obs, ak, bk, ps, phis, delp, ua, va, u_dt, v_dt, pt, q, nwat, zvir, gridstruct)
795 subroutine ps_nudging(dt, factor, npz, ak, bk, ps_obs, mask, tm, ps, phis, delp, ua, va, pt, nwat, q, bd, npx, npy, gridstruct, domain)
797 real,
intent(in):: dt, factor
798 integer,
intent(in):: npz, nwat, npx, npy
799 real,
intent(in),
dimension(npz+1):: ak, bk
800 type(fv_grid_bounds_type),
intent(IN) :: bd
801 real,
intent(in):: phis(isd:ied,jsd:jed)
802 real,
intent(in),
dimension(is:ie,js:je):: ps_obs, mask, tm
803 type(fv_grid_type),
intent(IN),
target :: gridstruct
804 type(domain2d),
intent(INOUT) :: domain
806 real,
intent(inout),
dimension(isd:ied,jsd:jed):: ps
807 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: delp, pt, ua, va
808 real,
intent(inout):: q(isd:ied,jsd:jed,npz,nwat)
810 real,
dimension(is:ie,js:je):: ps_dt
811 integer,
parameter:: kmax = 100
812 real:: pn0(kmax+1), pk0(kmax+1), pe0(kmax+1)
813 real,
dimension(is:ie,npz+1):: pe2, peln
814 real:: pst, dbk, pt0, rdt, bias
817 real,
pointer,
dimension(:,:) :: area
823 area => gridstruct%area
826 if ( kmax < km )
call mpp_error(fatal,
'==> KMAX must be larger than km')
836 pe0(k) =
ak0(k) +
bk0(k)*ps_obs(i,j)
837 pk0(k) = (
ak0(k) +
bk0(k)*ps_obs(i,j))**kappa
839 if( phis(i,j)>
gz0(i,j) )
then 841 if( phis(i,j)<
gz3(i,j,k) .and. phis(i,j) >=
gz3(i,j,k+1) )
then 842 pst = pk0(k) + (pk0(k+1)-pk0(k))*(
gz3(i,j,k)-phis(i,j))/(
gz3(i,j,k)-
gz3(i,j,k+1))
847 pn0(km ) = log(
ak0(km) +
bk0(km)*ps_obs(i,j))
848 pn0(km+1) = log(ps_obs(i,j))
851 pkx = (pk0(km+1)-pk0(km))/(kappa*(pn0(km+1)-pn0(km)))
852 pt0 = tm(i,j)/exp(kappax(km)*log(pkx))
853 pkx = exp((kappax(km)-1.0)*log(pkx))
854 pst = pk0(km+1) + (
gz0(i,j)-phis(i,j))/(cp_air*pt0*pkx)
856 pt0 = tm(i,j)/(pk0(km+1)-pk0(km))*(kappa*(pn0(km+1)-pn0(km)))
857 pst = pk0(km+1) + (
gz0(i,j)-phis(i,j))/(cp_air*pt0)
862 666 ps_dt(i,j) = pst**(1./(kappa*kappax(k))) - ps(i,j)
864 666 ps_dt(i,j) = pst**(1./kappa) - ps(i,j)
873 ps_dt(i,j) = sign( min(10.e2, abs(ps_dt(i,j))), ps_dt(i,j) )
874 ps_dt(i,j) = mask(i,j)*ps_dt(i,j) * max( 0., 1.-abs(
gz0(i,j)-phis(i,j))/(grav*500.) )
886 q(i,j,k,iq) = q(i,j,k,iq) * delp(i,j,k)
895 ua(i,j,k) = ua(i,j,k) * delp(i,j,k)
896 va(i,j,k) = va(i,j,k) * delp(i,j,k)
904 peln(i,1) = log(pe2(i,1))
908 pe2(i,k) = pe2(i,k-1) + delp(i,j,k-1)
909 peln(i,k) = log(pe2(i,k))
914 pt(i,j,k) = pt(i,j,k)*(peln(i,k+1)-peln(i,k))
927 rdt = dt / (
tau_ps/factor + dt)
929 dbk = rdt*(bk(k+1) - bk(k))
932 delp(i,j,k) = delp(i,j,k) + dbk*ps_dt(i,j)
933 ps(i,j) = delp(i,j,k) + ps(i,j)
942 pe2(i,k) = pe2(i,k-1) + delp(i,j,k-1)
943 peln(i,k) = log(pe2(i,k))
948 pt(i,j,k) = pt(i,j,k)/(peln(i,k+1)-peln(i,k))
958 q(i,j,k,iq) = q(i,j,k,iq) / delp(i,j,k)
967 ua(i,j,k) = ua(i,j,k) / delp(i,j,k)
968 va(i,j,k) = va(i,j,k) / delp(i,j,k)
977 integer,
intent(IN) :: is, ie, js, je
978 integer,
intent(IN) :: isd, ied, jsd,jed
979 real,
intent(inout):: ps_dt(is:ie,js:je)
980 real,
intent(IN),
dimension(isd:ied,jsd:jed) :: area
982 real:: esl, total_area
986 total_area = 4.*pi*
radius**2
989 bias =
g0_sum(ps_dt, is, ie, js, je, 1, .true., isd, ied, jsd, jed, area)
991 if ( abs(bias) < esl )
then 992 if(
master .and.
nudge_debug)
write(*,*)
'Very small PS bias=', -bias,
' No bais adjustment' 998 if ( bias > 0. )
then 1002 if ( ps_dt(i,j) > 0. )
then 1003 psum = psum + area(i,j)
1008 call mp_reduce_sum( psum )
1009 bias = bias * total_area / psum
1013 if ( ps_dt(i,j) > 0.0 )
then 1014 ps_dt(i,j) = max(0.0, ps_dt(i,j) - bias)
1022 if ( ps_dt(i,j) < 0. )
then 1023 psum = psum + area(i,j)
1028 call mp_reduce_sum( psum )
1029 bias = bias * total_area / psum
1033 if ( ps_dt(i,j) < 0.0 )
then 1034 ps_dt(i,j) = min(0.0, ps_dt(i,j) - bias)
1044 real function g0_sum(p, ifirst, ilast, jfirst, jlast, mode, reproduce, isd, ied, jsd, jed, area)
1046 integer,
intent(IN) :: ifirst, ilast
1047 integer,
intent(IN) :: jfirst, jlast
1048 integer,
intent(IN) :: mode
1049 logical,
intent(IN) :: reproduce
1050 real,
intent(IN) :: p(ifirst:ilast,jfirst:jlast)
1060 gsum = gsum + p(i,j)*area(i,j)
1064 call mp_reduce_sum(gsum)
1072 if ( reproduce )
g0_sum =
real(g0_sum, 4) 1077 subroutine compute_slp(isc, iec, jsc, jec, tm, ps, gz, slp)
1078 integer,
intent(in):: isc, iec, jsc, jec
1079 real,
intent(in),
dimension(isc:iec,jsc:jec):: tm, ps, gz
1081 real,
intent(out):: slp(isc:iec,jsc:jec)
1086 slp(i,j) = ps(i,j) * exp( gz(i,j)/(rdgas*(tm(i,j) + 3.25e-3*gz(i,j)/grav)) )
1093 subroutine get_obs(Time, dt, zvir, ak, bk, ps, ts, ps_obs, delp, pt, nwat, q, u_obs, v_obs, t_obs, q_obs, &
1094 phis, ua, va, u_dt, v_dt, npx, npy, npz, factor, mask, ptop, bd, gridstruct, domain)
1095 type(time_type),
intent(in):: Time
1096 integer,
intent(in):: npz, nwat, npx, npy
1097 real,
intent(in):: zvir, ptop
1098 real,
intent(in):: dt, factor
1099 real,
intent(in),
dimension(npz+1):: ak, bk
1100 type(fv_grid_bounds_type),
intent(IN) :: bd
1101 real,
intent(in),
dimension(isd:ied,jsd:jed):: phis
1102 real,
intent(in),
dimension(is:ie,js:je):: mask
1103 real,
intent(inout),
dimension(isd:ied,jsd:jed):: ps
1104 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: delp, pt, ua, va, u_dt, v_dt
1105 real,
intent(inout)::q(isd:ied,jsd:jed,npz,*)
1106 real,
intent(out),
dimension(is:ie,js:je):: ts, ps_obs
1107 real,
intent(out),
dimension(is:ie,js:je,npz):: u_obs, v_obs, t_obs, q_obs
1108 type(fv_grid_type),
intent(IN) :: gridstruct
1109 type(domain2d),
intent(INOUT) :: domain
1111 real:: tm(is:ie,js:je)
1112 real(KIND=4),
allocatable,
dimension(:,:,:):: ut, vt, wt
1113 real,
allocatable,
dimension(:,:,:):: uu, vv
1114 integer :: seconds, days
1118 call get_time (time, seconds, days)
1121 seconds = seconds - nint(dt)
1133 forecast_mode = .true.
1137 if (
master)
write(*,*)
'*** L-S nudging Ended at', days, seconds
1154 call get_ncep_analysis (
ps_dat(:,:,2),
u_dat(:,:,:,2),
v_dat(:,:,:,2), &
1168 if ( beta < 0. .or. beta > (1.+1.e-7) )
then 1169 if (
master)
write(*,*)
'Nudging: beta=', beta
1170 call mpp_error(fatal,
'==> Error from get_obs:data out of range')
1176 beta = 1.; alpha = 0.
1177 if(
nudge_debug .and.
master)
write(*,*)
'Doing final adiabatic initialization/nudging' 1191 allocate ( wt(is:ie,js:je,
km) )
1192 wt(:,:,:) = alpha*
t_dat(:,:,:,1) + beta*
t_dat(:,:,:,2)
1194 call get_int_hght(npz, ak, bk, ps(is:ie,js:je), delp, ps_obs(is:ie,js:je), wt)
1197 tm(i,j) = wt(i,j,
km)
1202 allocate ( uu(isd:ied,jsd:jed,npz) )
1203 allocate ( vv(isd:ied,jsd:jed,npz) )
1206 call ps_nudging(dt, factor, npz, ak, bk, ps_obs, mask, tm, ps, phis, delp, uu, vv, pt, nwat, q, bd, npx, npy, gridstruct, domain)
1210 u_dt(i,j,k) = u_dt(i,j,k) + (uu(i,j,k) - ua(i,j,k)) / dt
1211 v_dt(i,j,k) = v_dt(i,j,k) + (vv(i,j,k) - va(i,j,k)) / dt
1219 allocate ( ut(is:ie,js:je,npz) )
1220 allocate ( vt(is:ie,js:je,npz) )
1224 call remap_uv(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1225 km,
ps_dat(is:ie,js:je,1),
u_dat(:,:,:,1),
v_dat(:,:,:,1), ptop )
1227 u_obs(:,:,:) = alpha*ut(:,:,:)
1228 v_obs(:,:,:) = alpha*vt(:,:,:)
1230 call remap_uv(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1231 km,
ps_dat(is:ie,js:je,2),
u_dat(:,:,:,2),
v_dat(:,:,:,2), ptop )
1233 u_obs(:,:,:) = u_obs(:,:,:) + beta*ut(:,:,:)
1234 v_obs(:,:,:) = v_obs(:,:,:) + beta*vt(:,:,:)
1237 call remap_tq(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1238 km,
ps_dat(is:ie,js:je,1),
t_dat(:,:,:,1),
q_dat(:,:,:,1), zvir, ptop)
1240 t_obs(:,:,:) = alpha*ut(:,:,:)
1241 q_obs(:,:,:) = alpha*vt(:,:,:)
1243 call remap_tq(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1244 km,
ps_dat(is:ie,js:je,2),
t_dat(:,:,:,2),
q_dat(:,:,:,2), zvir, ptop)
1246 t_obs(:,:,:) = t_obs(:,:,:) + beta*ut(:,:,:)
1247 q_obs(:,:,:) = q_obs(:,:,:) + beta*vt(:,:,:)
1256 subroutine fv_nwp_nudge_init(time, axes, npz, zvir, ak, bk, ts, phis, gridstruct, ks, npx, neststruct, bd)
1257 character(len=17) :: mod_name =
'fv_nudge' 1258 type(time_type),
intent(in):: time
1259 integer,
intent(in):: axes(4)
1260 integer,
intent(in):: npz
1261 real,
intent(in):: zvir
1263 real,
intent(in),
dimension(bd%isd:bd%ied,bd%jsd:bd%jed):: phis
1264 real,
intent(in),
dimension(npz+1):: ak, bk
1265 real,
intent(out),
dimension(bd%is:bd%ie,bd%js:bd%je):: ts
1267 integer,
intent(in) :: ks, npx
1270 real :: missing_value = -1.e10
1273 integer :: i, j, j1, f_unit, unit, io, ierr, nt, k
1277 integer :: input_fname_unit
1278 character(len=128) :: fname_tmp
1281 real,
pointer,
dimension(:,:,:) :: agrid
1294 agrid => gridstruct%agrid
1301 if (neststruct%nested)
then 1303 grid_size = 1.e7/(neststruct%npx_global*neststruct%refinement_of_global)
1314 if( file_exist(
'input.nml' ) )
then 1315 unit = open_namelist_file()
1317 do while ( io .ne. 0 )
1318 read( unit, nml = fv_nwp_nudge_nml, iostat = io, end = 10 )
1319 ierr = check_nml_error(io,
'fv_nwp_nudge_nml')
1321 10
call close_file ( unit )
1323 call write_version_number (
'FV_NUDGE_MOD', version )
1326 write( f_unit, nml = fv_nwp_nudge_nml )
1327 write(*,*)
'NWP nudging initialized.' 1331 input_fname_unit = get_unit()
1336 do while ( io .eq. 0 )
1337 read ( input_fname_unit,
'(a)', iostat = io, end = 101 ) fname_tmp
1338 if( trim(fname_tmp) .ne.
"" )
then 1351 if( nt .eq. 1 )
then 1354 write(*,*)
'From NCEP file list: ', nt,
file_names(nt)
1362 if( nt .eq. 1 )
then 1365 write(*,*)
'From NCEP file list: ', nt,
file_names(nt)
1372 101
close( input_fname_unit )
1377 if (
file_names(nt) ==
"No_File_specified" )
then 1384 id_ht_err = register_diag_field( mod_name,
'ht_error', axes(1:2), time, &
1385 'height_error',
'DEG K', missing_value=missing_value )
1401 if(
master)
write(*,*)
'NCEP analysis dimensions:',
im,
jm,
km 1408 allocate (
lon(
im) )
1409 allocate (
lat(
jm) )
1411 call _get_var1 (ncid,
'lon',
im,
lon )
1412 call _get_var1 (ncid,
'lat',
jm,
lat )
1422 allocate (
ak0(
km+1) )
1423 allocate (
bk0(
km+1) )
1425 call _get_var1 (ncid,
'hyai',
km+1,
ak0, found )
1426 if ( .not. found )
ak0(:) = 0.
1428 call _get_var1 (ncid,
'hybi',
km+1,
bk0 )
1434 if (
bk0(1) < 1.e-9 )
ak0(1) = max(1.e-9,
ak0(1))
1438 write(*,*) k, 0.5*(ak(k)+ak(k+1))+0.5*(bk(k)+bk(k+1))*1.e5,
'del-B=', bk(k+1)-bk(k)
1473 call get_ncep_analysis (
ps_dat(:,:,nt),
u_dat(:,:,:,nt),
v_dat(:,:,:,nt), &
1486 real,
intent(in):: zvir
1487 character(len=128),
intent(in):: fname
1488 integer,
intent(inout):: nfile
1490 type(fv_grid_bounds_type),
intent(IN) :: bd
1491 real,
intent(out),
dimension(is:ie,js:je):: ts, ps
1492 real(KIND=4),
intent(out),
dimension(is:ie,js:je,km):: u, v, t, q
1494 real(kind=4),
allocatable:: wk0(:,:), wk1(:,:), wk2(:,:), wk3(:,:,:)
1496 integer:: i, j, k, npt
1497 integer:: i1, i2, j1, ncid
1499 logical:: read_ts = .true.
1500 logical:: land_ts = .false.
1502 integer:: status, var3id
1503 #include <netcdf.inc> 1506 if( .not. file_exist(fname) )
then 1507 call mpp_error(fatal,
'==> Error from get_ncep_analysis: file not found')
1510 if(
master)
write(*,*)
'Reading NCEP anlysis file:', fname
1514 allocate ( wk1(
im,
jm) )
1519 if ( .not. land_ts )
then 1520 allocate ( wk0(
im,
jm) )
1524 status = nf_inq_varid(ncid,
'ORO', var3id)
1525 if (status .eq. nf_noerr)
then 1529 status = nf_inq_varid(ncid,
'LAND', var3id)
1530 if (status .eq. nf_noerr)
then 1533 call mpp_error(fatal,
'Neither ORO nor LAND exists in re-analysis data')
1542 if( abs(wk0(i,j)-1.) > 0.99 )
then 1543 tmean = tmean + wk1(i,j)
1550 if ( npt /= 0 )
then 1551 tmean= tmean /
real(npt)
1553 if( abs(wk0(i,j)-1.) <= 0.99 )
then 1556 elseif ( i==
im )
then 1561 if ( abs(wk0(i2,j)-1.)>0.99 )
then 1562 wk1(i,j) = wk1(i2,j)
1563 elseif ( abs(wk0(i1,j)-1.)>0.99 )
then 1564 wk1(i,j) = wk1(i1,j)
1580 ts(i,j) =
s2c(i,j,1)*wk1(i1,j1 ) +
s2c(i,j,2)*wk1(i2,j1 ) + &
1581 s2c(i,j,3)*wk1(i2,j1+1) +
s2c(i,j,4)*wk1(i1,j1+1)
1589 if(
master)
call pmaxmin(
'SST_ncep', sst_ncep, i_sst, j_sst, 1.)
1593 if (
master)
write(*,*)
'Done processing NCEP SST' 1609 ps(i,j) =
s2c(i,j,1)*wk2(i1,j1 ) +
s2c(i,j,2)*wk2(i2,j1 ) + &
1610 s2c(i,j,3)*wk2(i2,j1+1) +
s2c(i,j,4)*wk2(i1,j1+1)
1616 status = nf_inq_varid(ncid,
'PHIS', var3id)
1617 if (status .eq. nf_noerr)
then 1621 status = nf_inq_varid(ncid,
'PHI', var3id)
1622 if (status .eq. nf_noerr)
then 1626 call mpp_error(fatal,
'Neither PHIS nor PHI exists in re-analysis data')
1637 gz0(i,j) =
s2c(i,j,1)*wk2(i1,j1 ) +
s2c(i,j,2)*wk2(i2,j1 ) + &
1638 s2c(i,j,3)*wk2(i2,j1+1) +
s2c(i,j,4)*wk2(i1,j1+1)
1658 u(i,j,k) =
s2c(i,j,1)*wk3(i1,j1 ,k) +
s2c(i,j,2)*wk3(i2,j1 ,k) + &
1659 s2c(i,j,3)*wk3(i2,j1+1,k) +
s2c(i,j,4)*wk3(i1,j1+1,k)
1672 v(i,j,k) =
s2c(i,j,1)*wk3(i1,j1 ,k) +
s2c(i,j,2)*wk3(i2,j1 ,k) + &
1673 s2c(i,j,3)*wk3(i2,j1+1,k) +
s2c(i,j,4)*wk3(i1,j1+1,k)
1689 q(i,j,k) =
s2c(i,j,1)*wk3(i1,j1 ,k) +
s2c(i,j,2)*wk3(i2,j1 ,k) + &
1690 s2c(i,j,3)*wk3(i2,j1+1,k) +
s2c(i,j,4)*wk3(i1,j1+1,k)
1704 t(i,j,k) =
s2c(i,j,1)*wk3(i1,j1 ,k) +
s2c(i,j,2)*wk3(i2,j1 ,k) + &
1705 s2c(i,j,3)*wk3(i2,j1+1,k) +
s2c(i,j,4)*wk3(i1,j1+1,k)
1717 t(i,j,k) = t(i,j,k)*(1.+zvir*q(i,j,k))
1734 real,
intent(IN),
dimension(isd:ied,jsd:jed,2) :: agrid
1740 integer i,j, i1, i2, jc, i0, j0
1743 rdlon(i) = 1. / (
lon(i+1) -
lon(i))
1745 rdlon(im) = 1. / (
lon(1) + 2.*pi -
lon(im))
1748 rdlat(j) = 1. / (
lat(j+1) -
lat(j))
1756 if ( agrid(i,j,1)>
lon(im) )
then 1758 a1 = (agrid(i,j,1)-
lon(im)) * rdlon(im)
1759 elseif ( agrid(i,j,1)<
lon(1) )
then 1761 a1 = (agrid(i,j,1)+2.*pi-
lon(im)) * rdlon(im)
1764 if ( agrid(i,j,1)>=
lon(i0) .and. agrid(i,j,1)<=
lon(i0+1) )
then 1766 a1 = (agrid(i,j,1)-
lon(i1)) * rdlon(i0)
1773 if ( agrid(i,j,2)<
lat(1) )
then 1776 elseif ( agrid(i,j,2)>
lat(jm) )
then 1781 if ( agrid(i,j,2)>=
lat(j0) .and. agrid(i,j,2)<=
lat(j0+1) )
then 1783 b1 = (agrid(i,j,2)-
lat(jc)) * rdlat(jc)
1790 if ( a1<0.0 .or. a1>1.0 .or. b1<0.0 .or. b1>1.0 )
then 1791 write(*,*)
'gid=', mpp_pe(), i,j,a1, b1
1794 s2c(i,j,1) = (1.-a1) * (1.-b1)
1795 s2c(i,j,2) = a1 * (1.-b1)
1796 s2c(i,j,3) = a1 * b1
1797 s2c(i,j,4) = (1.-a1) * b1
1809 real(kind=4),
intent(in):: sst(im,jm)
1816 real:: c1, c2, c3, c4
1817 integer i,j, i1, i2, jc, i0, j0, it, jt
1820 rdlon(i) = 1. / (
lon(i+1) -
lon(i))
1822 rdlon(im) = 1. / (
lon(1) + 2.*pi -
lon(im))
1825 rdlat(j) = 1. / (
lat(j+1) -
lat(j))
1832 delx = 360./
real(i_sst) 1833 dely = 180./
real(j_sst) 1838 yc = (-90. + dely * (0.5+
real(j-1))) * deg2rad
1839 if ( yc<
lat(1) )
then 1842 elseif ( yc>
lat(jm) )
then 1847 if ( yc>=
lat(j0) .and. yc<=
lat(j0+1) )
then 1850 b1 = (yc-
lat(jc)) * rdlat(jc)
1859 xc = delx * (0.5+
real(i-1)) * deg2rad
1860 if ( xc>
lon(im) )
then 1862 a1 = (xc-
lon(im)) * rdlon(im)
1863 elseif ( xc<
lon(1) )
then 1865 a1 = (xc+2.*pi-
lon(im)) * rdlon(im)
1868 if ( xc>=
lon(i0) .and. xc<=
lon(i0+1) )
then 1871 a1 = (xc-
lon(i1)) * rdlon(i0)
1881 c1 = (1.-a1) * (1.-b1)
1886 sst_ncep(i,j) = c1*sst(i1,jc ) + c2*sst(i2,jc ) + &
1887 c3*sst(i2,jc+1) + c4*sst(i1,jc+1)
1894 subroutine get_int_hght(npz, ak, bk, ps, delp, ps0, tv)
1895 integer,
intent(in):: npz
1896 real,
intent(in):: ak(npz+1), bk(npz+1)
1897 real,
intent(in),
dimension(is:ie,js:je):: ps, ps0
1898 real,
intent(in),
dimension(isd:ied,jsd:jed,npz):: delp
1899 real(KIND=4),
intent(in),
dimension(is:ie,js:je,km):: tv
1901 real,
dimension(is:ie,km+1):: pn0
1910 pn0(i,k) = log(
ak0(k) +
bk0(k)*ps0(i,j) )
1919 gz3(i,j,k) =
gz3(i,j,k+1) + rdgas*tv(i,j,k)*(pn0(i,k+1)-pn0(i,k))
1929 subroutine remap_tq( npz, ak, bk, ps, delp, t, q, &
1930 kmd, ps0, ta, qa, zvir, ptop)
1931 integer,
intent(in):: npz, kmd
1932 real,
intent(in):: zvir, ptop
1933 real,
intent(in):: ak(npz+1), bk(npz+1)
1934 real,
intent(in),
dimension(is:ie,js:je):: ps0
1935 real,
intent(inout),
dimension(is:ie,js:je):: ps
1936 real,
intent(in),
dimension(isd:ied,jsd:jed,npz):: delp
1937 real(KIND=4),
intent(in),
dimension(is:ie,js:je,kmd):: ta, qa
1938 real(KIND=4),
intent(out),
dimension(is:ie,js:je,npz):: t
1941 real(KIND=4),
intent(out),
dimension(is:ie,js:je,npz):: q
1943 real,
dimension(is:ie,kmd):: tp, qp
1944 real,
dimension(is:ie,kmd+1):: pe0, pn0
1945 real,
dimension(is:ie,npz):: qn1
1946 real,
dimension(is:ie,npz+1):: pe1, pn1
1953 pe0(i,k) =
ak0(k) +
bk0(k)*ps0(i,j)
1954 pn0(i,k) = log(pe0(i,k))
1965 pe1(i,k+1) = pe1(i,k) + delp(i,j,k)
1969 ps(i,j) = pe1(i,npz+1)
1973 pn1(i,k) = log(pe1(i,k))
1983 call mappm(kmd, pe0, qp, npz, pe1, qn1,
is,
ie, 0,
kord_data, ptop)
1996 call mappm(kmd, pn0, tp, npz, pn1, qn1,
is,
ie, 1,
kord_data, ptop)
2009 subroutine remap_uv(npz, ak, bk, ps, delp, u, v, kmd, ps0, u0, v0, ptop)
2010 integer,
intent(in):: npz
2011 real,
intent(IN):: ptop
2012 real,
intent(in):: ak(npz+1), bk(npz+1)
2013 real,
intent(inout):: ps(is:ie,js:je)
2014 real,
intent(in),
dimension(isd:ied,jsd:jed,npz):: delp
2015 real(KIND=4),
intent(inout),
dimension(is:ie,js:je,npz):: u, v
2017 integer,
intent(in):: kmd
2018 real,
intent(in):: ps0(is:ie,js:je)
2019 real(KIND=4),
intent(in),
dimension(is:ie,js:je,kmd):: u0, v0
2022 real,
dimension(is:ie,kmd+1):: pe0
2023 real,
dimension(is:ie,npz+1):: pe1
2024 real,
dimension(is:ie,kmd):: qt
2025 real,
dimension(is:ie,npz):: qn1
2034 pe0(i,k) =
ak0(k) +
bk0(k)*ps0(i,j)
2045 pe1(i,k+1) = pe1(i,k) + delp(i,j,k)
2049 ps(i,j) = pe1(i,npz+1)
2059 call mappm(kmd, pe0, qt, npz, pe1, qn1, is,ie, -1,
kord_data, ptop)
2073 call mappm(kmd, pe0, qt, npz, pe1, qn1, is,ie, -1,
kord_data, ptop)
2086 deallocate (
t_dat )
2087 deallocate (
q_dat )
2090 deallocate (
u_dat )
2091 deallocate (
v_dat )
2111 real :: slp_mask = 100900.
2113 type(time_type),
intent(in):: time
2114 real,
intent(inout):: mask(is:ie,js:je)
2115 real(kind=R_GRID),
intent(IN),
dimension(isd:ied,jsd:jed,2) :: agrid
2117 real(kind=R_GRID):: pos(2)
2128 call get_slp_obs(time,
nobs_tc(n),
x_obs(1,n),
y_obs(1,n),
wind_obs(1,n),
mslp_obs(1,n),
mslp_out(1,n),
rad_out(1,n), &
2129 time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_vor)
2131 if ( slp_o<880.e2 .or. slp_o>min(
slp_env,slp_mask) .or. abs(pos(2))*
rad2deg>40. )
goto 5000
2133 if ( r_vor < 30.e3 )
then 2139 dist = great_circle_dist(pos, agrid(i,j,1:2),
radius)
2140 if( dist < 6.*r_vor )
then 2141 mask(i,j) = mask(i,j) * ( 1. -
mask_fac*exp(-(0.5*dist/r_vor)**2)*min(1.,(
slp_env-slp_o)/10.e2) )
2155 subroutine breed_slp_inline(nstep, dt, npz, ak, bk, phis, pe, pk, peln, pkz, delp, u, v, pt, q, nwat, &
2156 zvir, gridstruct, ks, domain_local, bd, hydrostatic)
2158 integer,
intent(in):: nstep, npz, nwat, ks
2159 real,
intent(in):: dt
2160 real,
intent(in):: zvir
2161 real,
intent(in),
dimension(npz+1):: ak, bk
2162 logical,
intent(in):: hydrostatic
2165 type(domain2d),
intent(INOUT) :: domain_local
2169 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: delp, pt
2172 real,
intent(inout):: pk(
is:
ie,
js:
je, npz+1)
2173 real,
intent(inout):: pe(
is-1:
ie+1, npz+1,
js-1:
je+1)
2174 real,
intent(inout):: pkz(
is:
ie,
js:
je,npz)
2175 real,
intent(out):: peln(
is:
ie,npz+1,
js:
je)
2179 type(time_type):: time
2184 real(kind=R_GRID):: pos(2)
2188 real:: relx0, relx, f1, pbreed, pbtop, delp0, dp0
2189 real:: ratio, p_count, p_sum, a_sum, mass_sink, delps
2190 real:: p_lo, p_hi, tau_vt, mslp0
2191 real:: split_time, fac, pdep, r2, r3
2192 integer year, month, day, hour, minute, second
2193 integer n, i, j, k, iq, k0
2195 real,
pointer :: area(:,:)
2196 real(kind=R_GRID),
pointer :: agrid(:,:,:)
2201 if ( forecast_mode )
return 2203 agrid => gridstruct%agrid_64
2204 area => gridstruct%area
2207 if(
master)
write(*,*)
'NO TC data to process' 2215 call get_date(
fv_time, year, month, day, hour, minute, second)
2217 if (
master)
write(*,*)
'Warning: The year in storm track data is not the same as model year' 2221 split_time =
calday(year, month, day, hour, minute, second) + dt*
real(nstep)/86400.
2227 if (
master)
write(*,*)
'*** Vortext Breeding Ended at', day, hour, minute, second
2240 ps(i,j) = ps(i,j) + delp(i,j,k)
2245 tm(i,j) = pkz(i,j,npz)*pt(i,j,npz) / cp_air
2257 u(i,j,k) = u(i,j,k) * (delp(i,j-1,k)+delp(i,j,k))
2262 v(i,j,k) = v(i,j,k) * (delp(i-1,j,k)+delp(i,j,k))
2271 pt(i,j,k) = pt(i,j,k)*pkz(i,j,k)*delp(i,j,k)
2280 q(i,j,k,iq) = q(i,j,k,iq) * delp(i,j,k)
2294 mslp0 = min( mslp0,
mslp_obs(i,n) )
2301 call get_slp_obs(time,
nobs_tc(n),
x_obs(1,n),
y_obs(1,n),
wind_obs(1,n),
mslp_obs(1,n),
mslp_out(1,n),
rad_out(1,n), &
2302 time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_env, stime=split_time, fact=fac)
2304 if ( slp_o<87500. .or. slp_o>
slp_env .or. abs(pos(2))*
rad2deg>45. )
then 2310 write(*,*)
'Vortex breeding for TC:', n,
' slp=',slp_o/100.,pos(1)*
rad2deg,pos(2)*
rad2deg 2317 if ( slp_o > 1000.e2 )
then 2322 pbtop = max(100.e2, 850.e2-25.*(1000.e2-slp_o))
2324 if ( slp_o > 1000.e2 )
then 2327 pbtop = max(100.e2, 750.e2-20.*(1000.e2-slp_o))
2332 if ( slp_o > 1000.e2 )
then 2335 pbtop = max(500.e2, 900.e2-5.*(1000.e2-slp_o))
2340 pbreed = ak(k) + bk(k)*1.e5
2341 if ( pbreed>pbtop )
then 2350 dist(i,j) = great_circle_dist( pos, agrid(i,j,1:2),
radius)
2356 if ( r_vor < 30.e3 .or. p_env<900.e2 )
then 2367 if( dist(i,j)<(r_vor+
del_r) .and. dist(i,j)>r_vor .and. phis(i,j)<250.*grav )
then 2368 p_count = p_count + 1.
2369 p_sum = p_sum + slp(i,j)*area(i,j)
2370 a_sum = a_sum + area(i,j)
2375 call mp_reduce_sum(p_count)
2377 if ( p_count<32. )
then 2378 if(
nudge_debug .and.
master)
write(*,*) p_count,
'Skipping obs: too few p_count' 2382 call mp_reduce_sum(p_sum)
2383 call mp_reduce_sum(a_sum)
2384 p_env = p_sum / a_sum
2386 if(
nudge_debug .and.
master)
write(*,*)
'Environmental SLP=', p_env/100.,
' computed radius=', r_vor/1.e3
2388 if ( p_env>1015.e2 .or. p_env<920.e2 )
then 2390 if(
master)
write(*,*)
'Environmental SLP out of bound; skipping obs. p_count=', p_count, p_sum
2399 if ( p_env < max(
pre0_env, slp_o + 200.0) )
then 2401 r_vor = r_vor + 25.e3
2404 write(*,*)
'Computed environmental SLP too low' 2405 write(*,*)
' ', p_env/100., slp_o/100.,pos(1)*
rad2deg, pos(2)*
rad2deg 2408 if ( slp_o > 1003.e2 .and. r_vor > 1000.e3 )
then 2413 if ( r_vor < 1250.e3 )
goto 123
2420 tau_vt =
tau_vt_slp * (1. + (960.e2-slp_o)/100.e2 )
2421 tau_vt = max(abs(dt), tau_vt)
2424 relx0 = min(1., 2.*abs(dt)/tau_vt)
2426 relx0 = min(1., abs(dt)/tau_vt)
2436 if( dist(i,j) < r_vor .and. phis(i,j)<250.*grav )
then 2437 f1 = dist(i,j)/r_vor
2439 p_hi = p_env - (p_env-slp_o) * exp( -
r_hi*f1**2 )
2440 p_lo = p_env - (p_env-slp_o) * exp( -
r_lo*f1**2 )
2442 if ( ps(i,j) > p_hi .and. tm(i,j) <
tm_max )
then 2446 delps = relx*(ps(i,j) - p_hi) * min( (
tm_max-tm(i,j))/10., 1.)
2453 elseif ( slp(i,j) < p_lo )
then 2455 relx = max(0.5, relx0)
2456 delps = relx*(slp(i,j) - p_lo)
2464 pbreed = pbreed + delp(i,j,k)
2466 f1 = 1. - delps/(ps(i,j)-pbreed)
2468 delp(i,j,k) = delp(i,j,k)*f1
2470 mass_sink = mass_sink + delps*area(i,j)
2472 if ( delps > 0. )
then 2475 pbreed = pbreed + delp(i,j,k)
2477 f1 = 1. - delps/(ps(i,j)-pbreed)
2479 delp(i,j,k) = delp(i,j,k)*f1
2481 mass_sink = mass_sink + delps*area(i,j)
2485 if ( abs(delps) < 1. )
then 2486 delp(i,j,k) = delp(i,j,k) - delps
2487 mass_sink = mass_sink + delps*area(i,j)
2490 pdep = max(1.0, min(abs(0.4*delps), 0.2*dp0, 0.02*delp(i,j,k)))
2491 pdep = - min(pdep, abs(delps))
2492 delp(i,j,k) = delp(i,j,k) - pdep
2493 delps = delps - pdep
2494 mass_sink = mass_sink + pdep*area(i,j)
2505 call mp_reduce_sum(mass_sink)
2506 if ( abs(mass_sink)<1.e-40 )
goto 5000
2509 r3 = min( 4.*r_vor, max(2.*r_vor, 2500.e3) ) +
del_r 2514 if( dist(i,j)<r3 .and. dist(i,j)>r2 )
then 2515 p_sum = p_sum + area(i,j)
2520 call mp_reduce_sum(p_sum)
2521 mass_sink = mass_sink / p_sum
2522 if(
master .and.
nudge_debug)
write(*,*)
'TC#',n,
'Mass tele-ported (pa)=', mass_sink
2528 if( dist(i,j)<r3 .and. dist(i,j)>r2 )
then 2531 pbreed = pbreed + delp(i,j,k)
2533 f1 = 1. + mass_sink/(ps(i,j)-pbreed)
2535 delp(i,j,k) = delp(i,j,k)*f1
2548 ps(i,j) = ps(i,j) + delp(i,j,k)
2558 call mpp_update_domains(delp, domain_local, complete=.true.)
2567 pe(i,k,j) = pe(i,k-1,j) + delp(i,j,k-1)
2585 peln(i,k,j) = log(pe(i,k,j))
2586 pk(i,j,k) = pe(i,k,j)**kappa
2601 u(i,j,k) = u(i,j,k) / (delp(i,j-1,k)+delp(i,j,k))
2606 v(i,j,k) = v(i,j,k) / (delp(i-1,j,k)+delp(i,j,k))
2615 pkz(i,j,k) = (pk(i,j,k+1)-pk(i,j,k))/(kappa*(peln(i,k+1,j)-peln(i,k,j)))
2617 pkz(i,j,k) = exp(kappax(i,j,k)*log(pkz(i,j,k)))
2619 pt(i,j,k) = pt(i,j,k) / (pkz(i,j,k)*delp(i,j,k))
2633 q(i,j,k,iq) = q(i,j,k,iq) / delp(i,j,k)
2639 if(hydrostatic)
call mpp_update_domains(pt, domain_local, complete=.true.)
2648 subroutine breed_srf_w10(time, dt, npz, ak, bk, ps, phis, slp, delp, u, v, gridstruct)
2649 type(time_type),
intent(in):: time
2650 integer,
intent(in):: npz
2651 real,
intent(in):: dt
2652 real,
intent(in),
dimension(npz+1):: ak, bk
2653 real,
intent(in):: phis(isd:ied,jsd:jed)
2654 real,
intent(in):: ps(isd:ied,jsd:jed)
2655 real,
intent(in),
dimension(is:ie,js:je):: slp
2656 type(fv_grid_type),
intent(IN),
target :: gridstruct
2657 real,
intent(inout):: u(isd:ied,jsd:jed+1,npz)
2658 real,
intent(inout):: v(isd:ied+1,jsd:jed,npz)
2660 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: delp
2662 real,
dimension(is:ie,js:je):: us, vs
2663 real wu(is:ie, js:je+1)
2664 real wv(is:ie+1,js:je)
2665 real u1(is:ie), v1(is:ie)
2666 real:: dist(is:ie,js:je), wind(is:ie,js:je)
2667 real(kind=R_GRID):: pos(2)
2670 real:: r_vor, pc, p_count
2671 real:: r_max, speed, ut, vt, speed_local
2672 real:: u_bg, v_bg, mass, t_mass
2673 real:: relx0, relx, f1
2680 real,
pointer,
dimension(:,:) :: dx, dy, rdxa, rdya, a11, a21, a12, a22, area
2681 real(kind=R_GRID),
pointer,
dimension(:,:,:) :: agrid, vlon, vlat
2685 rdxa => gridstruct%rdxa
2686 rdya => gridstruct%rdya
2687 a11 => gridstruct%a11
2688 a21 => gridstruct%a21
2689 a12 => gridstruct%a12
2690 a22 => gridstruct%a22
2691 area => gridstruct%area
2692 agrid => gridstruct%agrid_64
2693 vlon => gridstruct%vlon
2694 vlat => gridstruct%vlat
2698 if(
master)
write(*,*)
'NO TC data to process' 2706 wu(i,j) = u(i,j,npz)*dx(i,j)
2712 wv(i,j) = v(i,j,npz)*dy(i,j)
2721 u1(i) = (wu(i,j) + wu(i,j+1)) * rdxa(i,j)
2722 v1(i) = (wv(i,j) + wv(i+1,j)) * rdya(i,j)
2724 us(i,j) = a11(i,j)*u1(i) + a12(i,j)*v1(i)
2725 vs(i,j) = a21(i,j)*u1(i) + a22(i,j)*v1(i)
2727 wind(i,j) = sqrt( us(i,j)**2 + vs(i,j)**2 )
2740 mslp0 = min( mslp0,
mslp_obs(i,n) )
2747 call get_slp_obs(time,
nobs_tc(n),
x_obs(1,n),
y_obs(1,n),
wind_obs(1,n),
mslp_obs(1,n),
mslp_out(1,n),
rad_out(1,n), &
2748 time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_env)
2750 if ( slp_o<90000. .or. slp_o>
slp_env .or. abs(pos(2))*
rad2deg>35. )
goto 3000
2755 dist(i,j) = great_circle_dist( pos, agrid(i,j,1:2),
radius )
2771 if( dist(i,j) < r_vor )
then 2774 if ( dist(i,j)<0.5*r_vor .and. abs(phis(i,j))>2.*grav ) p_count = p_count + 1.
2776 pc = min(pc, slp(i,j))
2778 if ( speed_local < wind(i,j) )
then 2779 speed_local = wind(i,j)
2787 call mp_reduce_sum(p_count)
2788 if ( p_count>32 )
goto 3000
2790 if ( w10_o < 0. )
then 2792 w10_o = 3.446778 * (1010.-slp_o/100.)**0.644
2796 call mp_reduce_max(speed)
2797 call mp_reduce_min(pc)
2799 if ( speed_local < speed )
then 2802 call mp_reduce_max(r_max)
2803 if( r_max<0. )
call mpp_error(fatal,
'==> Error in r_max')
2810 if ( w10_o > 12.5 )
then 2811 z0 = (0.085*w10_o - 0.58) * 1.e-3
2816 z0 = 0.0185/grav*(0.001*w10_o**2 + 0.028*w10_o)**2
2821 wind_fac = log(zz/z0) / log(10./z0)
2823 if( wind_fac<1. )
call mpp_error(fatal,
'==> Error in wind_fac')
2825 if ( pc < slp_o )
then 2833 speed = wind_fac*w10_o
2836 speed = max(wind_fac*w10_o, speed)
2837 if ( pc>1009.e2 ) r_max = 0.5 * r_vor
2844 r_max = min(0.75*r_vor, r_max)
2854 if( dist(i,j) <= min(r_vor,2.*
r_fac*r_max) .and. phis(i,j)<2.*grav )
then 2855 mass = area(i,j)*delp(i,j,npz)
2857 u_bg = u_bg + us(i,j)*mass
2858 v_bg = v_bg + vs(i,j)*mass
2860 t_mass = t_mass + mass
2861 p_count = p_count + 1.
2865 call mp_reduce_sum(p_count)
2866 if ( p_count<16. )
go to 3000
2868 call mp_reduce_sum(t_mass)
2869 call mp_reduce_sum(u_bg)
2870 call mp_reduce_sum(v_bg)
2871 u_bg = u_bg / t_mass
2872 v_bg = v_bg / t_mass
2881 if( dist(i,j) <= min(r_vor,
r_fac*r_max) .and. phis(i,j)<2.*grav )
then 2882 f1 = dist(i,j)/r_max
2884 if( dist(i,j)<=r_max )
then 2885 speed_local = speed * f1
2887 speed_local = speed / f1**0.75
2889 call tangent_wind(vlon(i,j,1:3), vlat(i,j,1:3), speed_local, pos, agrid(i,j,1:2), ut, vt)
2893 us(i,j) = relx*(ut-us(i,j))
2894 vs(i,j) = relx*(vt-vs(i,j))
2902 if( dist(i,j) <= min(r_vor,
r_fac*r_max) .and. phis(i,j)<2.*grav )
then 2903 f1 = dist(i,j)/r_max
2905 if( dist(i,j)<=r_max )
then 2906 speed_local = speed * f1
2908 speed_local = speed / f1**0.75
2910 call tangent_wind(vlon(i,j,1:3), vlat(i,j,1:3), speed_local, pos, agrid(i,j,1:2), ut, vt)
2914 us(i,j) = relx*(ut-us(i,j))
2915 vs(i,j) = relx*(vt-vs(i,j))
2927 subroutine breed_srf_winds(time, dt, npz, u_obs, v_obs, ak, bk, ps, phis, delp, ua, va, u_dt, v_dt, pt, q, nwat, zvir, gridstruct)
2929 type(time_type),
intent(in):: time
2930 integer,
intent(in):: npz, nwat
2931 real,
intent(in):: dt
2932 real,
intent(in):: zvir
2933 real,
intent(in),
dimension(npz+1):: ak, bk
2934 real,
intent(in),
dimension(isd:ied,jsd:jed):: phis, ps
2935 real,
intent(in),
dimension(is:ie,js:je,npz):: u_obs, v_obs
2936 type(fv_grid_type),
intent(IN),
target :: gridstruct
2938 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: delp, pt, ua, va, u_dt, v_dt
2939 real,
intent(inout)::q(isd:ied,jsd:jed,npz,nwat)
2941 real:: dist(is:ie,js:je), wind(is:ie,js:je)
2942 real:: slp(is:ie,js:je)
2943 real(kind=R_GRID):: pos(2)
2946 real:: r_vor, pc, p_count
2947 real:: r_max, speed, ut, vt, speed_local
2948 real:: u_bg, v_bg, mass, t_mass
2949 real:: relx0, relx, f1, rdt
2953 integer n, i, j, k, iq
2955 real,
pointer :: area(:,:)
2956 real(kind=R_GRID),
pointer :: vlon(:,:,:), vlat(:,:,:), agrid(:,:,:)
2958 area => gridstruct%area
2959 vlon => gridstruct%vlon
2960 vlat => gridstruct%vlat
2961 agrid => gridstruct%agrid_64
2964 if(
master)
write(*,*)
'NO TC data to process' 2974 slp(i,j) = ps(i,j)*exp(phis(i,j)/(rdgas*(pt(i,j,npz)+3.25e-3*phis(i,j)/grav)))
2975 wind(i,j) = sqrt( ua(i,j,npz)**2 + va(i,j,npz)**2 )
2987 mslp0 = min( mslp0,
mslp_obs(i,n) )
2994 call get_slp_obs(time,
nobs_tc(n),
x_obs(1,n),
y_obs(1,n),
wind_obs(1,n),
mslp_obs(1,n),
mslp_out(1,n),
rad_out(1,n), &
2995 time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_env)
2997 if ( slp_o<90000. .or. slp_o>
slp_env .or. abs(pos(2))*
rad2deg>35. )
goto 3000
3002 dist(i,j) = great_circle_dist( pos, agrid(i,j,1:2),
radius )
3019 if( dist(i,j) < r_vor )
then 3022 if ( dist(i,j)<0.5*r_vor .and. abs(phis(i,j))>2.*grav ) p_count = p_count + 1.
3024 pc = min(pc, slp(i,j))
3026 if ( speed_local < wind(i,j) )
then 3027 speed_local = wind(i,j)
3035 call mp_reduce_sum(p_count)
3036 if ( p_count>32 )
goto 3000
3038 if ( w10_o < 0. )
then 3040 w10_o = 3.446778 * (1010.-slp_o/100.)**0.644
3044 call mp_reduce_max(speed)
3045 call mp_reduce_min(pc)
3047 if ( speed_local < speed )
then 3050 call mp_reduce_max(r_max)
3051 if( r_max<0. )
call mpp_error(fatal,
'==> Error in r_max')
3058 if ( w10_o > 12.5 )
then 3059 z0 = (0.085*w10_o - 0.58) * 1.e-3
3064 z0 = 0.0185/grav*(0.001*w10_o**2 + 0.028*w10_o)**2
3069 wind_fac = log(zz/z0) / log(10./z0)
3071 if( wind_fac<1. )
call mpp_error(fatal,
'==> Error in wind_fac')
3073 if ( pc < slp_o )
then 3081 speed = wind_fac*w10_o
3084 speed = max(wind_fac*w10_o, speed)
3085 if ( pc>1009.e2 ) r_max = 0.5 * r_vor
3092 r_max = min(0.75*r_vor, r_max)
3102 if( dist(i,j) <= min(r_vor,2.*
r_fac*r_max) .and. phis(i,j)<2.*grav )
then 3103 mass = area(i,j)*delp(i,j,npz)
3109 u_bg = u_bg + u_obs(i,j,npz)*mass
3110 v_bg = v_bg + v_obs(i,j,npz)*mass
3111 t_mass = t_mass + mass
3112 p_count = p_count + 1.
3116 call mp_reduce_sum(p_count)
3117 if ( p_count<16. )
go to 3000
3119 call mp_reduce_sum(t_mass)
3120 call mp_reduce_sum(u_bg)
3121 call mp_reduce_sum(v_bg)
3122 u_bg = u_bg / t_mass
3123 v_bg = v_bg / t_mass
3132 if( dist(i,j) <= min(r_vor,
r_fac*r_max) .and. phis(i,j)<2.*grav )
then 3133 f1 = dist(i,j)/r_max
3135 if( dist(i,j)<=r_max )
then 3136 speed_local = speed * f1
3138 speed_local = speed / f1**0.75
3140 call tangent_wind(vlon(i,j,1:3), vlat(i,j,1:3), speed_local, pos, agrid(i,j,1:2), ut, vt)
3143 u_dt(i,j,k) = u_dt(i,j,k) + relx*(ut-ua(i,j,k)) * rdt
3144 v_dt(i,j,k) = v_dt(i,j,k) + relx*(vt-va(i,j,k)) * rdt
3146 ua(i,j,k) = ua(i,j,k) + relx*(ut-ua(i,j,k))
3147 va(i,j,k) = va(i,j,k) + relx*(vt-va(i,j,k))
3157 subroutine tangent_wind ( elon, elat, speed, po, pp, ut, vt )
3158 real,
intent(in):: speed
3159 real(kind=R_GRID),
intent(in):: po(2), pp(2)
3160 real(kind=R_GRID),
intent(in):: elon(3), elat(3)
3161 real,
intent(out):: ut, vt
3163 real(kind=R_GRID):: e1(3), eo(3), ep(3), op(3)
3168 op(:) = ep(:) - eo(:)
3173 ut = speed * (e1(1)*elon(1) + e1(2)*elon(2) + e1(3)*elon(3))
3174 vt = speed * (e1(1)*elat(1) + e1(2)*elat(2) + e1(3)*elat(3))
3177 if ( po(2) < 0. )
then 3185 subroutine get_slp_obs(time, nobs, lon_obs, lat_obs, w10, mslp, slp_out, r_out, time_obs, &
3186 x_o, y_o, w10_o, slp_o, r_vor, p_vor, stime, fact)
3188 type(time_type),
intent(in):: time
3189 integer,
intent(in):: nobs
3190 real(KIND=4),
intent(in):: lon_obs(nobs)
3191 real(KIND=4),
intent(in):: lat_obs(nobs)
3192 real(KIND=4),
intent(in):: w10(nobs)
3193 real(KIND=4),
intent(in):: mslp(nobs)
3194 real(KIND=4),
intent(in):: slp_out(nobs)
3195 real(KIND=4),
intent(in):: r_out(nobs)
3196 real(KIND=4),
intent(in):: time_obs(nobs)
3197 real,
optional,
intent(in):: stime
3198 real,
optional,
intent(out):: fact
3200 real(kind=R_GRID),
intent(out):: x_o , y_o
3201 real,
intent(out):: w10_o
3202 real,
intent(out):: slp_o
3203 real,
intent(out):: r_vor, p_vor
3206 real(kind=R_GRID):: p1(2), p2(2)
3208 real(kind=R_GRID) fac
3209 integer year, month, day, hour, minute, second, n
3211 t_thresh = 600./86400.
3220 if (
present(stime) )
then 3223 call get_date(time, year, month, day, hour, minute, second)
3226 if (
master)
write(*,*)
'Warning: The year in storm track data is not the same as model year' 3230 time_model =
calday(year, month, day, hour, minute, second)
3239 if ( time_model <= (time_obs(1)-t_thresh) .or. time_model >= time_obs(nobs) )
return 3241 if ( time_model <= time_obs(1) )
then 3249 if (
present(fact) ) fact = 1.25
3252 if( time_model >= time_obs(n) .and. time_model <= time_obs(n+1) )
then 3253 fac = (time_model-time_obs(n)) / (time_obs(n+1)-time_obs(n))
3254 w10_o = w10(n) + ( w10(n+1)- w10(n)) * fac
3255 slp_o = mslp(n) + ( mslp(n+1)- mslp(n)) * fac
3260 p1(1) = lon_obs(n); p1(2) = lat_obs(n)
3261 p2(1) = lon_obs(n+1); p2(2) = lat_obs(n+1)
3262 call intp_great_circle(fac, p1, p2, x_o, y_o)
3264 if (
present(fact) ) fact = 1. + 0.25*cos(fac*2.*pi)
3279 integer:: unit, n, nobs
3280 character(len=3):: GMT
3281 character(len=9):: ts_name
3282 character(len=19):: comment
3283 integer:: mmddhh, yr, year, month, day, hour, MPH, islp
3284 integer:: it, i1, i2, p_ring, d_ring
3285 real:: lon_deg, lat_deg, cald, slp, mps
3291 x_obs(:,:) = - 100.*pi
3292 y_obs(:,:) = - 100.*pi
3298 if(
master)
write(*,*)
'No TC track file specified' 3306 if(
master)
write(*,*)
'Reading TC track data for YEAR=', year
3320 read(unit, *) ts_name, nobs, yr, month, day, hour
3322 if ( yr /= year )
then 3323 if(
master)
write(*, *)
'Year inconsistency found !!!' 3324 call mpp_error(fatal,
'==> Error in reading best track data')
3327 do while ( ts_name==
'start' )
3334 read(unit, *) lon_deg, lat_deg, mps, slp, yr, month, day, hour
3338 cald =
calday(yr, month, day, hour, 0, 0)
3340 if(
nudge_debug .and.
master)
write(*, 100) cald, month, day, hour, lon_deg, lat_deg, mps, slp
3348 read(unit, *) ts_name, nobs, yr, month, day, hour
3350 100
format(1x, f9.2, 1x, i3, 1x, i2, 1x, i2, 1x, f6.1, 1x, f6.1, 1x, f4.1, 1x, f6.1)
3354 do while ( month /= 0 )
3356 read(unit, *) month, day, hour, gmt, lat_deg, lon_deg, mph, islp, comment
3363 if(
master)
write(*, *)
'Reading data for TC#',
nstorms, comment
3365 if(
master)
write(*, *)
'End of record reached' 3368 cald =
calday(year, month, day, hour, 0, 0)
3372 if(
master)
write(*, 200) nobs, cald, month, day, hour, gmt, lat_deg, lon_deg, mph, islp, comment
3375 if ( gmt ==
'GMT' )
then 3377 if ( lon_deg < 0 )
then 3382 elseif ( gmt ==
'PAC' )
then 3394 write(*,*)
'TC vortex breeding: total storms=',
nstorms 3397 write(*, *)
'TC#',n,
' contains ',
nobs_tc(n),
' observations' 3402 200
format(i3, 1x,f9.4, 1x, i2, 1x, i2, 1x, i2, 1x, a3, f5.1, 1x, f5.1, 1x, i3, 1x, i4, 1x, a19)
3408 real function calday(year, month, day, hour, minute, sec)
3410 integer,
intent(in):: year, month, day, hour
3411 integer,
intent(in):: minute, sec
3413 integer n, m, ds, nday
3416 data days /31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31/
3420 if( month /= 1 )
then 3436 calday =
real((year-year_track_data)*nday + ds) +
real(hour*3600 + minute*60 + sec)/86400.
3443 integer,
intent(in):: ny
3448 parameter( ny00 = 0000 )
3450 if( ny >= ny00 )
then 3451 if( mod(ny,100) == 0. .and. mod(ny,400) == 0. )
then 3453 elseif( mod(ny,4) == 0. .and. mod(ny,100) /= 0. )
then 3465 subroutine pmaxmin( qname, a, imax, jmax, fac )
3467 character(len=*) qname
3472 real qmin(jmax), qmax(jmax)
3480 pmax = max(pmax, a(i,j))
3481 pmin = min(pmin, a(i,j))
3492 pmax = max(pmax, qmax(j))
3493 pmin = min(pmin, qmin(j))
3496 write(*,*) qname,
' max = ', pmax*fac,
' min = ', pmin*fac
3501 subroutine del2_uv(du, dv, cd, kmd, ntimes, bd, npx, npy, gridstruct, domain)
3502 integer,
intent(in):: kmd
3503 integer,
intent(in):: ntimes
3504 real,
intent(in):: cd
3505 type(fv_grid_bounds_type),
intent(IN) :: bd
3506 real,
intent(inout),
dimension(is:ie,js:je,kmd):: du, dv
3507 integer,
intent(IN) :: npx, npy
3508 type(fv_grid_type),
intent(IN),
target :: gridstruct
3509 type(domain2d),
intent(INOUT) :: domain
3511 real(kind=R_GRID),
pointer,
dimension(:,:,:) :: vlon, vlat
3512 real,
dimension(is:ie,js:je,kmd):: v1, v2, v3
3515 vlon => gridstruct%vlon
3516 vlat => gridstruct%vlat
3523 v1(i,j,k) = du(i,j,k)*vlon(i,j,1) + dv(i,j,k)*vlat(i,j,1)
3524 v2(i,j,k) = du(i,j,k)*vlon(i,j,2) + dv(i,j,k)*vlat(i,j,2)
3525 v3(i,j,k) = du(i,j,k)*vlon(i,j,3) + dv(i,j,k)*vlat(i,j,3)
3531 call del2_scalar( v1(
is,
js,1), cd, kmd, ntimes, bd, npx, npy, gridstruct, domain )
3532 call del2_scalar( v2(
is,
js,1), cd, kmd, ntimes, bd, npx, npy, gridstruct, domain )
3533 call del2_scalar( v3(
is,
js,1), cd, kmd, ntimes, bd, npx, npy, gridstruct, domain )
3540 du(i,j,k) = v1(i,j,k)*vlon(i,j,1) + v2(i,j,k)*vlon(i,j,2) + v3(i,j,k)*vlon(i,j,3)
3541 dv(i,j,k) = v1(i,j,k)*vlat(i,j,1) + v2(i,j,k)*vlat(i,j,2) + v3(i,j,k)*vlat(i,j,3)
3549 subroutine del2_scalar(qdt, cd, kmd, nmax, bd, npx, npy, gridstruct, domain)
3550 integer,
intent(in):: kmd
3551 integer,
intent(in):: nmax
3552 real,
intent(in):: cd
3553 type(fv_grid_bounds_type),
intent(IN) :: bd
3554 real,
intent(inout):: qdt(is:ie,js:je,kmd)
3555 integer,
intent(IN) :: npx, npy
3556 type(fv_grid_type),
intent(IN),
target :: gridstruct
3557 type(domain2d),
intent(INOUT) :: domain
3559 real:: q(isd:ied,jsd:jed,kmd)
3560 real:: fx(isd:ied+1,jsd:jed), fy(isd:ied,jsd:jed+1)
3561 integer i,j,k, n, nt, ntimes
3564 real,
pointer,
dimension(:,:) :: rarea, area
3565 real,
pointer,
dimension(:,:) :: sina_u, sina_v
3566 real,
pointer,
dimension(:,:,:) :: sin_sg
3568 real,
pointer,
dimension(:,:) :: dx, dy, rdxc, rdyc
3570 real(kind=R_GRID),
pointer :: da_min
3572 logical,
pointer :: nested, sw_corner, se_corner, nw_corner, ne_corner
3574 area => gridstruct%area
3575 rarea => gridstruct%rarea
3577 sina_u => gridstruct%sina_u
3578 sina_v => gridstruct%sina_v
3579 sin_sg => gridstruct%sin_sg
3583 rdxc => gridstruct%rdxc
3584 rdyc => gridstruct%rdyc
3586 da_min => gridstruct%da_min
3588 nested => gridstruct%nested
3589 sw_corner => gridstruct%sw_corner
3590 se_corner => gridstruct%se_corner
3591 nw_corner => gridstruct%nw_corner
3592 ne_corner => gridstruct%ne_corner
3594 ntimes = min(3, nmax)
3602 q(i,j,k) = qdt(i,j,k)
3607 call mpp_update_domains(q, domain, complete=.true.)
3620 if(nt>0)
call copy_corners(q(isd,jsd,k), npx, npy, 1, nested, bd, &
3621 sw_corner, se_corner, nw_corner, ne_corner)
3624 fx(i,j) = dy(i,j)*sina_u(i,j)*(q(i-1,j,k)-q(i,j,k))*rdxc(i,j)
3626 if (is == 1) fx(i,j) = dy(is,j)*(q(is-1,j,k)-q(is,j,k))*rdxc(is,j)* &
3627 0.5*(sin_sg(1,j,1) + sin_sg(0,j,3))
3628 if (ie+1 == npx) fx(i,j) = dy(ie+1,j)*(q(ie,j,k)-q(ie+1,j,k))*rdxc(ie+1,j)* &
3629 0.5*(sin_sg(npx,j,1) + sin_sg(npx-1,j,3))
3632 if(nt>0)
call copy_corners(q(isd,jsd,k), npx, npy, 2, nested, bd, &
3633 sw_corner, se_corner, nw_corner, ne_corner)
3635 if (j == 1 .OR. j == npy)
then 3637 fy(i,j) = dx(i,j)*(q(i,j-1,k)-q(i,j,k))*rdyc(i,j) &
3638 *0.5*(sin_sg(i,j,2) + sin_sg(i,j-1,4) )
3642 fy(i,j) = dx(i,j)*sina_v(i,j)*(q(i,j-1,k)-q(i,j,k))*rdyc(i,j)
3650 qdt(i,j,k) = q(i,j,k) + damp*rarea(i,j)*(fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))
3656 q(i,j,k) = q(i,j,k) + damp*rarea(i,j)*(fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))
3666 subroutine rmse_bias(a, rms, bias, area)
3667 real,
intent(in):: a(is:ie,js:je)
3668 real,
intent(IN) :: area(isd:ied,jsd:jed)
3669 real,
intent(out):: rms, bias
3673 total_area = 4.*pi*
radius**2
3679 bias = bias + area(i,j) * a(i,j)
3680 rms = rms + area(i,j) * a(i,j)**2
3683 call mp_reduce_sum(bias)
3684 call mp_reduce_sum(rms)
3686 bias = bias / total_area
3687 rms = sqrt( rms / total_area )
3692 subroutine corr(a, b, co, area)
3693 real,
intent(in):: a(is:ie,js:je), b(is:ie,js:je)
3694 real,
intent(in):: area(isd:ied,jsd:jed)
3695 real,
intent(out):: co
3696 real:: m_a, m_b, std_a, std_b
3700 total_area = 4.*pi*
radius**2
3703 call std(a, m_a, std_a, area)
3704 call std(b, m_b, std_b, area)
3710 co = co + area(i,j) * (a(i,j)-m_a)*(b(i,j)-m_b)
3713 call mp_reduce_sum(co)
3714 co = co / (total_area*std_a*std_b )
3718 subroutine std(a, mean, stdv, area)
3719 real,
intent(in):: a(is:ie,js:je)
3720 real,
intent(IN) :: area(isd:ied,jsd:jed)
3721 real,
intent(out):: mean, stdv
3725 total_area = 4.*pi*
radius**2
3730 mean = mean + area(i,j) * a(i,j)
3733 call mp_reduce_sum(mean)
3734 mean = mean / total_area
3739 stdv = stdv + area(i,j) * (a(i,j)-mean)**2
3742 call mp_reduce_sum(stdv)
3743 stdv = sqrt( stdv / total_area )
real(kind=4), dimension(nobs_max, max_storm) mslp_out
outer ring SLP in mb
real, dimension(:,:), allocatable gz0
logical module_is_initialized
real, dimension(:,:,:), allocatable ps_dat
logical, public do_adiabatic_init
pure real function, public vicpqd(q)
real, dimension(:), allocatable lon
subroutine, public get_var3_r4(ncid, var3_name, is, ie, js, je, ks, ke, var3, time_slice)
subroutine breed_srf_winds(time, dt, npz, u_obs, v_obs, ak, bk, ps, phis, delp, ua, va, u_dt, v_dt, pt, q, nwat, zvir, gridstruct)
The subroutine 'breed_srf_winds' performs vortex breeding by nudging 1-m winds.
character(len=128) analysis_file_first
the first NCEP analysis file to be used for nudging, by default, the first file in the "input_fname_l...
subroutine compute_slp(isc, iec, jsc, jec, tm, ps, gz, slp)
real mask_fac
[0,1] 0: no mask; 1: full strength
subroutine remap_tq(npz, ak, bk, ps, delp, t, q, kmd, ps0, ta, qa, zvir, ptop)
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
subroutine timing_off(blk_name)
The subroutine 'timing_off' stops a timer.
The type 'fv_grid_type' is made up of grid-dependent information from fv_grid_tools and fv_grid_utils...
integer km
Data z-dimension.
subroutine del2_uv(du, dv, cd, kmd, ntimes, bd, npx, npy, gridstruct, domain)
The subroutine 'del2_uv' filters the wind tendency.
subroutine, public normalize_vect(e)
The subroutine 'normalize_vect' makes 'e' a unit vector.
subroutine get_ncep_analysis(ps, u, v, t, q, zvir, ts, nfile, fname, bd)
real(kind=4), dimension(:,:,:), allocatable gz3
work array
integer, dimension(:,:), allocatable jdc
subroutine, public mappm(km, pe1, q1, kn, pe2, q2, i1, i2, iv, kord, ptop)
The subroutine 'mappm' is a general-purpose routine for remapping one set of vertical levels to anoth...
subroutine, public intp_great_circle(beta, p1, p2, x_o, y_o)
subroutine tangent_wind(elon, elat, speed, po, pp, ut, vt)
The module 'multi_gases' peforms multi constitutents computations.
subroutine std(a, mean, stdv, area)
real, dimension(:,:,:), allocatable s2c
real(kind=4), dimension(nobs_max, max_storm) time_tc
start time of the track
real(kind=4), dimension(nobs_max, max_storm) mslp_obs
observed SLP in mb
real, dimension(:), allocatable bk0
real slp_env
storm environment pressure (pa)
real pre0_env
critical storm environment pressure (pa) for size computation
real(kind=r_grid), parameter radius
real p_norelax
from P_norelax upwards, no nudging
real(kind=4), dimension(:,:,:,:), allocatable q_dat
integer, dimension(max_storm) nobs_tc
integer, parameter, public r_grid
integer jm
Data y-dimension.
real(kind=4), dimension(nobs_max, max_storm) wind_obs
observed 10-m wind speed (m/s)
integer, dimension(:,:), allocatable id2
pure real function, public virqd(q)
The module fv_nwp_nudge contains routines for nudging to input analyses. note This module is currentl...
subroutine, public get_var1_real(ncid, var1_name, im, var1, var_exist)
integer im
Data x-dimension.
The module 'fv_timing' contains FV3 timers.
pure real function, public virq(q)
subroutine get_int_hght(npz, ak, bk, ps, delp, ps0, tv)
The module 'tp_core' is a collection of routines to support FV transport.
subroutine, public fv_nwp_nudge(Time, dt, npx, npy, npz, ps_dt, u_dt, v_dt, t_dt, q_dt, zvir, ptop, ak, bk, ts, ps, delp, ua, va, pt, nwat, q, phis, gridstruct, bd, domain)
Ths subroutine 'fv_nwp_nudge' computes and returns time tendencies for nudging to analysis...
real(kind=4), dimension(:,:,:,:), allocatable v_dat
subroutine, public get_var1_double(ncid, var1_name, im, var1, var_exist)
The 'get_var' subroutines read in variables from netcdf files.
character(len=128) track_file_name
subroutine remap_coef(agrid)
The module 'fv_mapz' contains the vertical mapping routines .
subroutine del2_scalar(qdt, cd, kmd, nmax, bd, npx, npy, gridstruct, domain)
The subroutine 'del2_scalar' filters the physics tendency.
subroutine, public copy_corners(q, npx, npy, dir, nested, bd, sw_corner, se_corner, nw_corner, ne_corner)
subroutine, public breed_slp_inline(nstep, dt, npz, ak, bk, phis, pe, pk, peln, pkz, delp, u, v, pt, q, nwat, zvir, gridstruct, ks, domain_local, bd, hydrostatic)
The subroutine 'breed_slp_inline' performs vortex breeding by nudging sea level pressure toward singl...
subroutine ps_bias_correction(ps_dt, is, ie, js, je, isd, ied, jsd, jed, area)
subroutine, public fv_nwp_nudge_end
The subroutine 'fv_nwp_nudge_end' terminates the nudging module.
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
real function, public great_circle_dist(q1, q2, radius)
subroutine pmaxmin(qname, a, imax, jmax, fac)
subroutine get_slp_obs(time, nobs, lon_obs, lat_obs, w10, mslp, slp_out, r_out, time_obs, x_o, y_o, w10_o, slp_o, r_vor, p_vor, stime, fact)
The module 'sim_nc' is a netcdf file reader.
real function g0_sum(p, ifirst, ilast, jfirst, jlast, mode, reproduce, isd, ied, jsd, jed, area)
Fast version of global sum that is reproduced if the result is rounded to 4-byte. ...
real dps_min
maximum PS increment (pa; each step) due to inline breeding
real, dimension(:), allocatable lat
real, dimension(:), allocatable ak0
integer, parameter nobs_max
Max # of observations per storm.
subroutine remap_uv(npz, ak, bk, ps, delp, u, v, kmd, ps0, u0, v0, ptop)
subroutine get_obs(Time, dt, zvir, ak, bk, ps, ts, ps_obs, delp, pt, nwat, q, u_obs, v_obs, t_obs, q_obs, phis, ua, va, u_dt, v_dt, npx, npy, npz, factor, mask, ptop, bd, gridstruct, domain)
type(time_type), public fv_time
subroutine, public open_ncfile(iflnm, ncid)
subroutine, public vect_cross(e, p1, p2)
The subroutine 'vect_cross performs cross products of 3D vectors: e = P1 X P2.
subroutine timing_on(blk_name)
The subroutine 'timing_on' starts a timer.
real p_wvp
cutoff level for specific humidity nudging
real(kind=4), dimension(nobs_max, max_storm) y_obs
latitude in degrees
integer, parameter nfile_max
maximum: 20-year analysis data, 4*366*20=29280
real p_relax
from P_relax upwards, nudging is reduced linearly proportional to pfull/P_relax
subroutine corr(a, b, co, area)
character(len=128) analysis_file_last
the last NCEP analysis file to be used for nudging by default, the last file in the "input_fname_list...
integer time_interval
dataset time interval (seconds)
@ The module 'fv_diagnostics' contains routines to compute diagnosic fields.
integer nfile_total
=5 for 1-day (if datasets are 6-hr apart)
The module 'fv_grid_utils' contains routines for setting up and computing grid-related quantities...
subroutine, public prt_maxmin(qname, q, is, ie, js, je, n_g, km, fac)
real function calday(year, month, day, hour, minute, sec)
The function 'calday' performs time interpolation:
real(kind=4), dimension(:,:,:,:), allocatable t_dat
subroutine, public fv_nwp_nudge_init(time, axes, npz, zvir, ak, bk, ts, phis, gridstruct, ks, npx, neststruct, bd)
The subroutine 'fv_nwp_nudge_init' initializes the nudging module.
integer, parameter max_storm
max # of storms to process
real(kind=4), dimension(nobs_max, max_storm) rad_out
outer ring radius in meters
subroutine get_tc_mask(time, mask, agrid)
real(kind=4), dimension(nobs_max, max_storm) x_obs
longitude in degrees
subroutine, public get_ncdim1(ncid, var1_name, im)
subroutine breed_srf_w10(time, dt, npz, ak, bk, ps, phis, slp, delp, u, v, gridstruct)
The subroutine 'breed_srf_w10' performs vortex breeding by nudging 10-m winds.
character(len=128), dimension(nfile_max) file_names
integer, dimension(:,:), allocatable id1
character(len=128) input_fname_list
a file lists the input NCEP analysis data
logical function leap_year(ny)
subroutine rmse_bias(a, rms, bias, area)
real(kind=4), dimension(:,:,:,:), allocatable u_dat
subroutine ps_nudging(dt, factor, npz, ak, bk, ps_obs, mask, tm, ps, phis, delp, ua, va, pt, nwat, q, bd, npx, npy, gridstruct, domain)
subroutine, public latlon2xyz(p, e, id)
The subroutine 'latlon2xyz' maps (lon, lat) to (x,y,z)
subroutine, public close_ncfile(ncid)