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 94 of file mod_cnv2d.F90.

94  use scale_prc, only: &
95  prc_abort
96  implicit none
97 
98  !---------------------------------------------------------------------------
99 
100  log_newline
101  log_info("CNV2D_setup",*) 'Setup'
102 
103  !--- read namelist
104 !!$ rewind(IO_FID_CONF)
105 !!$ read(IO_FID_CONF,nml=PARAM_CNV2D,iostat=ierr)
106 !!$ if( ierr < 0 ) then !--- missing
107 !!$ LOG_INFO("CNV2D_setup",*) 'Not found namelist. Default used.'
108 !!$ elseif( ierr > 0 ) then !--- fatal error
109 !!$ LOG_ERROR("CNV2D_setup",*) 'Not appropriate names in namelist PARAM_CNV2D. Check!'
110 !!$ call PRC_abort
111 !!$ endif
112 !!$ LOG_NML(PARAM_CNV2D)
113 
114 
115  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 126 of file mod_cnv2d.F90.

126  use scale_prc, only: &
127  prc_abort
128  use scale_const, only: &
129  undef => const_undef, &
130  d2r => const_d2r
131  use scale_file_tiledata, only: &
134  file_tiledata_get_data
135  use scale_atmos_grid_cartesc_real, only: &
138  implicit none
139  character(len=*), intent(in) :: dtype
140  real(RP), intent(in) :: dlat, dlon
141  character(len=*), intent(in) :: dir
142  character(len=*), intent(in) :: catalogue
143  character(len=*), intent(in) :: interp_type
144 
145  integer, intent(in), optional :: interp_level
146  integer, intent(in), optional :: nmax
147 
148  real(RP) :: DOMAIN_LATS, DOMAIN_LATE
149  real(RP) :: DOMAIN_LONS, DOMAIN_LONE
150 
151  character(len=H_LONG) :: fname
152 
153  integer :: ierr
154  integer :: i, j
155  !---------------------------------------------------------------------------
156 
157  if ( present(nmax) ) then
158  tile_nlim = nmax
159  else
160  tile_nlim = 1000
161  end if
162 
163  domain_lats = minval( latxv(:,:) )
164  domain_late = maxval( latxv(:,:) )
165  domain_lons = minval( lonuy(:,:) )
166  domain_lone = maxval( lonuy(:,:) )
167 
168  log_info("CNV2D_setup",*) 'Domain Information'
169  log_info_cont(*) 'Domain (LAT) :', domain_lats/d2r, domain_late/d2r
170  log_info_cont(*) ' (LON) :', domain_lons/d2r, domain_lone/d2r
171 
172  tile_dlat = dlat * d2r
173  tile_dlon = dlon * d2r
174 
175  if ( .not. first ) then
176  deallocate( tile_fname, tile_hit )
177  deallocate( tile_js, tile_je, tile_is, tile_ie )
178  deallocate( lat_1d, lon_1d )
179  end if
180 
181  allocate( tile_fname(tile_nlim), tile_hit(tile_nlim) )
182  allocate( tile_js(tile_nlim), tile_je(tile_nlim), tile_is(tile_nlim), tile_ie(tile_nlim) )
183 
184  ! catalogue file
185  fname = trim(dir)//'/'//trim(catalogue)
186 
187  call file_tiledata_get_info( tile_nlim, & ! [IN]
188  tile_dlat, tile_dlon, & ! [IN]
189  domain_lats, domain_late, domain_lons, domain_lone, & ! [IN]
190  fname, & ! [IN]
191  global_ia, & ! [OUT]
192  tile_nmax, & ! [OUT]
193  tile_fname(:), tile_hit(:), & ! [OUT]
194  tile_js(:), tile_je(:), tile_is(:), tile_ie(:), & ! [OUT]
195  nlat, nlon, dom_js, dom_je, dom_is, dom_ie, & ! [OUT]
196  zonal, pole ) ! [OUT]
197 
198  allocate( lat_1d(nlat) )
199  allocate( lon_1d(nlon) )
200 
201  call file_tiledata_get_latlon( nlat, nlon, & ! [IN]
202  dom_js, dom_is, & ! [IN]
203  tile_dlat, tile_dlon, & ! [IN]
204  lat_1d(:), lon_1d(:) ) ! [OUT]
205 
206  cnv2d_ftype = i_tile
207  cnv2d_tile_dtype = dtype
208  cnv2d_tile_dir = dir
209 
210  call cnv2d_init( interp_type, &
211  interp_level = interp_level, &
212  ll_struct = .true. )
213 
214  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_setup().

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), optional  interp_level,
real(rp), intent(in), optional  search_limit,
character(len=*), intent(in), optional  POSTFIX 
)

Definition at line 227 of file mod_cnv2d.F90.

