module Convert 2D data
More...
|
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) |
|
module Convert 2D data
- Description
- subroutines for preparing 2D data (convert from external file)
- Author
- Team SCALE
◆ cnv2d_setup()
subroutine, public mod_cnv2d::cnv2d_setup |
◆ 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.
132 file_tiledata_get_data
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
143 integer,
intent(in),
optional :: interp_level
144 integer,
intent(in),
optional :: nmax
146 real(RP) :: DOMAIN_LATS, DOMAIN_LATE
147 real(RP) :: DOMAIN_LONS, DOMAIN_LONE
149 character(len=H_LONG) :: fname
155 if (
present(nmax) )
then
161 domain_lats = minval( latxv(:,:) )
162 domain_late = maxval( latxv(:,:) )
163 domain_lons = minval( lonuy(:,:) )
164 domain_lone = maxval( lonuy(:,:) )
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
170 tile_dlat = dlat * d2r
171 tile_dlon = dlon * d2r
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)
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) )
186 fname = trim(dir)//
'/'//trim(catalogue)
189 tile_dlat, tile_dlon, &
190 domain_lats, domain_late, domain_lons, domain_lone, &
194 tile_fname(:), tile_hit(:), &
195 tile_js(:), tile_je(:), tile_is(:), tile_ie(:), &
196 nlat, nlon, dom_js, dom_je, dom_is, dom_ie, &
199 allocate( lat_1d(nlat) )
200 allocate( lon_1d(nlon) )
204 tile_dlat, tile_dlon, &
205 lat_1d(:), lon_1d(:) )
208 cnv2d_tile_dtype = dtype
211 call cnv2d_init( interp_type, &
212 interp_level = interp_level, &
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().
◆ 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.
230 file_grads_get_shape, &
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
246 integer :: file_id, var_id
254 call file_grads_get_shape( file_id, var_name, &
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) )
270 call file_grads_read( file_id, var_id, &
275 lat_org(i,j) = lat_1d(j) * d2r
279 call file_grads_read( file_id, var_id, &
285 lat_org(i,j) = lat_org(i,j) * d2r
294 call file_grads_read( file_id, var_id, &
299 lon_org(i,j) = lon_1d(i) * d2r
303 call file_grads_read( file_id, var_id, &
309 lon_org(i,j) = lon_org(i,j) * d2r
317 cnv2d_ftype = i_grads
318 cnv2d_grads_fid = file_id
319 cnv2d_grads_vid = var_id
321 call cnv2d_init( interp_type, &
322 interp_level = interp_level, &
323 search_limit = search_limit, &
324 ll_struct = .false. )
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().
◆ 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.
335 file_tiledata_get_data
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
346 select case ( cnv2d_ftype )
348 call file_tiledata_get_data( nlat, nlon, &
352 tile_dlat, tile_dlon, &
353 tile_fname(:), tile_hit(:), &
354 tile_js(:), tile_je(:), tile_is(:), tile_ie(:), &
355 dom_js, dom_je, dom_is, dom_ie, &
359 min_value = min_value, yrevers = yrevers )
361 call file_grads_read( cnv2d_grads_fid, cnv2d_grads_vid, &
370 idx_i(:,:,:), idx_j(:,:,:), &
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().
subroutine, public interp_interp2d(npoints, IA_ref, JA_ref, IA, JA, idx_i, idx_j, hfact, val_ref, val, threshold_undef, wsum, val2)
subroutine, public prc_abort
Abort Process.
logical function, public file_grads_isoned(file_id, var_id)
subroutine, public file_grads_open(file_name, file_id)
Open.
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
subroutine, public file_grads_varid(file_id, var_name, var_id)
module Atmosphere GRID CartesC Real(real space)
subroutine, public file_tiledata_get_latlon(nLAT, nLON, jsh, ish, TILE_DLAT, TILE_DLON, LAT, LON)
get tile data
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
real(rp), public const_d2r
degree to radian
real(rp), public const_undef
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]