FV3DYCORE  Version 2.0.0
fv_treat_da_inc.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 
35 !
36 ! REVISION HISTORY:
37 ! 02/12/2016 - Initial Version
38 !-------------------------------------------------------------------------------
39 
40 #ifdef OVERLOAD_R4
41 #define _GET_VAR1 get_var1_real
42 #else
43 #define _GET_VAR1 get_var1_double
44 #endif
45 
47 
48 ! <table>
49 ! <tr>
50 ! <th>Module Name</th>
51 ! <th>Functions Included</th>
52 ! </tr>
53 ! <tr>
54 ! <td>constants_mod</td>
55 ! <td>pi=>pi_8, omega, grav, kappa, rdgas, rvgas, cp_air</td>
56 ! </tr>
57 ! <tr>
58 ! <td>field_manager_mod</td>
59 ! <td>MODEL_ATMOS</td>
60 ! </tr>
61 ! <tr>
62 ! <td>fms_mod</td>
63 ! <td>file_exist, open_namelist_file,close_file, error_mesg, FATAL,
64 ! check_nml_error, stdlog,write_version_number,set_domain,
65 ! mpp_clock_id, mpp_clock_begin, mpp_clock_end, CLOCK_SUBCOMPONENT,
66 ! clock_flag_default, nullify_domain</td>
67 ! </tr>
68 ! <tr>
69 ! <td>fv_arrays_mod</td>
70 ! <td>fv_atmos_type, fv_grid_type, fv_grid_bounds_type,
71 ! R_GRID</td>
72 ! </tr>
73 ! <tr>
74 ! <td>fv_control_mod</td>
75 ! <td>fv_init, fv_end, ngrids</td>
76 ! </tr>
77 ! <tr>
78 ! <td>fv_grid_utils_mod</td>
79 ! <td>ptop_min, g_sum, mid_pt_sphere, get_unit_vect2,
80 ! get_latlon_vector, inner_prod, cubed_to_latlon</td>
81 ! </tr>
82 ! <tr>
83 ! <td>fms_mod</td>
84 ! <td>file_exist, read_data, field_exist, write_version_number</td>
85 ! </tr>
86 ! <tr>
87 ! <td>fv_mp_mod</td>
88 ! <td>ng,is_master,fill_corners,YDir,mp_reduce_min, mp_reduce_max</td>
89 ! </tr>
90 ! <tr>
91 ! <td>mpp_mod</td>
92 ! <td>mpp_error, FATAL, NOTE, mpp_pe</td>
93 ! </tr>
94 ! <tr>
95 ! <td>mpp_domains_mod</td>
96 ! <td>mpp_get_tile_id, domain2d, mpp_update_domains,NORTH, EAST</td>
97 ! </tr>
98 ! <tr>
99 ! <td>sim_nc_mod</td>
100 ! <td>open_ncfile, close_ncfile, get_ncdim1, get_var1_double,
101 ! get_var2_real, get_var3_r4, get_var1_real</td>
102 ! </tr>
103 ! <tr>
104 ! <td>tracer_manager_mod</td>
105 ! <td>get_tracer_names, get_number_tracers, get_tracer_index</td>
106 ! </tr>
107 ! </table>
108 
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, &
113  domain2d, &
114  mpp_update_domains, &
115  north, &
116  east
117  use tracer_manager_mod,only: get_tracer_names, &
118  get_number_tracers, &
119  get_tracer_index
120  use field_manager_mod, only: model_atmos
121 
122  use constants_mod, only: pi=>pi_8, omega, grav, kappa, &
123  rdgas, rvgas, cp_air
124  use fv_arrays_mod, only: fv_atmos_type, &
125  fv_grid_type, &
127  r_grid
128  use fv_grid_utils_mod, only: ptop_min, g_sum, &
132  use fv_mp_mod, only: is_master, &
133  fill_corners, &
134  ydir, &
135  mp_reduce_min, &
136  mp_reduce_max
137  use sim_nc_mod, only: open_ncfile, &
138  close_ncfile, &
139  get_ncdim1, &
140  get_var1_double, &
141  get_var2_real, &
142  get_var3_r4, &
143  get_var1_real, &
145  implicit none
146  private
147 
148  public :: read_da_inc,remap_coef
149 
150 contains
151  !=============================================================================
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 )
160  type(fv_atmos_type), intent(inout) :: Atm
161  type(domain2d), intent(inout) :: fv_domain
162  type(fv_grid_bounds_type), intent(IN) :: bd
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 ! D grid zonal wind (m/s)
165  real, intent(inout), dimension(is_in:ie_in+1,js_in:je_in ,npz_in):: v ! D grid meridional wind (m/s)
166  real, intent(inout) :: delp(is_in:ie_in ,js_in:je_in ,npz_in) ! pressure thickness (pascal)
167  real, intent(inout) :: pt( is_in:ie_in ,js_in:je_in ,npz_in) ! temperature (K)
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) !
170  ! local
171 
172  real :: deg2rad
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):: &
184  id1, id2, jdc
185  integer, dimension(Atm%bd%is:Atm%bd%ie+1,Atm%bd%js:Atm%bd%je)::&
186  id1_c, id2_c, jdc_c
187  integer, dimension(Atm%bd%is:Atm%bd%ie,Atm%bd%js:Atm%bd%je+1)::&
188  id1_d, id2_d, jdc_d
189 
190  integer:: i, j, k, im, jm, km, npt
191  integer:: i1, i2, j1, ncid
192  integer:: jbeg, jend
193  integer tsize(3)
194  real(kind=R_GRID), dimension(2):: p1, p2, p3
195  real(kind=R_GRID), dimension(3):: e1, e2, ex, ey
196 
197  logical:: found
198  integer :: is, ie, js, je
199  integer :: isd, ied, jsd, jed
200  integer :: sphum, liq_wat
201 #ifdef MULTI_GASES
202  integer :: spfo, spfo2, spfo3
203 #else
204  integer :: o3mr
205 #endif
206 
207  is = atm%bd%is
208  ie = atm%bd%ie
209  js = atm%bd%js
210  je = atm%bd%je
211  isd = atm%bd%isd
212  ied = atm%bd%ied
213  jsd = atm%bd%jsd
214  jed = atm%bd%jed
215 
216  deg2rad = pi/180.
217 
218  fname = 'INPUT/'//atm%flagstruct%res_latlon_dynamics
219 
220  if( file_exist(fname) ) then
221  call open_ncfile( fname, ncid ) ! open the file
222  call get_ncdim1( ncid, 'lon', tsize(1) )
223  call get_ncdim1( ncid, 'lat', tsize(2) )
224  call get_ncdim1( ncid, 'lev', tsize(3) )
225 
226  im = tsize(1); jm = tsize(2); km = tsize(3)
227 
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')
232  endif
233 
234  if(is_master()) write(*,*) fname, ' DA increment dimensions:', tsize
235 
236  allocate ( lon(im) )
237  allocate ( lat(jm) )
238 
239  call _get_var1 (ncid, 'lon', im, lon )
240  call _get_var1 (ncid, 'lat', jm, lat )
241 
242  ! Convert to radian
243  do i=1,im
244  lon(i) = lon(i) * deg2rad ! lon(1) = 0.
245  enddo
246  do j=1,jm
247  lat(j) = lat(j) * deg2rad
248  enddo
249 
250  else
251  call mpp_error(fatal,'==> Error in read_da_inc: Expected file '&
252  //trim(fname)//' for DA increment does not exist')
253  endif
254 
255  ! Initialize lat-lon to Cubed bi-linear interpolation coeff:
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)
259 
260  ! Find bounding latitudes:
261  jbeg = jm-1; jend = 2
262  do j=js,je
263  do i=is,ie
264  j1 = jdc(i,j)
265  jbeg = min(jbeg, j1)
266  jend = max(jend, j1+1)
267  enddo
268  enddo
269 
270  sphum = get_tracer_index(model_atmos, 'sphum')
271 #ifdef MULTI_GASES
272  spfo3 = get_tracer_index(model_atmos, 'spfo3')
273  spfo = get_tracer_index(model_atmos, 'spfo')
274  spfo2 = get_tracer_index(model_atmos, 'spfo2')
275 #else
276  o3mr = get_tracer_index(model_atmos, 'o3mr')
277 #endif
278  liq_wat = get_tracer_index(model_atmos, 'liq_wat')
279 
280  ! perform increments on scalars
281  allocate ( wk3(1:im,jbeg:jend, 1:km) )
282  allocate ( tp(is:ie,js:je,km) )
283 
284  call apply_inc_on_3d_scalar('T_inc',pt, is_in, js_in, ie_in, je_in)
285  call apply_inc_on_3d_scalar('delp_inc',delp, is_in, js_in, ie_in, je_in)
286  if ( .not. atm%flagstruct%hydrostatic ) then
287  call apply_inc_on_3d_scalar('delz_inc',delz, isc_in, jsc_in, iec_in, jec_in)
288  endif
289  call apply_inc_on_3d_scalar('sphum_inc',q(:,:,:,sphum), is_in, js_in, ie_in, je_in)
290  call apply_inc_on_3d_scalar('liq_wat_inc',q(:,:,:,liq_wat), is_in, js_in, ie_in, je_in)
291 #ifdef MULTI_GASES
292  call apply_inc_on_3d_scalar('spfo3_inc',q(:,:,:,spfo3), is_in, js_in, ie_in, je_in)
293  call apply_inc_on_3d_scalar('spfo_inc',q(:,:,:,spfo), is_in, js_in, ie_in, je_in)
294  call apply_inc_on_3d_scalar('spfo2_inc',q(:,:,:,spfo2), is_in, js_in, ie_in, je_in)
295 #else
296  call apply_inc_on_3d_scalar('o3mr_inc',q(:,:,:,o3mr), is_in, js_in, ie_in, je_in)
297 #endif
298 
299  deallocate ( tp )
300  deallocate ( wk3 )
301 
302  ! perform increments on winds
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))
307 
308  call get_staggered_grid( &
309  is, ie, js, je, &
310  isd, ied, jsd, jed, &
311  atm%gridstruct%grid, pt_c, pt_d)
312 
313  !------ pt_c part ------
314  ! Initialize lat-lon to Cubed bi-linear interpolation coeff:
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, &
317  pt_c)
318 
319  ! Find bounding latitudes:
320  jbeg = jm-1; jend = 2
321  do j=js,je
322  do i=is,ie+1
323  j1 = jdc_c(i,j)
324  jbeg = min(jbeg, j1)
325  jend = max(jend, j1+1)
326  enddo
327  enddo
328 
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) )
333 
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 )
336 
337  do k=1,km
338  do j=js,je
339  do i=is,ie+1
340  i1 = id1_c(i,j)
341  i2 = id2_c(i,j)
342  j1 = jdc_c(i,j)
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)
353  call mid_pt_sphere(p1, p2, p3)
354  call get_unit_vect2(p1, p2, e2)
355  call get_latlon_vector(p3, ex, ey)
356  vd_inc(i,j,k) = u_inc(i,j,k)*inner_prod(e2,ex) + &
357  v_inc(i,j,k)*inner_prod(e2,ey)
358  v(i,j,k) = v(i,j,k) + vd_inc(i,j,k)
359  enddo
360  enddo
361  enddo
362 
363  deallocate ( u_inc, v_inc )
364  deallocate ( wk3_u, wk3_v )
365 
366  !------ pt_d part ------
367  ! Initialize lat-lon to Cubed bi-linear interpolation coeff:
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, &
370  pt_d)
371 
372  ! Find bounding latitudes:
373  jbeg = jm-1; jend = 2
374  do j=js,je+1
375  do i=is,ie
376  j1 = jdc_d(i,j)
377  jbeg = min(jbeg, j1)
378  jend = max(jend, j1+1)
379  enddo
380  enddo
381 
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) )
386 
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 )
389 
390  do k=1,km
391  do j=js,je+1
392  do i=is,ie
393  i1 = id1_d(i,j)
394  i2 = id2_d(i,j)
395  j1 = jdc_d(i,j)
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)
406  call mid_pt_sphere(p1, p2, p3)
407  call get_unit_vect2(p1, p2, e1)
408  call get_latlon_vector(p3, ex, ey)
409  ud_inc(i,j,k) = u_inc(i,j,k)*inner_prod(e1,ex) + &
410  v_inc(i,j,k)*inner_prod(e1,ey)
411  u(i,j,k) = u(i,j,k) + ud_inc(i,j,k)
412  enddo
413  enddo
414  enddo
415 
416  deallocate ( u_inc, v_inc )
417  deallocate ( wk3_u, wk3_v )
418 
419 !rab The following is not necessary as ua/va will be re-calculated during model startup
420 !rab call cubed_to_latlon(Atm%u, Atm%v, Atm%ua, Atm%va, &
421 !rab Atm%gridstruct, Atm%flagstruct%npx, Atm%flagstruct%npy, &
422 !rab Atm%flagstruct%npz, 1, Atm%gridstruct%grid_type, &
423 !rab fv_domain, Atm%gridstruct%nested, &
424 !rab Atm%flagstruct%c2l_ord, Atm%bd)
425 
426  !------ winds clean up ------
427  deallocate ( pt_c, pt_d, ud_inc, vd_inc )
428  !------ all clean up ------
429  deallocate ( lon, lat )
430 
431  contains
432  !---------------------------------------------------------------------------
435  subroutine apply_inc_on_3d_scalar(field_name,var, is_in, js_in, ie_in, je_in)
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
439  integer :: ierr
440 
441  call check_var_exists(ncid, field_name, ierr)
442  if (ierr == 0) then
443  call get_var3_r4( ncid, field_name, 1,im, jbeg,jend, 1,km, wk3 )
444  else
445  if (is_master()) print *,'warning: no increment for ',trim(field_name),' found, assuming zero'
446  wk3 = 0.
447  endif
448 
449  do k=1,km
450  do j=js,je
451  do i=is,ie
452  i1 = id1(i,j)
453  i2 = id2(i,j)
454  j1 = jdc(i,j)
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)
458  enddo
459  enddo
460  enddo
461 
462  end subroutine apply_inc_on_3d_scalar
463  !---------------------------------------------------------------------------
464  end subroutine read_da_inc
465 
466  !=============================================================================
468 
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)
478  ! local:
479  real :: rdlon(im)
480  real :: rdlat(jm)
481  real:: a1, b1
482  integer i,j, i1, i2, jc, i0, j0
483  do i=1,im-1
484  rdlon(i) = 1. / (lon(i+1) - lon(i))
485  enddo
486  rdlon(im) = 1. / (lon(1) + 2.*pi - lon(im))
487 
488  do j=1,jm-1
489  rdlat(j) = 1. / (lat(j+1) - lat(j))
490  enddo
491 
492  ! * Interpolate to cubed sphere cell center
493  do 5000 j=js,je
494 
495  do i=is,ie
496 
497  if ( agrid(i,j,1)>lon(im) ) then
498  i1 = im; i2 = 1
499  a1 = (agrid(i,j,1)-lon(im)) * rdlon(im)
500  elseif ( agrid(i,j,1)<lon(1) ) then
501  i1 = im; i2 = 1
502  a1 = (agrid(i,j,1)+2.*pi-lon(im)) * rdlon(im)
503  else
504  do i0=1,im-1
505  if ( agrid(i,j,1)>=lon(i0) .and. agrid(i,j,1)<=lon(i0+1) ) then
506  i1 = i0; i2 = i0+1
507  a1 = (agrid(i,j,1)-lon(i1)) * rdlon(i0)
508  go to 111
509  endif
510  enddo
511  endif
512 111 continue
513 
514  if ( agrid(i,j,2)<lat(1) ) then
515  jc = 1
516  b1 = 0.
517  elseif ( agrid(i,j,2)>lat(jm) ) then
518  jc = jm-1
519  b1 = 1.
520  else
521  do j0=1,jm-1
522  if ( agrid(i,j,2)>=lat(j0) .and. agrid(i,j,2)<=lat(j0+1) ) then
523  jc = j0
524  b1 = (agrid(i,j,2)-lat(jc)) * rdlat(jc)
525  go to 222
526  endif
527  enddo
528  endif
529 222 continue
530 
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
533  endif
534 
535  s2c(i,j,1) = (1.-a1) * (1.-b1)
536  s2c(i,j,2) = a1 * (1.-b1)
537  s2c(i,j,3) = a1 * b1
538  s2c(i,j,4) = (1.-a1) * b1
539  id1(i,j) = i1
540  id2(i,j) = i2
541  jdc(i,j) = jc
542  enddo !i-loop
543 5000 continue ! j-loop
544 
545  end subroutine remap_coef
546 
547  !=============================================================================
550  subroutine get_staggered_grid( &
551  is, ie, js, je, &
552  isd, ied, jsd, jed, &
553  pt_b, pt_c, pt_d)
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
558  ! local
559  real(kind=R_GRID), dimension(2):: p1, p2, p3
560  integer :: i, j
561 
562  do j = js,je+1
563  do i = is,ie
564  p1(:) = pt_b(i, j,1:2)
565  p2(:) = pt_b(i+1,j,1:2)
566  call mid_pt_sphere(p1, p2, p3)
567  pt_d(i,j,1:2) = p3(:)
568  enddo
569  enddo
570  do j = js,je
571  do i = is,ie+1
572  p1(:) = pt_b(i,j ,1:2)
573  p2(:) = pt_b(i,j+1,1:2)
574  call mid_pt_sphere(p1, p2, p3)
575  pt_c(i,j,1:2) = p3(:)
576  enddo
577  enddo
578 
579  end subroutine get_staggered_grid
580  !=============================================================================
581 end module fv_treat_da_inc_mod
582 
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)
Definition: sim_nc_mod.F90:246
The module &#39;fv_mp_mod&#39; is a single program multiple data (SPMD) parallel decompostion/communication m...
Definition: fv_mp_mod.F90:24
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
real function, public g_sum(domain, p, ifirst, ilast, jfirst, jlast, ngc, area, mode, reproduce)
The function &#39;g_sum&#39; is the fast version of &#39;globalsum&#39;.
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
Definition: fv_arrays.F90:34
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 &#39;read_da_inc&#39; reads the increments of the diagnostic variables from the DA-generated f...
subroutine, public get_var1_real(ncid, var1_name, im, var1, var_exist)
Definition: sim_nc_mod.F90:116
subroutine apply_inc_on_3d_scalar(field_name, var, is_in, js_in, ie_in, je_in)
The subroutine &#39;apply_inc_on3d_scalar&#39; applies the input increments to the prognostic variables...
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
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 &#39;fv_arrays&#39; contains the &#39;fv_atmos_type&#39; and associated datatypes.
Definition: fv_arrays.F90:24
The module &#39;sim_nc&#39; is a netcdf file reader.
Definition: sim_nc_mod.F90:28
subroutine, public get_var2_real(ncid, var_name, im, jm, var2)
Definition: sim_nc_mod.F90:137
&#39;The module &#39;tread_da_increment&#39; contains routines for treating the increments of the prognostic vari...
subroutine, public open_ncfile(iflnm, ncid)
Definition: sim_nc_mod.F90:55
subroutine get_staggered_grid(is, ie, js, je, isd, ied, jsd, jed, pt_b, pt_c, pt_d)
The subroutine &#39;get_staggered_grid&#39; gets the lat-lon locations of the staggered points.
The module &#39;fv_grid_utils&#39; 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 &#39;remap_coef&#39; calculates the coefficients for horizonal regridding. ...
subroutine, public get_ncdim1(ncid, var1_name, im)
Definition: sim_nc_mod.F90:78
subroutine, public check_var_exists(ncid, var_name, status)
Definition: sim_nc_mod.F90:238
subroutine, public close_ncfile(ncid)
Definition: sim_nc_mod.F90:67