SCALE-RM
Functions/Subroutines
mod_cnv2d Module Reference

module Convert 2D data More...

Functions/Subroutines

subroutine, public cnv2d_setup
 Setup. More...
 
subroutine, public cnv2d_tile_init (dtype, dlat, dlon, dir, catalogue, interp_type, interp_level, nmax)
 
subroutine, public cnv2d_grads_init (FILE_NAME, VAR_NAME, LAT_NAME, LON_NAME, interp_type, interp_level, search_limit, POSTFIX)
 
subroutine, public cnv2d_exec (var, step, min_value, yrevers)
 

Detailed Description

module Convert 2D data

Description
subroutines for preparing 2D data (convert from external file)
Author
Team SCALE

Function/Subroutine Documentation

◆ cnv2d_setup()

subroutine, public mod_cnv2d::cnv2d_setup

Setup.

Definition at line 92 of file mod_cnv2d.F90.

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

References scale_prc::prc_abort().

Here is the call graph for this function:

◆ cnv2d_tile_init()

subroutine, public mod_cnv2d::cnv2d_tile_init ( character(len=*), intent(in)  dtype,
real(rp), intent(in)  dlat,
real(rp), intent(in)  dlon,
character(len=*), intent(in)  dir,
character(len=*), intent(in)  catalogue,
character(len=*), intent(in)  interp_type,
integer, intent(in), optional  interp_level,
integer, intent(in), optional  nmax 
)

Definition at line 124 of file mod_cnv2d.F90.

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

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy, scale_const::const_d2r, scale_const::const_undef, scale_file_tiledata::file_tiledata_get_info(), scale_file_tiledata::file_tiledata_get_latlon(), and scale_prc::prc_abort().

Referenced by mod_cnvtopo::cnvtopo(), and mod_cnvuser::cnvuser().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cnv2d_grads_init()

subroutine, public mod_cnv2d::cnv2d_grads_init ( character(len=*), intent(in)  FILE_NAME,
character(len=*), intent(in)  VAR_NAME,
character(len=*), intent(in)  LAT_NAME,
character(len=*), intent(in)  LON_NAME,
character(len=*), intent(in)  interp_type,
integer, intent(in)  interp_level,
real(rp), intent(in), optional  search_limit,
character(len=*), intent(in), optional  POSTFIX 
)

Definition at line 228 of file mod_cnv2d.F90.

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

References scale_const::const_d2r, scale_file_grads::file_grads_isoned(), scale_file_grads::file_grads_open(), and scale_file_grads::file_grads_varid().

Referenced by mod_cnvtopo::cnvtopo(), and mod_cnvuser::cnvuser().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cnv2d_exec()

subroutine, public mod_cnv2d::cnv2d_exec ( real(rp), dimension(ia,ja), intent(out)  var,
integer, intent(in), optional  step,
real(rp), intent(in), optional  min_value,
logical, intent(in), optional  yrevers 
)

Definition at line 334 of file mod_cnv2d.F90.

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

References scale_atmos_grid_cartesc::atmos_grid_cartesc_cx, scale_atmos_grid_cartesc::atmos_grid_cartesc_cy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon, scale_const::const_pi, scale_atmos_grid_cartesc_index::ia, scale_interp::interp_factor2d_linear_latlon(), scale_interp::interp_factor2d_linear_xy(), scale_interp::interp_factor2d_weight(), scale_interp::interp_interp2d(), scale_atmos_grid_cartesc_index::ja, and scale_prc::prc_abort().

Referenced by mod_cnvtopo::cnvtopo(), and mod_cnvuser::cnvuser().

Here is the call graph for this function:
Here is the caller graph for this function:
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
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
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_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_file_grads
module file_grads
Definition: scale_file_grads.F90:11
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_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_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
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