227  use scale_file_grads, only: &
228  file_grads_open, &
229  file_grads_get_shape, &
232  file_grads_read
233  use scale_const, only: &
234  d2r => const_d2r
235  implicit none
236  character(len=*), intent(in) :: FILE_NAME
237  character(len=*), intent(in) :: VAR_NAME
238  character(len=*), intent(in) :: LAT_NAME
239  character(len=*), intent(in) :: LON_NAME
240  character(len=*), intent(in) :: INTERP_TYPE
241  integer, intent(in), optional :: interp_level
242  real(RP), intent(in), optional :: search_limit
243  character(len=*), intent(in), optional :: POSTFIX
244 
245  integer :: file_id, var_id
246  integer :: shape(2)
247 
248  integer :: i, j
249 
250  call file_grads_open( file_name, & ! [IN]
251  file_id ) ! [OUT]
252 
253  call file_grads_get_shape( file_id, var_name, &
254  shape(:) )
255  nlon = shape(1)
256  nlat = shape(2)
257 
258  if ( .not. first ) deallocate( lat_org, lon_org, lat_1d, lon_1d )
259  allocate( lat_org(nlon,nlat), lon_org(nlon,nlat) )
260  allocate( lat_1d(nlat), lon_1d(nlon) )
261 
262  ! lat
263  call file_grads_varid( file_id, lat_name, & ! (in)
264  var_id ) ! (out)
265  if ( file_grads_isoned( file_id, var_id ) ) then
266  call file_grads_read( file_id, var_id, & ! (in)
267  lat_1d(:) ) ! (out)
268  !$omp parallel do collapse(2)
269  do j = 1, nlat
270  do i = 1, nlon
271  lat_org(i,j) = lat_1d(j) * d2r
272  end do
273  end do
274  else
275  call file_grads_read( file_id, var_id, & ! (in)
276  lat_org(:,:), & ! (out)
277  postfix = postfix ) ! (in)
278  !$omp parallel do collapse(2)
279  do j = 1, nlat
280  do i = 1, nlon
281  lat_org(i,j) = lat_org(i,j) * d2r
282  end do
283  end do
284  end if
285 
286  ! lon
287  call file_grads_varid( file_id, lon_name, & ! (in)
288  var_id ) ! (out)
289  if ( file_grads_isoned( file_id, var_id ) ) then
290  call file_grads_read( file_id, var_id, & ! (in)
291  lon_1d(:) ) ! (out)
292  !$omp parallel do collapse(2)
293  do j = 1, nlat
294  do i = 1, nlon
295  lon_org(i,j) = lon_1d(i) * d2r
296  end do
297  end do
298  else
299  call file_grads_read( file_id, var_id, & ! (in)
300  lon_org(:,:), & ! (out)
301  postfix = postfix ) ! (in)
302  !$omp parallel do collapse(2)
303  do j = 1, nlat
304  do i = 1, nlon
305  lon_org(i,j) = lon_org(i,j) * d2r
306  end do
307  end do
308  end if
309 
310  call file_grads_varid( file_id, var_name, & ! [IN]
311  var_id ) ! [OUT]
312 
313  cnv2d_ftype = i_grads
314  cnv2d_grads_fid = file_id
315  cnv2d_grads_vid = var_id
316 
317  call cnv2d_init( interp_type, &
318  interp_level = interp_level, &
319  search_limit = search_limit, &
320  ll_struct = .false. )
321 
322  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_setup().

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 330 of file mod_cnv2d.F90.

330  use scale_file_tiledata, only: &
331  file_tiledata_get_data
332  use scale_file_grads, only: &
333  file_grads_read
334  use scale_interp, only: &
336  implicit none
337  real(RP), intent(out) :: var(IA,JA)
338  integer, intent(in), optional :: step
339  real(RP), intent(in), optional :: min_value
340  logical, intent(in), optional :: yrevers
341 
342  select case ( cnv2d_ftype )
343  case ( i_tile )
344  call file_tiledata_get_data( nlat, nlon, & ! [IN]
345  cnv2d_tile_dir, & ! [IN]
346  global_ia, & ! [IN]
347  tile_nmax, & ! [IN]
348  tile_dlat, tile_dlon, & ! [IN]
349  tile_fname(:), tile_hit(:), & ! [IN]
350  tile_js(:), tile_je(:), tile_is(:), tile_ie(:), & ! [IN]
351  dom_js, dom_je, dom_is, dom_ie, & ! [IN]
352  cnv2d_tile_dtype, & ! [IN]
353  data_org(:,:), & ! [OUT]
354  step = step, & ! [IN]
355  min_value = min_value, yrevers = yrevers ) ! [IN]
356  case ( i_grads )
357  call file_grads_read( cnv2d_grads_fid, cnv2d_grads_vid, & ! [IN]
358  data_org(:,:), & ! [OUT]
359  step = step ) ! [IN]
360  end select
361 
362 
363  call interp_interp2d( itp_lev, & ! [IN]
364  nlon,nlat, & ! [IN]
365  ia, ja, & ! [IN]
366  idx_i(:,:,:), idx_j(:,:,:), & ! [IN]
367  hfact(:,:,:), & ! [IN]
368  data_org(:,:), & ! [IN]
369  var(:,:) ) ! [OUT]
370 
371  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:1232
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_file_grads::file_grads_isoned
logical function, public file_grads_isoned(file_id, var_id)
Definition: scale_file_grads.F90:332
scale_file_grads::file_grads_open
subroutine, public file_grads_open(file_name, file_id)
Open.
Definition: scale_file_grads.F90:102
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:300
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:175
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:54
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
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:41
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:49