41 #define _GET_VAR1 get_var1_real 43 #define _GET_VAR1 get_var1_double 109 use fms_mod
, only: file_exist, read_data, &
110 field_exist, write_version_number
111 use mpp_mod
, only: mpp_error, fatal, note, mpp_pe
112 use mpp_domains_mod
, only: mpp_get_tile_id, &
114 mpp_update_domains, &
117 use tracer_manager_mod
,only: get_tracer_names, &
118 get_number_tracers, &
120 use field_manager_mod
, only: model_atmos
122 use constants_mod
, only: pi=>pi_8, omega, grav, kappa, &
158 subroutine read_da_inc(Atm, fv_domain, bd, npz_in, nq, u, v, q, delp, pt, delz, is_in, js_in, ie_in, je_in, &
159 isc_in, jsc_in, iec_in, jec_in )
161 type(domain2d),
intent(inout) :: fv_domain
163 integer,
intent(IN) :: npz_in, nq, is_in, js_in, ie_in, je_in, isc_in, jsc_in, iec_in, jec_in
164 real,
intent(inout),
dimension(is_in:ie_in, js_in:je_in+1,npz_in):: u
165 real,
intent(inout),
dimension(is_in:ie_in+1,js_in:je_in ,npz_in):: v
166 real,
intent(inout) :: delp(is_in:ie_in ,js_in:je_in ,npz_in)
167 real,
intent(inout) :: pt( is_in:ie_in ,js_in:je_in ,npz_in)
168 real,
intent(inout) :: q( is_in:ie_in ,js_in:je_in ,npz_in, nq)
169 real,
intent(inout) :: delz(isc_in:iec_in ,jsc_in:jec_in ,npz_in)
173 character(len=128) :: fname
174 real(kind=4),
allocatable:: wk1(:), wk2(:,:), wk3(:,:,:)
175 real(kind=4),
allocatable:: wk3_u(:,:,:), wk3_v(:,:,:)
176 real,
allocatable:: tp(:,:,:), qp(:,:,:)
177 real,
dimension(:,:,:),
allocatable:: u_inc, v_inc, ud_inc, vd_inc
178 real,
allocatable:: lat(:), lon(:)
179 real,
allocatable:: pt_c(:,:,:), pt_d(:,:,:)
180 real:: s2c(atm%bd%is:atm%bd%ie,atm%bd%js:atm%bd%je,4)
181 real:: s2c_c(atm%bd%is:atm%bd%ie+1,atm%bd%js:atm%bd%je,4)
182 real:: s2c_d(atm%bd%is:atm%bd%ie,atm%bd%js:atm%bd%je+1,4)
183 integer,
dimension(Atm%bd%is:Atm%bd%ie,Atm%bd%js:Atm%bd%je):: &
185 integer,
dimension(Atm%bd%is:Atm%bd%ie+1,Atm%bd%js:Atm%bd%je)::&
187 integer,
dimension(Atm%bd%is:Atm%bd%ie,Atm%bd%js:Atm%bd%je+1)::&
190 integer:: i, j, k, im, jm, km, npt
191 integer:: i1, i2, j1, ncid
194 real(kind=R_GRID),
dimension(2):: p1, p2, p3
195 real(kind=R_GRID),
dimension(3):: e1, e2, ex, ey
198 integer :: is, ie, js, je
199 integer :: isd, ied, jsd, jed
200 integer :: sphum, liq_wat
202 integer :: spfo, spfo2, spfo3
218 fname =
'INPUT/'//atm%flagstruct%res_latlon_dynamics
220 if( file_exist(fname) )
then 226 im = tsize(1); jm = tsize(2); km = tsize(3)
228 if (km.ne.npz_in)
then 229 if (is_master()) print *,
'km = ', km
230 call mpp_error(fatal, &
231 '==> Error in read_da_inc: km is not equal to npz_in')
234 if(is_master())
write(*,*) fname,
' DA increment dimensions:', tsize
239 call _get_var1 (ncid,
'lon', im, lon )
240 call _get_var1 (ncid,
'lat', jm, lat )
244 lon(i) = lon(i) * deg2rad
247 lat(j) = lat(j) * deg2rad
251 call mpp_error(fatal,
'==> Error in read_da_inc: Expected file '&
252 //trim(fname)//
' for DA increment does not exist')
256 call remap_coef( is, ie, js, je, isd, ied, jsd, jed, &
257 im, jm, lon, lat, id1, id2, jdc, s2c, &
258 atm%gridstruct%agrid)
261 jbeg = jm-1; jend = 2
266 jend = max(jend, j1+1)
270 sphum = get_tracer_index(model_atmos,
'sphum')
272 spfo3 = get_tracer_index(model_atmos,
'spfo3')
273 spfo = get_tracer_index(model_atmos,
'spfo')
274 spfo2 = get_tracer_index(model_atmos,
'spfo2')
276 o3mr = get_tracer_index(model_atmos,
'o3mr')
278 liq_wat = get_tracer_index(model_atmos,
'liq_wat')
281 allocate ( wk3(1:im,jbeg:jend, 1:km) )
282 allocate ( tp(is:ie,js:je,km) )
286 if ( .not. atm%flagstruct%hydrostatic )
then 303 allocate (pt_c(isd:ied+1,jsd:jed ,2))
304 allocate (pt_d(isd:ied ,jsd:jed+1,2))
305 allocate (ud_inc(is:ie , js:je+1, km))
306 allocate (vd_inc(is:ie+1, js:je , km))
310 isd, ied, jsd, jed, &
311 atm%gridstruct%grid, pt_c, pt_d)
315 call remap_coef( is, ie+1, js, je, isd, ied+1, jsd, jed, &
316 im, jm, lon, lat, id1_c, id2_c, jdc_c, s2c_c, &
320 jbeg = jm-1; jend = 2
325 jend = max(jend, j1+1)
329 allocate ( wk3_u(1:im,jbeg:jend, 1:km) )
330 allocate ( wk3_v(1:im,jbeg:jend, 1:km) )
331 allocate ( u_inc(is:ie+1,js:je,km) )
332 allocate ( v_inc(is:ie+1,js:je,km) )
334 call get_var3_r4( ncid,
'u_inc', 1,im, jbeg,jend, 1,km, wk3_u )
335 call get_var3_r4( ncid,
'v_inc', 1,im, jbeg,jend, 1,km, wk3_v )
343 u_inc(i,j,k) = s2c_c(i,j,1)*wk3_u(i1,j1 ,k) + &
344 s2c_c(i,j,2)*wk3_u(i2,j1 ,k) + &
345 s2c_c(i,j,3)*wk3_u(i2,j1+1,k) + &
346 s2c_c(i,j,4)*wk3_u(i1,j1+1,k)
347 v_inc(i,j,k) = s2c_c(i,j,1)*wk3_v(i1,j1 ,k) + &
348 s2c_c(i,j,2)*wk3_v(i2,j1 ,k) + &
349 s2c_c(i,j,3)*wk3_v(i2,j1+1,k) + &
350 s2c_c(i,j,4)*wk3_v(i1,j1+1,k)
351 p1(:) = atm%gridstruct%grid(i,j ,1:2)
352 p2(:) = atm%gridstruct%grid(i,j+1,1:2)
356 vd_inc(i,j,k) = u_inc(i,j,k)*
inner_prod(e2,ex) + &
358 v(i,j,k) = v(i,j,k) + vd_inc(i,j,k)
363 deallocate ( u_inc, v_inc )
364 deallocate ( wk3_u, wk3_v )
368 call remap_coef( is, ie, js, je+1, isd, ied, jsd, jed+1, &
369 im, jm, lon, lat, id1_d, id2_d, jdc_d, s2c_d, &
373 jbeg = jm-1; jend = 2
378 jend = max(jend, j1+1)
382 allocate ( wk3_u(1:im,jbeg:jend, 1:km) )
383 allocate ( wk3_v(1:im,jbeg:jend, 1:km) )
384 allocate ( u_inc(is:ie,js:je+1,km) )
385 allocate ( v_inc(is:ie,js:je+1,km) )
387 call get_var3_r4( ncid,
'u_inc', 1,im, jbeg,jend, 1,km, wk3_u )
388 call get_var3_r4( ncid,
'v_inc', 1,im, jbeg,jend, 1,km, wk3_v )
396 u_inc(i,j,k) = s2c_d(i,j,1)*wk3_u(i1,j1 ,k) + &
397 s2c_d(i,j,2)*wk3_u(i2,j1 ,k) + &
398 s2c_d(i,j,3)*wk3_u(i2,j1+1,k) + &
399 s2c_d(i,j,4)*wk3_u(i1,j1+1,k)
400 v_inc(i,j,k) = s2c_d(i,j,1)*wk3_v(i1,j1 ,k) + &
401 s2c_d(i,j,2)*wk3_v(i2,j1 ,k) + &
402 s2c_d(i,j,3)*wk3_v(i2,j1+1,k) + &
403 s2c_d(i,j,4)*wk3_v(i1,j1+1,k)
404 p1(:) = atm%gridstruct%grid(i, j,1:2)
405 p2(:) = atm%gridstruct%grid(i+1,j,1:2)
409 ud_inc(i,j,k) = u_inc(i,j,k)*
inner_prod(e1,ex) + &
411 u(i,j,k) = u(i,j,k) + ud_inc(i,j,k)
416 deallocate ( u_inc, v_inc )
417 deallocate ( wk3_u, wk3_v )
427 deallocate ( pt_c, pt_d, ud_inc, vd_inc )
429 deallocate ( lon, lat )
436 character(len=*),
intent(in) :: field_name
437 integer,
intent(IN) :: is_in, js_in, ie_in, je_in
438 real,
dimension(is_in:ie_in,js_in:je_in,1:km),
intent(inout) :: var
443 call get_var3_r4( ncid, field_name, 1,im, jbeg,jend, 1,km, wk3 )
445 if (is_master()) print *,
'warning: no increment for ',trim(field_name),
' found, assuming zero' 455 tp(i,j,k) = s2c(i,j,1)*wk3(i1,j1 ,k) + s2c(i,j,2)*wk3(i2,j1 ,k)+&
456 s2c(i,j,3)*wk3(i2,j1+1,k) + s2c(i,j,4)*wk3(i1,j1+1,k)
457 var(i,j,k) = var(i,j,k)+tp(i,j,k)
469 subroutine remap_coef( is, ie, js, je, isd, ied, jsd, jed, &
470 im, jm, lon, lat, id1, id2, jdc, s2c, agrid )
472 integer,
intent(in):: is, ie, js, je, isd, ied, jsd, jed
473 integer,
intent(in):: im, jm
474 real,
intent(in):: lon(im), lat(jm)
475 real,
intent(out):: s2c(is:ie,js:je,4)
476 integer,
intent(out),
dimension(is:ie,js:je):: id1, id2, jdc
477 real,
intent(in):: agrid(isd:ied,jsd:jed,2)
482 integer i,j, i1, i2, jc, i0, j0
484 rdlon(i) = 1. / (lon(i+1) - lon(i))
486 rdlon(im) = 1. / (lon(1) + 2.*pi - lon(im))
489 rdlat(j) = 1. / (lat(j+1) - lat(j))
497 if ( agrid(i,j,1)>lon(im) )
then 499 a1 = (agrid(i,j,1)-lon(im)) * rdlon(im)
500 elseif ( agrid(i,j,1)<lon(1) )
then 502 a1 = (agrid(i,j,1)+2.*pi-lon(im)) * rdlon(im)
505 if ( agrid(i,j,1)>=lon(i0) .and. agrid(i,j,1)<=lon(i0+1) )
then 507 a1 = (agrid(i,j,1)-lon(i1)) * rdlon(i0)
514 if ( agrid(i,j,2)<lat(1) )
then 517 elseif ( agrid(i,j,2)>lat(jm) )
then 522 if ( agrid(i,j,2)>=lat(j0) .and. agrid(i,j,2)<=lat(j0+1) )
then 524 b1 = (agrid(i,j,2)-lat(jc)) * rdlat(jc)
531 if ( a1<0.0 .or. a1>1.0 .or. b1<0.0 .or. b1>1.0 )
then 532 write(*,*)
'gid=', mpp_pe(), i,j,a1, b1
535 s2c(i,j,1) = (1.-a1) * (1.-b1)
536 s2c(i,j,2) = a1 * (1.-b1)
538 s2c(i,j,4) = (1.-a1) * b1
552 isd, ied, jsd, jed, &
554 integer,
intent(in) :: is, ie, js, je, isd, ied, jsd, jed
555 real,
dimension(isd:ied+1,jsd:jed+1,2),
intent(in) :: pt_b
556 real,
dimension(isd:ied+1,jsd:jed ,2),
intent(out) :: pt_c
557 real,
dimension(isd:ied ,jsd:jed+1,2),
intent(out) :: pt_d
559 real(kind=R_GRID),
dimension(2):: p1, p2, p3
564 p1(:) = pt_b(i, j,1:2)
565 p2(:) = pt_b(i+1,j,1:2)
567 pt_d(i,j,1:2) = p3(:)
572 p1(:) = pt_b(i,j ,1:2)
573 p2(:) = pt_b(i,j+1,1:2)
575 pt_c(i,j,1:2) = p3(:)
subroutine, public mid_pt_sphere(p1, p2, pm)
subroutine, public get_var3_r4(ncid, var3_name, is, ie, js, je, ks, ke, var3, time_slice)
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
The type 'fv_grid_type' is made up of grid-dependent information from fv_grid_tools and fv_grid_utils...
real function, public g_sum(domain, p, ifirst, ilast, jfirst, jlast, ngc, area, mode, reproduce)
The function 'g_sum' is the fast version of 'globalsum'.
real function, public inner_prod(v1, v2)
subroutine, public get_latlon_vector(pp, elon, elat)
subroutine, public get_unit_vect2(e1, e2, uc)
integer, parameter, public r_grid
subroutine, public read_da_inc(Atm, fv_domain, bd, npz_in, nq, u, v, q, delp, pt, delz, is_in, js_in, ie_in, je_in, isc_in, jsc_in, iec_in, jec_in)
The subroutine 'read_da_inc' reads the increments of the diagnostic variables from the DA-generated f...
subroutine, public get_var1_real(ncid, var1_name, im, var1, var_exist)
subroutine apply_inc_on_3d_scalar(field_name, var, is_in, js_in, ie_in, je_in)
The subroutine 'apply_inc_on3d_scalar' applies the input increments to the prognostic variables...
subroutine, public get_var1_double(ncid, var1_name, im, var1, var_exist)
The 'get_var' subroutines read in variables from netcdf files.
real, parameter, public ptop_min
subroutine, public cubed_to_latlon(u, v, ua, va, gridstruct, npx, npy, km, mode, grid_type, domain, bounded_domain, c2l_ord, bd)
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
The module 'sim_nc' is a netcdf file reader.
subroutine, public get_var2_real(ncid, var_name, im, jm, var2)
'The module 'tread_da_increment' contains routines for treating the increments of the prognostic vari...
subroutine, public open_ncfile(iflnm, ncid)
subroutine get_staggered_grid(is, ie, js, je, isd, ied, jsd, jed, pt_b, pt_c, pt_d)
The subroutine 'get_staggered_grid' gets the lat-lon locations of the staggered points.
The module 'fv_grid_utils' contains routines for setting up and computing grid-related quantities...
subroutine, public remap_coef(is, ie, js, je, isd, ied, jsd, jed, im, jm, lon, lat, id1, id2, jdc, s2c, agrid)
The subroutine 'remap_coef' calculates the coefficients for horizonal regridding. ...
subroutine, public get_ncdim1(ncid, var1_name, im)
subroutine, public check_var_exists(ncid, var_name, status)
subroutine, public close_ncfile(ncid)