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
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35  logical, public :: cnv2d_donothing
36  logical, public :: cnv2d_usegrads = .true.
37 
38  !-----------------------------------------------------------------------------
39  !
40  !++ Private procedure
41  !
42  private :: cnv2d_grads_bilinear
43  private :: cnv2d_grads_areaweighted
44 
45  private :: cnv2d_write
46 
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private parameters & variables
50  !
51  real(RP), private :: cnv2d_unittile_ddeg = 0.0_rp ! dx for unit tile [deg]
52  real(RP), private :: cnv2d_oversampling_factor = 2.0_rp ! factor of min. dx against the unit tile
53 
54  character(len=H_SHORT), private :: cnv2d_interpolation_type = 'bilinear' ! interpolation type
55  ! 'bilinear' bi-linear interpolation (coarse input->fine model grid)
56  ! 'areaweighted' area-weighted average (fine input->coarse model grid)
57  ! 'nearestneighbor' nearest neighbor (both)
58 
59  character(len=H_LONG), private :: cnv2d_out_basename = '' ! basename of the output file
60  character(len=H_MID), private :: cnv2d_out_title = 'SCALE-RM 2D Boundary' ! title of the output file
61  character(len=H_SHORT), private :: cnv2d_out_varname = '' ! name of the variable
62  character(len=H_MID), private :: cnv2d_out_vardesc = '' ! title of the variable
63  character(len=H_SHORT), private :: cnv2d_out_varunit = '' ! units of the variable
64  character(len=H_SHORT), private :: cnv2d_out_dtype = 'DEFAULT' ! REAL4 or REAL8
65 
66  !-----------------------------------------------------------------------------
67 contains
68  !-----------------------------------------------------------------------------
70  subroutine cnv2d_setup
71  use scale_prc, only: &
72  prc_abort
73  use scale_const, only: &
74  d2r => const_d2r
75  use scale_statistics, only: &
76  statistics_horizontal_min
78  real_dlat => atmos_grid_cartesc_real_dlat, &
79  real_dlon => atmos_grid_cartesc_real_dlon
80  implicit none
81 
82  namelist / param_cnv2d / &
84  cnv2d_interpolation_type, &
85  cnv2d_unittile_ddeg, &
86  cnv2d_oversampling_factor, &
87  cnv2d_out_basename, &
88  cnv2d_out_title, &
89  cnv2d_out_varname, &
90  cnv2d_out_vardesc, &
91  cnv2d_out_varunit, &
92  cnv2d_out_dtype
93 
94  real(RP) :: drad(ia,ja)
95  real(RP) :: drad_min
96 
97  integer :: ierr
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  cnv2d_donothing = .true.
115 
116  if ( cnv2d_usegrads ) then
117  cnv2d_donothing = .false.
118  log_info("CNV2D_setup",*) 'Use GrADS format file for input'
119  endif
120 
121  if ( cnv2d_interpolation_type == 'bilinear' ) then
122  log_info("CNV2D_setup",*) 'Interpolation type : bi-linear interpolation (coarse input->fine model grid)'
123  elseif( cnv2d_interpolation_type == 'areaweighted' ) then
124  log_info("CNV2D_setup",*) 'Interpolation type : area-weighted average (fine input->coarse model grid)'
125  elseif( cnv2d_interpolation_type == 'nearestneighbor' ) then
126  log_info("CNV2D_setup",*) 'Interpolation type : nearest neighbor'
127  else
128  log_error("CNV2D_setup",*) 'Not appropriate interpolation type. Check!', trim(cnv2d_interpolation_type)
129  call prc_abort
130  endif
131 
132  if ( cnv2d_donothing ) then
133  log_info("CNV2D_setup",*) 'Do nothing for 2D data conversion'
134  else
135  drad(:,:) = min( real_dlat(:,:), real_dlon(:,:) )
136  call statistics_horizontal_min( ia, is, ie, ja, js, je, &
137  drad(:,:), drad_min )
138 
139  if ( cnv2d_unittile_ddeg > 0.0_rp ) then
140  cnv2d_oversampling_factor = ( drad_min / d2r ) / cnv2d_unittile_ddeg
141  endif
142  cnv2d_oversampling_factor = max( 1.0_rp, cnv2d_oversampling_factor )
143  cnv2d_unittile_ddeg = ( drad_min / d2r ) / cnv2d_oversampling_factor
144 
145  log_info("CNV2D_setup",*) 'The size of tile [deg] = ', cnv2d_unittile_ddeg
146  log_info("CNV2D_setup",*) 'oversampling factor = ', cnv2d_oversampling_factor
147  endif
148 
149  return
150  end subroutine cnv2d_setup
151 
152  !-----------------------------------------------------------------------------
154  subroutine cnv2d
155  implicit none
156  !---------------------------------------------------------------------------
157 
158  if ( cnv2d_donothing ) then
159  log_newline
160  log_progress(*) 'skip convert topography data'
161  else
162  log_newline
163  log_progress(*) 'start convert topography data'
164 
165  if ( cnv2d_usegrads ) then
166  call cnv2d_grads
167  endif
168 
169  log_progress(*) 'end convert topography data'
170  endif
171 
172  return
173  end subroutine cnv2d
174 
175  !-----------------------------------------------------------------------------
177  subroutine cnv2d_grads
178  use scale_prc, only: &
179  prc_abort
180  use scale_const, only: &
181  undef => const_undef, &
182  pi => const_pi, &
183  d2r => const_d2r
184  use scale_calendar, only: &
186  use scale_atmos_grid_cartesc_real, only: &
187  real_latxv => atmos_grid_cartesc_real_latxv, &
188  real_lonuy => atmos_grid_cartesc_real_lonuy
189  implicit none
190 
191  integer :: GrADS_NLAT = -1 ! number of latitude tile
192  integer :: GrADS_NLON = -1 ! number of longitude tile
193  real(RP) :: GrADS_DLAT = -1.0_rp ! width of latitude tile [deg.]
194  real(RP) :: GrADS_DLON = -1.0_rp ! width of longitude tile [deg.]
195  character(len=H_LONG) :: GrADS_IN_DIR = '.' ! directory contains data files (GrADS format)
196  character(len=H_LONG) :: GrADS_IN_CATALOGUE = '' ! catalogue file
197  character(len=H_LONG) :: GrADS_IN_FILENAME = '' ! single data file (GrADS format)
198  character(len=H_LONG) :: GrADS_IN_DATATYPE = 'REAL4' ! datatype (REAL4,REAL8,INT2)
199  logical :: GrADS_LATORDER_N2S = .false. ! data of the latitude direction is stored in ordar of North->South?
200  real(RP) :: GrADS_MISSINGVAL ! missing value
201  real(RP) :: GrADS_LAT_START = -90.0_rp ! (for single file) start latitude of domain in input data
202  real(RP) :: GrADS_LAT_END = 90.0_rp ! (for single file) end latitude of domain in input data
203  real(RP) :: GrADS_LON_START = 0.0_rp ! (for single file) start longitude of domain in input data
204  real(RP) :: GrADS_LON_END = 360.0_rp ! (for single file) end longitude of domain in input data
205  integer :: GrADS_NSTEP = 1 ! number of steps
206  real(DP) :: GrADS_DT = 0.0_dp ! time interval
207  character(len=H_SHORT) :: GrADS_DT_UNIT = "SEC" ! time unit for GrADS_DT
208 
209  namelist / param_cnv2d_grads / &
210  grads_nlat, &
211  grads_nlon, &
212  grads_dlat, &
213  grads_dlon, &
214  grads_in_catalogue, &
215  grads_in_dir, &
216  grads_in_filename, &
217  grads_in_datatype, &
218  grads_latorder_n2s, &
219  grads_missingval, &
220  grads_lat_start, &
221  grads_lat_end, &
222  grads_lon_start, &
223  grads_lon_end, &
224  grads_nstep, &
225  grads_dt, &
226  grads_dt_unit
227 
228  real(RP) :: VAR2D(ia,ja)
229  real(DP) :: VAR2D_DTSEC
230 
231  real(RP) :: REAL_LONUY_mod(0:ia,ja)
232  real(RP) :: DOMAIN_LATS, DOMAIN_LATE
233  real(RP) :: DOMAIN_LONS, DOMAIN_LONE
234  integer :: DOMAIN_LONSLOC(2), DOMAIN_LONELOC(2)
235  logical :: check_IDL
236 
237  integer :: TILE_NLAT, TILE_NLON
238  real(RP) :: TILE_DLAT, TILE_DLON
239 
240  ! data catalogue list
241  integer, parameter :: TILE_nlim = 1000
242  integer :: TILE_nmax
243  character(len=H_LONG) :: TILE_fname(tile_nlim)
244  real(RP) :: TILE_LATS (tile_nlim)
245  real(RP) :: TILE_LATE (tile_nlim)
246  real(RP) :: TILE_LONS (tile_nlim)
247  real(RP) :: TILE_LONE (tile_nlim)
248 
249  character(len=H_LONG) :: fname
250  integer :: t, index
251  integer :: nowstep
252  integer :: fid, ierr
253  !---------------------------------------------------------------------------
254 
255  log_newline
256  log_info("CNV2D_GrADS",*) 'Setup'
257 
258  grads_missingval = undef
259 
260  !--- read namelist
261  rewind(io_fid_conf)
262  read(io_fid_conf,nml=param_cnv2d_grads,iostat=ierr)
263  if( ierr < 0 ) then !--- missing
264  log_info("CNV2D_GrADS",*) 'Not found namelist. Default used.'
265  elseif( ierr > 0 ) then !--- fatal error
266  log_error("CNV2D_GrADS",*) 'Not appropriate names in namelist PARAM_CNV2D_GrADS. Check!'
267  call prc_abort
268  endif
269  log_nml(param_cnv2d_grads)
270 
271  if ( grads_nlat <= 0 ) then
272  log_error("CNV2D_GrADS",*) 'GrADS_NLAT (number of latitude tile) should be positive. Check!', grads_nlat
273  call prc_abort
274  endif
275 
276  if ( grads_nlon <= 0 ) then
277  log_error("CNV2D_GrADS",*) 'GrADS_NLON (number of longitude tile) should be positive. Check!', grads_nlon
278  call prc_abort
279  endif
280 
281  if ( grads_dlat <= 0.0_rp ) then
282  log_error("CNV2D_GrADS",*) 'GrADS_DLAT (width (deg.) of latitude tile) should be positive. Check!', grads_dlat
283  call prc_abort
284  endif
285 
286  if ( grads_dlon <= 0.0_rp ) then
287  log_error("CNV2D_GrADS",*) 'GrADS_DLON (width (deg.) of longitude tile) should be positive. Check!', grads_dlon
288  call prc_abort
289  endif
290 
291  if ( grads_in_catalogue == '' &
292  .AND. grads_in_filename == '' ) then
293  log_error("CNV2D_GrADS",*) 'Neither catalogue file nor single file do not specified. Check!'
294  call prc_abort
295  endif
296 
297  if ( grads_in_datatype == 'REAL8' ) then
298  log_info("CNV2D_GrADS",*) 'type of input data : REAL8'
299  elseif( grads_in_datatype == 'REAL4' ) then
300  log_info("CNV2D_GrADS",*) 'type of input data : REAL4'
301  elseif( grads_in_datatype == 'INT2' ) then
302  log_info("CNV2D_GrADS",*) 'type of input data : INT2'
303  else
304  log_error("CNV2D_GrADS",*) 'Not appropriate type for GrADS_IN_DATATYPE. Check!'
305  log_error_cont(*) 'REAL8, REAL4, INT2 are available. requested:', trim(grads_in_datatype)
306  call prc_abort
307  endif
308 
309  if ( grads_latorder_n2s ) then
310  log_info("CNV2D_GrADS",*) 'data ordar of the latitude direction : North -> South'
311  else
312  log_info("CNV2D_GrADS",*) 'data ordar of the latitude direction : South -> North'
313  endif
314 
315  log_info("CNV2D_GrADS",*) 'Number of steps : ', grads_nstep
316  if ( grads_nstep > 1 ) then
317  call calendar_unit2sec( var2d_dtsec, grads_dt, grads_dt_unit )
318  log_info("CNV2D_GrADS",*) 'Time interval : ', var2d_dtsec
319  else
320  var2d_dtsec = 0.0_dp
321  endif
322 
323 
324 
325  real_lonuy_mod(:,:) = mod( real_lonuy(:,:)+3.0_dp*pi, 2.0_dp*pi ) - pi ! [-> 0..180,-180..0]
326 
327  domain_lats = minval(real_latxv(:,:))
328  domain_late = maxval(real_latxv(:,:))
329  domain_lons = minval(real_lonuy_mod(:,:))
330  domain_lone = maxval(real_lonuy_mod(:,:))
331  domain_lonsloc = minloc(real_lonuy_mod(:,:))
332  domain_loneloc = maxloc(real_lonuy_mod(:,:))
333 
334  check_idl = .false.
335  if ( domain_lons < real_lonuy_mod( 0,domain_lonsloc(2)) &
336  .OR. domain_lone > real_lonuy_mod(ia,domain_loneloc(2)) ) then
337  check_idl = .true.
338  domain_lons = minval(real_lonuy_mod(:,:),mask=(real_lonuy_mod>0.0_rp))
339  domain_lone = maxval(real_lonuy_mod(:,:),mask=(real_lonuy_mod<0.0_rp))
340  endif
341 
342  tile_nlat = grads_nlat
343  tile_nlon = grads_nlon
344  log_info("CNV2D_GrADS",*) 'Size of data in each tile (j) = ', tile_nlat
345  log_info("CNV2D_GrADS",*) 'Size of data in each tile (i) = ', tile_nlon
346 
347  tile_dlat = grads_dlat * d2r
348  tile_dlon = grads_dlon * d2r
349  log_info("CNV2D_GrADS",*) 'TILE_DLAT :', tile_dlat/d2r
350  log_info("CNV2D_GrADS",*) 'TILE_DLON :', tile_dlon/d2r
351 
352  !---< READ from external files >---
353 
354  if ( grads_in_catalogue /= '' ) then
355 
356  ! input from catalogue file
357  fname = trim(grads_in_dir)//'/'//trim(grads_in_catalogue)
358 
359  log_newline
360  log_info("CNV2D_GrADS",*) '+ Input catalogue file:', trim(fname)
361 
362  fid = io_get_available_fid()
363  open( fid, &
364  file = trim(fname), &
365  form = 'formatted', &
366  status = 'old', &
367  iostat = ierr )
368 
369  if ( ierr /= 0 ) then
370  log_error("CNV2D_GrADS",*) 'catalogue file not found!', trim(fname)
371  call prc_abort
372  endif
373 
374  do t = 1, tile_nlim
375  read(fid,*,iostat=ierr) index, tile_lats(t), tile_late(t), & ! South->North
376  tile_lons(t), tile_lone(t), & ! WEST->EAST
377  tile_fname(t)
378  if ( ierr /= 0 ) exit
379  enddo
380 
381  tile_nmax = t - 1
382  close(fid)
383 
384  elseif( grads_in_filename /= '' ) then
385 
386  ! input from single file
387  tile_nmax = 1
388  tile_fname(1) = grads_in_filename
389  tile_lats(1) = grads_lat_start
390  tile_late(1) = grads_lat_end
391  tile_lons(1) = grads_lon_start
392  tile_lone(1) = grads_lon_end
393 
394  endif
395 
396  do t = 1, tile_nmax
397  if ( tile_lons(t) >= 180.0_rp ) then
398  tile_lons(t) = tile_lons(t) - 360.0_rp
399  tile_lone(t) = tile_lone(t) - 360.0_rp
400  endif
401  if ( tile_lons(t) < -180.0_rp ) tile_lons(t) = tile_lons(t) + 360.0_rp
402  if ( tile_lone(t) < -180.0_rp ) tile_lone(t) = tile_lone(t) + 360.0_rp
403  enddo
404 
405  do nowstep = 1, grads_nstep
406 
407  log_info("CNV2D_GrADS",*) 'step = ', nowstep
408 
409  var2d(:,:) = 0.0_rp
410 
411  if ( cnv2d_interpolation_type == 'bilinear' &
412  .OR. cnv2d_interpolation_type == 'nearestneighbor' ) then
413 
414  call cnv2d_grads_bilinear ( domain_lats, & ! [IN]
415  domain_late, & ! [IN]
416  domain_lons, & ! [IN]
417  domain_lone, & ! [IN]
418  tile_nmax, & ! [IN]
419  grads_in_dir, & ! [IN]
420  tile_fname(1:tile_nmax), & ! [IN]
421  tile_lats(1:tile_nmax), & ! [IN]
422  tile_late(1:tile_nmax), & ! [IN]
423  tile_lons(1:tile_nmax), & ! [IN]
424  tile_lone(1:tile_nmax), & ! [IN]
425  tile_nlat, & ! [IN]
426  tile_nlon, & ! [IN]
427  tile_dlat, & ! [IN]
428  tile_dlon, & ! [IN]
429  grads_in_datatype, & ! [IN]
430  check_idl, & ! [IN]
431  grads_latorder_n2s, & ! [IN]
432  grads_missingval, & ! [IN]
433  nowstep, & ! [IN]
434  var2d(:,:) ) ! [OUT]
435 
436  elseif( cnv2d_interpolation_type == 'areaweighted' ) then
437 
438  call cnv2d_grads_areaweighted( domain_lats, & ! [IN]
439  domain_late, & ! [IN]
440  domain_lons, & ! [IN]
441  domain_lone, & ! [IN]
442  tile_nmax, & ! [IN]
443  grads_in_dir, & ! [IN]
444  tile_fname(1:tile_nmax), & ! [IN]
445  tile_lats(1:tile_nmax), & ! [IN]
446  tile_late(1:tile_nmax), & ! [IN]
447  tile_lons(1:tile_nmax), & ! [IN]
448  tile_lone(1:tile_nmax), & ! [IN]
449  tile_nlat, & ! [IN]
450  tile_nlon, & ! [IN]
451  tile_dlat, & ! [IN]
452  tile_dlon, & ! [IN]
453  grads_in_datatype, & ! [IN]
454  check_idl, & ! [IN]
455  grads_latorder_n2s, & ! [IN]
456  grads_missingval, & ! [IN]
457  nowstep, & ! [IN]
458  var2d(:,:) ) ! [OUT]
459 
460  endif
461 
462  ! output 2D file
463  call cnv2d_write( var2d(:,:), & ! [IN]
464  var2d_dtsec, & ! [IN]
465  nowstep ) ! [IN]
466  enddo
467 
468  return
469  end subroutine cnv2d_grads
470 
471  !-----------------------------------------------------------------------------
473  subroutine cnv2d_grads_bilinear( &
474  DOMAIN_LATS, &
475  DOMAIN_LATE, &
476  DOMAIN_LONS, &
477  DOMAIN_LONE, &
478  TILE_nmax, &
479  TILE_dir, &
480  TILE_fname, &
481  TILE_LATS, &
482  TILE_LATE, &
483  TILE_LONS, &
484  TILE_LONE, &
485  TILE_NLAT, &
486  TILE_NLON, &
487  TILE_DLAT, &
488  TILE_DLON, &
489  TILE_DATATYPE, &
490  TILE_check_IDL, &
491  TILE_LATORDER_N2S, &
492  TILE_MISSINGVAL, &
493  nowstep, &
494  VAR2D )
495  use scale_prc, only: &
496  prc_abort
497  use scale_const, only: &
498  radius => const_radius, &
499  pi => const_pi, &
500  d2r => const_d2r
501  use scale_atmos_grid_cartesc_real, only: &
502  real_lat => atmos_grid_cartesc_real_lat, &
503  real_lon => atmos_grid_cartesc_real_lon
504  implicit none
505 
506  real(RP), intent(in) :: DOMAIN_LATS
507  real(RP), intent(in) :: DOMAIN_LATE
508  real(RP), intent(in) :: DOMAIN_LONS
509  real(RP), intent(in) :: DOMAIN_LONE
510  integer, intent(in) :: TILE_nmax
511  character(len=H_LONG), intent(in) :: TILE_dir
512  character(len=H_LONG), intent(in) :: TILE_fname(tile_nmax)
513  real(RP), intent(in) :: TILE_LATS (tile_nmax)
514  real(RP), intent(in) :: TILE_LATE (tile_nmax)
515  real(RP), intent(in) :: TILE_LONS (tile_nmax)
516  real(RP), intent(in) :: TILE_LONE (tile_nmax)
517  integer, intent(in) :: TILE_NLAT
518  integer, intent(in) :: TILE_NLON
519  real(RP), intent(in) :: TILE_DLAT
520  real(RP), intent(in) :: TILE_DLON
521  character(len=H_LONG), intent(in) :: TILE_DATATYPE
522  logical, intent(in) :: TILE_check_IDL
523  logical, intent(in) :: TILE_LATORDER_N2S
524  real(RP), intent(in) :: TILE_MISSINGVAL
525  integer, intent(in) :: nowstep
526  real(RP), intent(out) :: VAR2D(ia,ja)
527 
528  real(RP) :: REAL_LON_mod(ia,ja)
529 
530  real(RP) :: TILE_LAT (0:tile_nlat+1)
531  real(RP) :: TILE_LON (0:tile_nlon+1)
532  real(RP) :: TILE_VALUE (tile_nlon,tile_nlat)
533  real(DP) :: TILE_VALUE_DP(tile_nlon,tile_nlat)
534  real(SP) :: TILE_VALUE_SP(tile_nlon,tile_nlat)
535  integer(2) :: TILE_VALUE_I2(tile_nlon,tile_nlat)
536 
537  integer :: jloc_b
538  integer :: jloc_t
539  integer :: iloc_l
540  integer :: iloc_r
541  real(RP) :: jfrac_t ! fraction for jloc_t
542  real(RP) :: ifrac_r ! fraction for iloc_r
543 
544  character(len=H_LONG) :: fname
545 
546  logical :: hit_lat, hit_lon
547  integer :: fid, ierr
548  integer :: i, j, ii, jj, t
549  !---------------------------------------------------------------------------
550 
551  real_lon_mod(:,:) = mod( real_lon(:,:)+3.0_dp*pi, 2.0_dp*pi ) - pi ! [-> 0..180,-180..0]
552 
553  ! data file
554  do t = 1, tile_nmax
555  hit_lat = .false.
556  hit_lon = .false.
557 
558  if ( ( tile_lats(t)*d2r >= domain_lats .AND. tile_lats(t)*d2r < domain_late ) &
559  .OR. ( tile_late(t)*d2r >= domain_lats .AND. tile_late(t)*d2r < domain_late ) ) then
560  hit_lat = .true.
561  endif
562 
563  if ( ( domain_lats >= tile_lats(t)*d2r .AND. domain_lats < tile_late(t)*d2r ) &
564  .OR. ( domain_late >= tile_lats(t)*d2r .AND. domain_late < tile_late(t)*d2r ) ) then
565  hit_lat = .true.
566  endif
567 
568  if ( tile_check_idl ) then
569  if ( ( tile_lons(t)*d2r >= domain_lons .AND. tile_lons(t)*d2r < pi ) &
570  .OR. ( tile_lons(t)*d2r >= -pi .AND. tile_lons(t)*d2r < domain_lone ) &
571  .OR. ( tile_lone(t)*d2r >= domain_lons .AND. tile_lone(t)*d2r < pi ) &
572  .OR. ( tile_lone(t)*d2r >= -pi .AND. tile_lone(t)*d2r < domain_lone ) ) then
573  hit_lon = .true.
574  endif
575  else
576  if ( ( tile_lons(t)*d2r >= domain_lons .AND. tile_lons(t)*d2r < domain_lone ) &
577  .OR. ( tile_lone(t)*d2r >= domain_lons .AND. tile_lone(t)*d2r < domain_lone ) ) then
578  hit_lon = .true.
579  endif
580  endif
581 
582  if ( ( domain_lons >= tile_lons(t)*d2r .AND. domain_lons < tile_lone(t)*d2r ) &
583  .OR. ( domain_lone >= tile_lons(t)*d2r .AND. domain_lone < tile_lone(t)*d2r ) ) then
584  hit_lon = .true.
585  endif
586 
587  if ( ( domain_lons+2.0_rp*pi >= tile_lons(t)*d2r .AND. domain_lons+2.0_rp*pi < tile_lone(t)*d2r ) &
588  .OR. ( domain_lone+2.0_rp*pi >= tile_lons(t)*d2r .AND. domain_lone+2.0_rp*pi < tile_lone(t)*d2r ) ) then
589  hit_lon = .true.
590  endif
591 
592  if ( hit_lat .AND. hit_lon ) then
593  fname = trim(tile_dir)//'/'//trim(tile_fname(t))
594 
595  log_newline
596  log_info("CNV2D_GrADS_bilinear",*) '+ Input data file :', trim(fname)
597  log_info_cont(*) 'Domain (LAT) :', domain_lats/d2r, domain_late/d2r
598  log_info_cont(*) ' (LON) :', domain_lons/d2r, domain_lone/d2r
599  if ( tile_check_idl ) then
600  log_info_cont(*) '(Date line exists within the domain)'
601  endif
602  log_info_cont(*) 'Tile (LAT) :', tile_lats(t), tile_late(t)
603  log_info_cont(*) ' (LON) :', tile_lons(t), tile_lone(t)
604 
605  if ( tile_datatype == 'REAL8' ) then
606 
607  fid = io_get_available_fid()
608  open( fid, &
609  file = trim(fname), &
610  form = 'unformatted', &
611  access = 'direct', &
612  status = 'old', &
613  recl = tile_nlon*tile_nlat*8, &
614  iostat = ierr )
615 
616  if ( ierr /= 0 ) then
617  log_error("CNV2D_GrADS_bilinear",*) 'data file not found!'
618  call prc_abort
619  endif
620 
621  read(fid,rec=nowstep) tile_value_dp(:,:)
622  close(fid)
623 
624  if ( tile_latorder_n2s ) then
625  do j = 1, tile_nlat
626  do i = 1, tile_nlon
627  tile_value(i,j) = real( TILE_VALUE_DP(i,TILE_NLAT-j+1), kind=RP ) ! reverse order
628  enddo
629  enddo
630  else
631  do j = 1, tile_nlat
632  do i = 1, tile_nlon
633  tile_value(i,j) = real( TILE_VALUE_DP(i,j), kind=rp )
634  enddo
635  enddo
636  endif
637 
638  elseif( tile_datatype == 'REAL4' ) then
639 
640  fid = io_get_available_fid()
641  open( fid, &
642  file = trim(fname), &
643  form = 'unformatted', &
644  access = 'direct', &
645  status = 'old', &
646  recl = tile_nlon*tile_nlat*4, &
647  iostat = ierr )
648 
649  if ( ierr /= 0 ) then
650  log_error("CNV2D_GrADS_bilinear",*) 'data file not found!'
651  call prc_abort
652  endif
653 
654  read(fid,rec=nowstep) tile_value_sp(:,:)
655  close(fid)
656 
657  if ( tile_latorder_n2s ) then
658  do j = 1, tile_nlat
659  do i = 1, tile_nlon
660  tile_value(i,j) = real( TILE_VALUE_SP(i,TILE_NLAT-j+1), kind=RP ) ! reverse order
661  enddo
662  enddo
663  else
664  do j = 1, tile_nlat
665  do i = 1, tile_nlon
666  tile_value(i,j) = real( TILE_VALUE_SP(i,j), kind=rp )
667  enddo
668  enddo
669  endif
670 
671  elseif( tile_datatype == 'INT2' ) then
672 
673  fid = io_get_available_fid()
674  open( fid, &
675  file = trim(fname), &
676  form = 'unformatted', &
677  access = 'direct', &
678  status = 'old', &
679  recl = tile_nlon*tile_nlat*2, &
680  iostat = ierr )
681 
682  if ( ierr /= 0 ) then
683  log_error("CNV2D_GrADS_bilinear",*) 'data file not found!'
684  call prc_abort
685  endif
686 
687  read(fid,rec=nowstep) tile_value_i2(:,:)
688  close(fid)
689 
690  if ( tile_latorder_n2s ) then
691  do j = 1, tile_nlat
692  do i = 1, tile_nlon
693  tile_value(i,j) = real( TILE_VALUE_I2(i,TILE_NLAT-j+1), kind=RP ) ! reverse order
694  enddo
695  enddo
696  else
697  do j = 1, tile_nlat
698  do i = 1, tile_nlon
699  tile_value(i,j) = real( TILE_VALUE_I2(i,j), kind=rp )
700  enddo
701  enddo
702  endif
703 
704  endif
705 
706  call prof_valcheck('CNV2D','VAR',tile_value(:,:))
707 
708  tile_lat(0) = tile_lats(t) * d2r - tile_dlat
709  do jj = 1, tile_nlat+1
710  tile_lat(jj) = tile_lat(jj-1) + tile_dlat
711  !LOG_INFO("CNV2D_GrADS_bilinear",*) jj, TILE_LAT(jj) / D2R
712  enddo
713 
714  tile_lon(0) = tile_lons(t) * d2r - tile_dlon
715  do ii = 1, tile_nlon+1
716  tile_lon(ii) = tile_lon(ii-1) + tile_dlon
717  if( tile_lon(ii) > pi ) tile_lon(ii) = tile_lon(ii) - 2.0_rp*pi
718  !LOG_INFO("CNV2D_GrADS_bilinear",*) ii, TILE_LON(ii) / D2R
719  enddo
720 
721  ! match and calc fraction
722  do jj = 1, tile_nlat+1
723  do ii = 1, tile_nlon+1
724  if ( tile_lat(jj ) < domain_lats &
725  .OR. tile_lat(jj-1) > domain_late ) then
726  cycle
727  endif
728 
729  if ( tile_check_idl .OR. tile_lon(ii-1) >= tile_lon(ii) ) then
730  ! do nothing
731  else
732  if ( tile_lon(ii ) < domain_lons &
733  .OR. tile_lon(ii-1) > domain_lone ) then
734  cycle
735  endif
736  endif
737 
738  do j = 1, ja
739  do i = 1, ia
740 
741  iloc_l = -1
742  iloc_r = -1
743  jloc_b = -1
744  jloc_t = -1
745 
746  ! normal
747  if ( real_lon_mod(i,j) >= tile_lon(ii-1) &
748  .AND. real_lon_mod(i,j) < tile_lon(ii ) &
749  .AND. real_lat(i,j) >= tile_lat(jj-1) &
750  .AND. real_lat(i,j) < tile_lat(jj ) ) then
751 
752  iloc_l = ii-1
753  iloc_r = ii
754  ifrac_r = min( real_lon_mod(i,j)-tile_lon(ii-1), tile_dlon ) / tile_dlon
755 
756  jloc_b = jj-1
757  jloc_t = jj
758  jfrac_t = min( real_lat(i,j)-tile_lat(jj-1), tile_dlat ) / tile_dlat
759 
760  endif
761 
762  ! across IDL
763  if ( tile_lon(ii-1) >= tile_lon(ii) &
764  .AND. real_lat(i,j) >= tile_lat(jj-1) &
765  .AND. real_lat(i,j) < tile_lat(jj ) ) then
766 
767  if ( real_lon_mod(i,j) >= tile_lon(ii-1) &
768  .AND. real_lon_mod(i,j) < pi ) then
769 
770  iloc_l = ii-1
771  iloc_r = ii
772  ifrac_r = min( real_lon_mod(i,j)-tile_lon(ii-1), tile_dlon ) / tile_dlon
773 
774  jloc_b = jj-1
775  jloc_t = jj
776  jfrac_t = min( real_lat(i,j)-tile_lat(jj-1), tile_dlat ) / tile_dlat
777 
778  elseif( real_lon_mod(i,j) >= -pi &
779  .AND. real_lon_mod(i,j) < tile_lon(ii ) ) then
780 
781  iloc_l = ii-1
782  iloc_r = ii
783  ifrac_r = min( real_lon_mod(i,j)+2.0_rp*pi-tile_lon(ii-1), tile_dlon ) / tile_dlon
784 
785  jloc_b = jj-1
786  jloc_t = jj
787  jfrac_t = min( real_lat(i,j)-tile_lat(jj-1), tile_dlat ) / tile_dlat
788 
789  endif
790 
791  endif
792 
793  if( iloc_r == -1 .AND. jloc_t == -1 ) cycle
794 
795  if( tile_lon(0) <= 0.0_rp .AND. iloc_l == 0 ) iloc_l = tile_nlon ! around prime meridian
796  if( tile_lon(tile_nlon+1) >= 0.0_rp .AND. iloc_r == tile_nlon+1 ) iloc_r = 1 ! around prime meridian
797 
798  if( tile_lat(0) <= -0.5_rp*pi .AND. jloc_b == 0 ) jloc_b = 1 ! around south pole
799  if( tile_lat(tile_nlat+1) >= 0.5_rp*pi .AND. jloc_t == tile_nlat+1 ) jloc_t = tile_nlat ! around north pole
800 
801  if ( cnv2d_interpolation_type == 'nearestneighbor' ) then
802  if ( ifrac_r >= 0.5d0 ) then
803  iloc_l = iloc_r
804  else
805  iloc_r = iloc_l
806  endif
807 
808  if ( jfrac_t >= 0.5d0 ) then
809  jloc_b = jloc_t
810  else
811  jloc_t = jloc_b
812  endif
813  endif
814 
815  var2d(i,j) = ( ifrac_r) * ( jfrac_t) * tile_value(iloc_r,jloc_t) &
816  + (1.0_rp-ifrac_r) * ( jfrac_t) * tile_value(iloc_l,jloc_t) &
817  + ( ifrac_r) * (1.0_rp-jfrac_t) * tile_value(iloc_r,jloc_b) &
818  + (1.0_rp-ifrac_r) * (1.0_rp-jfrac_t) * tile_value(iloc_l,jloc_b) + var2d(i,j)
819 
820  enddo
821  enddo
822 
823  enddo
824  enddo
825 
826  endif
827  enddo ! tile loop
828 
829  return
830  end subroutine cnv2d_grads_bilinear
831 
832  !-----------------------------------------------------------------------------
834  subroutine cnv2d_grads_areaweighted( &
835  DOMAIN_LATS, &
836  DOMAIN_LATE, &
837  DOMAIN_LONS, &
838  DOMAIN_LONE, &
839  TILE_nmax, &
840  TILE_dir, &
841  TILE_fname, &
842  TILE_LATS, &
843  TILE_LATE, &
844  TILE_LONS, &
845  TILE_LONE, &
846  TILE_NLAT, &
847  TILE_NLON, &
848  TILE_DLAT, &
849  TILE_DLON, &
850  TILE_DATATYPE, &
851  TILE_check_IDL, &
852  TILE_LATORDER_N2S, &
853  TILE_MISSINGVAL, &
854  nowstep, &
855  VAR2D )
856  use scale_prc, only: &
857  prc_abort
858  use scale_const, only: &
859  radius => const_radius, &
860  pi => const_pi, &
861  eps => const_eps, &
862  d2r => const_d2r
863  use scale_atmos_grid_cartesc_real, only: &
864  real_latxv => atmos_grid_cartesc_real_latxv, &
865  real_lonuy => atmos_grid_cartesc_real_lonuy
866  implicit none
867 
868  real(RP), intent(in) :: DOMAIN_LATS
869  real(RP), intent(in) :: DOMAIN_LATE
870  real(RP), intent(in) :: DOMAIN_LONS
871  real(RP), intent(in) :: DOMAIN_LONE
872  integer, intent(in) :: TILE_nmax
873  character(len=H_LONG), intent(in) :: TILE_dir
874  character(len=H_LONG), intent(in) :: TILE_fname(tile_nmax)
875  real(RP), intent(in) :: TILE_LATS (tile_nmax)
876  real(RP), intent(in) :: TILE_LATE (tile_nmax)
877  real(RP), intent(in) :: TILE_LONS (tile_nmax)
878  real(RP), intent(in) :: TILE_LONE (tile_nmax)
879  integer, intent(in) :: TILE_NLAT
880  integer, intent(in) :: TILE_NLON
881  real(RP), intent(in) :: TILE_DLAT
882  real(RP), intent(in) :: TILE_DLON
883  character(len=H_LONG), intent(in) :: TILE_DATATYPE
884  logical, intent(in) :: TILE_check_IDL
885  logical, intent(in) :: TILE_LATORDER_N2S
886  real(RP), intent(in) :: TILE_MISSINGVAL
887  integer, intent(in) :: nowstep
888  real(RP), intent(out) :: VAR2D(ia,ja)
889 
890  real(RP) :: REAL_LONUY_mod(0:ia,ja)
891 
892  real(RP) :: TILE_LATH (0:tile_nlat)
893  real(RP) :: TILE_LONH (0:tile_nlon)
894  real(RP) :: TILE_VALUE (tile_nlon,tile_nlat)
895  real(DP) :: TILE_VALUE_DP(tile_nlon,tile_nlat)
896  real(SP) :: TILE_VALUE_SP(tile_nlon,tile_nlat)
897  integer(2) :: TILE_VALUE_I2(tile_nlon,tile_nlat)
898 
899  integer :: jloc, iloc
900  real(RP) :: jfrac_b ! fraction for jloc
901  real(RP) :: ifrac_l ! fraction for iloc
902 
903  real(RP) :: val_sum (ia,ja)
904  real(RP) :: area_sum(ia,ja)
905  real(RP) :: val, mask
906  real(RP) :: area, area_fraction
907 
908  character(len=H_LONG) :: fname
909 
910  real(RP) :: zerosw
911  logical :: hit_lat, hit_lon
912  integer :: fid, ierr
913  integer :: i, j, ii, jj, t
914  !---------------------------------------------------------------------------
915 
916  real_lonuy_mod(:,:) = mod( real_lonuy(:,:)+3.0_dp*pi, 2.0_dp*pi ) - pi ! [-> 0..180,-180..0]
917 
918  ! data file
919  do t = 1, tile_nmax
920  hit_lat = .false.
921  hit_lon = .false.
922 
923  if ( ( tile_lats(t)*d2r >= domain_lats .AND. tile_lats(t)*d2r < domain_late ) &
924  .OR. ( tile_late(t)*d2r >= domain_lats .AND. tile_late(t)*d2r < domain_late ) ) then
925  hit_lat = .true.
926  endif
927 
928  if ( ( domain_lats >= tile_lats(t)*d2r .AND. domain_lats < tile_late(t)*d2r ) &
929  .OR. ( domain_late >= tile_lats(t)*d2r .AND. domain_late < tile_late(t)*d2r ) ) then
930  hit_lat = .true.
931  endif
932 
933  if ( tile_check_idl ) then
934  if ( ( tile_lons(t)*d2r >= domain_lons .AND. tile_lons(t)*d2r < pi ) &
935  .OR. ( tile_lons(t)*d2r >= -pi .AND. tile_lons(t)*d2r < domain_lone ) &
936  .OR. ( tile_lone(t)*d2r >= domain_lons .AND. tile_lone(t)*d2r < pi ) &
937  .OR. ( tile_lone(t)*d2r >= -pi .AND. tile_lone(t)*d2r < domain_lone ) ) then
938  hit_lon = .true.
939  endif
940  else
941  if ( ( tile_lons(t)*d2r >= domain_lons .AND. tile_lons(t)*d2r < domain_lone ) &
942  .OR. ( tile_lone(t)*d2r >= domain_lons .AND. tile_lone(t)*d2r < domain_lone ) ) then
943  hit_lon = .true.
944  endif
945  endif
946 
947  if ( ( domain_lons >= tile_lons(t)*d2r .AND. domain_lons < tile_lone(t)*d2r ) &
948  .OR. ( domain_lone >= tile_lons(t)*d2r .AND. domain_lone < tile_lone(t)*d2r ) ) then
949  hit_lon = .true.
950  endif
951 
952  if ( hit_lat .AND. hit_lon ) then
953  fname = trim(tile_dir)//'/'//trim(tile_fname(t))
954 
955  log_newline
956  log_info("CNV2D_GrADS_areaweighted",*) '+ Input data file :', trim(fname)
957  log_info_cont(*) 'Domain (LAT) :', domain_lats/d2r, domain_late/d2r
958  log_info_cont(*) ' (LON) :', domain_lons/d2r, domain_lone/d2r
959  if ( tile_check_idl ) then
960  log_info_cont(*) '(Date line exists within the domain)'
961  endif
962  log_info_cont(*) 'Tile (LAT) :', tile_lats(t), tile_late(t)
963  log_info_cont(*) ' (LON) :', tile_lons(t), tile_lone(t)
964 
965  if ( tile_datatype == 'REAL8' ) then
966 
967  fid = io_get_available_fid()
968  open( fid, &
969  file = trim(fname), &
970  form = 'unformatted', &
971  access = 'direct', &
972  status = 'old', &
973  recl = tile_nlon*tile_nlat*8, &
974  iostat = ierr )
975 
976  if ( ierr /= 0 ) then
977  log_error("CNV2D_GrADS_areaweighted",*) 'data file not found!'
978  call prc_abort
979  endif
980 
981  read(fid,rec=nowstep) tile_value_dp(:,:)
982  close(fid)
983 
984  if ( tile_latorder_n2s ) then
985  do j = 1, tile_nlat
986  do i = 1, tile_nlon
987  tile_value(i,j) = real( TILE_VALUE_DP(i,TILE_NLAT-j+1), kind=RP ) ! reverse order
988  enddo
989  enddo
990  else
991  do j = 1, tile_nlat
992  do i = 1, tile_nlon
993  tile_value(i,j) = real( TILE_VALUE_DP(i,j), kind=rp )
994  enddo
995  enddo
996  endif
997 
998  elseif( tile_datatype == 'REAL4' ) then
999 
1000  fid = io_get_available_fid()
1001  open( fid, &
1002  file = trim(fname), &
1003  form = 'unformatted', &
1004  access = 'direct', &
1005  status = 'old', &
1006  recl = tile_nlon*tile_nlat*4, &
1007  iostat = ierr )
1008 
1009  if ( ierr /= 0 ) then
1010  log_error("CNV2D_GrADS_areaweighted",*) 'data file not found!'
1011  call prc_abort
1012  endif
1013 
1014  read(fid,rec=nowstep) tile_value_sp(:,:)
1015  close(fid)
1016 
1017  if ( tile_latorder_n2s ) then
1018  do j = 1, tile_nlat
1019  do i = 1, tile_nlon
1020  tile_value(i,j) = real( TILE_VALUE_SP(i,TILE_NLAT-j+1), kind=RP ) ! reverse order
1021  enddo
1022  enddo
1023  else
1024  do j = 1, tile_nlat
1025  do i = 1, tile_nlon
1026  tile_value(i,j) = real( TILE_VALUE_SP(i,j), kind=rp )
1027  enddo
1028  enddo
1029  endif
1030 
1031  elseif( tile_datatype == 'INT2' ) then
1032 
1033  fid = io_get_available_fid()
1034  open( fid, &
1035  file = trim(fname), &
1036  form = 'unformatted', &
1037  access = 'direct', &
1038  status = 'old', &
1039  recl = tile_nlon*tile_nlat*2, &
1040  iostat = ierr )
1041 
1042  if ( ierr /= 0 ) then
1043  log_error("CNV2D_GrADS_areaweighted",*) 'data file not found!'
1044  call prc_abort
1045  endif
1046 
1047  read(fid,rec=nowstep) tile_value_i2(:,:)
1048  close(fid)
1049 
1050  if ( tile_latorder_n2s ) then
1051  do j = 1, tile_nlat
1052  do i = 1, tile_nlon
1053  tile_value(i,j) = real( TILE_VALUE_I2(i,TILE_NLAT-j+1), kind=RP ) ! reverse order
1054  enddo
1055  enddo
1056  else
1057  do j = 1, tile_nlat
1058  do i = 1, tile_nlon
1059  tile_value(i,j) = real( TILE_VALUE_I2(i,j), kind=rp )
1060  enddo
1061  enddo
1062  endif
1063 
1064  endif
1065 
1066  call prof_valcheck('CNV2D','VAR',tile_value(:,:))
1067 
1068  tile_lath(0) = tile_lats(t) * d2r
1069  do jj = 1, tile_nlat
1070  tile_lath(jj) = tile_lath(jj-1) + tile_dlat
1071  !LOG_INFO("CNV2D_GrADS_areaweighted",*) jj, TILE_LATH(jj) / D2R
1072  enddo
1073 
1074  tile_lonh(0) = tile_lons(t) * d2r
1075  do ii = 1, tile_nlon
1076  tile_lonh(ii) = tile_lonh(ii-1) + tile_dlon
1077  if( tile_lonh(ii) > pi ) tile_lonh(ii) = tile_lonh(ii) - 2.0_rp*pi
1078  !LOG_INFO("CNV2D_GrADS_areaweighted",*) ii, TILE_LONH(ii) / D2R
1079  enddo
1080 
1081  ! match and calc fraction
1082  do jj = 1, tile_nlat
1083  do ii = 1, tile_nlon
1084 
1085  iloc = 1 ! Z_sfc(1,1) is used for dummy grid
1086  ifrac_l = 1.0_rp
1087 
1088  jloc = 1 ! Z_sfc(1,1) is used for dummy grid
1089  jfrac_b = 1.0_rp
1090 
1091  if ( tile_lath(jj ) < domain_lats &
1092  .OR. tile_lath(jj-1) > domain_late ) then
1093  cycle
1094  endif
1095 
1096  if ( tile_check_idl ) then
1097  if ( tile_lonh(ii ) < domain_lons &
1098  .AND. tile_lonh(ii-1) > domain_lone ) then
1099  cycle
1100  endif
1101  else
1102  if ( tile_lonh(ii ) < domain_lons &
1103  .OR. tile_lonh(ii-1) > domain_lone ) then
1104  cycle
1105  endif
1106  endif
1107 
1108  jloop: do j = js-1, je+1
1109  iloop: do i = is-1, ie+1
1110  if ( tile_lonh(ii-1) >= real_lonuy_mod(i-1,j ) &
1111  .AND. tile_lonh(ii-1) < real_lonuy_mod(i ,j ) &
1112  .AND. tile_lath(jj-1) >= real_latxv(i ,j-1) &
1113  .AND. tile_lath(jj-1) < real_latxv(i ,j ) ) then
1114 
1115  iloc = i
1116  ifrac_l = min( real_lonuy_mod(i,j)-tile_lonh(ii-1), tile_dlon ) / tile_dlon
1117 
1118  jloc = j
1119  jfrac_b = min( real_latxv(i,j)-tile_lath(jj-1), tile_dlat ) / tile_dlat
1120  exit jloop
1121 
1122  endif
1123 
1124  if ( real_lonuy_mod(i-1,j) >= real_lonuy_mod(i ,j ) &
1125  .AND. tile_lath(jj-1) >= real_latxv(i ,j-1) &
1126  .AND. tile_lath(jj-1) < real_latxv(i ,j ) ) then ! across the IDL
1127 
1128  if ( tile_lonh(ii-1) >= real_lonuy_mod(i-1,j) &
1129  .AND. tile_lonh(ii-1) < pi ) then
1130 
1131  iloc = i
1132  ifrac_l = min( real_lonuy_mod(i,j)-tile_lonh(ii-1)+2.0_rp*pi, tile_dlon ) / tile_dlon
1133 
1134  jloc = j
1135  jfrac_b = min( real_latxv(i,j)-tile_lath(jj-1), tile_dlat ) / tile_dlat
1136  exit jloop
1137 
1138  elseif( tile_lonh(ii-1) >= -pi &
1139  .AND. tile_lonh(ii-1) < real_lonuy_mod(i ,j) ) then
1140 
1141  iloc = i
1142  ifrac_l = min( real_lonuy_mod(i,j)-tile_lonh(ii-1), tile_dlon ) / tile_dlon
1143 
1144  jloc = j
1145  jfrac_b = min( real_latxv(i,j)-tile_lath(jj-1), tile_dlat ) / tile_dlat
1146  exit jloop
1147 
1148  endif
1149 
1150  endif
1151  enddo iloop
1152  enddo jloop
1153 
1154  if( iloc == 1 .AND. jloc == 1 ) cycle
1155 
1156  val = tile_value(ii,jj)
1157  mask = 0.0_rp
1158 
1159  if ( abs(val-tile_missingval) < eps ) mask = 1.0_rp! if missing value, mask = 1
1160 
1161  area = radius * radius * tile_dlon * ( sin(tile_lath(jj))-sin(tile_lath(jj-1)) ) * ( 1.0_rp - mask )
1162 
1163 ! LOG_INFO("CNV2D_GrADS_areaweighted",*) ii, jj, area, iloc, jloc, ifrac_l, jfrac_b, TILE_VALUE(ii,jj)
1164 
1165  area_fraction = ( ifrac_l) * ( jfrac_b) * area
1166  area_sum(iloc ,jloc ) = area_sum(iloc ,jloc ) + area_fraction
1167  val_sum(iloc ,jloc ) = val_sum(iloc ,jloc ) + area_fraction * val
1168 
1169  area_fraction = (1.0_rp-ifrac_l) * ( jfrac_b) * area
1170  area_sum(iloc+1,jloc ) = area_sum(iloc+1,jloc ) + area_fraction
1171  val_sum(iloc+1,jloc ) = val_sum(iloc+1,jloc ) + area_fraction * val
1172 
1173  area_fraction = ( ifrac_l) * (1.0_rp-jfrac_b) * area
1174  area_sum(iloc ,jloc+1) = area_sum(iloc ,jloc+1) + area_fraction
1175  val_sum(iloc ,jloc+1) = val_sum(iloc ,jloc+1) + area_fraction * val
1176 
1177  area_fraction = (1.0_rp-ifrac_l) * (1.0_rp-jfrac_b) * area
1178  area_sum(iloc+1,jloc+1) = area_sum(iloc+1,jloc+1) + area_fraction
1179  val_sum(iloc+1,jloc+1) = val_sum(iloc+1,jloc+1) + area_fraction * val
1180  enddo
1181  enddo
1182 
1183  endif
1184  enddo ! tile loop
1185 
1186  do j = js, je
1187  do i = is, ie
1188  zerosw = 0.5_rp - sign( 0.5_rp, area_sum(i,j)-eps )
1189  var2d(i,j) = val_sum(i,j) * ( 1.0_rp-zerosw ) / ( area_sum(i,j)-zerosw )
1190  enddo
1191  enddo
1192 
1193  return
1194  end subroutine cnv2d_grads_areaweighted
1195 
1196  !-----------------------------------------------------------------------------
1198  subroutine cnv2d_fillhalo( VAR, FILL_BND )
1199  use scale_comm_cartesc, only: &
1200  comm_vars8, &
1201  comm_wait
1202  implicit none
1203 
1204  real(RP), intent(inout) :: VAR(ia,ja)
1205 
1206  logical, intent(in), optional :: FILL_BND
1207 
1208  logical :: FILL_BND_
1209  !---------------------------------------------------------------------------
1210 
1211  fill_bnd_ = .false.
1212  if ( present(fill_bnd) ) fill_bnd_ = fill_bnd
1213 
1214  call comm_vars8( var(:,:), 1 )
1215  call comm_wait ( var(:,:), 1, fill_bnd_ )
1216 
1217  return
1218  end subroutine cnv2d_fillhalo
1219 
1220  !-----------------------------------------------------------------------------
1222  subroutine cnv2d_write( &
1223  VAR2D, &
1224  timeintv, &
1225  istep )
1226  use scale_time, only: &
1227  time_nowdate
1228  use scale_file_cartesc, only: &
1229  file_cartesc_write
1230  implicit none
1231 
1232  real(RP), intent(in) :: VAR2D(ia,ja)
1233  real(DP), intent(in) :: timeintv
1234  integer, intent(in) :: istep
1235 
1236  real(RP) :: work(ia,ja,1)
1237  real(DP) :: timeofs
1238  !---------------------------------------------------------------------------
1239 
1240  if ( cnv2d_out_basename /= '' ) then
1241  log_newline
1242  log_info("CNV2D_write",*) 'Output converted 2D file '
1243 
1244  work(:,:,1) = var2d(:,:)
1245 
1246  call cnv2d_fillhalo( work(:,:,1), fill_bnd=.false. )
1247 
1248  timeofs = real(istep-1,kind=DP) * timeintv
1249 
1250  call file_cartesc_write( work(:,:,:), & ! [IN]
1251  cnv2d_out_basename, & ! [IN]
1252  cnv2d_out_title, & ! [IN]
1253  cnv2d_out_varname, & ! [IN]
1254  cnv2d_out_vardesc, & ! [IN]
1255  cnv2d_out_varunit, & ! [IN]
1256  'XYT', & ! [IN]
1257  cnv2d_out_dtype, & ! [IN]
1258  timeintv, & ! [IN]
1259  time_nowdate, & ! [IN]
1260  timeofs=timeofs ) ! [IN]
1261  endif
1262 
1263  return
1264  end subroutine cnv2d_write
1265 
1266 end module mod_cnv2d
subroutine, public calendar_unit2sec(second, value, unit)
Convert several units to second.
integer, public ia
of whole cells: x, local, with HALO
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
real(rp), public const_radius
radius of the planet [m]
Definition: scale_const.F90:44
integer, public ja
of whole cells: y, local, with HALO
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
real(rp), public const_undef
Definition: scale_const.F90:41
module COMMUNICATION
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module TRACER
logical, public cnv2d_donothing
Definition: mod_cnv2d.F90:35
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_dlon
delta longitude
module atmosphere / grid / cartesC index
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:313
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
module TIME
Definition: scale_time.F90:16
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]
module Convert 2D data
Definition: mod_cnv2d.F90:11
logical, public cnv2d_usegrads
Definition: mod_cnv2d.F90:36
subroutine, public cnv2d_setup
Setup.
Definition: mod_cnv2d.F90:71
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_dlat
delta latitude
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
module profiler
Definition: scale_prof.F90:11
subroutine cnv2d_fillhalo(VAR, FILL_BND)
HALO Communication.
Definition: mod_cnv2d.F90:1199
real(rp), public const_eps
small number
Definition: scale_const.F90:33
module Atmosphere GRID CartesC Real(real space)
module PRECISION
module file / cartesianC
subroutine, public cnv2d
Driver.
Definition: mod_cnv2d.F90:155
real(rp), public const_pi
pi
Definition: scale_const.F90:31
module CALENDAR
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
module Statistics
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:69
module STDIO
Definition: scale_io.F90:10
integer, parameter, public rp