FV3DYCORE  Version 2.0.0
fv_nudge.F90
Go to the documentation of this file.
1 
2 !***********************************************************************
3 !* GNU Lesser General Public License
4 !*
5 !* This file is part of the FV3 dynamical core.
6 !*
7 !* The FV3 dynamical core is free software: you can redistribute it
8 !* and/or modify it under the terms of the
9 !* GNU Lesser General Public License as published by the
10 !* Free Software Foundation, either version 3 of the License, or
11 !* (at your option) any later version.
12 !*
13 !* The FV3 dynamical core is distributed in the hope that it will be
14 !* useful, but WITHOUT ANYWARRANTY; without even the implied warranty
15 !* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
16 !* See the GNU General Public License for more details.
17 !*
18 !* You should have received a copy of the GNU Lesser General Public
19 !* License along with the FV3 dynamical core.
20 !* If not, see <http://www.gnu.org/licenses/>.
21 !***********************************************************************
22 #ifdef OVERLOAD_R4
23 #define _GET_VAR1 get_var1_real
24 #else
25 #define _GET_VAR1 get_var1_double
26 #endif
27 
31 
33 
34 ! <table>
35 ! <tr>
36 ! <th>Module Name</th>
37 ! <th>Functions Included</th>
38 ! </tr>
39 ! <tr>
40 ! <td>constants_mod</td>
41 ! <td>pi=>pi_8, grav, rdgas, cp_air, kappa, cnst_radius =>radius</td>
42 ! </tr>
43 ! <tr>
44 ! <td>external_sst_mod</td>
45 ! <td>i_sst, j_sst, sst_ncep, sst_anom, forecast_mode</td>
46 ! </tr>
47 ! <tr>
48 ! <td>diag_manager_mod</td>
49 ! <td>register_diag_field, send_data</td>
50 ! </tr>
51 ! <tr>
52 ! <td>fms_mod</td>
53 ! <td>write_version_number, open_namelist_file, check_nml_error,
54 ! file_exist, close_file</td>
55 ! </tr>
56 ! <tr>
57 ! <td>fv_arrays_mod</td>
58 ! <td>fv_grid_type, fv_grid_bounds_type, fv_nest_type, R_GRID</td>
59 ! </tr>
60 ! <tr>
61 ! <td>fv_diagnostics_mod</td>
62 ! <td>prt_maxmin, fv_time</td>
63 ! </tr>
64 ! <tr>
65 ! <td>fv_grid_utils_mod</td>
66 ! <td>great_circle_dist, intp_great_circle,
67 ! latlon2xyz, vect_cross, normalize_vect</td>
68 ! </tr>
69 ! <tr>
70 ! <td>fv_mp_mod</td>
71 ! <td>nmp_reduce_sum, mp_reduce_min, mp_reduce_max, is_master</td>
72 ! </tr>
73 ! <tr>
74 ! <td>fv_mapz_mod</td>
75 ! <td>mappm</td>
76 ! </tr>
77 ! <tr>
78 ! <td>fv_timing_mod</td>
79 ! <td>timing_on, timing_off</td>
80 ! </tr>
81 ! <tr>
82 ! <td>mpp_mod</td>
83 ! <td>mpp_error, FATAL, stdlog, get_unit, mpp_pe</td>
84 ! </tr>
85 ! <tr>
86 ! <td>sim_nc_mod</td>
87 ! <td>open_ncfile, close_ncfile, get_ncdim1, get_var1_double,
88 ! get_var3_r4, get_var1_real</td>
89 ! </tr>
90 ! <tr>
91 ! <td>mpp_domains_mod</td>
92 ! <td>mpp_update_domains, domain2d</td>
93 ! </tr>
94 ! <tr>
95 ! <td>time_manager_mod</td>
96 ! <td>time_type, get_time, get_date</td>
97 ! </tr>
98 ! <tr>
99 ! <td>tp_core_mod</td>
100 ! <td>copy_corners</td>
101 ! </tr>
102 ! </table>
103 
104  use external_sst_mod, only: i_sst, j_sst, sst_ncep, sst_anom, forecast_mode
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
109 !use fms_io_mod, only: field_size
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
113 
117  use tp_core_mod, only: copy_corners
118  use fv_mapz_mod, only: mappm
119  use fv_mp_mod, only: mp_reduce_sum, mp_reduce_min, mp_reduce_max, is_master
121 
125 #ifdef MULTI_GASES
126  use multi_gases_mod, only: virq, virqd, vicpqd, num_gas
127 #endif
128 
129  implicit none
130  private
131 
132  real(kind=R_GRID), parameter :: radius = cnst_radius
133 
134 ! version number of this module
135 ! Include variable "version" to be written to log file.
136 #include<file_version.h>
137 
138  logical :: do_adiabatic_init
139 
141  public do_adiabatic_init
142  integer im
143  integer jm
144  integer km
145  real, allocatable:: ak0(:), bk0(:)
146  real, allocatable:: lat(:), lon(:)
147 
148  logical :: module_is_initialized = .false.
149  logical :: master
150  logical :: no_obs
151  real :: deg2rad, rad2deg
152  real :: time_nudge = 0.
153  integer :: time_interval = 6*3600
154 ! ---> h1g, enhance the max. analysis data files, 2012-10-22
155 ! integer, parameter :: nfile_max = 125
156  integer, parameter :: nfile_max = 29280
157 ! <--- h1g, 2012-10-22
158  integer :: nfile
159 
160  integer :: k_breed = 0
161  integer :: k_trop = 0
162  real :: p_trop = 950.e2
163  real :: dps_min = 50.
164  real :: del2_cd = 0.16
165 
166  real, allocatable:: s2c(:,:,:)
167  integer, allocatable:: id1(:,:), id2(:,:), jdc(:,:)
168  real, allocatable :: ps_dat(:,:,:)
169  real(KIND=4), allocatable, dimension(:,:,:,:):: u_dat, v_dat, t_dat, q_dat
170  real(KIND=4), allocatable, dimension(:,:,:):: gz3
171  real, allocatable:: gz0(:,:)
172 
173 ! Namelist variables:
174 ! ---> h1g, add the list of input NCEP analysis data files, 2012-10-22
175  character(len=128):: input_fname_list =""
176  character(len=128):: analysis_file_first =""
178  character(len=128):: analysis_file_last=""
180 
181  real :: p_relax = 30.e2
183 
184  real :: p_norelax = 0.0
185 ! <--- h1g, 2012-10-22
186 
187  character(len=128):: file_names(nfile_max)
188  character(len=128):: track_file_name
189  integer :: nfile_total = 0
190  real :: p_wvp = 100.e2
191  integer :: kord_data = 8
192 
193  real :: mask_fac = 0.25
194 
195  logical :: t_is_tv = .false.
196  logical :: use_pt_inc = .false.
197  logical :: use_high_top = .false.
198  logical :: add_bg_wind = .true.
199  logical :: conserve_mom = .true.
200  logical :: conserve_hgt = .true.
201  logical :: tc_mask = .false.
202  logical :: strong_mask = .false.
203  logical :: ibtrack = .true.
204  logical :: nudge_debug = .false.
205  logical :: do_ps_bias = .false.
206  logical :: nudge_ps = .false.
207  logical :: nudge_q = .false.
208  logical :: nudge_winds = .true.
209  logical :: nudge_virt = .true.
210  logical :: nudge_hght = .true.
211  logical :: time_varying = .true.
212  logical :: print_end_breed = .true.
213  logical :: print_end_nudge = .true.
214 
215 !Nudging time-scales (seconds):
216 !note, however, the effective time-scale is 2X smaller (stronger) due
217 ! to the use of the time-varying weighting factor
218  real :: tau_ps = 21600.
219  real :: tau_q = 86400.
220  real :: tau_winds = 21600.
221  real :: tau_virt = 43200.
222  real :: tau_hght = 43200.
223 
224  real :: q_min = 1.e-8
225 
226  integer :: jbeg, jend
227  integer :: nf_uv = 0
228  integer :: nf_ps = 0
229  integer :: nf_t = 2
230  integer :: nf_ht = 1
231 
232 ! starting layer (top layer is sponge layer and is skipped)
233  integer :: kstart = 2
234 
235 ! skip "kbot" layers
236  integer :: kbot_winds = 0
237  integer :: kbot_t = 0
238  integer :: kbot_q = 0
239  logical :: analysis_time
240 
241 !-- Tropical cyclones --------------------------------------------------------------------
242 
243 ! track dataset: 'INPUT/tropical_cyclones.txt'
244 
245  logical :: breed_srf_w = .false.
246  real :: grid_size = 28.e3
247  real :: tau_vt_slp = 1200.
248  real :: tau_vt_wind = 1200.
249  real :: tau_vt_rad = 4.0
250 
251  real :: pt_lim = 0.2
252  real :: slp_env = 101010.
253  real :: pre0_env = 100000.
254  real, parameter:: tm_max = 315.
255 !------------------
256  real:: r_lo = 2.0
257  real:: r_hi = 5.0 ! try 4.0?
258 !------------------
259  real:: r_fac = 1.2
260  real :: r_min = 200.e3
261  real :: r_inc = 25.e3
262  real, parameter:: del_r = 50.e3
263  real:: elapsed_time = 0.0
264  real:: nudged_time = 1.e12
265  ! usage example: set to 43200. to do inline vortex breeding
266  ! for only the first 12 hours
267  ! In addition, specify only 3 analysis files (12 hours)
268  integer:: year_track_data
269  integer, parameter:: max_storm = 140
270  integer, parameter:: nobs_max = 125
271 
272  integer :: nstorms = 0
273  integer :: nobs_tc(max_storm)
274  integer :: min_nobs = 16
275  real :: min_mslp = 1009.e2
276  real(KIND=4):: x_obs(nobs_max,max_storm)
277  real(KIND=4):: y_obs(nobs_max,max_storm)
278  real(KIND=4):: wind_obs(nobs_max,max_storm)
279  real(KIND=4):: mslp_obs(nobs_max,max_storm)
280  real(KIND=4):: mslp_out(nobs_max,max_storm)
281  real(KIND=4):: rad_out(nobs_max,max_storm)
282  real(KIND=4):: time_tc(nobs_max,max_storm)
283 !------------------------------------------------------------------------------------------
284  integer :: id_ht_err
285 
286  !!! CLEANUP Module pointers
287  integer :: is, ie, js, je
288  integer :: isd, ied, jsd, jed
289 
290  namelist /fv_nwp_nudge_nml/t_is_tv, nudge_ps, nudge_virt, nudge_hght, nudge_q, nudge_winds, &
297  input_fname_list, analysis_file_first, analysis_file_last, p_relax, p_norelax !h1g, add 3 namelist variables, 2012-20-22
298 
299  contains
300 
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, &
305  bd, domain )
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
314  type(fv_grid_bounds_type), intent(IN) :: bd
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
318 ! pt as input is true tempeature
319 #ifdef MULTI_GASES
320  real, intent(inout):: q(isd:ied,jsd:jed,npz,*)
321 #else
322  real, intent(inout):: q(isd:ied,jsd:jed,npz,nwat)
323 #endif
324  real, intent(inout), dimension(isd:ied,jsd:jed):: ps
325 ! Accumulated tendencies
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
329 
330  type(fv_grid_type), intent(INOUT), target :: gridstruct
331 ! local:
332  real:: h2(is:ie,js:je)
333  real:: m_err(is:ie,js:je)
334  real:: slp_n(is:ie,js:je)
335  real:: slp_m(is:ie,js:je)
336  real:: mask(is:ie,js:je)
337  real:: tv(is:ie,js:je)
338  real:: peln(is:ie,npz+1)
339  real:: pe2(is:ie, npz+1)
340  real:: ptmp
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
344  real:: ps_fac(is:ie,js:je)
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
349  logical used
350 
351 
352  real(kind=R_GRID), pointer, dimension(:,:,:) :: agrid
353  real, pointer, dimension(:,:) :: rarea, area
354 
355  real, pointer, dimension(:,:) :: sina_u, sina_v
356  real, pointer, dimension(:,:,:) :: sin_sg
357  real(kind=R_GRID), pointer, dimension(:,:,:) :: vlon, vlat
358 
359  real, pointer, dimension(:,:) :: dx, dy, rdxc, rdyc
360 
361  real(kind=R_GRID), pointer :: da_min
362 
363  logical, pointer :: bounded_domain, sw_corner, se_corner, nw_corner, ne_corner
364 
365  if ( .not. module_is_initialized ) then
366  call mpp_error(fatal,'==> Error from fv_nwp_nudge: module not initialized')
367  endif
368  agrid => gridstruct%agrid_64
369  area => gridstruct%area
370  rarea => gridstruct%rarea
371 
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
377 
378  dx => gridstruct%dx
379  dy => gridstruct%dy
380  rdxc => gridstruct%rdxc
381  rdyc => gridstruct%rdyc
382 
383  da_min => gridstruct%da_min
384 
385 
386  sw_corner => gridstruct%sw_corner
387  se_corner => gridstruct%se_corner
388  nw_corner => gridstruct%nw_corner
389  ne_corner => gridstruct%ne_corner
390 
391  if ( no_obs ) then
392 #ifndef DYCORE_SOLO
393  forecast_mode = .true.
394 #endif
395  return
396  endif
397 
398  call get_time (time, seconds, days)
399 
400  do j=js,je
401  do i=is,ie
402  mask(i,j) = 1.
403  enddo
404  enddo
405  if ( tc_mask ) call get_tc_mask(time, mask, agrid)
406 
407 ! The following profile is suitable only for nwp purposes; if the analysis has a good representation
408 ! of the strat-meso-sphere the profile for upper layers should be changed.
409 
410  profile(:) = 1.
411 
412 !$OMP parallel do default(none) shared(npz,press,ak,bk,P_relax,P_norelax,profile)
413  do k=1,npz
414  press(k) = 0.5*(ak(k) + ak(k+1)) + 0.5*(bk(k)+bk(k+1))*1.e5
415  if ( press(k) < p_relax ) then
416  profile(k) = max(0.01, press(k)/p_relax)
417  endif
418 
419  ! above P_norelax, no nudging. added by h1g
420  if( press(k) < p_norelax ) profile(k) = 0.0
421  enddo
422  profile(1) = 0.
423 
424 ! Thermodynamics:
425  prof_t(:) = 1.
426 !$OMP parallel do default(none) shared(npz,press,prof_t)
427  do k=1,npz
428  if ( press(k) < 10.e2 ) then
429  prof_t(k) = max(0.01, press(k)/10.e2)
430  endif
431  enddo
432  prof_t(1) = 0.
433 
434 ! Water vapor:
435  prof_q(:) = 1.
436 !$OMP parallel do default(none) shared(npz,press,prof_q)
437  do k=1,npz
438  if ( press(k) < 300.e2 ) then
439  prof_q(k) = max(0., press(k)/300.e2)
440  endif
441  enddo
442  prof_q(1) = 0.
443 
444 ! Height
445  if ( k_trop == 0 ) then
446  k_trop = 2
447  do k=2,npz-1
448  ptmp = ak(k+1) + bk(k+1)*1.e5
449  if ( ptmp > p_trop ) then
450  k_trop = k
451  exit
452  endif
453  enddo
454  endif
455 
456  if ( nudge_virt .and. nudge_hght ) then
457  kht = k_trop
458  else
459  kht = npz-kbot_t
460  endif
461 
462  if ( time_varying ) then
463  factor = 1. + cos(real(mod(seconds,time_interval))/real(time_interval)*2.*pi)
464  factor = max(1.e-5, factor)
465  else
466  factor = 1.
467  endif
468 
469  if ( do_adiabatic_init ) factor = 2.*factor
470 
471  allocate (ps_obs(is:ie,js:je) )
472  allocate ( t_obs(is:ie,js:je,npz) )
473  allocate ( q_obs(is:ie,js:je,npz) )
474 
475  if ( nudge_winds ) then
476  allocate (u_obs(is:ie,js:je,npz) )
477  allocate (v_obs(is:ie,js:je,npz) )
478  endif
479 
480 
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)
483 ! *t_obs* is virtual temperature
484 
485  if ( no_obs ) then
486  deallocate (ps_obs)
487  deallocate (t_obs)
488  deallocate (q_obs)
489  if ( nudge_winds ) then
490  deallocate (u_obs)
491  deallocate (v_obs)
492  endif
493 #ifndef DYCORE_SOLO
494  forecast_mode = .true.
495 #endif
496  return
497  endif
498 
499  do j=js,je
500  do i=is,ie
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))
503  else
504  ps_fac(i,j) = 1.
505  endif
506  enddo
507  enddo
508 
509  if( analysis_time ) then
510 !-------------------------------------------
511 ! Compute RMSE, bias, and correlation of SLP
512 !-------------------------------------------
513  do j=js,je
514  do i=is,ie
515 #ifdef MULTI_GASES
516  tv(i,j) = pt(i,j,npz)*virq(q(i,j,npz,1:num_gas))
517 #else
518  tv(i,j) = pt(i,j,npz)*(1.+zvir*q(i,j,npz,1))
519 #endif
520  enddo
521  enddo
522  call compute_slp(is, ie, js, je, tv, ps(is:ie,js:je), phis(is:ie,js:je), slp_m)
523  call compute_slp(is, ie, js, je, t_obs(is,js,npz), ps_obs, gz0, slp_n)
524 
525  if ( nudge_debug ) then
526  if(master) write(*,*) 'kht=', kht
527  call prt_maxmin('PS_o', ps_obs, is, ie, js, je, 0, 1, 0.01)
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
530  call prt_maxmin('SLP_m', slp_m, is, ie, js, je, 0, 1, 0.01)
531  call prt_maxmin('SLP_o', slp_n, is, ie, js, je, 0, 1, 0.01)
532  endif
533 
534 !$OMP parallel do default(none) shared(is,ie,js,je,phis,gz0,m_err,mask,slp_m,slp_n)
535  do j=js,je
536  do i=is,ie
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))
539  else
540  m_err(i,j) = mask(i,j)*(slp_m(i,j) - slp_n(i,j))
541  endif
542  enddo
543  enddo
544 
545  call rmse_bias(m_err, rms, bias, area)
546  call corr(slp_m, slp_n, co, area)
547 
548  call prt_maxmin('SLP Error (mb)=', m_err, is, ie, js, je, 0, 1, 0.01)
549  if(master) write(*,*) 'SLP (Pa): RMS=', rms, ' Bias=', bias
550  if(master) write(*,*) 'SLP correlation=',co
551  endif
552 
553  if ( nudge_winds ) then
554 
555  allocate (du_obs(is:ie,js:je,npz) )
556  allocate (dv_obs(is:ie,js:je,npz) )
557 
558  du_obs = 0.
559  dv_obs = 0.
560 
561 ! Compute tendencies:
562  rdt = 1. / (tau_winds/factor + dt)
563 !$OMP parallel do default(none) shared(kstart,npz,kbot_winds,is,ie,js,je,du_obs,dv_obs, &
564 !$OMP profile,u_obs,v_obs,ua,va,rdt)
565  do k=kstart, npz - kbot_winds
566  do j=js,je
567  do i=is,ie
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
570  enddo
571  enddo
572  enddo
573 
574  if ( nf_uv>0 ) call del2_uv(du_obs, dv_obs, del2_cd, npz, nf_uv, bd, npx, npy, gridstruct, domain)
575 
576 !$OMP parallel do default(none) shared(kstart,kbot_winds,npz,is,ie,js,je,du_obs,dv_obs, &
577 !$OMP mask,ps_fac,u_dt,v_dt,ua,va,dt)
578  do k=kstart, npz - kbot_winds
579  if ( k==npz ) then
580  do j=js,je
581  do i=is,ie
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)
584 !
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
589  enddo
590  enddo
591  else
592  do j=js,je
593  do i=is,ie
594 ! Apply TC mask
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)
597 !
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
602  enddo
603  enddo
604  endif
605  enddo
606  endif
607 
608 !$OMP parallel do default(none) shared(is,ie,js,je,npz,t_dt)
609  do k=1,npz
610  do j=js,je
611  do i=is,ie
612  t_dt(i,j,k) = 0.
613  enddo
614  enddo
615  enddo
616 
617  if ( nudge_virt ) then
618  rdt = 1./(tau_virt/factor + dt)
619 !$OMP parallel do default(none) shared(is,ie,js,je,npz,kstart,kht,t_dt,prof_t,t_obs,zvir, &
620 #ifdef MULTI_GASES
621 !$OMP num_gas, &
622 #endif
623 !$OMP q,pt,rdt,ps_fac)
624  do k=kstart, kht
625  if ( k==npz ) then
626  do j=js,je
627  do i=is,ie
628 #ifdef MULTI_GASES
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)
630 #else
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)
632 #endif
633  enddo
634  enddo
635  else
636  do j=js,je
637  do i=is,ie
638 #ifdef MULTI_GASES
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
640 #else
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
642 #endif
643  enddo
644  enddo
645  endif
646  enddo
647  endif
648 
649  if ( nudge_hght .and. kht<npz ) then ! averaged (in log-p) temperature
650  rdt = 1. / (tau_hght/factor + dt)
651 !$OMP parallel do default(none) shared(is,ie,js,je,npz,ak,h2,delp,kht,t_obs,pt,zvir, &
652 #ifdef MULTI_GASES
653 !$OMP num_gas, &
654 #endif
655 !$OMP q,rdt,ps_fac,mask,t_dt) &
656 !$OMP private(pe2, peln )
657  do j=js,je
658  do i=is,ie
659  pe2(i,1) = ak(1)
660  h2(i,j) = 0.
661  enddo
662  do k=1, npz
663  do i=is,ie
664  pe2(i,k+1) = pe2(i,k) + delp(i,j,k)
665  enddo
666  enddo
667  do k=kht+1, npz+1
668  do i=is,ie
669  peln(i,k) = log(pe2(i,k))
670  enddo
671  enddo
672 
673  do k=kht+1, npz
674  do i=is,ie
675 ! Difference between "Mean virtual tempearture" (pseudo height)
676 #ifdef MULTI_GASES
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))
678 #else
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))
680 #endif
681  enddo
682  enddo
683  do i=is,ie
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)
686  enddo
687 
688  do k=kht+1, npz
689  do i=is,ie
690 #ifdef MULTI_GASES
691  t_dt(i,j,k) = h2(i,j) / virq(q(i,j,k,1:num_gas))
692 #else
693  t_dt(i,j,k) = h2(i,j) / (1.+zvir*q(i,j,k,1))
694 #endif
695  enddo
696  enddo
697  enddo ! j-loop
698  if(nudge_debug) call prt_maxmin('H2 increment=', h2, is, ie, js, je, 0, 1, dt)
699  endif
700 
701  if ( nudge_virt .or. nudge_hght ) then
702  if ( nf_t>0 ) call del2_scalar(t_dt, del2_cd, npz, nf_t, bd, npx, npy, gridstruct, domain)
703 !$OMP parallel do default(none) shared(kstart,npz,is,ie,js,je,pt,t_dt,dt,mask)
704  do k=kstart, npz
705  do j=js,je
706  do i=is,ie
707  pt(i,j,k) = pt(i,j,k) + t_dt(i,j,k)*dt*mask(i,j)
708  enddo
709  enddo
710  enddo
711  endif
712 
713  q_dt(:,:,:) = 0.
714  if ( nudge_q ) then
715  rdt = 1./(tau_q/factor + dt)
716 !$OMP parallel do default(none) shared(kstart,npz,kbot_q,is,ie,js,je,press,p_wvp,nwat, &
717 !$OMP q,delp,q_dt,prof_q,q_min,q_obs,rdt,mask,dt)
718  do k=kstart, npz - kbot_q
719  if ( press(k) > p_wvp ) then
720  do iq=2,nwat
721  do j=js,je
722  do i=is,ie
723  q(i,j,k,iq) = q(i,j,k,iq)*delp(i,j,k)
724  enddo
725  enddo
726  enddo
727 ! Specific humidity:
728  do j=js,je
729  do i=is,ie
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))
734  enddo
735  enddo
736  do iq=2,nwat
737  do j=js,je
738  do i=is,ie
739  q(i,j,k,iq) = q(i,j,k,iq)/delp(i,j,k)
740  enddo
741  enddo
742  enddo
743  endif
744  enddo
745  endif
746 
747  ps_dt(:,:) = 0.
748 
749  deallocate ( t_obs )
750  deallocate ( q_obs )
751  deallocate ( ps_obs )
752 
753  if ( breed_srf_w .and. nudge_winds ) &
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)
755 
756  if ( nudge_debug) then
757  call prt_maxmin('T increment=', t_dt, is, ie, js, je, 0, npz, dt)
758  call prt_maxmin('U increment=', du_obs, is, ie, js, je, 0, npz, dt)
759  call prt_maxmin('V increment=', dv_obs, is, ie, js, je, 0, npz, dt)
760  endif
761 
762  if ( nudge_winds ) then
763  deallocate ( u_obs )
764  deallocate ( v_obs )
765  deallocate (du_obs)
766  deallocate (dv_obs)
767  endif
768 
769  nullify(agrid)
770  nullify(area)
771  nullify(rarea)
772 
773  nullify(vlon)
774  nullify(vlat)
775  nullify(sina_u)
776  nullify(sina_v)
777  nullify(sin_sg)
778 
779  nullify(dx)
780  nullify(dy)
781  nullify(rdxc)
782  nullify(rdyc)
783 
784  nullify(da_min)
785 
786  nullify(sw_corner)
787  nullify(se_corner)
788  nullify(nw_corner)
789  nullify(ne_corner)
790 
791  end subroutine fv_nwp_nudge
792 
793 
794  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)
795 ! Input
796  real, intent(in):: dt, factor
797  integer, intent(in):: npz, nwat, npx, npy
798  real, intent(in), dimension(npz+1):: ak, bk
799  type(fv_grid_bounds_type), intent(IN) :: bd
800  real, intent(in):: phis(isd:ied,jsd:jed)
801  real, intent(in), dimension(is:ie,js:je):: ps_obs, mask, tm
802  type(fv_grid_type), intent(IN), target :: gridstruct
803  type(domain2d), intent(INOUT) :: domain
804 ! Input/Output
805  real, intent(inout), dimension(isd:ied,jsd:jed):: ps
806  real, intent(inout), dimension(isd:ied,jsd:jed,npz):: delp, pt, ua, va
807  real, intent(inout):: q(isd:ied,jsd:jed,npz,nwat)
808 ! local
809  real, dimension(is:ie,js:je):: ps_dt
810  integer, parameter:: kmax = 100
811  real:: pn0(kmax+1), pk0(kmax+1), pe0(kmax+1)
812  real, dimension(is:ie,npz+1):: pe2, peln
813  real:: pst, dbk, pt0, rdt, bias
814  integer i, j, k, iq
815 
816  real, pointer, dimension(:,:) :: area
817 #ifdef MULTI_GASES
818  real :: kappax(km)
819  real :: pkx
820 #endif
821 
822  area => gridstruct%area
823 
824 ! Adjust interpolated ps to model terrain
825  if ( kmax < km ) call mpp_error(fatal,'==> KMAX must be larger than km')
826 
827  do j=js,je
828  do 666 i=is,ie
829 #ifdef MULTI_GASES
830  do k=1,km
831  kappax(k)= virqd(q(i,j,k,1:num_gas))/vicpqd(q(i,j,k,1:num_gas))
832  enddo
833 #endif
834  do k=1, km+1
835  pe0(k) = ak0(k) + bk0(k)*ps_obs(i,j)
836  pk0(k) = (ak0(k) + bk0(k)*ps_obs(i,j))**kappa
837  enddo
838  if( phis(i,j)>gz0(i,j) ) then
839  do k=km,1,-1
840  if( phis(i,j)<gz3(i,j,k) .and. phis(i,j) >= gz3(i,j,k+1) ) then
841  pst = pk0(k) + (pk0(k+1)-pk0(k))*(gz3(i,j,k)-phis(i,j))/(gz3(i,j,k)-gz3(i,j,k+1))
842  go to 666
843  endif
844  enddo
845  else
846  pn0(km ) = log(ak0(km) + bk0(km)*ps_obs(i,j))
847  pn0(km+1) = log(ps_obs(i,j))
848 ! Extrapolation into the ground using only the lowest layer potential temp
849 #ifdef MULTI_GASES
850  pkx = (pk0(km+1)-pk0(km))/(kappa*(pn0(km+1)-pn0(km)))
851  pt0 = tm(i,j)/exp(kappax(km)*log(pkx))
852  pkx = exp((kappax(km)-1.0)*log(pkx))
853  pst = pk0(km+1) + (gz0(i,j)-phis(i,j))/(cp_air*pt0*pkx)
854 #else
855  pt0 = tm(i,j)/(pk0(km+1)-pk0(km))*(kappa*(pn0(km+1)-pn0(km)))
856  pst = pk0(km+1) + (gz0(i,j)-phis(i,j))/(cp_air*pt0)
857 #endif
858  endif
859 #ifdef MULTI_GASES
860  k=km
861 666 ps_dt(i,j) = pst**(1./(kappa*kappax(k))) - ps(i,j)
862 #else
863 666 ps_dt(i,j) = pst**(1./kappa) - ps(i,j)
864 #endif
865  enddo ! j-loop
866 
867  if( nf_ps>0 ) call del2_scalar(ps_dt, del2_cd, 1, nf_ps, bd, npx, npy, gridstruct, domain)
868 
869  do j=js,je
870  do i=is,ie
871 ! Cap the errors:
872  ps_dt(i,j) = sign( min(10.e2, abs(ps_dt(i,j))), ps_dt(i,j) )
873  ps_dt(i,j) = mask(i,j)*ps_dt(i,j) * max( 0., 1.-abs(gz0(i,j)-phis(i,j))/(grav*500.) )
874  enddo
875  enddo
876 
877  if( do_ps_bias ) call ps_bias_correction( ps_dt, is, ie, js, je, isd, ied, jsd, jed, area)
878 
879 #ifdef CONSV_HGHT
880 ! Convert tracer moist mixing ratio to mass
881  do iq=1,nwat
882  do k=1,npz
883  do j=js,je
884  do i=is,ie
885  q(i,j,k,iq) = q(i,j,k,iq) * delp(i,j,k)
886  enddo
887  enddo
888  enddo
889  enddo
890 
891  do k=1,npz
892  do j=js,je
893  do i=is,ie
894  ua(i,j,k) = ua(i,j,k) * delp(i,j,k)
895  va(i,j,k) = va(i,j,k) * delp(i,j,k)
896  enddo
897  enddo
898  enddo
899 
900  do j=js,je
901  do i=is,ie
902  pe2(i,1) = ak(1)
903  peln(i,1) = log(pe2(i,1))
904  enddo
905  do k=2,npz+1
906  do i=is,ie
907  pe2(i,k) = pe2(i,k-1) + delp(i,j,k-1)
908  peln(i,k) = log(pe2(i,k))
909  enddo
910  enddo
911  do k=1,npz
912  do i=is,ie
913  pt(i,j,k) = pt(i,j,k)*(peln(i,k+1)-peln(i,k))
914  enddo
915  enddo
916  enddo
917 #endif
918 
919 ! Update ps:
920  do j=js,je
921  do i=is,ie
922  ps(i,j) = ak(1)
923  enddo
924  enddo
925 
926  rdt = dt / (tau_ps/factor + dt)
927  do k=1,npz
928  dbk = rdt*(bk(k+1) - bk(k))
929  do j=js,je
930  do i=is,ie
931  delp(i,j,k) = delp(i,j,k) + dbk*ps_dt(i,j)
932  ps(i,j) = delp(i,j,k) + ps(i,j)
933  enddo
934  enddo
935  enddo
936 
937 #ifdef CONSV_HGHT
938  do j=js,je
939  do k=2,npz+1
940  do i=is,ie
941  pe2(i,k) = pe2(i,k-1) + delp(i,j,k-1)
942  peln(i,k) = log(pe2(i,k))
943  enddo
944  enddo
945  do k=1,npz
946  do i=is,ie
947  pt(i,j,k) = pt(i,j,k)/(peln(i,k+1)-peln(i,k))
948  enddo
949  enddo
950  enddo
951 
952 ! Convert tracer density back to moist mixing ratios
953  do iq=1,nwat
954  do k=1,npz
955  do j=js,je
956  do i=is,ie
957  q(i,j,k,iq) = q(i,j,k,iq) / delp(i,j,k)
958  enddo
959  enddo
960  enddo
961  enddo
962 
963  do k=1,npz
964  do j=js,je
965  do i=is,ie
966  ua(i,j,k) = ua(i,j,k) / delp(i,j,k)
967  va(i,j,k) = va(i,j,k) / delp(i,j,k)
968  enddo
969  enddo
970  enddo
971 #endif
972 
973  end subroutine ps_nudging
974 
975  subroutine ps_bias_correction ( ps_dt, is, ie, js, je, isd, ied, jsd, jed, area )
976  integer, intent(IN) :: is, ie, js, je
977  integer, intent(IN) :: isd, ied, jsd,jed
978  real, intent(inout):: ps_dt(is:ie,js:je)
979  real, intent(IN), dimension(isd:ied,jsd:jed) :: area
980 !
981  real:: esl, total_area
982  real:: bias, psum
983  integer:: i, j
984 
985  total_area = 4.*pi*radius**2
986  esl = 0.01 ! Pascal
987 
988  bias = g0_sum(ps_dt, is, ie, js, je, 1, .true., isd, ied, jsd, jed, area)
989 
990  if ( abs(bias) < esl ) then
991  if(master .and. nudge_debug) write(*,*) 'Very small PS bias=', -bias, ' No bais adjustment'
992  return
993  else
994  if(master .and. nudge_debug) write(*,*) 'Significant PS bias=', -bias
995  endif
996 
997  if ( bias > 0. ) then
998  psum = 0.
999  do j=js,je
1000  do i=is,ie
1001  if ( ps_dt(i,j) > 0. ) then
1002  psum = psum + area(i,j)
1003  endif
1004  enddo
1005  enddo
1006 
1007  call mp_reduce_sum( psum )
1008  bias = bias * total_area / psum
1009 
1010  do j=js,je
1011  do i=is,ie
1012  if ( ps_dt(i,j) > 0.0 ) then
1013  ps_dt(i,j) = max(0.0, ps_dt(i,j) - bias)
1014  endif
1015  enddo
1016  enddo
1017  else
1018  psum = 0.
1019  do j=js,je
1020  do i=is,ie
1021  if ( ps_dt(i,j) < 0. ) then
1022  psum = psum + area(i,j)
1023  endif
1024  enddo
1025  enddo
1026 
1027  call mp_reduce_sum( psum )
1028  bias = bias * total_area / psum
1029 
1030  do j=js,je
1031  do i=is,ie
1032  if ( ps_dt(i,j) < 0.0 ) then
1033  ps_dt(i,j) = min(0.0, ps_dt(i,j) - bias)
1034  endif
1035  enddo
1036  enddo
1037  endif
1038 
1039  end subroutine ps_bias_correction
1040 
1043  real function g0_sum(p, ifirst, ilast, jfirst, jlast, mode, reproduce, isd, ied, jsd, jed, area)
1044 ! Fast version of global sum; reproduced if result rounded to 4-byte
1045  integer, intent(IN) :: ifirst, ilast
1046  integer, intent(IN) :: jfirst, jlast
1047  integer, intent(IN) :: mode
1048  logical, intent(IN) :: reproduce
1049  real, intent(IN) :: p(ifirst:ilast,jfirst:jlast)
1050  integer, intent(IN) :: isd, ied, jsd, jed
1051  real, intent(IN) :: area(isd:ied,jsd:jed)
1052 
1053  integer :: i,j
1054  real gsum
1055 
1056  gsum = 0.
1057  do j=jfirst,jlast
1058  do i=ifirst,ilast
1059  gsum = gsum + p(i,j)*area(i,j)
1060  enddo
1061  enddo
1062 ! call mpp_sum(gsum) ! does this work?
1063  call mp_reduce_sum(gsum)
1064 
1065  if ( mode==1 ) then
1066  g0_sum = gsum / (4.*pi*radius**2)
1067  else
1068  g0_sum = gsum
1069  endif
1070 
1071  if ( reproduce ) g0_sum = real(g0_sum, 4) ! convert to 4-byte real
1072 
1073  end function g0_sum
1074 
1075 
1076  subroutine compute_slp(isc, iec, jsc, jec, tm, ps, gz, slp)
1077  integer, intent(in):: isc, iec, jsc, jec
1078  real, intent(in), dimension(isc:iec,jsc:jec):: tm, ps, gz
1079 ! tm: virtual temperature required as input
1080  real, intent(out):: slp(isc:iec,jsc:jec)
1081  integer:: i,j
1082 
1083  do j=jsc,jec
1084  do i=isc,iec
1085  slp(i,j) = ps(i,j) * exp( gz(i,j)/(rdgas*(tm(i,j) + 3.25e-3*gz(i,j)/grav)) )
1086  enddo
1087  enddo
1088 
1089  end subroutine compute_slp
1090 
1091 
1092  subroutine get_obs(Time, dt, zvir, ak, bk, ps, ts, ps_obs, delp, pt, nwat, q, u_obs, v_obs, t_obs, q_obs, &
1093  phis, ua, va, u_dt, v_dt, npx, npy, npz, factor, mask, ptop, bd, gridstruct, domain)
1094  type(time_type), intent(in):: Time
1095  integer, intent(in):: npz, nwat, npx, npy
1096  real, intent(in):: zvir, ptop
1097  real, intent(in):: dt, factor
1098  real, intent(in), dimension(npz+1):: ak, bk
1099  type(fv_grid_bounds_type), intent(IN) :: bd
1100  real, intent(in), dimension(isd:ied,jsd:jed):: phis
1101  real, intent(in), dimension(is:ie,js:je):: mask
1102  real, intent(inout), dimension(isd:ied,jsd:jed):: ps
1103  real, intent(inout), dimension(isd:ied,jsd:jed,npz):: delp, pt, ua, va, u_dt, v_dt
1104  real, intent(inout)::q(isd:ied,jsd:jed,npz,*)
1105  real, intent(out), dimension(is:ie,js:je):: ts, ps_obs
1106  real, intent(out), dimension(is:ie,js:je,npz):: u_obs, v_obs, t_obs, q_obs
1107  type(fv_grid_type), intent(IN) :: gridstruct
1108  type(domain2d), intent(INOUT) :: domain
1109 ! local:
1110  real:: tm(is:ie,js:je)
1111  real(KIND=4), allocatable,dimension(:,:,:):: ut, vt, wt
1112  real, allocatable,dimension(:,:,:):: uu, vv
1113  integer :: seconds, days
1114  integer :: i,j,k
1115  real :: alpha, beta
1116 
1117  call get_time (time, seconds, days)
1118 
1119  if ( do_adiabatic_init ) goto 333
1120  seconds = seconds - nint(dt)
1121 
1122 ! Data must be "time_interval" (hr) apart; keep two time levels in memory
1123 
1124  no_obs = .false.
1125  analysis_time = .false.
1126 
1127  if ( mod(seconds, time_interval) == 0 ) then
1128 
1129  if ( nfile == nfile_total ) then
1130  no_obs = .true.
1131 #ifndef DYCORE_SOLO
1132  forecast_mode = .true.
1133 #endif
1134  if(print_end_nudge) then
1135  print_end_nudge = .false.
1136  if (master) write(*,*) '*** L-S nudging Ended at', days, seconds
1137  endif
1138  return ! free-running mode
1139  endif
1140 
1141  ps_dat(:,:,1) = ps_dat(:,:,2)
1142  if ( nudge_winds ) then
1143  u_dat(:,:,:,1) = u_dat(:,:,:,2)
1144  v_dat(:,:,:,1) = v_dat(:,:,:,2)
1145  endif
1146  t_dat(:,:,:,1) = t_dat(:,:,:,2)
1147  q_dat(:,:,:,1) = q_dat(:,:,:,2)
1148 
1149 !---------------
1150 ! Read next data
1151 !---------------
1152  nfile = nfile + 1
1153  call get_ncep_analysis ( ps_dat(:,:,2), u_dat(:,:,:,2), v_dat(:,:,:,2), &
1154  t_dat(:,:,:,2), q_dat(:,:,:,2), zvir, &
1155  ts, nfile, file_names(nfile), bd )
1156  time_nudge = dt
1157  else
1158  time_nudge = time_nudge + dt
1159  endif
1160 
1161 !--------------------
1162 ! Time interpolation:
1163 !--------------------
1164 
1165  beta = time_nudge / real(time_interval)
1166 
1167  if ( beta < 0. .or. beta > (1.+1.e-7) ) then
1168  if (master) write(*,*) 'Nudging: beta=', beta
1169  call mpp_error(fatal,'==> Error from get_obs:data out of range')
1170  endif
1171 
1172 333 continue
1173 
1174  if ( do_adiabatic_init ) then
1175  beta = 1.; alpha = 0.
1176  if( nudge_debug .and. master) write(*,*) 'Doing final adiabatic initialization/nudging'
1177  else
1178  alpha = 1. - beta
1179  if( abs(alpha)<1.e-7 ) analysis_time = .true.
1180  endif
1181 
1182 ! Warning: ps_data are not adjusted for the differences in terrain yet
1183  ps_obs(:,:) = alpha*ps_dat(:,:,1) + beta*ps_dat(:,:,2)
1184 
1185 !---------------------------------
1186 !*** nudge & update ps & delp here
1187 !---------------------------------
1188  if (nudge_ps) then
1189 
1190  allocate ( wt(is:ie,js:je,km) )
1191  wt(:,:,:) = alpha*t_dat(:,:,:,1) + beta*t_dat(:,:,:,2)
1192 ! Needs gz3 for ps_nudging
1193  call get_int_hght(npz, ak, bk, ps(is:ie,js:je), delp, ps_obs(is:ie,js:je), wt)
1194  do j=js,je
1195  do i=is,ie
1196  tm(i,j) = wt(i,j,km)
1197  enddo
1198  enddo
1199  deallocate ( wt )
1200 
1201  allocate ( uu(isd:ied,jsd:jed,npz) )
1202  allocate ( vv(isd:ied,jsd:jed,npz) )
1203  uu = ua
1204  vv = va
1205  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)
1206  do k=1,npz
1207  do j=js,je
1208  do i=is,ie
1209  u_dt(i,j,k) = u_dt(i,j,k) + (uu(i,j,k) - ua(i,j,k)) / dt
1210  v_dt(i,j,k) = v_dt(i,j,k) + (vv(i,j,k) - va(i,j,k)) / dt
1211  enddo
1212  enddo
1213  enddo
1214  deallocate (uu )
1215  deallocate (vv )
1216  endif
1217 
1218  allocate ( ut(is:ie,js:je,npz) )
1219  allocate ( vt(is:ie,js:je,npz) )
1220 
1221  if ( nudge_winds ) then
1222 
1223  call remap_uv(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1224  km, ps_dat(is:ie,js:je,1), u_dat(:,:,:,1), v_dat(:,:,:,1), ptop )
1225 
1226  u_obs(:,:,:) = alpha*ut(:,:,:)
1227  v_obs(:,:,:) = alpha*vt(:,:,:)
1228 
1229  call remap_uv(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1230  km, ps_dat(is:ie,js:je,2), u_dat(:,:,:,2), v_dat(:,:,:,2), ptop )
1231 
1232  u_obs(:,:,:) = u_obs(:,:,:) + beta*ut(:,:,:)
1233  v_obs(:,:,:) = v_obs(:,:,:) + beta*vt(:,:,:)
1234  endif
1235 
1236  call remap_tq(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1237  km, ps_dat(is:ie,js:je,1), t_dat(:,:,:,1), q_dat(:,:,:,1), zvir, ptop)
1238 
1239  t_obs(:,:,:) = alpha*ut(:,:,:)
1240  q_obs(:,:,:) = alpha*vt(:,:,:)
1241 
1242  call remap_tq(npz, ak, bk, ps(is:ie,js:je), delp, ut, vt, &
1243  km, ps_dat(is:ie,js:je,2), t_dat(:,:,:,2), q_dat(:,:,:,2), zvir, ptop)
1244 
1245  t_obs(:,:,:) = t_obs(:,:,:) + beta*ut(:,:,:)
1246  q_obs(:,:,:) = q_obs(:,:,:) + beta*vt(:,:,:)
1247 
1248  deallocate ( ut )
1249  deallocate ( vt )
1250 
1251  end subroutine get_obs
1252 
1255  subroutine fv_nwp_nudge_init(time, axes, npz, zvir, ak, bk, ts, phis, gridstruct, ks, npx, neststruct, bd)
1256  character(len=17) :: mod_name = 'fv_nudge'
1257  type(time_type), intent(in):: time
1258  integer, intent(in):: axes(4)
1259  integer, intent(in):: npz
1260  real, intent(in):: zvir
1261  type(fv_grid_bounds_type), intent(IN) :: bd
1262  real, intent(in), dimension(bd%isd:bd%ied,bd%jsd:bd%jed):: phis
1263  real, intent(in), dimension(npz+1):: ak, bk
1264  real, intent(out), dimension(bd%is:bd%ie,bd%js:bd%je):: ts
1265  type(fv_grid_type), target :: gridstruct
1266  integer, intent(in) :: ks, npx
1267  type(fv_nest_type) :: neststruct
1268 
1269  real :: missing_value = -1.e10
1270  logical found
1271  integer tsize(4)
1272  integer :: i, j, j1, f_unit, unit, io, ierr, nt, k
1273  integer :: ncid
1274 
1275 ! ---> h1g, read NCEP analysis data list, 2012-10-22
1276  integer :: input_fname_unit
1277  character(len=128) :: fname_tmp
1278 ! <--- h1g, 2012-10-22
1279 
1280  real, pointer, dimension(:,:,:) :: agrid
1281 
1282  is = bd%is
1283  ie = bd%ie
1284  js = bd%js
1285  je = bd%je
1286 
1287  isd = bd%isd
1288  ied = bd%ied
1289  jsd = bd%jsd
1290  jed = bd%jed
1291 
1292 
1293  agrid => gridstruct%agrid
1294 
1295  master = is_master()
1296  do_adiabatic_init = .false.
1297  deg2rad = pi/180.
1298  rad2deg = 180./pi
1299 
1300  if (neststruct%nested) then
1301 !!! Assumes no grid stretching
1302  grid_size = 1.e7/(neststruct%npx_global*neststruct%refinement_of_global)
1303  else
1304  grid_size = 1.e7/real(npx-1) ! mean grid size
1305  endif
1306 
1307  do nt=1,nfile_max
1308  file_names(nt) = "No_File_specified"
1309  enddo
1310 
1311  track_file_name = "No_File_specified"
1312 
1313  if( file_exist( 'input.nml' ) ) then
1314  unit = open_namelist_file()
1315  io = 1
1316  do while ( io .ne. 0 )
1317  read( unit, nml = fv_nwp_nudge_nml, iostat = io, end = 10 )
1318  ierr = check_nml_error(io,'fv_nwp_nudge_nml')
1319  end do
1320 10 call close_file ( unit )
1321  end if
1322  call write_version_number ( 'FV_NUDGE_MOD', version )
1323  if ( master ) then
1324  f_unit=stdlog()
1325  write( f_unit, nml = fv_nwp_nudge_nml )
1326  write(*,*) 'NWP nudging initialized.'
1327  endif
1328 ! ---> h1g, specify the NCEP analysis data for nudging, 2012-10-22
1329  if ( trim(input_fname_list) .ne. "" ) then
1330  input_fname_unit = get_unit()
1331  open( input_fname_unit, file = input_fname_list )
1332  nt = 0
1333  io = 0
1334 
1335  do while ( io .eq. 0 )
1336  read ( input_fname_unit, '(a)', iostat = io, end = 101 ) fname_tmp
1337  if( trim(fname_tmp) .ne. "" ) then ! escape any empty record
1338  if ( trim(fname_tmp) == trim(analysis_file_last) ) then
1339  nt = nt + 1
1340  file_names(nt) = 'INPUT/'//trim(fname_tmp)
1341  if(master .and. nudge_debug) write(*,*) 'From NCEP file list, last file: ', nt, file_names(nt)
1342  nt = 0
1343  goto 101 ! read last analysis data and then close file
1344  endif ! trim(fname_tmp) == trim(analysis_file_last)
1345 
1346  if ( trim(analysis_file_first) == "" ) then
1347  nt = nt + 1
1348  file_names(nt) = 'INPUT/'//trim(fname_tmp)
1349  if(master .and. nudge_debug) then
1350  if( nt .eq. 1 ) then
1351  write(*,*) 'From NCEP file list, first file: ', nt, file_names(nt),trim(analysis_file_first)
1352  else
1353  write(*,*) 'From NCEP file list: ', nt, file_names(nt)
1354  endif ! nt .eq. 1
1355  endif ! master .and. nudge_debug
1356  else
1357  if ( trim(fname_tmp) == trim(analysis_file_first) .or. nt > 0 ) then
1358  nt = nt + 1
1359  file_names(nt) = 'INPUT/'//trim(fname_tmp)
1360  if(master .and. nudge_debug) then
1361  if( nt .eq. 1 ) then
1362  write(*,*) 'From NCEP file list, first file: ', nt, file_names(nt),trim(analysis_file_first)
1363  else
1364  write(*,*) 'From NCEP file list: ', nt, file_names(nt)
1365  endif ! nt .eq. 1
1366  endif ! master .and. nudge_debug
1367  endif ! trim(fname_tmp) == trim(analysis_file_first) .or. nt > 0
1368  endif ! trim(analysis_file_first) == ""
1369  endif ! trim(fname_tmp) .ne. ""
1370  end do ! io .eq. 0
1371 101 close( input_fname_unit )
1372  endif
1373 ! <--- h1g, 2012-10-22
1374 
1375  do nt=1,nfile_max
1376  if ( file_names(nt) == "No_File_specified" ) then
1377  nfile_total = nt - 1
1378  if(master) write(*,*) 'Total of NCEP files specified=', nfile_total
1379  exit
1380  endif
1381  enddo
1382 
1383  id_ht_err = register_diag_field( mod_name, 'ht_error', axes(1:2), time, &
1384  'height_error', 'DEG K', missing_value=missing_value )
1385 
1386 
1387 ! Initialize remapping coefficients:
1388 
1389 ! call field_size(file_names(1), 'T', tsize, field_found=found)
1390 ! if ( found ) then
1391 ! im = tsize(1); jm = tsize(2); km = tsize(3)
1392 ! if(master) write(*,*) 'NCEP analysis dimensions:', tsize
1393 ! else
1394 ! call mpp_error(FATAL,'==> Error from get_ncep_analysis: T field not found')
1395 ! endif
1396  call open_ncfile( file_names(1), ncid ) ! open the file
1397  call get_ncdim1( ncid, 'lon', im )
1398  call get_ncdim1( ncid, 'lat', jm )
1399  call get_ncdim1( ncid, 'lev', km )
1400  if(master) write(*,*) 'NCEP analysis dimensions:', im, jm, km
1401 
1402  allocate ( s2c(is:ie,js:je,4) )
1403  allocate ( id1(is:ie,js:je) )
1404  allocate ( id2(is:ie,js:je) )
1405  allocate ( jdc(is:ie,js:je) )
1406 
1407  allocate ( lon(im) )
1408  allocate ( lat(jm) )
1409 
1410  call _get_var1 (ncid, 'lon', im, lon )
1411  call _get_var1 (ncid, 'lat', jm, lat )
1412 
1413 ! Convert to radian
1414  do i=1,im
1415  lon(i) = lon(i) * deg2rad ! lon(1) = 0.
1416  enddo
1417  do j=1,jm
1418  lat(j) = lat(j) * deg2rad
1419  enddo
1420 
1421  allocate ( ak0(km+1) )
1422  allocate ( bk0(km+1) )
1423 
1424  call _get_var1 (ncid, 'hyai', km+1, ak0, found )
1425  if ( .not. found ) ak0(:) = 0.
1426 
1427  call _get_var1 (ncid, 'hybi', km+1, bk0 )
1428  call close_ncfile( ncid )
1429 
1430 ! Note: definition of NCEP hybrid is p(k) = a(k)*1.E5 + b(k)*ps
1431  ak0(:) = ak0(:) * 1.e5
1432 ! Limiter to prevent NAN at top during remapping
1433  if ( bk0(1) < 1.e-9 ) ak0(1) = max(1.e-9, ak0(1))
1434 
1435  if ( master ) then
1436  do k=1,npz
1437  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)
1438  enddo
1439  endif
1440 
1441  if ( k_breed==0 ) k_breed = max(1, ks)
1442 
1443  call slp_obs_init
1444 
1445 !-----------------------------------------------------------
1446 ! Initialize lat-lon to Cubed bi-linear interpolation coeff:
1447 !-----------------------------------------------------------
1448  call remap_coef(agrid)
1449 
1450 ! Find bounding latitudes:
1451  jbeg = jm-1; jend = 2
1452  do j=js,je
1453  do i=is,ie
1454  j1 = jdc(i,j)
1455  jbeg = min(jbeg, j1)
1456  jend = max(jend, j1+1)
1457  enddo
1458  enddo
1459 
1460  allocate ( gz0(is:ie,js:je) )
1461  allocate ( gz3(is:ie,js:je,km+1) )
1462  allocate (ps_dat(is:ie,js:je,2) )
1463  allocate ( u_dat(is:ie,js:je,km,2) )
1464  allocate ( v_dat(is:ie,js:je,km,2) )
1465  allocate ( t_dat(is:ie,js:je,km,2) )
1466  allocate ( q_dat(is:ie,js:je,km,2) )
1467 
1468 
1469 ! Get first dataset
1470  nt = 2
1471  nfile = 1
1472  call get_ncep_analysis ( ps_dat(:,:,nt), u_dat(:,:,:,nt), v_dat(:,:,:,nt), &
1473  t_dat(:,:,:,nt), q_dat(:,:,:,nt), zvir, &
1474  ts, nfile, file_names(nfile), bd )
1475 
1476 
1477  module_is_initialized = .true.
1478 
1479  nullify(agrid)
1480 
1481  end subroutine fv_nwp_nudge_init
1482 
1483 
1484  subroutine get_ncep_analysis ( ps, u, v, t, q, zvir, ts, nfile, fname, bd )
1485  real, intent(in):: zvir
1486  character(len=128), intent(in):: fname
1487  integer, intent(inout):: nfile
1488 !
1489  type(fv_grid_bounds_type), intent(IN) :: bd
1490  real, intent(out), dimension(is:ie,js:je):: ts, ps
1491  real(KIND=4), intent(out), dimension(is:ie,js:je,km):: u, v, t, q
1492 ! local:
1493  real(kind=4), allocatable:: wk0(:,:), wk1(:,:), wk2(:,:), wk3(:,:,:)
1494  real tmean
1495  integer:: i, j, k, npt
1496  integer:: i1, i2, j1, ncid
1497  logical found
1498  logical:: read_ts = .true.
1499  logical:: land_ts = .false.
1500 
1501  integer:: status, var3id ! h1g, 2016-08-10
1502 #include <netcdf.inc>
1503 
1504 
1505  if( .not. file_exist(fname) ) then
1506  call mpp_error(fatal,'==> Error from get_ncep_analysis: file not found: '//fname)
1507  else
1508  call open_ncfile( fname, ncid ) ! open the file
1509  if(master) write(*,*) 'Reading NCEP anlysis file:', fname
1510  endif
1511 
1512  if ( read_ts ) then ! read skin temperature; could be used for SST
1513  allocate ( wk1(im,jm) )
1514 
1515  call get_var3_r4( ncid, 'TS', 1,im, 1,jm, 1,1, wk1 )
1516 ! if ( master ) write(*,*) 'Done reading NCEP TS data'
1517 
1518  if ( .not. land_ts ) then
1519  allocate ( wk0(im,jm) )
1520 ! Read NCEP ORO (1; land; 0: ocean; 2: sea_ice)
1521 
1522 ! ---> h1g, read either 'ORO' or 'LAND', 2016-08-10
1523  status = nf_inq_varid(ncid, 'ORO', var3id)
1524  if (status .eq. nf_noerr) then
1525  call get_var3_r4( ncid, 'ORO', 1,im, 1,jm, 1,1, wk0 )
1526 
1527  else !there is no 'ORO'
1528  status = nf_inq_varid(ncid, 'LAND', var3id)
1529  if (status .eq. nf_noerr) then
1530  call get_var3_r4( ncid, 'LAND', 1,im, 1,jm, 1,1, wk0 )
1531  else
1532  call mpp_error(fatal,'Neither ORO nor LAND exists in re-analysis data')
1533  endif
1534  endif
1535 ! <--- h1g, 2016-08-10
1536 
1537  do j=1,jm
1538  tmean = 0.
1539  npt = 0
1540  do i=1,im
1541  if( abs(wk0(i,j)-1.) > 0.99 ) then ! non-land points
1542  tmean = tmean + wk1(i,j)
1543  npt = npt + 1
1544  endif
1545  enddo
1546 !-------------------------------------------------------
1547 ! Replace TS over interior land with zonal mean SST/Ice
1548 !-------------------------------------------------------
1549  if ( npt /= 0 ) then
1550  tmean= tmean / real(npt)
1551  do i=1,im
1552  if( abs(wk0(i,j)-1.) <= 0.99 ) then ! land points
1553  if ( i==1 ) then
1554  i1 = im; i2 = 2
1555  elseif ( i==im ) then
1556  i1 = im-1; i2 = 1
1557  else
1558  i1 = i-1; i2 = i+1
1559  endif
1560  if ( abs(wk0(i2,j)-1.)>0.99 ) then ! east side has priority
1561  wk1(i,j) = wk1(i2,j)
1562  elseif ( abs(wk0(i1,j)-1.)>0.99 ) then ! west side
1563  wk1(i,j) = wk1(i1,j)
1564  else
1565  wk1(i,j) = tmean
1566  endif
1567  endif
1568  enddo
1569  endif
1570  enddo
1571  deallocate ( wk0 )
1572  endif ! land_ts
1573 
1574  do j=js,je
1575  do i=is,ie
1576  i1 = id1(i,j)
1577  i2 = id2(i,j)
1578  j1 = jdc(i,j)
1579  ts(i,j) = s2c(i,j,1)*wk1(i1,j1 ) + s2c(i,j,2)*wk1(i2,j1 ) + &
1580  s2c(i,j,3)*wk1(i2,j1+1) + s2c(i,j,4)*wk1(i1,j1+1)
1581  enddo
1582  enddo
1583  call prt_maxmin('SST_model', ts, is, ie, js, je, 0, 1, 1.)
1584 
1585 #ifndef DYCORE_SOLO
1586 ! Perform interp to FMS SST format/grid
1587  call ncep2fms( wk1 )
1588  if(master) call pmaxmin( 'SST_ncep', sst_ncep, i_sst, j_sst, 1.)
1589 ! if(nfile/=1 .and. master) call pmaxmin( 'SST_anom', sst_anom, i_sst, j_sst, 1.)
1590 #endif
1591  deallocate ( wk1 )
1592  if (master) write(*,*) 'Done processing NCEP SST'
1593 
1594  endif ! read_ts
1595 
1596 !----------------------------------
1597 ! remap surface pressure and height:
1598 !----------------------------------
1599 
1600  allocate ( wk2(im,jbeg:jend) )
1601  call get_var3_r4( ncid, 'PS', 1,im, jbeg,jend, 1,1, wk2 )
1602 
1603  do j=js,je
1604  do i=is,ie
1605  i1 = id1(i,j)
1606  i2 = id2(i,j)
1607  j1 = jdc(i,j)
1608  ps(i,j) = s2c(i,j,1)*wk2(i1,j1 ) + s2c(i,j,2)*wk2(i2,j1 ) + &
1609  s2c(i,j,3)*wk2(i2,j1+1) + s2c(i,j,4)*wk2(i1,j1+1)
1610  enddo
1611  enddo
1612 
1613 
1614 ! ---> h1g, read either 'PHIS' or 'PHI', 2016-08-10
1615  status = nf_inq_varid(ncid, 'PHIS', var3id)
1616  if (status .eq. nf_noerr) then
1617  call get_var3_r4( ncid, 'PHIS', 1,im, jbeg,jend, 1,1, wk2 )
1618 
1619  else !there is no 'PHIS'
1620  status = nf_inq_varid(ncid, 'PHI', var3id)
1621  if (status .eq. nf_noerr) then
1622  call get_var3_r4( ncid, 'PHI', 1,im, jbeg,jend, 1,1, wk2 )
1623  wk2 = wk2 * grav ! convert unit from geopotential meter (m) to geopotential height (m2/s2)
1624  else
1625  call mpp_error(fatal,'Neither PHIS nor PHI exists in re-analysis data')
1626  endif
1627  endif
1628 ! <--- h1g, 2016-08-10
1629 
1630 
1631  do j=js,je
1632  do i=is,ie
1633  i1 = id1(i,j)
1634  i2 = id2(i,j)
1635  j1 = jdc(i,j)
1636  gz0(i,j) = s2c(i,j,1)*wk2(i1,j1 ) + s2c(i,j,2)*wk2(i2,j1 ) + &
1637  s2c(i,j,3)*wk2(i2,j1+1) + s2c(i,j,4)*wk2(i1,j1+1)
1638  enddo
1639  enddo
1640  call prt_maxmin('ZS_ncep', gz0, is, ie, js, je, 0, 1, 1./grav)
1641  deallocate ( wk2 )
1642 
1643 
1644  allocate ( wk3(1:im, jbeg:jend, 1:km) )
1645 
1646 ! Winds:
1647  if ( nudge_winds ) then
1648 
1649  call get_var3_r4( ncid, 'U', 1,im, jbeg,jend, 1,km, wk3 )
1650 
1651  do k=1,km
1652  do j=js,je
1653  do i=is,ie
1654  i1 = id1(i,j)
1655  i2 = id2(i,j)
1656  j1 = jdc(i,j)
1657  u(i,j,k) = s2c(i,j,1)*wk3(i1,j1 ,k) + s2c(i,j,2)*wk3(i2,j1 ,k) + &
1658  s2c(i,j,3)*wk3(i2,j1+1,k) + s2c(i,j,4)*wk3(i1,j1+1,k)
1659  enddo
1660  enddo
1661  enddo
1662 
1663  call get_var3_r4( ncid, 'V', 1,im, jbeg,jend, 1,km, wk3 )
1664 
1665  do k=1,km
1666  do j=js,je
1667  do i=is,ie
1668  i1 = id1(i,j)
1669  i2 = id2(i,j)
1670  j1 = jdc(i,j)
1671  v(i,j,k) = s2c(i,j,1)*wk3(i1,j1 ,k) + s2c(i,j,2)*wk3(i2,j1 ,k) + &
1672  s2c(i,j,3)*wk3(i2,j1+1,k) + s2c(i,j,4)*wk3(i1,j1+1,k)
1673  enddo
1674  enddo
1675  enddo
1676 
1677  endif
1678 
1679 ! Read in tracers: only sphum at this point
1680  call get_var3_r4( ncid, 'Q', 1,im, jbeg,jend, 1,km , wk3 )
1681 
1682  do k=1,km
1683  do j=js,je
1684  do i=is,ie
1685  i1 = id1(i,j)
1686  i2 = id2(i,j)
1687  j1 = jdc(i,j)
1688  q(i,j,k) = s2c(i,j,1)*wk3(i1,j1 ,k) + s2c(i,j,2)*wk3(i2,j1 ,k) + &
1689  s2c(i,j,3)*wk3(i2,j1+1,k) + s2c(i,j,4)*wk3(i1,j1+1,k)
1690  enddo
1691  enddo
1692  enddo
1693 
1694  call get_var3_r4( ncid, 'T', 1,im, jbeg,jend, 1,km , wk3 )
1695  call close_ncfile ( ncid )
1696 
1697  do k=1,km
1698  do j=js,je
1699  do i=is,ie
1700  i1 = id1(i,j)
1701  i2 = id2(i,j)
1702  j1 = jdc(i,j)
1703  t(i,j,k) = s2c(i,j,1)*wk3(i1,j1 ,k) + s2c(i,j,2)*wk3(i2,j1 ,k) + &
1704  s2c(i,j,3)*wk3(i2,j1+1,k) + s2c(i,j,4)*wk3(i1,j1+1,k)
1705  enddo
1706  enddo
1707  enddo
1708 
1709  if ( .not. t_is_tv ) then
1710  do k=1,km
1711  do j=js,je
1712  do i=is,ie
1713 ! The field T in Larry H.'s post processing of NCEP analysis is actually virtual temperature
1714 ! before Dec 1, 2005
1715 ! Convert t to virtual temperature:
1716  t(i,j,k) = t(i,j,k)*(1.+zvir*q(i,j,k))
1717  enddo
1718  enddo
1719  enddo
1720  endif
1721 
1722 ! endif
1723 
1724  deallocate ( wk3 )
1725 
1726 ! nfile = nfile + 1
1727 
1728  end subroutine get_ncep_analysis
1729 
1730 
1731  subroutine remap_coef(agrid)
1733  real, intent(IN), dimension(isd:ied,jsd:jed,2) :: agrid
1734 
1735 ! local:
1736  real :: rdlon(im)
1737  real :: rdlat(jm)
1738  real:: a1, b1
1739  integer i,j, i1, i2, jc, i0, j0
1740 
1741  do i=1,im-1
1742  rdlon(i) = 1. / (lon(i+1) - lon(i))
1743  enddo
1744  rdlon(im) = 1. / (lon(1) + 2.*pi - lon(im))
1745 
1746  do j=1,jm-1
1747  rdlat(j) = 1. / (lat(j+1) - lat(j))
1748  enddo
1749 
1750 ! * Interpolate to cubed sphere cell center
1751  do 5000 j=js,je
1752 
1753  do i=is,ie
1754 
1755  if ( agrid(i,j,1)>lon(im) ) then
1756  i1 = im; i2 = 1
1757  a1 = (agrid(i,j,1)-lon(im)) * rdlon(im)
1758  elseif ( agrid(i,j,1)<lon(1) ) then
1759  i1 = im; i2 = 1
1760  a1 = (agrid(i,j,1)+2.*pi-lon(im)) * rdlon(im)
1761  else
1762  do i0=1,im-1
1763  if ( agrid(i,j,1)>=lon(i0) .and. agrid(i,j,1)<=lon(i0+1) ) then
1764  i1 = i0; i2 = i0+1
1765  a1 = (agrid(i,j,1)-lon(i1)) * rdlon(i0)
1766  go to 111
1767  endif
1768  enddo
1769  endif
1770 111 continue
1771 
1772  if ( agrid(i,j,2)<lat(1) ) then
1773  jc = 1
1774  b1 = 0.
1775  elseif ( agrid(i,j,2)>lat(jm) ) then
1776  jc = jm-1
1777  b1 = 1.
1778  else
1779  do j0=1,jm-1
1780  if ( agrid(i,j,2)>=lat(j0) .and. agrid(i,j,2)<=lat(j0+1) ) then
1781  jc = j0
1782  b1 = (agrid(i,j,2)-lat(jc)) * rdlat(jc)
1783  go to 222
1784  endif
1785  enddo
1786  endif
1787 222 continue
1788 
1789  if ( a1<0.0 .or. a1>1.0 .or. b1<0.0 .or. b1>1.0 ) then
1790  write(*,*) 'gid=', mpp_pe(), i,j,a1, b1
1791  endif
1792 
1793  s2c(i,j,1) = (1.-a1) * (1.-b1)
1794  s2c(i,j,2) = a1 * (1.-b1)
1795  s2c(i,j,3) = a1 * b1
1796  s2c(i,j,4) = (1.-a1) * b1
1797  id1(i,j) = i1
1798  id2(i,j) = i2
1799  jdc(i,j) = jc
1800  enddo !i-loop
1801 5000 continue ! j-loop
1802 
1803  end subroutine remap_coef
1804 
1805 
1806 #ifndef DYCORE_SOLO
1807  subroutine ncep2fms( sst )
1808  real(kind=4), intent(in):: sst(im,jm)
1809 ! local:
1810  real :: rdlon(im)
1811  real :: rdlat(jm)
1812  real:: a1, b1
1813  real:: delx, dely
1814  real:: xc, yc ! "data" location
1815  real:: c1, c2, c3, c4
1816  integer i,j, i1, i2, jc, i0, j0, it, jt
1817 
1818  do i=1,im-1
1819  rdlon(i) = 1. / (lon(i+1) - lon(i))
1820  enddo
1821  rdlon(im) = 1. / (lon(1) + 2.*pi - lon(im))
1822 
1823  do j=1,jm-1
1824  rdlat(j) = 1. / (lat(j+1) - lat(j))
1825  enddo
1826 
1827 ! * Interpolate to "FMS" 1x1 SST data grid
1828 ! lon: 0.5, 1.5, ..., 359.5
1829 ! lat: -89.5, -88.5, ... , 88.5, 89.5
1830 
1831  delx = 360./real(i_sst)
1832  dely = 180./real(j_sst)
1833 
1834  jt = 1
1835  do 5000 j=1,j_sst
1836 
1837  yc = (-90. + dely * (0.5+real(j-1))) * deg2rad
1838  if ( yc<lat(1) ) then
1839  jc = 1
1840  b1 = 0.
1841  elseif ( yc>lat(jm) ) then
1842  jc = jm-1
1843  b1 = 1.
1844  else
1845  do j0=jt,jm-1
1846  if ( yc>=lat(j0) .and. yc<=lat(j0+1) ) then
1847  jc = j0
1848  jt = j0
1849  b1 = (yc-lat(jc)) * rdlat(jc)
1850  go to 222
1851  endif
1852  enddo
1853  endif
1854 222 continue
1855  it = 1
1856 
1857  do i=1,i_sst
1858  xc = delx * (0.5+real(i-1)) * deg2rad
1859  if ( xc>lon(im) ) then
1860  i1 = im; i2 = 1
1861  a1 = (xc-lon(im)) * rdlon(im)
1862  elseif ( xc<lon(1) ) then
1863  i1 = im; i2 = 1
1864  a1 = (xc+2.*pi-lon(im)) * rdlon(im)
1865  else
1866  do i0=it,im-1
1867  if ( xc>=lon(i0) .and. xc<=lon(i0+1) ) then
1868  i1 = i0; i2 = i0+1
1869  it = i0
1870  a1 = (xc-lon(i1)) * rdlon(i0)
1871  go to 111
1872  endif
1873  enddo
1874  endif
1875 111 continue
1876 
1877 ! if ( a1<0.0 .or. a1>1.0 .or. b1<0.0 .or. b1>1.0 ) then
1878 ! write(*,*) 'gid=', mpp_pe(), i,j,a1, b1
1879 ! endif
1880  c1 = (1.-a1) * (1.-b1)
1881  c2 = a1 * (1.-b1)
1882  c3 = a1 * b1
1883  c4 = (1.-a1) * b1
1884 ! Interpolated surface pressure
1885  sst_ncep(i,j) = c1*sst(i1,jc ) + c2*sst(i2,jc ) + &
1886  c3*sst(i2,jc+1) + c4*sst(i1,jc+1)
1887  enddo !i-loop
1888 5000 continue ! j-loop
1889 
1890  end subroutine ncep2fms
1891 #endif
1892 
1893  subroutine get_int_hght(npz, ak, bk, ps, delp, ps0, tv)
1894  integer, intent(in):: npz
1895  real, intent(in):: ak(npz+1), bk(npz+1)
1896  real, intent(in), dimension(is:ie,js:je):: ps, ps0
1897  real, intent(in), dimension(isd:ied,jsd:jed,npz):: delp
1898  real(KIND=4), intent(in), dimension(is:ie,js:je,km):: tv ! virtual temperature
1899 ! local:
1900  real, dimension(is:ie,km+1):: pn0
1901  integer i,j,k
1902 
1903 !$OMP parallel do default(none) shared(is,ie,js,je,km,ak0,bk0,ps0,gz3,gz0,tv) &
1904 !$OMP private(pn0)
1905  do 5000 j=js,je
1906 
1907  do k=1,km+1
1908  do i=is,ie
1909  pn0(i,k) = log( ak0(k) + bk0(k)*ps0(i,j) )
1910  enddo
1911  enddo
1912  do i=is,ie
1913  gz3(i,j,km+1) = gz0(i,j) ! Data Surface geopotential
1914  enddo
1915 
1916  do k=km,1,-1
1917  do i=is,ie
1918  gz3(i,j,k) = gz3(i,j,k+1) + rdgas*tv(i,j,k)*(pn0(i,k+1)-pn0(i,k))
1919  enddo
1920  enddo
1921 
1922 5000 continue
1923 
1924  end subroutine get_int_hght
1925 
1926 
1927 
1928  subroutine remap_tq( npz, ak, bk, ps, delp, t, q, &
1929  kmd, ps0, ta, qa, zvir, ptop)
1930  integer, intent(in):: npz, kmd
1931  real, intent(in):: zvir, ptop
1932  real, intent(in):: ak(npz+1), bk(npz+1)
1933  real, intent(in), dimension(is:ie,js:je):: ps0
1934  real, intent(inout), dimension(is:ie,js:je):: ps
1935  real, intent(in), dimension(isd:ied,jsd:jed,npz):: delp
1936  real(KIND=4), intent(in), dimension(is:ie,js:je,kmd):: ta, qa
1937  real(KIND=4), intent(out), dimension(is:ie,js:je,npz):: t
1938 ! on input "ta" is virtual temperature
1939 ! on output "t" is virtual temperature
1940  real(KIND=4), intent(out), dimension(is:ie,js:je,npz):: q
1941 ! local:
1942  real, dimension(is:ie,kmd):: tp, qp
1943  real, dimension(is:ie,kmd+1):: pe0, pn0
1944  real, dimension(is:ie,npz):: qn1
1945  real, dimension(is:ie,npz+1):: pe1, pn1
1946  integer i,j,k
1947 
1948  do 5000 j=js,je
1949 
1950  do k=1,kmd+1
1951  do i=is,ie
1952  pe0(i,k) = ak0(k) + bk0(k)*ps0(i,j)
1953  pn0(i,k) = log(pe0(i,k))
1954  enddo
1955  enddo
1956 !------
1957 ! Model
1958 !------
1959  do i=is,ie
1960  pe1(i,1) = ak(1)
1961  enddo
1962  do k=1,npz
1963  do i=is,ie
1964  pe1(i,k+1) = pe1(i,k) + delp(i,j,k)
1965  enddo
1966  enddo
1967  do i=is,ie
1968  ps(i,j) = pe1(i,npz+1)
1969  enddo
1970  do k=1,npz+1
1971  do i=is,ie
1972  pn1(i,k) = log(pe1(i,k))
1973  enddo
1974  enddo
1975 
1976  if ( nudge_q ) then
1977  do k=1,kmd
1978  do i=is,ie
1979  qp(i,k) = qa(i,j,k)
1980  enddo
1981  enddo
1982  call mappm(kmd, pe0, qp, npz, pe1, qn1, is,ie, 0, kord_data, ptop)
1983  do k=1,npz
1984  do i=is,ie
1985  q(i,j,k) = qn1(i,k)
1986  enddo
1987  enddo
1988  endif
1989 
1990  do k=1,kmd
1991  do i=is,ie
1992  tp(i,k) = ta(i,j,k)
1993  enddo
1994  enddo
1995  call mappm(kmd, pn0, tp, npz, pn1, qn1, is,ie, 1, kord_data, ptop)
1996 
1997  do k=1,npz
1998  do i=is,ie
1999  t(i,j,k) = qn1(i,k)
2000  enddo
2001  enddo
2002 
2003 5000 continue
2004 
2005  end subroutine remap_tq
2006 
2007 
2008  subroutine remap_uv(npz, ak, bk, ps, delp, u, v, kmd, ps0, u0, v0, ptop)
2009  integer, intent(in):: npz
2010  real, intent(IN):: ptop
2011  real, intent(in):: ak(npz+1), bk(npz+1)
2012  real, intent(inout):: ps(is:ie,js:je)
2013  real, intent(in), dimension(isd:ied,jsd:jed,npz):: delp
2014  real(KIND=4), intent(inout), dimension(is:ie,js:je,npz):: u, v
2015 !
2016  integer, intent(in):: kmd
2017  real, intent(in):: ps0(is:ie,js:je)
2018  real(KIND=4), intent(in), dimension(is:ie,js:je,kmd):: u0, v0
2019 !
2020 ! local:
2021  real, dimension(is:ie,kmd+1):: pe0
2022  real, dimension(is:ie,npz+1):: pe1
2023  real, dimension(is:ie,kmd):: qt
2024  real, dimension(is:ie,npz):: qn1
2025  integer i,j,k
2026 
2027  do 5000 j=js,je
2028 !------
2029 ! Data
2030 !------
2031  do k=1,kmd+1
2032  do i=is,ie
2033  pe0(i,k) = ak0(k) + bk0(k)*ps0(i,j)
2034  enddo
2035  enddo
2036 !------
2037 ! Model
2038 !------
2039  do i=is,ie
2040  pe1(i,1) = ak(1)
2041  enddo
2042  do k=1,npz
2043  do i=is,ie
2044  pe1(i,k+1) = pe1(i,k) + delp(i,j,k)
2045  enddo
2046  enddo
2047  do i=is,ie
2048  ps(i,j) = pe1(i,npz+1)
2049  enddo
2050 !------
2051 ! map u
2052 !------
2053  do k=1,kmd
2054  do i=is,ie
2055  qt(i,k) = u0(i,j,k)
2056  enddo
2057  enddo
2058  call mappm(kmd, pe0, qt, npz, pe1, qn1, is,ie, -1, kord_data, ptop)
2059  do k=1,npz
2060  do i=is,ie
2061  u(i,j,k) = qn1(i,k)
2062  enddo
2063  enddo
2064 !------
2065 ! map v
2066 !------
2067  do k=1,kmd
2068  do i=is,ie
2069  qt(i,k) = v0(i,j,k)
2070  enddo
2071  enddo
2072  call mappm(kmd, pe0, qt, npz, pe1, qn1, is,ie, -1, kord_data, ptop)
2073  do k=1,npz
2074  do i=is,ie
2075  v(i,j,k) = qn1(i,k)
2076  enddo
2077  enddo
2078 5000 continue
2079 
2080  end subroutine remap_uv
2081 
2083  subroutine fv_nwp_nudge_end
2084  deallocate ( ps_dat )
2085  deallocate ( t_dat )
2086  deallocate ( q_dat )
2087 
2088  if ( nudge_winds ) then
2089  deallocate ( u_dat )
2090  deallocate ( v_dat )
2091  endif
2092 
2093  deallocate ( s2c )
2094  deallocate ( id1 )
2095  deallocate ( id2 )
2096  deallocate ( jdc )
2097 
2098  deallocate ( ak0 )
2099  deallocate ( bk0 )
2100  deallocate ( lat )
2101  deallocate ( lon )
2102 
2103  deallocate ( gz3 )
2104  deallocate ( gz0 )
2105 
2106  end subroutine fv_nwp_nudge_end
2107 
2108 
2109  subroutine get_tc_mask(time, mask, agrid)
2110  real :: slp_mask = 100900.
2111 ! Input
2112  type(time_type), intent(in):: time
2113  real, intent(inout):: mask(is:ie,js:je)
2114  real(kind=R_GRID), intent(IN), dimension(isd:ied,jsd:jed,2) :: agrid
2115 ! local
2116  real(kind=R_GRID):: pos(2)
2117  real:: slp_o
2118  real:: w10_o
2119  real:: r_vor, p_vor
2120  real:: dist
2121  integer n, i, j
2122 
2123  do 5000 n=1,nstorms ! loop through all storms
2124 !----------------------------------------
2125 ! Obtain slp observation
2126 !----------------------------------------
2127  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), &
2128  time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_vor)
2129 
2130  if ( slp_o<880.e2 .or. slp_o>min(slp_env,slp_mask) .or. abs(pos(2))*rad2deg>40. ) goto 5000 ! next storm
2131 
2132  if ( r_vor < 30.e3 ) then
2133  r_vor = r_min + (slp_env-slp_o)/20.e2*r_inc ! radius of influence
2134  endif
2135 
2136  do j=js, je
2137  do i=is, ie
2138  dist = great_circle_dist(pos, agrid(i,j,1:2), radius)
2139  if( dist < 6.*r_vor ) then
2140  mask(i,j) = mask(i,j) * ( 1. - mask_fac*exp(-(0.5*dist/r_vor)**2)*min(1.,(slp_env-slp_o)/10.e2) )
2141  endif
2142  enddo ! i-loop
2143  enddo ! end j-loop
2144 
2145 5000 continue
2146 
2147  end subroutine get_tc_mask
2148 
2154  subroutine breed_slp_inline(nstep, dt, npz, ak, bk, phis, pe, pk, peln, pkz, delp, u, v, pt, q, nwat, &
2155  zvir, gridstruct, ks, domain_local, bd, hydrostatic)
2156 ! Input
2157  integer, intent(in):: nstep, npz, nwat, ks
2158  real, intent(in):: dt
2159  real, intent(in):: zvir
2160  real, intent(in), dimension(npz+1):: ak, bk
2161  logical, intent(in):: hydrostatic
2162  type(fv_grid_bounds_type), intent(IN) :: bd
2163  real, intent(in):: phis(isd:ied,jsd:jed)
2164  type(domain2d), intent(INOUT) :: domain_local
2165 ! Input/Output
2166  real, intent(inout):: u(isd:ied,jsd:jed+1,npz)
2167  real, intent(inout):: v(isd:ied+1,jsd:jed,npz)
2168  real, intent(inout), dimension(isd:ied,jsd:jed,npz):: delp, pt
2169  real, intent(inout)::q(isd:ied,jsd:jed,npz,*)
2170 
2171  real, intent(inout):: pk(is:ie,js:je, npz+1)
2172  real, intent(inout):: pe(is-1:ie+1, npz+1,js-1:je+1)
2173  real, intent(inout):: pkz(is:ie,js:je,npz)
2174  real, intent(out):: peln(is:ie,npz+1,js:je)
2175 
2176  type(fv_grid_type), target :: gridstruct
2177 ! local
2178  type(time_type):: time
2179  real:: ps(is:ie,js:je)
2180  real:: dist(is:ie,js:je)
2181  real:: tm(is:ie,js:je)
2182  real:: slp(is:ie,js:je)
2183  real(kind=R_GRID):: pos(2)
2184  real:: slp_o ! sea-level pressure (Pa)
2185  real:: w10_o, p_env
2186  real:: r_vor
2187  real:: relx0, relx, f1, pbreed, pbtop, delp0, dp0
2188  real:: ratio, p_count, p_sum, a_sum, mass_sink, delps
2189  real:: p_lo, p_hi, tau_vt, mslp0
2190  real:: split_time, fac, pdep, r2, r3
2191  integer year, month, day, hour, minute, second
2192  integer n, i, j, k, iq, k0
2193 
2194  real, pointer :: area(:,:)
2195  real(kind=R_GRID), pointer :: agrid(:,:,:)
2196 #ifdef MULTI_GASES
2197  real :: kappax(is:ie,js:je,npz)
2198 #endif
2199 
2200  if ( forecast_mode ) return
2201 
2202  agrid => gridstruct%agrid_64
2203  area => gridstruct%area
2204 
2205  if ( nstorms==0 ) then
2206  if(master) write(*,*) 'NO TC data to process'
2207  return
2208  endif
2209 
2210  if ( k_breed==0 ) k_breed = max(1, ks)
2211 
2212 
2213 ! Advance (local) time
2214  call get_date(fv_time, year, month, day, hour, minute, second)
2215  if ( year /= year_track_data ) then
2216  if (master) write(*,*) 'Warning: The year in storm track data is not the same as model year'
2217  return
2218  endif
2219  time = fv_time ! fv_time is the time at past time step (set in fv_diag)
2220  split_time = calday(year, month, day, hour, minute, second) + dt*real(nstep)/86400.
2221 
2222  elapsed_time = elapsed_time + dt
2223  if ( elapsed_time > nudged_time + 0.1 ) then
2224  if(print_end_breed) then
2225  print_end_breed = .false.
2226  if (master) write(*,*) '*** Vortext Breeding Ended at', day, hour, minute, second
2227  endif
2228  return ! time to return to forecast mode
2229  endif
2230 
2231 !$OMP parallel do default(none) shared(is,ie,js,je,npz,ps,ak,delp,tm,pkz,pt)
2232  do j=js,je
2233 ! ---- Compute ps
2234  do i=is,ie
2235  ps(i,j) = ak(1)
2236  enddo
2237  do k=1,npz
2238  do i=is,ie
2239  ps(i,j) = ps(i,j) + delp(i,j,k)
2240  enddo
2241  enddo
2242 ! Compute lowest layer air temperature:
2243  do i=is,ie
2244  tm(i,j) = pkz(i,j,npz)*pt(i,j,npz) / cp_air ! virtual temp
2245  enddo
2246  enddo
2247 ! call prt_maxmin('TM', tm, is, ie, js, je, 0, 1, 1.)
2248 
2249 !$OMP parallel do default(none) shared(k_breed,npz,conserve_mom,is,ie,js,je,u,v,delp, &
2250 !$OMP conserve_hgt,hydrostatic,pt,pkz,q,nwat)
2251  do k=k_breed+1,npz
2252 
2253  if ( conserve_mom ) then
2254  do j=js,je+1
2255  do i=is,ie
2256  u(i,j,k) = u(i,j,k) * (delp(i,j-1,k)+delp(i,j,k))
2257  enddo
2258  enddo
2259  do j=js,je
2260  do i=is,ie+1
2261  v(i,j,k) = v(i,j,k) * (delp(i-1,j,k)+delp(i,j,k))
2262  enddo
2263  enddo
2264  endif
2265 
2266  if ( conserve_hgt .and. hydrostatic ) then
2267  do j=js,je
2268  do i=is,ie
2269 ! Conserve total enthalpy
2270  pt(i,j,k) = pt(i,j,k)*pkz(i,j,k)*delp(i,j,k)
2271  enddo
2272  enddo
2273  endif
2274 
2275 ! Convert tracer moist mixing ratio to mass
2276  do iq=1,nwat
2277  do j=js,je
2278  do i=is,ie
2279  q(i,j,k,iq) = q(i,j,k,iq) * delp(i,j,k)
2280  enddo
2281  enddo
2282  enddo
2283 
2284  enddo
2285 
2286  do 5000 n=1,nstorms ! loop through all storms
2287 
2288  if ( nobs_tc(n) < min_nobs ) goto 5000
2289 
2290 ! Check min MSLP
2291  mslp0 = 1013.e2
2292  do i=1,nobs_tc(n)
2293  mslp0 = min( mslp0, mslp_obs(i,n) )
2294  enddo
2295  if ( mslp0 > min_mslp ) goto 5000
2296 
2297 !----------------------------------------
2298 ! Obtain slp observation
2299 !----------------------------------------
2300  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), &
2301  time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_env, stime=split_time, fact=fac)
2302 
2303  if ( slp_o<87500. .or. slp_o>slp_env .or. abs(pos(2))*rad2deg>45. ) then
2304  goto 5000 ! next storm
2305  endif
2306 
2307 
2308  if(nudge_debug .and. master) &
2309  write(*,*) 'Vortex breeding for TC:', n, ' slp=',slp_o/100.,pos(1)*rad2deg,pos(2)*rad2deg
2310 
2311 ! Determine pbtop (top pressure of vortex breeding)
2312  k0 = k_breed
2313 
2314  if ( use_high_top ) then
2315 #ifdef HIGH_TEST
2316  if ( slp_o > 1000.e2 ) then
2317  pbtop = 850.e2
2318  else
2319 ! pbtop = max(200.E2, 850.E2-20.*(1000.E2-slp_o))
2320 ! mp87:
2321  pbtop = max(100.e2, 850.e2-25.*(1000.e2-slp_o))
2322 #else
2323  if ( slp_o > 1000.e2 ) then
2324  pbtop = 750.e2
2325  else
2326  pbtop = max(100.e2, 750.e2-20.*(1000.e2-slp_o))
2327 #endif
2328  endif
2329  else
2330 ! Lower top for vrotex breeding
2331  if ( slp_o > 1000.e2 ) then
2332  pbtop = 900.e2
2333  else
2334  pbtop = max(500.e2, 900.e2-5.*(1000.e2-slp_o)) ! mp48
2335  endif
2336  endif
2337 
2338  do k=1,npz
2339  pbreed = ak(k) + bk(k)*1.e5
2340  if ( pbreed>pbtop ) then
2341  k0 = k
2342  exit
2343  endif
2344  enddo
2345  k0 = max(k0, k_breed)
2346 
2347  do j=js, je
2348  do i=is, ie
2349  dist(i,j) = great_circle_dist( pos, agrid(i,j,1:2), radius)
2350  enddo
2351  enddo
2352 
2353  call compute_slp(is, ie, js, je, tm, ps, phis(is:ie,js:je), slp)
2354 
2355  if ( r_vor < 30.e3 .or. p_env<900.e2 ) then
2356 
2357 ! Compute r_vor & p_env
2358  r_vor = r_min + (slp_env-slp_o)/25.e2*r_inc
2359 
2360 123 continue
2361  p_count = 0.
2362  p_sum = 0.
2363  a_sum = 0.
2364  do j=js, je
2365  do i=is, ie
2366  if( dist(i,j)<(r_vor+del_r) .and. dist(i,j)>r_vor .and. phis(i,j)<250.*grav ) then
2367  p_count = p_count + 1.
2368  p_sum = p_sum + slp(i,j)*area(i,j)
2369  a_sum = a_sum + area(i,j)
2370  endif
2371  enddo
2372  enddo
2373 
2374  call mp_reduce_sum(p_count)
2375 
2376  if ( p_count<32. ) then
2377  if(nudge_debug .and. master) write(*,*) p_count, 'Skipping obs: too few p_count'
2378  goto 5000
2379  endif
2380 
2381  call mp_reduce_sum(p_sum)
2382  call mp_reduce_sum(a_sum)
2383  p_env = p_sum / a_sum
2384 
2385  if(nudge_debug .and. master) write(*,*) 'Environmental SLP=', p_env/100., ' computed radius=', r_vor/1.e3
2386 
2387  if ( p_env>1015.e2 .or. p_env<920.e2 ) then
2388  if( nudge_debug ) then
2389  if(master) write(*,*) 'Environmental SLP out of bound; skipping obs. p_count=', p_count, p_sum
2390  call prt_maxmin('SLP_breeding', slp, is, ie, js, je, 0, 1, 0.01)
2391  endif
2392  goto 5000
2393  endif
2394 
2395  endif
2396 
2397 
2398  if ( p_env < max(pre0_env, slp_o + 200.0) ) then
2399 
2400  r_vor = r_vor + 25.e3
2401 
2402  if(nudge_debug .and. master) then
2403  write(*,*) 'Computed environmental SLP too low'
2404  write(*,*) ' ', p_env/100., slp_o/100.,pos(1)*rad2deg, pos(2)*rad2deg
2405  endif
2406 
2407  if ( slp_o > 1003.e2 .and. r_vor > 1000.e3 ) then
2408 ! if(master) write(*,*) 'Failed to size the Vortex for the weak storm'
2409  goto 5000
2410  endif
2411 
2412  if ( r_vor < 1250.e3 ) goto 123
2413 
2414 ! if(master) write(*,*) 'Failed to size the Vortex; skipping this storm'
2415  goto 5000
2416 
2417  endif
2418 
2419  tau_vt = tau_vt_slp * (1. + (960.e2-slp_o)/100.e2 )
2420  tau_vt = max(abs(dt), tau_vt)
2421 
2422  if ( do_adiabatic_init ) then
2423  relx0 = min(1., 2.*abs(dt)/tau_vt)
2424  else
2425  relx0 = min(1., abs(dt)/tau_vt)
2426  endif
2427 
2428  mass_sink = 0.
2429 
2430 !$OMP parallel do default(none) shared(is,ie,js,je,dist,r_vor,phis,p_env,slp_o,r_hi,r_lo, &
2431 !$OMP ps,tm,relx0,tau_vt_rad,dps_min,slp,ak,k0,delp,npz,area) &
2432 !$OMP private(f1, p_hi, p_lo, relx, delps, pbreed, mass_sink, dp0, pdep)
2433  do j=js, je
2434  do i=is, ie
2435  if( dist(i,j) < r_vor .and. phis(i,j)<250.*grav ) then
2436  f1 = dist(i,j)/r_vor
2437 ! Compute p_obs: assuming local radial distributions of slp are Gaussian
2438  p_hi = p_env - (p_env-slp_o) * exp( -r_hi*f1**2 ) ! upper bound
2439  p_lo = p_env - (p_env-slp_o) * exp( -r_lo*f1**2 ) ! lower bound
2440 
2441  if ( ps(i,j) > p_hi .and. tm(i,j) < tm_max ) then
2442 ! do nothing if lowest layer is too hot
2443 ! Under-development:
2444  relx = relx0*exp( -tau_vt_rad*f1**2 )
2445  delps = relx*(ps(i,j) - p_hi) * min( (tm_max-tm(i,j))/10., 1.)
2446 ! After mp115
2447 ! delps = relx*(ps(i,j) - p_hi) * min( (tm_max-tm(i,j))/5., 1.)
2448  ! Cap the increment to prevent overheating
2449  ! Note: ps is used here to prevent
2450  ! over deepening over terrain
2451  delps = min(delps, dps_min)
2452  elseif ( slp(i,j) < p_lo ) then
2453 ! Over-development:
2454  relx = max(0.5, relx0)
2455  delps = relx*(slp(i,j) - p_lo) ! Note: slp is used here
2456  else
2457  goto 400 ! do nothing; proceed to next storm
2458  endif
2459 
2460 #ifdef SIM_TEST
2461  pbreed = ak(1)
2462  do k=1,k0
2463  pbreed = pbreed + delp(i,j,k)
2464  enddo
2465  f1 = 1. - delps/(ps(i,j)-pbreed)
2466  do k=k0+1,npz
2467  delp(i,j,k) = delp(i,j,k)*f1
2468  enddo
2469  mass_sink = mass_sink + delps*area(i,j)
2470 #else
2471  if ( delps > 0. ) then
2472  pbreed = ak(1)
2473  do k=1,k0
2474  pbreed = pbreed + delp(i,j,k)
2475  enddo
2476  f1 = 1. - delps/(ps(i,j)-pbreed)
2477  do k=k0+1,npz
2478  delp(i,j,k) = delp(i,j,k)*f1
2479  enddo
2480  mass_sink = mass_sink + delps*area(i,j)
2481  else
2482  dp0 = abs(delps)
2483  do k=npz,k0+1,-1
2484  if ( abs(delps) < 1. ) then
2485  delp(i,j,k) = delp(i,j,k) - delps
2486  mass_sink = mass_sink + delps*area(i,j)
2487  go to 400
2488  else
2489  pdep = max(1.0, min(abs(0.4*delps), 0.2*dp0, 0.02*delp(i,j,k)))
2490  pdep = - min(pdep, abs(delps))
2491  delp(i,j,k) = delp(i,j,k) - pdep
2492  delps = delps - pdep
2493  mass_sink = mass_sink + pdep*area(i,j)
2494  endif
2495  enddo
2496  endif
2497 #endif
2498 
2499  endif
2500 400 continue
2501  enddo ! end i-loop
2502  enddo ! end j-loop
2503 
2504  call mp_reduce_sum(mass_sink)
2505  if ( abs(mass_sink)<1.e-40 ) goto 5000
2506 
2507  r2 = r_vor + del_r
2508  r3 = min( 4.*r_vor, max(2.*r_vor, 2500.e3) ) + del_r
2509 
2510  p_sum = 0.
2511  do j=js, je
2512  do i=is, ie
2513  if( dist(i,j)<r3 .and. dist(i,j)>r2 ) then
2514  p_sum = p_sum + area(i,j)
2515  endif
2516  enddo
2517  enddo
2518 
2519  call mp_reduce_sum(p_sum)
2520  mass_sink = mass_sink / p_sum ! mean delta pressure to be added back to the environment to conserve mass
2521  if(master .and. nudge_debug) write(*,*) 'TC#',n, 'Mass tele-ported (pa)=', mass_sink
2522 
2523 !$OMP parallel do default(none) shared(is,ie,js,je,dist,r3,r2,ak,k_breed,delp,ps,mass_sink,npz) &
2524 !$OMP private(pbreed, f1)
2525  do j=js, je
2526  do i=is, ie
2527  if( dist(i,j)<r3 .and. dist(i,j)>r2 ) then
2528  pbreed = ak(1)
2529  do k=1,k_breed
2530  pbreed = pbreed + delp(i,j,k)
2531  enddo
2532  f1 = 1. + mass_sink/(ps(i,j)-pbreed)
2533  do k=k_breed+1,npz
2534  delp(i,j,k) = delp(i,j,k)*f1
2535  enddo
2536  endif
2537  enddo
2538  enddo
2539 
2540 ! ---- re-compute ps
2541  do j=js,je
2542  do i=is,ie
2543  ps(i,j) = ak(1)
2544  enddo
2545  do k=1,npz
2546  do i=is,ie
2547  ps(i,j) = ps(i,j) + delp(i,j,k)
2548  enddo
2549  enddo
2550  enddo
2551 
2552 5000 continue
2553 
2554 !--------------------------
2555 ! Update delp halo regions:
2556 !--------------------------
2557  call mpp_update_domains(delp, domain_local, complete=.true.)
2558 
2559 !$OMP parallel do default(none) shared(is,ie,js,je,npz,pe,ak,delp,k_breed,peln,pk)
2560  do j=js-1,je+1
2561  do i=is-1,ie+1
2562  pe(i,1,j) = ak(1)
2563  enddo
2564  do k=2,npz+1
2565  do i=is-1,ie+1
2566  pe(i,k,j) = pe(i,k-1,j) + delp(i,j,k-1)
2567  enddo
2568  enddo
2569  enddo
2570 
2571 #ifdef MULTI_GASES
2572 !$OMP parallel do default(none) shared(is,ie,js,je,npz,k_breed,kappax,num_gas)
2573  do k=k_breed,npz
2574  do j=js,je
2575  do i=is,ie
2576  kappax(i,j,k)= virqd(q(i,j,k,1:num_gas))/vicpqd(q(i,j,k,1:num_gas))
2577  enddo
2578  enddo
2579  enddo
2580 #endif
2581  do k=k_breed+1,npz+1
2582  do j=js,je
2583  do i=is,ie
2584  peln(i,k,j) = log(pe(i,k,j))
2585  pk(i,j,k) = pe(i,k,j)**kappa
2586  enddo
2587  enddo
2588  enddo
2589 
2590 !$OMP parallel do default(none) shared(k_breed,npz,is,ie,js,je,conserve_mom,u,v,delp, &
2591 #ifdef MULTI_GASES
2592 !$OMP kappax, &
2593 #endif
2594 !$OMP conserve_hgt,hydrostatic,pkz,pk,pt,peln)
2595  do k=k_breed+1,npz
2596 
2597  if ( conserve_mom ) then
2598  do j=js,je+1
2599  do i=is,ie
2600  u(i,j,k) = u(i,j,k) / (delp(i,j-1,k)+delp(i,j,k))
2601  enddo
2602  enddo
2603  do j=js,je
2604  do i=is,ie+1
2605  v(i,j,k) = v(i,j,k) / (delp(i-1,j,k)+delp(i,j,k))
2606  enddo
2607  enddo
2608  endif
2609 
2610  if ( conserve_hgt .and. hydrostatic ) then
2611  do j=js,je
2612  do i=is,ie
2613 ! Conserve total enthalpy (static energy)
2614  pkz(i,j,k) = (pk(i,j,k+1)-pk(i,j,k))/(kappa*(peln(i,k+1,j)-peln(i,k,j)))
2615 #ifdef MULTI_GASES
2616  pkz(i,j,k) = exp(kappax(i,j,k)*log(pkz(i,j,k)))
2617 #endif
2618  pt(i,j,k) = pt(i,j,k) / (pkz(i,j,k)*delp(i,j,k))
2619  enddo
2620  enddo
2621 ! pkz last used
2622  endif
2623  enddo
2624 
2625 ! Convert tracer mass back to moist mixing ratio
2626 
2627 !$OMP parallel do default(none) shared(nwat,k_breed,npz,is,ie,js,je,q,delp)
2628  do iq=1,nwat
2629  do k=k_breed+1,npz
2630  do j=js,je
2631  do i=is,ie
2632  q(i,j,k,iq) = q(i,j,k,iq) / delp(i,j,k)
2633  enddo
2634  enddo
2635  enddo
2636  enddo
2637 
2638  if(hydrostatic) call mpp_update_domains(pt, domain_local, complete=.true.)
2639 
2640  nullify(agrid)
2641  nullify(area)
2642 
2643  end subroutine breed_slp_inline
2644 
2647  subroutine breed_srf_w10(time, dt, npz, ak, bk, ps, phis, slp, delp, u, v, gridstruct)
2648  type(time_type), intent(in):: time
2649  integer, intent(in):: npz
2650  real, intent(in):: dt
2651  real, intent(in), dimension(npz+1):: ak, bk
2652  real, intent(in):: phis(isd:ied,jsd:jed)
2653  real, intent(in):: ps(isd:ied,jsd:jed)
2654  real, intent(in), dimension(is:ie,js:je):: slp
2655  type(fv_grid_type), intent(IN), target :: gridstruct
2656  real, intent(inout):: u(isd:ied,jsd:jed+1,npz)
2657  real, intent(inout):: v(isd:ied+1,jsd:jed,npz)
2658 ! Input/Output
2659  real, intent(inout), dimension(isd:ied,jsd:jed,npz):: delp
2660 ! local
2661  real, dimension(is:ie,js:je):: us, vs
2662  real wu(is:ie, js:je+1)
2663  real wv(is:ie+1,js:je)
2664  real u1(is:ie), v1(is:ie)
2665  real:: dist(is:ie,js:je), wind(is:ie,js:je)
2666  real(kind=R_GRID):: pos(2)
2667  real:: slp_o
2668  real:: w10_o, p_env
2669  real:: r_vor, pc, p_count
2670  real:: r_max, speed, ut, vt, speed_local ! tangent wind speed
2671  real:: u_bg, v_bg, mass, t_mass
2672  real:: relx0, relx, f1
2673  real:: z0, mslp0
2674  real:: zz = 35.
2675  real:: wind_fac
2676  integer n, i, j
2677 
2678  ! Pointers
2679  real, pointer, dimension(:,:) :: dx, dy, rdxa, rdya, a11, a21, a12, a22, area
2680  real(kind=R_GRID), pointer, dimension(:,:,:) :: agrid, vlon, vlat
2681 
2682  dx => gridstruct%dx
2683  dy => gridstruct%dy
2684  rdxa => gridstruct%rdxa
2685  rdya => gridstruct%rdya
2686  a11 => gridstruct%a11
2687  a21 => gridstruct%a21
2688  a12 => gridstruct%a12
2689  a22 => gridstruct%a22
2690  area => gridstruct%area
2691  agrid => gridstruct%agrid_64
2692  vlon => gridstruct%vlon
2693  vlat => gridstruct%vlat
2694 
2695 
2696  if ( nstorms==0 ) then
2697  if(master) write(*,*) 'NO TC data to process'
2698  return
2699  endif
2700 
2701 ! Compute lat-lon winds on A grid
2702 !$OMP parallel do default(none) shared(is,ie,js,je,wu,u,npz,dx)
2703  do j=js,je+1
2704  do i=is,ie
2705  wu(i,j) = u(i,j,npz)*dx(i,j)
2706  enddo
2707  enddo
2708 !$OMP parallel do default(none) shared(is,ie,js,je,wv,v,npz,dy)
2709  do j=js,je
2710  do i=is,ie+1
2711  wv(i,j) = v(i,j,npz)*dy(i,j)
2712  enddo
2713  enddo
2714 
2715 !$OMP parallel do default(none) shared(is,ie,js,je,wu,rdxa,wv,rdya,us,vs,a11,a12,a21,a22,wind) &
2716 !$OMP private(u1,v1)
2717  do j=js, je
2718  do i=is, ie
2719 ! Co-variant to Co-variant "vorticity-conserving" interpolation
2720  u1(i) = (wu(i,j) + wu(i,j+1)) * rdxa(i,j)
2721  v1(i) = (wv(i,j) + wv(i+1,j)) * rdya(i,j)
2722 ! Cubed (cell center co-variant winds) to lat-lon:
2723  us(i,j) = a11(i,j)*u1(i) + a12(i,j)*v1(i)
2724  vs(i,j) = a21(i,j)*u1(i) + a22(i,j)*v1(i)
2725 ! Wind speed
2726  wind(i,j) = sqrt( us(i,j)**2 + vs(i,j)**2 )
2727  enddo
2728  enddo
2729 
2730  relx0 = min(1., dt/tau_vt_wind)
2731 
2732  do 3000 n=1,nstorms ! loop through all storms
2733 
2734  if ( nobs_tc(n) < min_nobs ) goto 3000
2735 
2736 ! Check min MSLP
2737  mslp0 = 1013.e2
2738  do i=1,nobs_tc(n)
2739  mslp0 = min( mslp0, mslp_obs(i,n) )
2740  enddo
2741  if ( mslp0 > min_mslp ) goto 3000
2742 
2743 !----------------------------------------
2744 ! Obtain slp observation
2745 !----------------------------------------
2746  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), &
2747  time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_env)
2748 
2749  if ( slp_o<90000. .or. slp_o>slp_env .or. abs(pos(2))*rad2deg>35. ) goto 3000 ! next storm
2750 
2751 
2752  do j=js, je
2753  do i=is, ie
2754  dist(i,j) = great_circle_dist( pos, agrid(i,j,1:2), radius )
2755  enddo
2756  enddo
2757 
2758  r_vor = r_min + (slp_env-slp_o)/25.e2*r_inc
2759 
2760 !----------------------------------------------------
2761 ! * Find model's SLP center nearest to the observation
2762 ! * Find maximum wind speed at the lowest model level
2763 
2764  speed_local = 0.
2765  r_max = -999.
2766  pc = 1013.e2
2767  p_count = 0.
2768  do j=js, je
2769  do i=is, ie
2770  if( dist(i,j) < r_vor ) then
2771 
2772 ! Counting the "land" points:
2773  if ( dist(i,j)<0.5*r_vor .and. abs(phis(i,j))>2.*grav ) p_count = p_count + 1.
2774 
2775  pc = min(pc, slp(i,j))
2776 
2777  if ( speed_local < wind(i,j) ) then
2778  speed_local = wind(i,j)
2779  r_max = dist(i,j)
2780  endif
2781 
2782  endif
2783  enddo
2784  enddo
2785 
2786  call mp_reduce_sum(p_count)
2787  if ( p_count>32 ) goto 3000 ! over/near rough land
2788 
2789  if ( w10_o < 0. ) then ! 10-m wind obs is not available
2790 ! Uses Atkinson_Holliday wind-pressure correlation
2791  w10_o = 3.446778 * (1010.-slp_o/100.)**0.644
2792  endif
2793 
2794  speed = speed_local
2795  call mp_reduce_max(speed) ! global max wind (near storm)
2796  call mp_reduce_min(pc)
2797 
2798  if ( speed_local < speed ) then
2799  r_max = -999.
2800  endif
2801  call mp_reduce_max(r_max)
2802  if( r_max<0. ) call mpp_error(fatal,'==> Error in r_max')
2803 
2804 ! ---------------------------------------------------
2805 ! Determine surface wind speed and radius for nudging
2806 ! ---------------------------------------------------
2807 
2808 ! Compute surface roughness z0 from w10, based on Eq (4) & (5) from Moon et al. 2007
2809  if ( w10_o > 12.5 ) then
2810  z0 = (0.085*w10_o - 0.58) * 1.e-3
2811 ! z0 (w10=40) = 2.82E-3
2812 ! z0 = min( z0, 2.82E-3 )
2813 
2814  else
2815  z0 = 0.0185/grav*(0.001*w10_o**2 + 0.028*w10_o)**2
2816  endif
2817 
2818 ! lowest layer height: zz
2819 
2820  wind_fac = log(zz/z0) / log(10./z0)
2821  if( nudge_debug .and. master ) write(*,*) 'Wind adjustment factor=', wind_fac
2822  if( wind_fac<1. ) call mpp_error(fatal,'==> Error in wind_fac')
2823 
2824  if ( pc < slp_o ) then
2825 !--
2826 ! The storm in the model is over developed
2827 ! if ( (pc+3.0E2)>slp_o .or. speed <= w10_o ) go to 3000 ! next storm
2828 !--
2829 ! using radius (r_max) as dtermined above;
2830 ! What if model's pressure center is very far from the observed?
2831 ! using obs wind
2832  speed = wind_fac*w10_o
2833  else
2834 ! The storm in the model is under developed; using max wind
2835  speed = max(wind_fac*w10_o, speed)
2836  if ( pc>1009.e2 ) r_max = 0.5 * r_vor
2837  endif
2838 
2839 ! More adjustment
2840 
2841 ! Some bounds on the radius of maximum wind:
2842  r_max = max(2.5*grid_size, r_max) ! at least 2.5X the grid size
2843  r_max = min(0.75*r_vor, r_max)
2844 
2845  t_mass = 0.
2846  u_bg = 0.
2847  v_bg = 0.
2848 
2849  if ( add_bg_wind ) then
2850  p_count = 0.
2851  do j=js, je
2852  do i=is, ie
2853  if( dist(i,j) <= min(r_vor,2.*r_fac*r_max) .and. phis(i,j)<2.*grav ) then
2854  mass = area(i,j)*delp(i,j,npz)
2855 !-- using model winds ----------------------------------
2856  u_bg = u_bg + us(i,j)*mass
2857  v_bg = v_bg + vs(i,j)*mass
2858 !-------------------------------------------------------
2859  t_mass = t_mass + mass
2860  p_count = p_count + 1.
2861  endif
2862  enddo
2863  enddo
2864  call mp_reduce_sum(p_count)
2865  if ( p_count<16. ) go to 3000
2866 
2867  call mp_reduce_sum(t_mass)
2868  call mp_reduce_sum(u_bg)
2869  call mp_reduce_sum(v_bg)
2870  u_bg = u_bg / t_mass
2871  v_bg = v_bg / t_mass
2872  if ( nudge_debug .and. master ) write(*,*) pos(2)*rad2deg, 'vortex bg wind=', u_bg, v_bg
2873  endif
2874 
2875  relx = relx0
2876 ! Nudge wind in the "inner core":
2877 #ifdef TTT
2878  do j=js, je
2879  do i=is, ie
2880  if( dist(i,j) <= min(r_vor, r_fac*r_max) .and. phis(i,j)<2.*grav ) then
2881  f1 = dist(i,j)/r_max
2882  relx = relx0*exp( -tau_vt_rad*f1**2 )
2883  if( dist(i,j)<=r_max ) then
2884  speed_local = speed * f1
2885  else
2886  speed_local = speed / f1**0.75
2887  endif
2888  call tangent_wind(vlon(i,j,1:3), vlat(i,j,1:3), speed_local, pos, agrid(i,j,1:2), ut, vt)
2889  ut = ut + u_bg
2890  vt = vt + v_bg
2891 ! Update:
2892  us(i,j) = relx*(ut-us(i,j))
2893  vs(i,j) = relx*(vt-vs(i,j))
2894  endif
2895 400 continue
2896  enddo ! end i-loop
2897  enddo ! end j-loop
2898 #else
2899  do j=js, je
2900  do i=is, ie
2901  if( dist(i,j) <= min(r_vor, r_fac*r_max) .and. phis(i,j)<2.*grav ) then
2902  f1 = dist(i,j)/r_max
2903  relx = relx0*exp( -tau_vt_rad*f1**2 )
2904  if( dist(i,j)<=r_max ) then
2905  speed_local = speed * f1
2906  else
2907  speed_local = speed / f1**0.75
2908  endif
2909  call tangent_wind(vlon(i,j,1:3), vlat(i,j,1:3), speed_local, pos, agrid(i,j,1:2), ut, vt)
2910  ut = ut + u_bg
2911  vt = vt + v_bg
2912 ! Update:
2913  us(i,j) = relx*(ut-us(i,j))
2914  vs(i,j) = relx*(vt-vs(i,j))
2915  endif
2916 400 continue
2917  enddo ! end i-loop
2918  enddo ! end j-loop
2919 #endif
2920 
2921 3000 continue
2922 
2923  end subroutine breed_srf_w10
2924 
2926  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)
2927 ! Input
2928  type(time_type), intent(in):: time
2929  integer, intent(in):: npz, nwat
2930  real, intent(in):: dt
2931  real, intent(in):: zvir
2932  real, intent(in), dimension(npz+1):: ak, bk
2933  real, intent(in), dimension(isd:ied,jsd:jed):: phis, ps
2934  real, intent(in), dimension(is:ie,js:je,npz):: u_obs, v_obs
2935  type(fv_grid_type), intent(IN), target :: gridstruct
2936 ! Input/Output
2937  real, intent(inout), dimension(isd:ied,jsd:jed,npz):: delp, pt, ua, va, u_dt, v_dt
2938  real, intent(inout)::q(isd:ied,jsd:jed,npz,nwat)
2939 ! local
2940  real:: dist(is:ie,js:je), wind(is:ie,js:je)
2941  real:: slp(is:ie,js:je)
2942  real(kind=R_GRID):: pos(2)
2943  real:: slp_o
2944  real:: w10_o, p_env
2945  real:: r_vor, pc, p_count
2946  real:: r_max, speed, ut, vt, speed_local ! tangent wind speed
2947  real:: u_bg, v_bg, mass, t_mass
2948  real:: relx0, relx, f1, rdt
2949  real:: z0, mslp0
2950  real:: zz = 35.
2951  real:: wind_fac
2952  integer n, i, j, k, iq
2953 
2954  real, pointer :: area(:,:)
2955  real(kind=R_GRID), pointer :: vlon(:,:,:), vlat(:,:,:), agrid(:,:,:)
2956 
2957  area => gridstruct%area
2958  vlon => gridstruct%vlon
2959  vlat => gridstruct%vlat
2960  agrid => gridstruct%agrid_64
2961 
2962  if ( nstorms==0 ) then
2963  if(master) write(*,*) 'NO TC data to process'
2964  return
2965  endif
2966 
2967  rdt = 1./dt
2968  relx0 = min(1., dt/tau_vt_wind)
2969 
2970 !$OMP parallel do default(none) shared(is,ie,js,je,slp,ps,phis,pt,npz,wind,ua,va)
2971  do j=js, je
2972  do i=is, ie
2973  slp(i,j) = ps(i,j)*exp(phis(i,j)/(rdgas*(pt(i,j,npz)+3.25e-3*phis(i,j)/grav)))
2974  wind(i,j) = sqrt( ua(i,j,npz)**2 + va(i,j,npz)**2 )
2975  enddo
2976  enddo
2977 
2978 !!!!!$OMP parallel do default(none) private(pos, w10_o, slp_o, r_vor, p_env)
2979  do 3000 n=1,nstorms ! loop through all storms
2980 
2981  if ( nobs_tc(n) < min_nobs ) goto 3000
2982 
2983 ! Check min MSLP
2984  mslp0 = 1013.e2
2985  do i=1,nobs_tc(n)
2986  mslp0 = min( mslp0, mslp_obs(i,n) )
2987  enddo
2988  if ( mslp0 > min_mslp ) goto 3000
2989 
2990 !----------------------------------------
2991 ! Obtain slp observation
2992 !----------------------------------------
2993  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), &
2994  time_tc(1,n), pos(1), pos(2), w10_o, slp_o, r_vor, p_env)
2995 
2996  if ( slp_o<90000. .or. slp_o>slp_env .or. abs(pos(2))*rad2deg>35. ) goto 3000 ! next storm
2997 
2998 
2999  do j=js, je
3000  do i=is, ie
3001  dist(i,j) = great_circle_dist( pos, agrid(i,j,1:2), radius )
3002  enddo
3003  enddo
3004 
3005  r_vor = r_min + (slp_env-slp_o)/25.e2*r_inc
3006 
3007 !----------------------------------------------------
3008 
3009 ! * Find model's SLP center nearest to the observation
3010 ! * Find maximum wind speed at the lowest model level
3011 
3012  speed_local = 0.
3013  r_max = -999.
3014  pc = 1013.e2
3015  p_count = 0.
3016  do j=js, je
3017  do i=is, ie
3018  if( dist(i,j) < r_vor ) then
3019 
3020 ! Counting the "land" points:
3021  if ( dist(i,j)<0.5*r_vor .and. abs(phis(i,j))>2.*grav ) p_count = p_count + 1.
3022 
3023  pc = min(pc, slp(i,j))
3024 
3025  if ( speed_local < wind(i,j) ) then
3026  speed_local = wind(i,j)
3027  r_max = dist(i,j)
3028  endif
3029 
3030  endif
3031  enddo
3032  enddo
3033 
3034  call mp_reduce_sum(p_count)
3035  if ( p_count>32 ) goto 3000 ! over/near rough land
3036 
3037  if ( w10_o < 0. ) then ! 10-m wind obs is not available
3038 ! Uses Atkinson_Holliday wind-pressure correlation
3039  w10_o = 3.446778 * (1010.-slp_o/100.)**0.644
3040  endif
3041 
3042  speed = speed_local
3043  call mp_reduce_max(speed) ! global max wind (near storm)
3044  call mp_reduce_min(pc)
3045 
3046  if ( speed_local < speed ) then
3047  r_max = -999.
3048  endif
3049  call mp_reduce_max(r_max)
3050  if( r_max<0. ) call mpp_error(fatal,'==> Error in r_max')
3051 
3052 ! ---------------------------------------------------
3053 ! Determine surface wind speed and radius for nudging
3054 ! ---------------------------------------------------
3055 
3056 ! Compute surface roughness z0 from w10, based on Eq (4) & (5) from Moon et al. 2007
3057  if ( w10_o > 12.5 ) then
3058  z0 = (0.085*w10_o - 0.58) * 1.e-3
3059 ! z0 (w10=40) = 2.82E-3
3060 ! z0 = min( z0, 2.82E-3 )
3061 
3062  else
3063  z0 = 0.0185/grav*(0.001*w10_o**2 + 0.028*w10_o)**2
3064  endif
3065 
3066 ! lowest layer height: zz
3067 
3068  wind_fac = log(zz/z0) / log(10./z0)
3069  if( nudge_debug .and. master ) write(*,*) 'Wind adjustment factor=', wind_fac
3070  if( wind_fac<1. ) call mpp_error(fatal,'==> Error in wind_fac')
3071 
3072  if ( pc < slp_o ) then
3073 !--
3074 ! The storm in the model is over developed
3075 ! if ( (pc+3.0E2)>slp_o .or. speed <= w10_o ) go to 3000 ! next storm
3076 !--
3077 ! using radius (r_max) as dtermined above;
3078 ! What if model's pressure center is very far from the observed?
3079 ! using obs wind
3080  speed = wind_fac*w10_o
3081  else
3082 ! The storm in the model is under developed; using max wind
3083  speed = max(wind_fac*w10_o, speed)
3084  if ( pc>1009.e2 ) r_max = 0.5 * r_vor
3085  endif
3086 
3087 ! More adjustment
3088 
3089 ! Some bounds on the radius of maximum wind:
3090  r_max = max(2.5*grid_size, r_max) ! at least 2.5X the grid size
3091  r_max = min(0.75*r_vor, r_max)
3092 
3093  t_mass = 0.
3094  u_bg = 0.
3095  v_bg = 0.
3096 
3097  if ( add_bg_wind ) then
3098  p_count = 0.
3099  do j=js, je
3100  do i=is, ie
3101  if( dist(i,j) <= min(r_vor,2.*r_fac*r_max) .and. phis(i,j)<2.*grav ) then
3102  mass = area(i,j)*delp(i,j,npz)
3103 !-- using model winds ----------------------------------
3104 ! u_bg = u_bg + ua(i,j,npz)*mass
3105 ! v_bg = v_bg + va(i,j,npz)*mass
3106 !-------------------------------------------------------
3107 ! Using analysis winds
3108  u_bg = u_bg + u_obs(i,j,npz)*mass
3109  v_bg = v_bg + v_obs(i,j,npz)*mass
3110  t_mass = t_mass + mass
3111  p_count = p_count + 1.
3112  endif
3113  enddo
3114  enddo
3115  call mp_reduce_sum(p_count)
3116  if ( p_count<16. ) go to 3000
3117 
3118  call mp_reduce_sum(t_mass)
3119  call mp_reduce_sum(u_bg)
3120  call mp_reduce_sum(v_bg)
3121  u_bg = u_bg / t_mass
3122  v_bg = v_bg / t_mass
3123 ! if ( master ) write(*,*) pos(2)*rad2deg, 'vortex bg wind=', u_bg, v_bg
3124  endif
3125 
3126  relx = relx0
3127  k = npz ! lowest layer only
3128 ! Nudge wind in the "inner core":
3129  do j=js, je
3130  do i=is, ie
3131  if( dist(i,j) <= min(r_vor, r_fac*r_max) .and. phis(i,j)<2.*grav ) then
3132  f1 = dist(i,j)/r_max
3133  relx = relx0*exp( -tau_vt_rad*f1**2 )
3134  if( dist(i,j)<=r_max ) then
3135  speed_local = speed * f1
3136  else
3137  speed_local = speed / f1**0.75
3138  endif
3139  call tangent_wind(vlon(i,j,1:3), vlat(i,j,1:3), speed_local, pos, agrid(i,j,1:2), ut, vt)
3140  ut = ut + u_bg
3141  vt = vt + v_bg
3142  u_dt(i,j,k) = u_dt(i,j,k) + relx*(ut-ua(i,j,k)) * rdt
3143  v_dt(i,j,k) = v_dt(i,j,k) + relx*(vt-va(i,j,k)) * rdt
3144 ! Update:
3145  ua(i,j,k) = ua(i,j,k) + relx*(ut-ua(i,j,k))
3146  va(i,j,k) = va(i,j,k) + relx*(vt-va(i,j,k))
3147  endif
3148 400 continue
3149  enddo ! end i-loop
3150  enddo ! end j-loop
3151 
3152 3000 continue
3153 
3154  end subroutine breed_srf_winds
3155 
3156  subroutine tangent_wind ( elon, elat, speed, po, pp, ut, vt )
3157  real, intent(in):: speed
3158  real(kind=R_GRID), intent(in):: po(2), pp(2)
3159  real(kind=R_GRID), intent(in):: elon(3), elat(3)
3160  real, intent(out):: ut, vt
3161 ! local
3162  real(kind=R_GRID):: e1(3), eo(3), ep(3), op(3)
3163 
3164  call latlon2xyz(po, eo)
3165  call latlon2xyz(pp, ep)
3166 
3167  op(:) = ep(:) - eo(:)
3168  call normalize_vect( op )
3169 
3170  call vect_cross(e1, ep, eo)
3171 
3172  ut = speed * (e1(1)*elon(1) + e1(2)*elon(2) + e1(3)*elon(3))
3173  vt = speed * (e1(1)*elat(1) + e1(2)*elat(2) + e1(3)*elat(3))
3174 
3175 ! SH:
3176  if ( po(2) < 0. ) then
3177  ut = -ut
3178  vt = -vt
3179  endif
3180 
3181  end subroutine tangent_wind
3182 
3183 
3184  subroutine get_slp_obs(time, nobs, lon_obs, lat_obs, w10, mslp, slp_out, r_out, time_obs, &
3185  x_o, y_o, w10_o, slp_o, r_vor, p_vor, stime, fact)
3186 ! Input
3187  type(time_type), intent(in):: time
3188  integer, intent(in):: nobs
3189  real(KIND=4), intent(in):: lon_obs(nobs)
3190  real(KIND=4), intent(in):: lat_obs(nobs)
3191  real(KIND=4), intent(in):: w10(nobs)
3192  real(KIND=4), intent(in):: mslp(nobs)
3193  real(KIND=4), intent(in):: slp_out(nobs)
3194  real(KIND=4), intent(in):: r_out(nobs)
3195  real(KIND=4), intent(in):: time_obs(nobs)
3196  real, optional, intent(in):: stime
3197  real, optional, intent(out):: fact
3198 ! Output
3199  real(kind=R_GRID), intent(out):: x_o , y_o
3200  real, intent(out):: w10_o
3201  real, intent(out):: slp_o
3202  real, intent(out):: r_vor, p_vor
3203 ! Internal:
3204  real:: t_thresh
3205  real(kind=R_GRID):: p1(2), p2(2)
3206  real time_model
3207  real(kind=R_GRID) fac
3208  integer year, month, day, hour, minute, second, n
3209 
3210  t_thresh = 600./86400. ! unit: days
3211 
3212  w10_o = -100000.
3213  slp_o = -100000.
3214  x_o = -100.*pi
3215  y_o = -100.*pi
3216  p_vor = -1.e10
3217  r_vor = -1.e10
3218 
3219  if ( present(stime) ) then
3220  time_model = stime
3221  else
3222  call get_date(time, year, month, day, hour, minute, second)
3223 
3224  if ( year /= year_track_data ) then
3225  if (master) write(*,*) 'Warning: The year in storm track data is not the same as model year'
3226  return
3227  endif
3228 
3229  time_model = calday(year, month, day, hour, minute, second)
3230 ! if(nudge_debug .and. master) write(*,*) 'Model:', time_model, year, month, day, hour, minute, second
3231  endif
3232 
3233 !-------------------------------------------------------------------------------------------
3234 ! if ( time_model <= time_obs(1) .or. time_model >= time_obs(nobs) ) then
3235 ! return
3236 !-------------------------------------------------------------------------------------------
3237 
3238  if ( time_model <= (time_obs(1)-t_thresh) .or. time_model >= time_obs(nobs) ) return
3239 
3240  if ( time_model <= time_obs(1) ) then
3241 !--
3242 ! This is an attempt to perform vortex breeding several minutes before the first available observation
3243 !--
3244  w10_o = w10(1)
3245  slp_o = mslp(1)
3246  x_o = lon_obs(1)
3247  y_o = lat_obs(1)
3248  if ( present(fact) ) fact = 1.25
3249  else
3250  do n=1,nobs-1
3251  if( time_model >= time_obs(n) .and. time_model <= time_obs(n+1) ) then
3252  fac = (time_model-time_obs(n)) / (time_obs(n+1)-time_obs(n))
3253  w10_o = w10(n) + ( w10(n+1)- w10(n)) * fac
3254  slp_o = mslp(n) + ( mslp(n+1)- mslp(n)) * fac
3255 ! Trajectory interpolation:
3256 ! Linear in (lon,lat) space
3257 ! x_o = lon_obs(n) + (lon_obs(n+1)-lon_obs(n)) * fac
3258 ! y_o = lat_obs(n) + (lat_obs(n+1)-lat_obs(n)) * fac
3259  p1(1) = lon_obs(n); p1(2) = lat_obs(n)
3260  p2(1) = lon_obs(n+1); p2(2) = lat_obs(n+1)
3261  call intp_great_circle(fac, p1, p2, x_o, y_o)
3262 !----------------------------------------------------------------------
3263  if ( present(fact) ) fact = 1. + 0.25*cos(fac*2.*pi)
3264 ! Additional data from the extended best track
3265 ! if ( slp_out(n)>0. .and. slp_out(n+1)>0. .and. r_out(n)>0. .and. r_out(n+1)>0. ) then
3266 ! p_vor = slp_out(n) + ( slp_out(n+1) - slp_out(n)) * fac
3267 ! r_vor = r_out(n) + ( r_out(n+1) - r_out(n)) * fac
3268 ! endif
3269  return
3270  endif
3271  enddo
3272  endif
3273 
3274  end subroutine get_slp_obs
3275 
3276 
3277  subroutine slp_obs_init
3278  integer:: unit, n, nobs
3279  character(len=3):: GMT
3280  character(len=9):: ts_name
3281  character(len=19):: comment
3282  integer:: mmddhh, yr, year, month, day, hour, MPH, islp
3283  integer:: it, i1, i2, p_ring, d_ring
3284  real:: lon_deg, lat_deg, cald, slp, mps
3285 
3286  nobs_tc(:) = 0
3287  time_tc(:,:) = 0.
3288  wind_obs(:,:) = -100000.
3289  mslp_obs(:,:) = -100000.
3290  x_obs(:,:) = - 100.*pi
3291  y_obs(:,:) = - 100.*pi
3292 
3293  mslp_out(:,:) = -1.e10
3294  rad_out(:,:) = -1.e10
3295 
3296  if( track_file_name == "No_File_specified" ) then
3297  if(master) write(*,*) 'No TC track file specified'
3298  return
3299  else
3300  unit = get_unit()
3301  open( unit, file=track_file_name)
3302  endif
3303 
3304  read(unit, *) year
3305  if(master) write(*,*) 'Reading TC track data for YEAR=', year
3306 
3307  year_track_data = year
3308 
3309  nstorms = 0
3310  nobs = 0
3311  month = 99
3312 
3313  if ( ibtrack ) then
3314 
3315 !---------------------------------------------------------------
3316 ! The data format is from Ming Zhoa's processed ibTrack datasets
3317 !---------------------------------------------------------------
3318 
3319  read(unit, *) ts_name, nobs, yr, month, day, hour
3320 
3321  if ( yr /= year ) then
3322  if(master) write(*, *) 'Year inconsistency found !!!'
3323  call mpp_error(fatal,'==> Error in reading best track data')
3324  endif
3325 
3326  do while ( ts_name=='start' )
3327 
3328  nstorms = nstorms + 1
3329  nobs_tc(nstorms) = nobs ! observation count for this storm
3330  if(nudge_debug .and. master) write(*, *) 'Read Data for TC#', nstorms, nobs
3331 
3332  do it=1, nobs
3333  read(unit, *) lon_deg, lat_deg, mps, slp, yr, month, day, hour
3334 ! if ( yr /= year ) then
3335 ! if(master) write(*, *) 'Extended to year + 1', yr
3336 ! endif
3337  cald = calday(yr, month, day, hour, 0, 0)
3338  time_tc(it,nstorms) = cald
3339  if(nudge_debug .and. master) write(*, 100) cald, month, day, hour, lon_deg, lat_deg, mps, slp
3340 
3341  wind_obs(it,nstorms) = mps ! m/s
3342  mslp_obs(it,nstorms) = 100.*slp
3343  y_obs(it,nstorms) = lat_deg * deg2rad
3344  x_obs(it,nstorms) = lon_deg * deg2rad
3345  enddo
3346 
3347  read(unit, *) ts_name, nobs, yr, month, day, hour
3348  enddo
3349 100 format(1x, f9.2, 1x, i3, 1x, i2, 1x, i2, 1x, f6.1, 1x, f6.1, 1x, f4.1, 1x, f6.1)
3350 
3351  else
3352 
3353  do while ( month /= 0 )
3354 
3355  read(unit, *) month, day, hour, gmt, lat_deg, lon_deg, mph, islp, comment
3356 
3357  select case (month)
3358 
3359  case (99) ! New storm record to start
3360  nstorms = nstorms + 1
3361  nobs = 0
3362  if(master) write(*, *) 'Reading data for TC#', nstorms, comment
3363  case ( 0) ! end of record
3364  if(master) write(*, *) 'End of record reached'
3365  case default
3366  nobs = nobs + 1
3367  cald = calday(year, month, day, hour, 0, 0)
3368  time_tc(nobs,nstorms) = cald
3369  nobs_tc(nstorms) = nobs ! observation count for this storm
3370 
3371  if(master) write(*, 200) nobs, cald, month, day, hour, gmt, lat_deg, lon_deg, mph, islp, comment
3372  mslp_obs(nobs,nstorms) = 100. * real(islp)
3373  y_obs(nobs,nstorms) = lat_deg * deg2rad
3374  if ( gmt == 'GMT' ) then
3375 ! Transfrom x from (-180 , 180) to (0, 360) then to radian
3376  if ( lon_deg < 0 ) then
3377  x_obs(nobs,nstorms) = (360.+lon_deg) * deg2rad
3378  else
3379  x_obs(nobs,nstorms) = (360.-lon_deg) * deg2rad
3380  endif
3381  elseif ( gmt == 'PAC' ) then ! Pacific storms
3382  x_obs(nobs,nstorms) = lon_deg * deg2rad
3383  endif
3384  end select
3385 
3386  enddo
3387 
3388  endif
3389 
3390  close(unit)
3391 
3392  if(master) then
3393  write(*,*) 'TC vortex breeding: total storms=', nstorms
3394  if ( nstorms/=0 ) then
3395  do n=1,nstorms
3396  write(*, *) 'TC#',n, ' contains ', nobs_tc(n),' observations'
3397  enddo
3398  endif
3399  endif
3400 
3401 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)
3402 
3403  end subroutine slp_obs_init
3404 
3406  ! Julian day (0 to 365 for non-leap year)
3407  real function calday(year, month, day, hour, minute, sec)
3408 ! input:
3409  integer, intent(in):: year, month, day, hour
3410  integer, intent(in):: minute, sec
3411 ! Local:
3412  integer n, m, ds, nday
3413  real tsec
3414  integer days(12)
3415  data days /31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31/
3416 
3417  ds = day - 1
3418 
3419  if( month /= 1 ) then
3420  do m=1, month-1
3421  if( m==2 .and. leap_year(year) ) then
3422  ds = ds + 29
3423  else
3424  ds = ds + days(m)
3425  endif
3426  enddo
3427  endif
3428 
3429  if ( leap_year(year_track_data) ) then
3430  nday = 366
3431  else
3432  nday = 365
3433  endif
3434 
3435  calday = real((year-year_track_data)*nday + ds) + real(hour*3600 + minute*60 + sec)/86400.
3436 
3437  end function calday
3438 
3441  logical function leap_year(ny)
3442  integer, intent(in):: ny
3443  integer ny00
3444 !
3445 ! No leap years prior to 0000
3446 !
3447  parameter( ny00 = 0000 )
3448 
3449  if( ny >= ny00 ) then
3450  if( mod(ny,100) == 0. .and. mod(ny,400) == 0. ) then
3451  leap_year = .true.
3452  elseif( mod(ny,4) == 0. .and. mod(ny,100) /= 0. ) then
3453  leap_year = .true.
3454  else
3455  leap_year = .false.
3456  endif
3457  else
3458  leap_year = .false.
3459  endif
3460 
3461  end function leap_year
3462 
3463 
3464  subroutine pmaxmin( qname, a, imax, jmax, fac )
3466  character(len=*) qname
3467  integer imax, jmax
3468  integer i, j
3469  real a(imax,jmax)
3470 
3471  real qmin(jmax), qmax(jmax)
3472  real pmax, pmin
3473  real fac ! multiplication factor
3474 
3475  do j=1,jmax
3476  pmax = a(1,j)
3477  pmin = a(1,j)
3478  do i=2,imax
3479  pmax = max(pmax, a(i,j))
3480  pmin = min(pmin, a(i,j))
3481  enddo
3482  qmax(j) = pmax
3483  qmin(j) = pmin
3484  enddo
3485 !
3486 ! Now find max/min of amax/amin
3487 !
3488  pmax = qmax(1)
3489  pmin = qmin(1)
3490  do j=2,jmax
3491  pmax = max(pmax, qmax(j))
3492  pmin = min(pmin, qmin(j))
3493  enddo
3494 
3495  write(*,*) qname, ' max = ', pmax*fac, ' min = ', pmin*fac
3496 
3497  end subroutine pmaxmin
3498 
3500  subroutine del2_uv(du, dv, cd, kmd, ntimes, bd, npx, npy, gridstruct, domain)
3501  integer, intent(in):: kmd
3502  integer, intent(in):: ntimes
3503  real, intent(in):: cd
3504  type(fv_grid_bounds_type), intent(IN) :: bd
3505  real, intent(inout), dimension(is:ie,js:je,kmd):: du, dv
3506  integer, intent(IN) :: npx, npy
3507  type(fv_grid_type), intent(IN), target :: gridstruct
3508  type(domain2d), intent(INOUT) :: domain
3509 ! local:
3510  real(kind=R_GRID), pointer, dimension(:,:,:) :: vlon, vlat
3511  real, dimension(is:ie,js:je,kmd):: v1, v2, v3
3512  integer i,j,k
3513 
3514  vlon => gridstruct%vlon
3515  vlat => gridstruct%vlat
3516 
3517 ! transform to 3D Cartesian:
3518 !$OMP parallel do default(none) shared(kmd,is,ie,js,je,v1,v2,v3,du,vlon,dv,vlat)
3519  do k=1,kmd
3520  do j=js,je
3521  do i=is,ie
3522  v1(i,j,k) = du(i,j,k)*vlon(i,j,1) + dv(i,j,k)*vlat(i,j,1)
3523  v2(i,j,k) = du(i,j,k)*vlon(i,j,2) + dv(i,j,k)*vlat(i,j,2)
3524  v3(i,j,k) = du(i,j,k)*vlon(i,j,3) + dv(i,j,k)*vlat(i,j,3)
3525  enddo
3526  enddo
3527  enddo
3528 
3529 ! Filter individual components as scalar:
3530  call del2_scalar( v1(is,js,1), cd, kmd, ntimes, bd, npx, npy, gridstruct, domain )
3531  call del2_scalar( v2(is,js,1), cd, kmd, ntimes, bd, npx, npy, gridstruct, domain )
3532  call del2_scalar( v3(is,js,1), cd, kmd, ntimes, bd, npx, npy, gridstruct, domain )
3533 
3534 ! Convert back to lat-lon components:
3535 !$OMP parallel do default(none) shared(kmd,is,ie,js,je,du,dv,v1,v2,v3,vlon,vlat)
3536  do k=1,kmd
3537  do j=js,je
3538  do i=is,ie
3539  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)
3540  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)
3541  enddo
3542  enddo
3543  enddo
3544 
3545  end subroutine del2_uv
3546 
3548  subroutine del2_scalar(qdt, cd, kmd, nmax, bd, npx, npy, gridstruct, domain)
3549  integer, intent(in):: kmd
3550  integer, intent(in):: nmax
3551  real, intent(in):: cd
3552  type(fv_grid_bounds_type), intent(IN) :: bd
3553  real, intent(inout):: qdt(is:ie,js:je,kmd)
3554  integer, intent(IN) :: npx, npy
3555  type(fv_grid_type), intent(IN), target :: gridstruct
3556  type(domain2d), intent(INOUT) :: domain
3557 ! local:
3558  real:: q(isd:ied,jsd:jed,kmd)
3559  real:: fx(isd:ied+1,jsd:jed), fy(isd:ied,jsd:jed+1)
3560  integer i,j,k, n, nt, ntimes
3561  real :: damp
3562 
3563  real, pointer, dimension(:,:) :: rarea, area
3564  real, pointer, dimension(:,:) :: sina_u, sina_v
3565  real, pointer, dimension(:,:,:) :: sin_sg
3566 
3567  real, pointer, dimension(:,:) :: dx, dy, rdxc, rdyc
3568 
3569  real(kind=R_GRID), pointer :: da_min
3570 
3571  logical, pointer :: bounded_domain, sw_corner, se_corner, nw_corner, ne_corner
3572 
3573  area => gridstruct%area
3574  rarea => gridstruct%rarea
3575 
3576  sina_u => gridstruct%sina_u
3577  sina_v => gridstruct%sina_v
3578  sin_sg => gridstruct%sin_sg
3579 
3580  dx => gridstruct%dx
3581  dy => gridstruct%dy
3582  rdxc => gridstruct%rdxc
3583  rdyc => gridstruct%rdyc
3584 
3585  da_min => gridstruct%da_min
3586 
3587  bounded_domain => gridstruct%bounded_domain
3588  sw_corner => gridstruct%sw_corner
3589  se_corner => gridstruct%se_corner
3590  nw_corner => gridstruct%nw_corner
3591  ne_corner => gridstruct%ne_corner
3592 
3593  ntimes = min(3, nmax)
3594 
3595  damp = cd * da_min
3596 
3597 !$OMP parallel do default(none) shared(is,ie,js,je,kmd,q,qdt)
3598  do k=1,kmd
3599  do j=js,je
3600  do i=is,ie
3601  q(i,j,k) = qdt(i,j,k)
3602  enddo
3603  enddo
3604  enddo
3605  call timing_on('COMM_TOTAL')
3606  call mpp_update_domains(q, domain, complete=.true.)
3607  call timing_off('COMM_TOTAL')
3608 
3609  do n=1,ntimes
3610 
3611  nt = ntimes - n
3612 
3613 !$OMP parallel do default(none) shared(is,ie,js,je,kmd,nt,dy,q,isd,jsd,npx,npy,bounded_domain, &
3614 !$OMP bd,sw_corner,se_corner,nw_corner,ne_corner, &
3615 !$OMP sina_u,rdxc,sin_sg,dx,rdyc,sina_v,qdt,damp,rarea) &
3616 !$OMP private(fx, fy)
3617  do k=1,kmd
3618 
3619  if(nt>0) call copy_corners(q(isd,jsd,k), npx, npy, 1, bounded_domain, bd, &
3620  sw_corner, se_corner, nw_corner, ne_corner)
3621  do j=js-nt,je+nt
3622  do i=is-nt,ie+1+nt
3623  fx(i,j) = dy(i,j)*sina_u(i,j)*(q(i-1,j,k)-q(i,j,k))*rdxc(i,j)
3624  enddo
3625  if (is == 1) fx(i,j) = dy(is,j)*(q(is-1,j,k)-q(is,j,k))*rdxc(is,j)* &
3626  0.5*(sin_sg(1,j,1) + sin_sg(0,j,3))
3627  if (ie+1 == npx) fx(i,j) = dy(ie+1,j)*(q(ie,j,k)-q(ie+1,j,k))*rdxc(ie+1,j)* &
3628  0.5*(sin_sg(npx,j,1) + sin_sg(npx-1,j,3))
3629  enddo
3630 
3631  if(nt>0) call copy_corners(q(isd,jsd,k), npx, npy, 2, bounded_domain, bd, &
3632  sw_corner, se_corner, nw_corner, ne_corner)
3633  do j=js-nt,je+1+nt
3634  if (j == 1 .OR. j == npy) then
3635  do i=is-nt,ie+nt
3636  fy(i,j) = dx(i,j)*(q(i,j-1,k)-q(i,j,k))*rdyc(i,j) &
3637  *0.5*(sin_sg(i,j,2) + sin_sg(i,j-1,4) )
3638  enddo
3639  else
3640  do i=is-nt,ie+nt
3641  fy(i,j) = dx(i,j)*sina_v(i,j)*(q(i,j-1,k)-q(i,j,k))*rdyc(i,j)
3642  enddo
3643  end if
3644  enddo
3645 
3646  if ( nt==0 ) then
3647  do j=js,je
3648  do i=is,ie
3649  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))
3650  enddo
3651  enddo
3652  else
3653  do j=js-nt,je+nt
3654  do i=is-nt,ie+nt
3655  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))
3656  enddo
3657  enddo
3658  endif
3659  enddo
3660 
3661  enddo
3662 
3663  end subroutine del2_scalar
3664 
3665  subroutine rmse_bias(a, rms, bias, area)
3666  real, intent(in):: a(is:ie,js:je)
3667  real, intent(IN) :: area(isd:ied,jsd:jed)
3668  real, intent(out):: rms, bias
3669  integer:: i,j
3670  real:: total_area
3671 
3672  total_area = 4.*pi*radius**2
3673 
3674  rms = 0.
3675  bias = 0.
3676  do j=js,je
3677  do i=is,ie
3678  bias = bias + area(i,j) * a(i,j)
3679  rms = rms + area(i,j) * a(i,j)**2
3680  enddo
3681  enddo
3682  call mp_reduce_sum(bias)
3683  call mp_reduce_sum(rms)
3684 
3685  bias = bias / total_area
3686  rms = sqrt( rms / total_area )
3687 
3688  end subroutine rmse_bias
3689 
3690 
3691  subroutine corr(a, b, co, area)
3692  real, intent(in):: a(is:ie,js:je), b(is:ie,js:je)
3693  real, intent(in):: area(isd:ied,jsd:jed)
3694  real, intent(out):: co
3695  real:: m_a, m_b, std_a, std_b
3696  integer:: i,j
3697  real:: total_area
3698 
3699  total_area = 4.*pi*radius**2
3700 
3701 ! Compute standard deviation:
3702  call std(a, m_a, std_a, area)
3703  call std(b, m_b, std_b, area)
3704 
3705 ! Compute correlation:
3706  co = 0.
3707  do j=js,je
3708  do i=is,ie
3709  co = co + area(i,j) * (a(i,j)-m_a)*(b(i,j)-m_b)
3710  enddo
3711  enddo
3712  call mp_reduce_sum(co)
3713  co = co / (total_area*std_a*std_b )
3714 
3715  end subroutine corr
3716 
3717  subroutine std(a, mean, stdv, area)
3718  real, intent(in):: a(is:ie,js:je)
3719  real, intent(IN) :: area(isd:ied,jsd:jed)
3720  real, intent(out):: mean, stdv
3721  integer:: i,j
3722  real:: total_area
3723 
3724  total_area = 4.*pi*radius**2
3725 
3726  mean = 0.
3727  do j=js,je
3728  do i=is,ie
3729  mean = mean + area(i,j) * a(i,j)
3730  enddo
3731  enddo
3732  call mp_reduce_sum(mean)
3733  mean = mean / total_area
3734 
3735  stdv = 0.
3736  do j=js,je
3737  do i=is,ie
3738  stdv = stdv + area(i,j) * (a(i,j)-mean)**2
3739  enddo
3740  enddo
3741  call mp_reduce_sum(stdv)
3742  stdv = sqrt( stdv / total_area )
3743 
3744  end subroutine std
3745 
3746 
3747 end module fv_nwp_nudge_mod
real(kind=4), dimension(nobs_max, max_storm) mslp_out
outer ring SLP in mb
Definition: fv_nudge.F90:280
logical do_ps_bias
Definition: fv_nudge.F90:205
real, dimension(:,:), allocatable gz0
Definition: fv_nudge.F90:171
logical module_is_initialized
Definition: fv_nudge.F90:148
real, dimension(:,:,:), allocatable ps_dat
Definition: fv_nudge.F90:168
logical, public do_adiabatic_init
Definition: fv_nudge.F90:138
pure real function, public vicpqd(q)
real, dimension(:), allocatable lon
Definition: fv_nudge.F90:146
subroutine, public get_var3_r4(ncid, var3_name, is, ie, js, je, ks, ke, var3, time_slice)
Definition: sim_nc_mod.F90:246
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 &#39;breed_srf_winds&#39; performs vortex breeding by nudging 1-m winds.
Definition: fv_nudge.F90:2927
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...
Definition: fv_nudge.F90:176
subroutine compute_slp(isc, iec, jsc, jec, tm, ps, gz, slp)
Definition: fv_nudge.F90:1077
real mask_fac
[0,1] 0: no mask; 1: full strength
Definition: fv_nudge.F90:193
subroutine remap_tq(npz, ak, bk, ps, delp, t, q, kmd, ps0, ta, qa, zvir, ptop)
Definition: fv_nudge.F90:1930
The module &#39;fv_mp_mod&#39; is a single program multiple data (SPMD) parallel decompostion/communication m...
Definition: fv_mp_mod.F90:24
subroutine timing_off(blk_name)
The subroutine &#39;timing_off&#39; stops a timer.
Definition: fv_timing.F90:180
The type &#39;fv_grid_type&#39; is made up of grid-dependent information from fv_grid_tools and fv_grid_utils...
Definition: fv_arrays.F90:123
integer km
Data z-dimension.
Definition: fv_nudge.F90:144
subroutine del2_uv(du, dv, cd, kmd, ntimes, bd, npx, npy, gridstruct, domain)
The subroutine &#39;del2_uv&#39; filters the wind tendency.
Definition: fv_nudge.F90:3501
integer id_ht_err
Definition: fv_nudge.F90:284
subroutine, public normalize_vect(e)
The subroutine &#39;normalize_vect&#39; makes &#39;e&#39; a unit vector.
subroutine get_ncep_analysis(ps, u, v, t, q, zvir, ts, nfile, fname, bd)
Definition: fv_nudge.F90:1485
real(kind=4), dimension(:,:,:), allocatable gz3
work array
Definition: fv_nudge.F90:170
integer, dimension(:,:), allocatable jdc
Definition: fv_nudge.F90:167
subroutine, public mappm(km, pe1, q1, kn, pe2, q2, i1, i2, iv, kord, ptop)
The subroutine &#39;mappm&#39; is a general-purpose routine for remapping one set of vertical levels to anoth...
Definition: fv_mapz.F90:3281
integer, public num_gas
Definition: multi_gases.F90:46
subroutine, public intp_great_circle(beta, p1, p2, x_o, y_o)
real tau_ps
1-day
Definition: fv_nudge.F90:218
subroutine tangent_wind(elon, elat, speed, po, pp, ut, vt)
Definition: fv_nudge.F90:3157
The module &#39;multi_gases&#39; peforms multi constitutents computations.
Definition: multi_gases.F90:25
logical nudge_hght
Definition: fv_nudge.F90:210
subroutine std(a, mean, stdv, area)
Definition: fv_nudge.F90:3718
real, dimension(:,:,:), allocatable s2c
Definition: fv_nudge.F90:166
real(kind=4), dimension(nobs_max, max_storm) time_tc
start time of the track
Definition: fv_nudge.F90:282
real(kind=4), dimension(nobs_max, max_storm) mslp_obs
observed SLP in mb
Definition: fv_nudge.F90:279
real, dimension(:), allocatable bk0
Definition: fv_nudge.F90:145
logical time_varying
Definition: fv_nudge.F90:211
real slp_env
storm environment pressure (pa)
Definition: fv_nudge.F90:252
real pre0_env
critical storm environment pressure (pa) for size computation
Definition: fv_nudge.F90:253
logical strong_mask
Definition: fv_nudge.F90:202
real(kind=r_grid), parameter radius
Definition: fv_nudge.F90:132
integer min_nobs
Definition: fv_nudge.F90:274
real p_norelax
from P_norelax upwards, no nudging
Definition: fv_nudge.F90:184
real(kind=4), dimension(:,:,:,:), allocatable q_dat
Definition: fv_nudge.F90:169
integer, dimension(max_storm) nobs_tc
Definition: fv_nudge.F90:273
integer, parameter, public r_grid
Definition: fv_arrays.F90:34
integer kord_data
Definition: fv_nudge.F90:191
integer jm
Data y-dimension.
Definition: fv_nudge.F90:143
real(kind=4), dimension(nobs_max, max_storm) wind_obs
observed 10-m wind speed (m/s)
Definition: fv_nudge.F90:278
real, parameter del_r
Definition: fv_nudge.F90:262
integer, dimension(:,:), allocatable id2
Definition: fv_nudge.F90:167
integer year_track_data
Definition: fv_nudge.F90:268
pure real function, public virqd(q)
real tau_q
1-day
Definition: fv_nudge.F90:219
The module fv_nwp_nudge contains routines for nudging to input analyses. note This module is currentl...
Definition: fv_nudge.F90:32
subroutine, public get_var1_real(ncid, var1_name, im, var1, var_exist)
Definition: sim_nc_mod.F90:116
integer im
Data x-dimension.
Definition: fv_nudge.F90:142
The module &#39;fv_timing&#39; contains FV3 timers.
Definition: fv_timing.F90:24
logical nudge_winds
Definition: fv_nudge.F90:208
pure real function, public virq(q)
logical, public t_is_tv
Definition: fv_nudge.F90:195
subroutine get_int_hght(npz, ak, bk, ps, delp, ps0, tv)
Definition: fv_nudge.F90:1894
The module &#39;tp_core&#39; is a collection of routines to support FV transport.
Definition: tp_core.F90:24
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 &#39;fv_nwp_nudge&#39; computes and returns time tendencies for nudging to analysis...
Definition: fv_nudge.F90:306
real(kind=4), dimension(:,:,:,:), allocatable v_dat
Definition: fv_nudge.F90:169
subroutine, public get_var1_double(ncid, var1_name, im, var1, var_exist)
The &#39;get_var&#39; subroutines read in variables from netcdf files.
Definition: sim_nc_mod.F90:93
character(len=128) track_file_name
Definition: fv_nudge.F90:188
subroutine remap_coef(agrid)
Definition: fv_nudge.F90:1732
logical use_high_top
Definition: fv_nudge.F90:197
The module &#39;fv_mapz&#39; contains the vertical mapping routines .
Definition: fv_mapz.F90:27
subroutine del2_scalar(qdt, cd, kmd, nmax, bd, npx, npy, gridstruct, domain)
The subroutine &#39;del2_scalar&#39; filters the physics tendency.
Definition: fv_nudge.F90:3549
logical breed_srf_w
Definition: fv_nudge.F90:245
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 &#39;breed_slp_inline&#39; performs vortex breeding by nudging sea level pressure toward singl...
Definition: fv_nudge.F90:2156
subroutine ps_bias_correction(ps_dt, is, ie, js, je, isd, ied, jsd, jed, area)
Definition: fv_nudge.F90:976
subroutine, public fv_nwp_nudge_end
The subroutine &#39;fv_nwp_nudge_end&#39; terminates the nudging module.
Definition: fv_nudge.F90:2084
The module &#39;fv_arrays&#39; contains the &#39;fv_atmos_type&#39; and associated datatypes.
Definition: fv_arrays.F90:24
real function, public great_circle_dist(q1, q2, radius)
subroutine pmaxmin(qname, a, imax, jmax, fac)
Definition: fv_nudge.F90:3465
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)
Definition: fv_nudge.F90:3186
The module &#39;sim_nc&#39; is a netcdf file reader.
Definition: sim_nc_mod.F90:28
logical print_end_breed
Definition: fv_nudge.F90:212
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. ...
Definition: fv_nudge.F90:1044
real dps_min
maximum PS increment (pa; each step) due to inline breeding
Definition: fv_nudge.F90:163
logical use_pt_inc
Definition: fv_nudge.F90:196
real, dimension(:), allocatable lat
Definition: fv_nudge.F90:146
real, dimension(:), allocatable ak0
Definition: fv_nudge.F90:145
integer, parameter nobs_max
Max # of observations per storm.
Definition: fv_nudge.F90:270
logical nudge_virt
Definition: fv_nudge.F90:209
subroutine remap_uv(npz, ak, bk, ps, delp, u, v, kmd, ps0, u0, v0, ptop)
Definition: fv_nudge.F90:2009
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)
Definition: fv_nudge.F90:1094
type(time_type), public fv_time
subroutine, public open_ncfile(iflnm, ncid)
Definition: sim_nc_mod.F90:55
subroutine ncep2fms(sst)
Definition: fv_nudge.F90:1808
subroutine, public vect_cross(e, p1, p2)
The subroutine &#39;vect_cross performs cross products of 3D vectors: e = P1 X P2.
subroutine timing_on(blk_name)
The subroutine &#39;timing_on&#39; starts a timer.
Definition: fv_timing.F90:116
logical nudge_debug
Definition: fv_nudge.F90:204
real p_wvp
cutoff level for specific humidity nudging
Definition: fv_nudge.F90:190
real(kind=4), dimension(nobs_max, max_storm) y_obs
latitude in degrees
Definition: fv_nudge.F90:277
integer, parameter nfile_max
maximum: 20-year analysis data, 4*366*20=29280
Definition: fv_nudge.F90:156
real p_relax
from P_relax upwards, nudging is reduced linearly proportional to pfull/P_relax
Definition: fv_nudge.F90:181
logical print_end_nudge
Definition: fv_nudge.F90:213
subroutine corr(a, b, co, area)
Definition: fv_nudge.F90:3692
real, parameter tm_max
Definition: fv_nudge.F90:254
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...
Definition: fv_nudge.F90:178
integer time_interval
dataset time interval (seconds)
Definition: fv_nudge.F90:153
@ The module &#39;fv_diagnostics&#39; contains routines to compute diagnosic fields.
integer nfile_total
=5 for 1-day (if datasets are 6-hr apart)
Definition: fv_nudge.F90:189
The module &#39;fv_grid_utils&#39; contains routines for setting up and computing grid-related quantities...
logical conserve_mom
Definition: fv_nudge.F90:199
subroutine, public prt_maxmin(qname, q, is, ie, js, je, n_g, km, fac)
logical conserve_hgt
Definition: fv_nudge.F90:200
real nudged_time
seconds
Definition: fv_nudge.F90:264
real function calday(year, month, day, hour, minute, sec)
The function &#39;calday&#39; performs time interpolation:
Definition: fv_nudge.F90:3408
subroutine, public copy_corners(q, npx, npy, dir, bounded_domain, bd, sw_corner, se_corner, nw_corner, ne_corner)
Definition: tp_core.F90:248
real(kind=4), dimension(:,:,:,:), allocatable t_dat
Definition: fv_nudge.F90:169
subroutine, public fv_nwp_nudge_init(time, axes, npz, zvir, ak, bk, ts, phis, gridstruct, ks, npx, neststruct, bd)
The subroutine &#39;fv_nwp_nudge_init&#39; initializes the nudging module.
Definition: fv_nudge.F90:1256
integer, parameter max_storm
max # of storms to process
Definition: fv_nudge.F90:269
logical nudge_ps
Definition: fv_nudge.F90:206
real(kind=4), dimension(nobs_max, max_storm) rad_out
outer ring radius in meters
Definition: fv_nudge.F90:281
subroutine get_tc_mask(time, mask, agrid)
Definition: fv_nudge.F90:2110
real(kind=4), dimension(nobs_max, max_storm) x_obs
longitude in degrees
Definition: fv_nudge.F90:276
subroutine, public get_ncdim1(ncid, var1_name, im)
Definition: sim_nc_mod.F90:78
subroutine slp_obs_init
Definition: fv_nudge.F90:3278
subroutine breed_srf_w10(time, dt, npz, ak, bk, ps, phis, slp, delp, u, v, gridstruct)
The subroutine &#39;breed_srf_w10&#39; performs vortex breeding by nudging 10-m winds.
Definition: fv_nudge.F90:2648
logical add_bg_wind
Definition: fv_nudge.F90:198
integer kbot_winds
Definition: fv_nudge.F90:236
real tau_winds
6-hr
Definition: fv_nudge.F90:220
character(len=128), dimension(nfile_max) file_names
Definition: fv_nudge.F90:187
integer, dimension(:,:), allocatable id1
Definition: fv_nudge.F90:167
character(len=128) input_fname_list
a file lists the input NCEP analysis data
Definition: fv_nudge.F90:175
logical analysis_time
Definition: fv_nudge.F90:239
logical function leap_year(ny)
Definition: fv_nudge.F90:3442
subroutine rmse_bias(a, rms, bias, area)
Definition: fv_nudge.F90:3666
real(kind=4), dimension(:,:,:,:), allocatable u_dat
Definition: fv_nudge.F90:169
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)
Definition: fv_nudge.F90:795
subroutine, public latlon2xyz(p, e, id)
The subroutine &#39;latlon2xyz&#39; maps (lon, lat) to (x,y,z)
subroutine, public close_ncfile(ncid)
Definition: sim_nc_mod.F90:67