SCALE-RM
mod_cnv2d.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
11 module mod_cnv2d
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
20  use scale_tracer
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: cnv2d_setup
29  public :: cnv2d_tile_init
30  public :: cnv2d_grads_init
31  public :: cnv2d_exec
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private procedure
40  !
41  private :: cnv2d_init
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  integer, parameter :: I_TILE = 1
48  integer, parameter :: I_GrADS = 2
49  integer :: CNV2D_ftype
50 
51  integer :: GLOBAL_IA
52 
53  integer :: nLON, nLAT
54  real(RP), allocatable :: data_org(:,:)
55  real(RP), allocatable :: LAT_org (:,:)
56  real(RP), allocatable :: LON_org (:,:)
57  real(RP), allocatable :: LAT_1d (:)
58  real(RP), allocatable :: LON_1d (:)
59 
60  ! interpolation
61  integer :: itp_lev
62  integer, allocatable :: idx_i(:,:,:)
63  integer, allocatable :: idx_j(:,:,:)
64  real(RP), allocatable :: hfact(:,:,:)
65  logical :: zonal, pole
66 
67  ! TILE data
68  character(len=H_SHORT) :: CNV2D_tile_dtype
69  character(len=H_LONG) :: CNV2D_tile_dir
70  real(RP) :: TILE_DLAT, TILE_DLON
71  integer :: TILE_nlim
72  integer :: TILE_nmax
73  character(len=H_LONG), allocatable :: TILE_fname(:)
74  logical, allocatable :: TILE_hit (:)
75  integer, allocatable :: TILE_JS (:)
76  integer, allocatable :: TILE_JE (:)
77  integer, allocatable :: TILE_IS (:)
78  integer, allocatable :: TILE_IE (:)
79 
80  integer :: dom_is, dom_ie, dom_js, dom_je
81 
82 
83  ! GrADS data
84  integer :: CNV2D_grads_fid
85  integer :: CNV2D_grads_vid
86 
87  !-----------------------------------------------------------------------------
88 contains
89  !-----------------------------------------------------------------------------
91  subroutine cnv2d_setup
92  use scale_prc, only: &
93  prc_abort
94  implicit none
95 
96  !---------------------------------------------------------------------------
97 
98  log_newline
99  log_info("CNV2D_setup",*) 'Setup'
100 
101  !--- read namelist
102 !!$ rewind(IO_FID_CONF)
103 !!$ read(IO_FID_CONF,nml=PARAM_CNV2D,iostat=ierr)
104 !!$ if( ierr < 0 ) then !--- missing
105 !!$ LOG_INFO("CNV2D_setup",*) 'Not found namelist. Default used.'
106 !!$ elseif( ierr > 0 ) then !--- fatal error
107 !!$ LOG_ERROR("CNV2D_setup",*) 'Not appropriate names in namelist PARAM_CNV2D. Check!'
108 !!$ call PRC_abort
109 !!$ endif
110 !!$ LOG_NML(PARAM_CNV2D)
111 
112 
113  return
114  end subroutine cnv2d_setup
115 
116  subroutine cnv2d_tile_init( &
117  dtype, &
118  dlat, dlon, &
119  dir, &
120  catalogue, &
121  interp_type, &
122  interp_level, &
123  nmax )
124  use scale_prc, only: &
125  prc_abort
126  use scale_const, only: &
127  undef => const_undef, &
128  d2r => const_d2r
129  use scale_file_tiledata, only: &
132  file_tiledata_get_data
133  use scale_atmos_grid_cartesc_real, only: &
136  implicit none
137  character(len=*), intent(in) :: dtype
138  real(rp), intent(in) :: dlat, dlon
139  character(len=*), intent(in) :: dir
140  character(len=*), intent(in) :: catalogue
141  character(len=*), intent(in) :: interp_type
142 
143  integer, intent(in), optional :: interp_level
144  integer, intent(in), optional :: nmax
145 
146  real(rp) :: domain_lats, domain_late
147  real(rp) :: domain_lons, domain_lone
148 
149  character(len=H_LONG) :: fname
150 
151  integer :: ierr
152  integer :: i, j
153  !---------------------------------------------------------------------------
154 
155  if ( present(nmax) ) then
156  tile_nlim = nmax
157  else
158  tile_nlim = 1000
159  end if
160 
161  domain_lats = minval( latxv(:,:) )
162  domain_late = maxval( latxv(:,:) )
163  domain_lons = minval( lonuy(:,:) )
164  domain_lone = maxval( lonuy(:,:) )
165 
166  log_info("CNV2D_setup",*) 'Domain Information'
167  log_info_cont(*) 'Domain (LAT) :', domain_lats/d2r, domain_late/d2r
168  log_info_cont(*) ' (LON) :', domain_lons/d2r, domain_lone/d2r
169 
170  tile_dlat = dlat * d2r
171  tile_dlon = dlon * d2r
172 
173  if (allocated(tile_fname)) deallocate(tile_fname)
174  if (allocated(tile_hit)) deallocate(tile_hit)
175  if (allocated(tile_js)) deallocate(tile_js)
176  if (allocated(tile_je)) deallocate(tile_je)
177  if (allocated(tile_is)) deallocate(tile_is)
178  if (allocated(tile_ie)) deallocate(tile_ie)
179  if (allocated(lat_1d)) deallocate(lat_1d)
180  if (allocated(lon_1d)) deallocate(lon_1d)
181 
182  allocate( tile_fname(tile_nlim), tile_hit(tile_nlim) )
183  allocate( tile_js(tile_nlim), tile_je(tile_nlim), tile_is(tile_nlim), tile_ie(tile_nlim) )
184 
185  ! catalogue file
186  fname = trim(dir)//'/'//trim(catalogue)
187 
188  call file_tiledata_get_info( tile_nlim, & ! [IN]
189  tile_dlat, tile_dlon, & ! [IN]
190  domain_lats, domain_late, domain_lons, domain_lone, & ! [IN]
191  fname, & ! [IN]
192  global_ia, & ! [OUT]
193  tile_nmax, & ! [OUT]
194  tile_fname(:), tile_hit(:), & ! [OUT]
195  tile_js(:), tile_je(:), tile_is(:), tile_ie(:), & ! [OUT]
196  nlat, nlon, dom_js, dom_je, dom_is, dom_ie, & ! [OUT]
197  zonal, pole ) ! [OUT]
198 
199  allocate( lat_1d(nlat) )
200  allocate( lon_1d(nlon) )
201 
202  call file_tiledata_get_latlon( nlat, nlon, & ! [IN]
203  dom_js, dom_is, & ! [IN]
204  tile_dlat, tile_dlon, & ! [IN]
205  lat_1d(:), lon_1d(:) ) ! [OUT]
206 
207  cnv2d_ftype = i_tile
208  cnv2d_tile_dtype = dtype
209  cnv2d_tile_dir = dir
210 
211  call cnv2d_init( interp_type, &
212  interp_level = interp_level, &
213  ll_struct = .true. )
214 
215  return
216  end subroutine cnv2d_tile_init
217 
218 
219  subroutine cnv2d_grads_init( &
220  FILE_NAME, &
221  VAR_NAME, &
222  LAT_NAME, &
223  LON_NAME, &
224  interp_type, &
225  interp_level, &
226  search_limit, &
227  POSTFIX )
228  use scale_file_grads, only: &
229  file_grads_open, &
230  file_grads_get_shape, &
233  file_grads_read
234  use scale_const, only: &
235  d2r => const_d2r
236  implicit none
237  character(len=*), intent(in) :: file_name
238  character(len=*), intent(in) :: var_name
239  character(len=*), intent(in) :: lat_name
240  character(len=*), intent(in) :: lon_name
241  character(len=*), intent(in) :: interp_type
242  integer, intent(in) :: interp_level
243  real(rp), intent(in), optional :: search_limit
244  character(len=*), intent(in), optional :: postfix
245 
246  integer :: file_id, var_id
247  integer :: shape(2)
248 
249  integer :: i, j
250 
251  call file_grads_open( file_name, & ! [IN]
252  file_id ) ! [OUT]
253 
254  call file_grads_get_shape( file_id, var_name, &
255  shape(:) )
256  nlon = shape(1)
257  nlat = shape(2)
258 
259  if (allocated(lat_org)) deallocate(lat_org)
260  if (allocated(lon_org)) deallocate(lon_org)
261  if (allocated(lat_1d)) deallocate(lat_1d)
262  if (allocated(lon_1d)) deallocate(lon_1d)
263  allocate( lat_org(nlon,nlat), lon_org(nlon,nlat) )
264  allocate( lat_1d(nlat), lon_1d(nlon) )
265 
266  ! lat
267  call file_grads_varid( file_id, lat_name, & ! (in)
268  var_id ) ! (out)
269  if ( file_grads_isoned( file_id, var_id ) ) then
270  call file_grads_read( file_id, var_id, & ! (in)
271  lat_1d(:) ) ! (out)
272  !$omp parallel do collapse(2)
273  do j = 1, nlat
274  do i = 1, nlon
275  lat_org(i,j) = lat_1d(j) * d2r
276  end do
277  end do
278  else
279  call file_grads_read( file_id, var_id, & ! (in)
280  lat_org(:,:), & ! (out)
281  postfix = postfix ) ! (in)
282  !$omp parallel do collapse(2)
283  do j = 1, nlat
284  do i = 1, nlon
285  lat_org(i,j) = lat_org(i,j) * d2r
286  end do
287  end do
288  end if
289 
290  ! lon
291  call file_grads_varid( file_id, lon_name, & ! (in)
292  var_id ) ! (out)
293  if ( file_grads_isoned( file_id, var_id ) ) then
294  call file_grads_read( file_id, var_id, & ! (in)
295  lon_1d(:) ) ! (out)
296  !$omp parallel do collapse(2)
297  do j = 1, nlat
298  do i = 1, nlon
299  lon_org(i,j) = lon_1d(i) * d2r
300  end do
301  end do
302  else
303  call file_grads_read( file_id, var_id, & ! (in)
304  lon_org(:,:), & ! (out)
305  postfix = postfix ) ! (in)
306  !$omp parallel do collapse(2)
307  do j = 1, nlat
308  do i = 1, nlon
309  lon_org(i,j) = lon_org(i,j) * d2r
310  end do
311  end do
312  end if
313 
314  call file_grads_varid( file_id, var_name, & ! [IN]
315  var_id ) ! [OUT]
316 
317  cnv2d_ftype = i_grads
318  cnv2d_grads_fid = file_id
319  cnv2d_grads_vid = var_id
320 
321  call cnv2d_init( interp_type, &
322  interp_level = interp_level, &
323  search_limit = search_limit, &
324  ll_struct = .false. )
325 
326  return
327  end subroutine cnv2d_grads_init
328 
329  subroutine cnv2d_exec( &
330  var, &
331  step, &
332  min_value, &
333  yrevers )
334  use scale_file_tiledata, only: &
335  file_tiledata_get_data
336  use scale_file_grads, only: &
337  file_grads_read
338  use scale_interp, only: &
340  implicit none
341  real(rp), intent(out) :: var(ia,ja)
342  integer, intent(in), optional :: step
343  real(rp), intent(in), optional :: min_value
344  logical, intent(in), optional :: yrevers
345 
346  select case ( cnv2d_ftype )
347  case ( i_tile )
348  call file_tiledata_get_data( nlat, nlon, & ! [IN]
349  cnv2d_tile_dir, & ! [IN]
350  global_ia, & ! [IN]
351  tile_nmax, & ! [IN]
352  tile_dlat, tile_dlon, & ! [IN]
353  tile_fname(:), tile_hit(:), & ! [IN]
354  tile_js(:), tile_je(:), tile_is(:), tile_ie(:), & ! [IN]
355  dom_js, dom_je, dom_is, dom_ie, & ! [IN]
356  cnv2d_tile_dtype, & ! [IN]
357  data_org(:,:), & ! [OUT]
358  step = step, & ! [IN]
359  min_value = min_value, yrevers = yrevers ) ! [IN]
360  case ( i_grads )
361  call file_grads_read( cnv2d_grads_fid, cnv2d_grads_vid, & ! [IN]
362  data_org(:,:), & ! [OUT]
363  step = step ) ! [IN]
364  end select
365 
366 
367  call interp_interp2d( itp_lev, & ! [IN]
368  nlon,nlat, & ! [IN]
369  ia, ja, & ! [IN]
370  idx_i(:,:,:), idx_j(:,:,:), & ! [IN]
371  hfact(:,:,:), & ! [IN]
372  data_org(:,:), & ! [IN]
373  var(:,:) ) ! [OUT]
374 
375  return
376  end subroutine cnv2d_exec
377 
378 
379  ! private
380 
381  subroutine cnv2d_init( &
382  interp_type, &
383  interp_level, &
384  search_limit, &
385  ll_struct )
386  use scale_prc, only: &
387  prc_abort
388  use scale_const, only: &
389  pi => const_pi
390  use scale_interp, only: &
394  use scale_mapprojection, only: &
395  mapprojection_lonlat2xy
396  use scale_atmos_grid_cartesc, only: &
397  cx => atmos_grid_cartesc_cx, &
399  use scale_atmos_grid_cartesc_real, only: &
402  implicit none
403 
404  character(len=*), intent(in) :: interp_type
405  integer, intent(in), optional :: interp_level
406  real(rp), intent(in), optional :: search_limit
407  logical, intent(in), optional :: ll_struct
408 
409  real(rp), allocatable :: x_org(:,:), y_org(:,:)
410  !---------------------------------------------------------------------------
411 
412  if (allocated(idx_i)) deallocate(idx_i)
413  if (allocated(idx_j)) deallocate(idx_j)
414  if (allocated(hfact)) deallocate(hfact)
415  if (allocated(data_org)) deallocate(data_org)
416 
417  select case ( interp_type )
418  case ( 'LINEAR' )
419  itp_lev = 4
420  case ( 'DIST-WEIGHT' )
421  if ( present(interp_level) ) then
422  itp_lev = interp_level
423  else
424  log_error('CNV2D_init',*) 'INTERP_LEVEL is required for the DIST-WEIGHT interpolation'
425  call prc_abort
426  end if
427  case default
428  log_error('CNV2D_init',*) 'INTERP_TYPE is invalid: ', trim(interp_type)
429  call prc_abort
430  end select
431 
432  ! interporation
433  allocate( idx_i(ia,ja,itp_lev) )
434  allocate( idx_j(ia,ja,itp_lev) )
435  allocate( hfact(ia,ja,itp_lev) )
436 
437  select case ( interp_type )
438  case ( 'LINEAR' )
439  if ( ll_struct ) then
440  call interp_factor2d_linear_latlon( nlon, nlat, & ! [IN]
441  ia, ja, & ! [IN]
442  lon_1d(:), lat_1d(:), & ! [IN]
443  lon(:,:), lat(:,:), & ! [IN]
444  idx_i(:,:,:), idx_j(:,:,:), & ! [OUT]
445  hfact(:,:,:) ) ! [OUT]
446  else
447  allocate( x_org(nlon,nlat), y_org(nlon,nlat) )
448  call mapprojection_lonlat2xy( nlon, 1, nlon, nlat, 1, nlat, & ! [IN]
449  lon_org(:,:), lat_org(:,:), & ! [IN]
450  x_org(:,:), y_org(:,:) ) ! [OUT]
451 
452  zonal = ( maxval(lon_org) - minval(lat_org) ) > 2.0_rp * pi * 0.9_rp
453  pole = ( maxval(lat_org) > pi * 0.5_rp * 0.9_rp ) .or. ( minval(lat_org) < - pi * 0.5_rp * 0.9_rp )
454  call interp_factor2d_linear_xy( nlon, nlat, & ! [IN]
455  ia, ja, & ! [IN]
456  x_org(:,:), y_org(:,:), & ! [IN]
457  cx(:), cy(:), & ! [IN]
458  idx_i(:,:,:), idx_j(:,:,:), & ! [OUT]
459  hfact(:,:,:), & ! [OUT]
460  zonal = zonal, pole = pole, & ! [IN]
461  missing = .true. ) ! [IN, optional]
462  deallocate( x_org, y_org )
463  end if
464  case ( 'DIST-WEIGHT' )
465  call interp_factor2d_weight( itp_lev, & ! [IN]
466  nlon, nlat, & ! [IN]
467  ia, ja, & ! [IN]
468  lon_org(:,:), lat_org(:,:), & ! [IN]
469  lon(:,:), lat(:,:), & ! [IN]
470  idx_i(:,:,:), idx_j(:,:,:), & ! [OUT]
471  hfact(:,:,:), & ! [OUT]
472  latlon_structure = ll_struct, & ! [IN]
473  lon_1d = lon_1d(:), & ! [IN]
474  lat_1d = lat_1d(:), & ! [IN]
475  search_limit = search_limit ) ! [IN]
476  end select
477 
478  allocate( data_org(nlon,nlat) )
479 
480  return
481  end subroutine cnv2d_init
482 
483  ! private
484 
485 end module mod_cnv2d
scale_interp::interp_interp2d
subroutine, public interp_interp2d(npoints, IA_ref, JA_ref, IA, JA, idx_i, idx_j, hfact, val_ref, val, threshold_undef, wsum, val2)
Definition: scale_interp.F90:1376
mod_cnv2d::cnv2d_grads_init
subroutine, public cnv2d_grads_init(FILE_NAME, VAR_NAME, LAT_NAME, LON_NAME, interp_type, interp_level, search_limit, POSTFIX)
Definition: mod_cnv2d.F90:228
mod_cnv2d
module Convert 2D data
Definition: mod_cnv2d.F90:11
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_file_grads::file_grads_isoned
logical function, public file_grads_isoned(file_id, var_id)
Definition: scale_file_grads.F90:373
scale_file_grads::file_grads_open
subroutine, public file_grads_open(file_name, file_id)
Open.
Definition: scale_file_grads.F90:104
scale_interp::interp_factor2d_weight
subroutine, public interp_factor2d_weight(npoints, IA_ref, JA_ref, IA, JA, lon_ref, lat_ref, lon, lat, idx_i, idx_j, hfact, search_limit, latlon_structure, lon_1d, lat_1d, weight_order)
Definition: scale_interp.F90:748
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_interp
module INTERPOLATION
Definition: scale_interp.F90:12
scale_file_tiledata::file_tiledata_get_info
subroutine, public file_tiledata_get_info(TILE_nlim, TILE_DLAT, TILE_DLON, DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE, catalog_fname, GLOBAL_IA, TILE_nmax, TILE_fname, TILE_hit, TILE_JS, TILE_JE, TILE_IS, TILE_IE, nLATH, nLONH, jsh, jeh, ish, ieh, zonal, pole, single_fname, LATS, LATE, LONS, LONE)
get tile information
Definition: scale_file_tiledata.F90:62
scale_file_grads::file_grads_varid
subroutine, public file_grads_varid(file_id, var_name, var_id)
Definition: scale_file_grads.F90:312
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:49
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_file_grads
module file_grads
Definition: scale_file_grads.F90:11
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_cnv2d::cnv2d_exec
subroutine, public cnv2d_exec(var, step, min_value, yrevers)
Definition: mod_cnv2d.F90:334
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_file_tiledata::file_tiledata_get_latlon
subroutine, public file_tiledata_get_latlon(nLAT, nLON, jsh, ish, TILE_DLAT, TILE_DLON, LAT, LON)
get tile data
Definition: scale_file_tiledata.F90:177
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:55
scale_interp::interp_factor2d_linear_xy
subroutine, public interp_factor2d_linear_xy(IA_ref, JA_ref, IA, JA, x_ref, y_ref, x, y, idx_i, idx_j, hfact, zonal, pole, missing)
Definition: scale_interp.F90:483
mod_cnv2d::cnv2d_setup
subroutine, public cnv2d_setup
Setup.
Definition: mod_cnv2d.F90:92
mod_cnv2d::cnv2d_tile_init
subroutine, public cnv2d_tile_init(dtype, dlat, dlon, dir, catalogue, interp_type, interp_level, nmax)
Definition: mod_cnv2d.F90:124
scale_atmos_grid_cartesc::atmos_grid_cartesc_cy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cy
center coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:57
scale_interp::interp_factor2d_linear_latlon
subroutine, public interp_factor2d_linear_latlon(IA_ref, JA_ref, IA, JA, lon_ref, lat_ref, lon, lat, idx_i, idx_j, hfact)
Definition: scale_interp.F90:357
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_file_tiledata
module file_tiledata
Definition: scale_file_tiledata.F90:12
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:53
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_grid_cartesc::atmos_grid_cartesc_cx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cx
center coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:56
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:50