FV3DYCORE  Version1.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: ng, &
133  is_master, &
134  fill_corners, &
135  ydir, &
136  mp_reduce_min, &
137  mp_reduce_max
138  use sim_nc_mod, only: open_ncfile, &
139  close_ncfile, &
140  get_ncdim1, &
141  get_var1_double, &
142  get_var2_real, &
143  get_var3_r4, &
144  get_var1_real, &
146  implicit none
147  private
148 
149  public :: read_da_inc,remap_coef
150 
151 contains
152  !=============================================================================
159  subroutine read_da_inc(Atm, fv_domain)
160  type(fv_atmos_type), intent(inout) :: atm(:)
161  type(domain2d), intent(inout) :: fv_domain
162  ! local
163  integer :: nq
164 
165  real :: deg2rad
166  character(len=128) :: fname
167  real(kind=4), allocatable:: wk1(:), wk2(:,:), wk3(:,:,:)
168  real(kind=4), allocatable:: wk3_u(:,:,:), wk3_v(:,:,:)
169  real, allocatable:: tp(:,:,:), qp(:,:,:)
170  real, dimension(:,:,:), allocatable:: u_inc, v_inc, ud_inc, vd_inc
171  real, allocatable:: lat(:), lon(:)
172  real, allocatable:: pt_c(:,:,:), pt_d(:,:,:)
173  real:: s2c(atm(1)%bd%is:atm(1)%bd%ie,atm(1)%bd%js:atm(1)%bd%je,4)
174  real:: s2c_c(atm(1)%bd%is:atm(1)%bd%ie+1,atm(1)%bd%js:atm(1)%bd%je,4)
175  real:: s2c_d(atm(1)%bd%is:atm(1)%bd%ie,atm(1)%bd%js:atm(1)%bd%je+1,4)
176  integer, dimension(Atm(1)%bd%is:Atm(1)%bd%ie,Atm(1)%bd%js:Atm(1)%bd%je):: &
177  id1, id2, jdc
178  integer, dimension(Atm(1)%bd%is:Atm(1)%bd%ie+1,Atm(1)%bd%js:Atm(1)%bd%je)::&
179  id1_c, id2_c, jdc_c
180  integer, dimension(Atm(1)%bd%is:Atm(1)%bd%ie,Atm(1)%bd%js:Atm(1)%bd%je+1)::&
181  id1_d, id2_d, jdc_d
182 
183  integer:: i, j, k, im, jm, km, npz, npt
184  integer:: i1, i2, j1, ncid
185  integer:: jbeg, jend
186  integer tsize(3)
187  real(kind=R_GRID), dimension(2):: p1, p2, p3
188  real(kind=R_GRID), dimension(3):: e1, e2, ex, ey
189 
190  logical:: found
191  integer :: is, ie, js, je
192  integer :: isd, ied, jsd, jed
193  integer :: sphum, liq_wat
194 #ifdef MULTI_GASES
195  integer :: spfo, spfo2, spfo3
196 #else
197  integer :: o3mr
198 #endif
199 
200  is = atm(1)%bd%is
201  ie = atm(1)%bd%ie
202  js = atm(1)%bd%js
203  je = atm(1)%bd%je
204  isd = atm(1)%bd%isd
205  ied = atm(1)%bd%ied
206  jsd = atm(1)%bd%jsd
207  jed = atm(1)%bd%jed
208 
209  deg2rad = pi/180.
210 
211  npz = atm(1)%npz
212 
213  fname = 'INPUT/'//atm(1)%flagstruct%res_latlon_dynamics
214 
215  if( file_exist(fname) ) then
216  call open_ncfile( fname, ncid ) ! open the file
217  call get_ncdim1( ncid, 'lon', tsize(1) )
218  call get_ncdim1( ncid, 'lat', tsize(2) )
219  call get_ncdim1( ncid, 'lev', tsize(3) )
220 
221  im = tsize(1); jm = tsize(2); km = tsize(3)
222 
223  if (km.ne.npz) then
224  if (is_master()) print *, 'km = ', km
225  call mpp_error(fatal, &
226  '==> Error in read_da_inc: km is not equal to npz')
227  endif
228 
229  if(is_master()) write(*,*) fname, ' DA increment dimensions:', tsize
230 
231  allocate ( lon(im) )
232  allocate ( lat(jm) )
233 
234  call _get_var1 (ncid, 'lon', im, lon )
235  call _get_var1 (ncid, 'lat', jm, lat )
236 
237  ! Convert to radian
238  do i=1,im
239  lon(i) = lon(i) * deg2rad ! lon(1) = 0.
240  enddo
241  do j=1,jm
242  lat(j) = lat(j) * deg2rad
243  enddo
244 
245  else
246  call mpp_error(fatal,'==> Error in read_da_inc: Expected file '&
247  //trim(fname)//' for DA increment does not exist')
248  endif
249 
250  ! Initialize lat-lon to Cubed bi-linear interpolation coeff:
251  call remap_coef( is, ie, js, je, &
252  im, jm, lon, lat, id1, id2, jdc, s2c, &
253  atm(1)%gridstruct%agrid(is:ie,js:je,:))
254 
255  ! Find bounding latitudes:
256  jbeg = jm-1; jend = 2
257  do j=js,je
258  do i=is,ie
259  j1 = jdc(i,j)
260  jbeg = min(jbeg, j1)
261  jend = max(jend, j1+1)
262  enddo
263  enddo
264 
265  sphum = get_tracer_index(model_atmos, 'sphum')
266 #ifdef MULTI_GASES
267  spfo3 = get_tracer_index(model_atmos, 'spfo3')
268  spfo = get_tracer_index(model_atmos, 'spfo')
269  spfo2 = get_tracer_index(model_atmos, 'spfo2')
270 #else
271  o3mr = get_tracer_index(model_atmos, 'o3mr')
272 #endif
273  liq_wat = get_tracer_index(model_atmos, 'liq_wat')
274 
275  ! perform increments on scalars
276  allocate ( wk3(1:im,jbeg:jend, 1:km) )
277  allocate ( tp(is:ie,js:je,km) )
278 
279  call apply_inc_on_3d_scalar('T_inc',atm(1)%pt)
280  call apply_inc_on_3d_scalar('delp_inc',atm(1)%delp)
281  if (.not. atm(1)%flagstruct%hydrostatic) then
282  call apply_inc_on_3d_scalar('delz_inc',atm(1)%delz)
283  endif
284  call apply_inc_on_3d_scalar('sphum_inc',atm(1)%q(:,:,:,sphum))
285  call apply_inc_on_3d_scalar('liq_wat_inc',atm(1)%q(:,:,:,liq_wat))
286 #ifdef MULTI_GASES
287  call apply_inc_on_3d_scalar('spfo3_inc',atm(1)%q(:,:,:,spfo3))
288  call apply_inc_on_3d_scalar('spfo_inc',atm(1)%q(:,:,:,spfo))
289  call apply_inc_on_3d_scalar('spfo2_inc',atm(1)%q(:,:,:,spfo2))
290 #else
291  call apply_inc_on_3d_scalar('o3mr_inc',atm(1)%q(:,:,:,o3mr))
292 #endif
293 
294  deallocate ( tp )
295  deallocate ( wk3 )
296 
297  ! perform increments on winds
298  allocate (pt_c(is:ie+1,js:je ,2))
299  allocate (pt_d(is:ie ,js:je+1,2))
300  allocate (ud_inc(is:ie , js:je+1, km))
301  allocate (vd_inc(is:ie+1, js:je , km))
302 
303  call get_staggered_grid( &
304  is, ie, js, je, &
305  isd, ied, jsd, jed, &
306  atm(1)%gridstruct%grid, pt_c, pt_d)
307 
308  !------ pt_c part ------
309  ! Initialize lat-lon to Cubed bi-linear interpolation coeff:
310  call remap_coef( is, ie+1, js, je, &
311  im, jm, lon, lat, id1_c, id2_c, jdc_c, s2c_c, &
312  pt_c)
313 
314  ! Find bounding latitudes:
315  jbeg = jm-1; jend = 2
316  do j=js,je
317  do i=is,ie+1
318  j1 = jdc_c(i,j)
319  jbeg = min(jbeg, j1)
320  jend = max(jend, j1+1)
321  enddo
322  enddo
323 
324  allocate ( wk3_u(1:im,jbeg:jend, 1:km) )
325  allocate ( wk3_v(1:im,jbeg:jend, 1:km) )
326  allocate ( u_inc(is:ie+1,js:je,km) )
327  allocate ( v_inc(is:ie+1,js:je,km) )
328 
329  call get_var3_r4( ncid, 'u_inc', 1,im, jbeg,jend, 1,km, wk3_u )
330  call get_var3_r4( ncid, 'v_inc', 1,im, jbeg,jend, 1,km, wk3_v )
331 
332  do k=1,km
333  do j=js,je
334  do i=is,ie+1
335  i1 = id1_c(i,j)
336  i2 = id2_c(i,j)
337  j1 = jdc_c(i,j)
338  u_inc(i,j,k) = s2c_c(i,j,1)*wk3_u(i1,j1 ,k) + &
339  s2c_c(i,j,2)*wk3_u(i2,j1 ,k) + &
340  s2c_c(i,j,3)*wk3_u(i2,j1+1,k) + &
341  s2c_c(i,j,4)*wk3_u(i1,j1+1,k)
342  v_inc(i,j,k) = s2c_c(i,j,1)*wk3_v(i1,j1 ,k) + &
343  s2c_c(i,j,2)*wk3_v(i2,j1 ,k) + &
344  s2c_c(i,j,3)*wk3_v(i2,j1+1,k) + &
345  s2c_c(i,j,4)*wk3_v(i1,j1+1,k)
346  p1(:) = atm(1)%gridstruct%grid(i,j ,1:2)
347  p2(:) = atm(1)%gridstruct%grid(i,j+1,1:2)
348  call mid_pt_sphere(p1, p2, p3)
349  call get_unit_vect2(p1, p2, e2)
350  call get_latlon_vector(p3, ex, ey)
351  vd_inc(i,j,k) = u_inc(i,j,k)*inner_prod(e2,ex) + &
352  v_inc(i,j,k)*inner_prod(e2,ey)
353  atm(1)%v(i,j,k) = atm(1)%v(i,j,k) + vd_inc(i,j,k)
354  enddo
355  enddo
356  enddo
357 
358  deallocate ( u_inc, v_inc )
359  deallocate ( wk3_u, wk3_v )
360 
361  !------ pt_d part ------
362  ! Initialize lat-lon to Cubed bi-linear interpolation coeff:
363  call remap_coef( is, ie, js, je+1,&
364  im, jm, lon, lat, id1_d, id2_d, jdc_d, s2c_d, &
365  pt_d)
366 
367  ! Find bounding latitudes:
368  jbeg = jm-1; jend = 2
369  do j=js,je+1
370  do i=is,ie
371  j1 = jdc_d(i,j)
372  jbeg = min(jbeg, j1)
373  jend = max(jend, j1+1)
374  enddo
375  enddo
376 
377  allocate ( wk3_u(1:im,jbeg:jend, 1:km) )
378  allocate ( wk3_v(1:im,jbeg:jend, 1:km) )
379  allocate ( u_inc(is:ie,js:je+1,km) )
380  allocate ( v_inc(is:ie,js:je+1,km) )
381 
382  call get_var3_r4( ncid, 'u_inc', 1,im, jbeg,jend, 1,km, wk3_u )
383  call get_var3_r4( ncid, 'v_inc', 1,im, jbeg,jend, 1,km, wk3_v )
384 
385  do k=1,km
386  do j=js,je+1
387  do i=is,ie
388  i1 = id1_d(i,j)
389  i2 = id2_d(i,j)
390  j1 = jdc_d(i,j)
391  u_inc(i,j,k) = s2c_d(i,j,1)*wk3_u(i1,j1 ,k) + &
392  s2c_d(i,j,2)*wk3_u(i2,j1 ,k) + &
393  s2c_d(i,j,3)*wk3_u(i2,j1+1,k) + &
394  s2c_d(i,j,4)*wk3_u(i1,j1+1,k)
395  v_inc(i,j,k) = s2c_d(i,j,1)*wk3_v(i1,j1 ,k) + &
396  s2c_d(i,j,2)*wk3_v(i2,j1 ,k) + &
397  s2c_d(i,j,3)*wk3_v(i2,j1+1,k) + &
398  s2c_d(i,j,4)*wk3_v(i1,j1+1,k)
399  p1(:) = atm(1)%gridstruct%grid(i, j,1:2)
400  p2(:) = atm(1)%gridstruct%grid(i+1,j,1:2)
401  call mid_pt_sphere(p1, p2, p3)
402  call get_unit_vect2(p1, p2, e1)
403  call get_latlon_vector(p3, ex, ey)
404  ud_inc(i,j,k) = u_inc(i,j,k)*inner_prod(e1,ex) + &
405  v_inc(i,j,k)*inner_prod(e1,ey)
406  atm(1)%u(i,j,k) = atm(1)%u(i,j,k) + ud_inc(i,j,k)
407  enddo
408  enddo
409  enddo
410 
411  deallocate ( u_inc, v_inc )
412  deallocate ( wk3_u, wk3_v )
413 
414 !rab The following is not necessary as ua/va will be re-calculated during model startup
415 !rab call cubed_to_latlon(Atm(1)%u, Atm(1)%v, Atm(1)%ua, Atm(1)%va, &
416 !rab Atm(1)%gridstruct, Atm(1)%flagstruct%npx, Atm(1)%flagstruct%npy, &
417 !rab Atm(1)%flagstruct%npz, 1, Atm(1)%gridstruct%grid_type, &
418 !rab fv_domain, Atm(1)%gridstruct%nested, &
419 !rab Atm(1)%flagstruct%c2l_ord, Atm(1)%bd)
420 
421  !------ winds clean up ------
422  deallocate ( pt_c, pt_d, ud_inc, vd_inc )
423  !------ all clean up ------
424  deallocate ( lon, lat )
425 
426  contains
427  !---------------------------------------------------------------------------
430  subroutine apply_inc_on_3d_scalar(field_name,var)
431  character(len=*), intent(in) :: field_name
432  real, dimension(isd:ied,jsd:jed,1:km), intent(inout) :: var
433  integer :: ierr
434 
435  call check_var_exists(ncid, field_name, ierr)
436  if (ierr == 0) then
437  call get_var3_r4( ncid, field_name, 1,im, jbeg,jend, 1,km, wk3 )
438  else
439  if (is_master()) print *,'warning: no increment for ',trim(field_name),' found, assuming zero'
440  wk3 = 0.
441  endif
442 
443  do k=1,km
444  do j=js,je
445  do i=is,ie
446  i1 = id1(i,j)
447  i2 = id2(i,j)
448  j1 = jdc(i,j)
449  tp(i,j,k) = s2c(i,j,1)*wk3(i1,j1 ,k) + s2c(i,j,2)*wk3(i2,j1 ,k)+&
450  s2c(i,j,3)*wk3(i2,j1+1,k) + s2c(i,j,4)*wk3(i1,j1+1,k)
451  var(i,j,k) = var(i,j,k)+tp(i,j,k)
452  enddo
453  enddo
454  enddo
455 
456  end subroutine apply_inc_on_3d_scalar
457  !---------------------------------------------------------------------------
458  end subroutine read_da_inc
459 
460  !=============================================================================
462 
463  subroutine remap_coef( is, ie, js, je,&
464  im, jm, lon, lat, id1, id2, jdc, s2c, agrid )
466  integer, intent(in):: is, ie, js, je
467  integer, intent(in):: im, jm
468  real, intent(in):: lon(im), lat(jm)
469  real, intent(out):: s2c(is:ie,js:je,4)
470  integer, intent(out), dimension(is:ie,js:je):: id1, id2, jdc
471  real, intent(in):: agrid(is:ie,js:je,2)
472  ! local:
473  real :: rdlon(im)
474  real :: rdlat(jm)
475  real:: a1, b1
476  integer i,j, i1, i2, jc, i0, j0
477  do i=1,im-1
478  rdlon(i) = 1. / (lon(i+1) - lon(i))
479  enddo
480  rdlon(im) = 1. / (lon(1) + 2.*pi - lon(im))
481 
482  do j=1,jm-1
483  rdlat(j) = 1. / (lat(j+1) - lat(j))
484  enddo
485 
486  ! * Interpolate to cubed sphere cell center
487  do 5000 j=js,je
488 
489  do i=is,ie
490 
491  if ( agrid(i,j,1)>lon(im) ) then
492  i1 = im; i2 = 1
493  a1 = (agrid(i,j,1)-lon(im)) * rdlon(im)
494  elseif ( agrid(i,j,1)<lon(1) ) then
495  i1 = im; i2 = 1
496  a1 = (agrid(i,j,1)+2.*pi-lon(im)) * rdlon(im)
497  else
498  do i0=1,im-1
499  if ( agrid(i,j,1)>=lon(i0) .and. agrid(i,j,1)<=lon(i0+1) ) then
500  i1 = i0; i2 = i0+1
501  a1 = (agrid(i,j,1)-lon(i1)) * rdlon(i0)
502  go to 111
503  endif
504  enddo
505  endif
506 111 continue
507 
508  if ( agrid(i,j,2)<lat(1) ) then
509  jc = 1
510  b1 = 0.
511  elseif ( agrid(i,j,2)>lat(jm) ) then
512  jc = jm-1
513  b1 = 1.
514  else
515  do j0=1,jm-1
516  if ( agrid(i,j,2)>=lat(j0) .and. agrid(i,j,2)<=lat(j0+1) ) then
517  jc = j0
518  b1 = (agrid(i,j,2)-lat(jc)) * rdlat(jc)
519  go to 222
520  endif
521  enddo
522  endif
523 222 continue
524 
525  if ( a1<0.0 .or. a1>1.0 .or. b1<0.0 .or. b1>1.0 ) then
526  write(*,*) 'gid=', mpp_pe(), i,j,a1, b1
527  endif
528 
529  s2c(i,j,1) = (1.-a1) * (1.-b1)
530  s2c(i,j,2) = a1 * (1.-b1)
531  s2c(i,j,3) = a1 * b1
532  s2c(i,j,4) = (1.-a1) * b1
533  id1(i,j) = i1
534  id2(i,j) = i2
535  jdc(i,j) = jc
536  enddo !i-loop
537 5000 continue ! j-loop
538 
539  end subroutine remap_coef
540 
541  !=============================================================================
544  subroutine get_staggered_grid( &
545  is, ie, js, je, &
546  isd, ied, jsd, jed, &
547  pt_b, pt_c, pt_d)
548  integer, intent(in) :: is, ie, js, je
549  integer, intent(in) :: isd, ied, jsd, jed
550  real, dimension(isd:ied+1,jsd:jed+1,2), intent(in) :: pt_b
551  real, dimension(is:ie+1,js:je ,2), intent(out) :: pt_c
552  real, dimension(is:ie ,js:je+1,2), intent(out) :: pt_d
553  ! local
554  real(kind=R_GRID), dimension(2):: p1, p2, p3
555  integer :: i, j
556 
557  do j = js,je+1
558  do i = is,ie
559  p1(:) = pt_b(i, j,1:2)
560  p2(:) = pt_b(i+1,j,1:2)
561  call mid_pt_sphere(p1, p2, p3)
562  pt_d(i,j,1:2) = p3(:)
563  enddo
564  enddo
565  do j = js,je
566  do i = is,ie+1
567  p1(:) = pt_b(i,j ,1:2)
568  p2(:) = pt_b(i,j+1,1:2)
569  call mid_pt_sphere(p1, p2, p3)
570  pt_c(i,j,1:2) = p3(:)
571  enddo
572  enddo
573 
574  end subroutine get_staggered_grid
575  !=============================================================================
576 end module fv_treat_da_inc_mod
577 
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:120
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:35
subroutine, public get_var1_real(ncid, var1_name, im, var1, var_exist)
Definition: sim_nc_mod.F90:116
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
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 apply_inc_on_3d_scalar(field_name, var)
The subroutine &#39;apply_inc_on3d_scalar&#39; applies the input increments to the prognostic variables...
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
integer, parameter, public ng
Definition: fv_mp_mod.F90:2716
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 read_da_inc(Atm, fv_domain)
The subroutine &#39;read_da_inc&#39; reads the increments of the diagnostic variables from the DA-generated f...
subroutine, public remap_coef(is, ie, js, je, 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 cubed_to_latlon(u, v, ua, va, gridstruct, npx, npy, km, mode, grid_type, domain, nested, c2l_ord, bd)
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