38 private :: copytopo_transgrid
39 private :: copytopo_setalpha
40 private :: copytopo_input_data
41 private :: copytopo_mix_data
47 integer,
private,
parameter :: handle = 1
49 character(len=H_SHORT),
private :: COPYTOPO_IN_FILETYPE =
''
50 character(len=H_LONG),
private :: COPYTOPO_IN_BASENAME =
''
51 integer,
private :: COPYTOPO_IN_PRC_NUM_X = 0
52 integer,
private :: COPYTOPO_IN_PRC_NUM_Y = 0
53 character(len=H_LONG),
private :: COPYTOPO_LATLON_CATALOGUE =
''
54 integer,
private :: COPYTOPO_DOMID
56 character(len=H_LONG),
private :: COPYTOPO_IN_POSTFIX =
''
57 character(len=H_SHORT),
private :: COPYTOPO_IN_VARNAME =
'topo'
58 character(len=H_SHORT),
private :: COPYTOPO_IN_LONNAME =
'lon'
59 character(len=H_SHORT),
private :: COPYTOPO_IN_LATNAME =
'lat'
60 real(RP),
private :: COPYTOPO_TRANSITION_DX = -1.0_rp
61 real(RP),
private :: COPYTOPO_TRANSITION_DY = -1.0_rp
62 real(RP),
private :: COPYTOPO_FRACX = 1.0_rp
63 real(RP),
private :: COPYTOPO_FRACY = 1.0_rp
64 real(RP),
private :: COPYTOPO_taux = 1.0_rp
65 real(RP),
private :: COPYTOPO_tauy = 1.0_rp
66 logical,
private :: COPYTOPO_ENTIRE_REGION = .false.
67 logical,
private :: COPYTOPO_LINEAR_H = .true.
68 real(RP),
private :: COPYTOPO_EXP_H = 2.0_rp
69 character(len=H_SHORT),
private :: COPYTOPO_INTRP_TYPE =
"LINEAR"
70 integer,
private :: COPYTOPO_hypdiff_ORDER = 4
71 integer,
private :: COPYTOPO_hypdiff_NITER = 20
83 real(rp),
intent(inout) :: topo_child(:,:)
87 real(rp) :: alpha (
ia,
ja)
88 real(rp) :: topo_parent(
ia,
ja)
90 namelist / param_copytopo / &
91 copytopo_in_filetype, &
92 copytopo_in_basename, &
93 copytopo_in_prc_num_x, &
94 copytopo_in_prc_num_y, &
95 copytopo_latlon_catalogue, &
96 copytopo_in_postfix, &
97 copytopo_in_varname, &
98 copytopo_in_lonname, &
99 copytopo_in_latname, &
100 copytopo_transition_dx, &
101 copytopo_transition_dy, &
106 copytopo_entire_region, &
109 copytopo_hypdiff_order, &
110 copytopo_hypdiff_niter
116 log_info(
"COPYTOPO",*)
'Setup'
122 log_info(
"COPYTOPO",*)
'Not found namelist. Default used.'
123 elseif( ierr > 0 )
then
124 log_error(
"COPYTOPO",*)
'Not appropriate names in namelist PARAM_COPYTOPO. Check!'
127 log_nml(param_copytopo)
131 call copytopo_transgrid( ctrx(:), ctry(:) )
133 call copytopo_setalpha ( ctrx(:), ctry(:), &
136 call copytopo_input_data( topo_parent(:,:) )
138 call copytopo_mix_data( topo_parent(:,:), &
147 subroutine copytopo_transgrid( &
166 real(rp),
intent(out) :: ctrx(
ia)
167 real(rp),
intent(out) :: ctry(
ja)
169 real(rp) :: ctrxg (
iag)
170 real(rp) :: ctryg (
jag)
172 real(rp) :: bufftotx, transtotx
173 real(rp) :: bufftoty, transtoty
174 integer :: imain, ibuff, itrans
175 integer :: jmain, jbuff, jtrans
176 integer :: copy_is, copy_ie
177 integer :: copy_js, copy_je
179 integer :: i, j, ii, jj
186 if( abs(cbfxg(i)) < eps )
exit
188 ibuff = i - 1 -
ihalo
189 if ( ibuff == 0 )
then
195 if ( copytopo_transition_dx < 0.0_rp )
then
196 copytopo_transition_dx = bufftotx
199 if( cxg(i) - bufftotx - fxg(
ihalo) >= copytopo_transition_dx )
exit
201 itrans = i - 1 -
ihalo - ibuff
202 if ( itrans == 0 )
then
205 transtotx = fxg(itrans+ibuff+
ihalo) - fxg(
ihalo) - bufftotx
208 imain =
iag - 2*ibuff - 2*itrans - 2*
ihalo
210 if ( imain < 1 )
then
211 log_error(
"COPYTOPO_transgrid",*)
'Not appropriate transition width for global domain(X).', copytopo_transition_dx
212 log_error_cont(*)
'# of buffer region (one side) = ', ibuff
213 log_error_cont(*)
'# of transion region (one side) = ', itrans
221 copy_ie =
ihalo+ibuff
223 do i = copy_is, copy_ie
227 if ( itrans > 0 )
then
228 copy_is =
ihalo+ibuff+1
229 copy_ie =
ihalo+ibuff+itrans
231 do i = copy_is, copy_ie
232 ctrxg(i) = (transtotx+bufftotx+fxg(
ihalo )-cxg(i)) / transtotx
234 copy_is =
iag -
ihalo - ibuff - itrans + 1
237 do i = copy_is, copy_ie
238 ctrxg(i) = (transtotx+bufftotx-fxg(
iag-
ihalo)+cxg(i)) / transtotx
245 do i = copy_is, copy_ie
249 ctrxg(:) = max( min( ctrxg(:), 1.0_rp ), 0.0_rp )
254 if( abs(cbfyg(j)) < eps )
exit
256 jbuff = j - 1 -
jhalo
257 if ( jbuff == 0 )
then
263 if ( copytopo_transition_dy < 0.0_rp )
then
264 copytopo_transition_dy = bufftoty
267 if( cyg(j) - bufftoty - fyg(
jhalo) >= copytopo_transition_dy )
exit
269 jtrans = j - 1 -
jhalo - jbuff
270 if ( jtrans == 0 )
then
273 transtoty = fyg(jtrans+jbuff+
jhalo) - fyg(
jhalo) - bufftoty
276 jmain =
jag - 2*jbuff - 2*jtrans - 2*
jhalo
278 if ( jmain < 1 )
then
279 log_error(
"COPYTOPO_transgrid",*)
'Not appropriate transition width for global domain(Y).', copytopo_transition_dy
280 log_error_cont(*)
'# of buffer region (one side)', jbuff
281 log_error_cont(*)
'# of transion region (one side)', jtrans
289 copy_je =
jhalo+jbuff
291 do j = copy_js, copy_je
295 if ( jtrans > 0 )
then
296 copy_js =
jhalo+jbuff+1
297 copy_je =
jhalo+jbuff+jtrans
299 do j = copy_js, copy_je
300 ctryg(j) = (transtoty+bufftoty+fyg(
jhalo )-cyg(j)) / transtoty
302 copy_js =
jag -
jhalo - jbuff - jtrans + 1
305 do j = copy_js, copy_je
306 ctryg(j) = (transtoty+bufftoty-fyg(
jag-
jhalo)+cyg(j)) / transtoty
313 do j = copy_js, copy_je
317 ctryg(:) = max( min( ctryg(:), 1.0_rp ), 0.0_rp )
336 log_info(
"COPYTOPO_transgrid",*)
'COPYTOPO grid Information'
337 log_info_cont(*)
'size of buffer region (x and y; one side) = ', bufftotx, bufftoty
338 log_info_cont(*)
'# of buffer region (x and y; one side) = ', ibuff, jbuff
339 log_info_cont(*)
'size of transion region (x and y; one side) = ', transtotx, transtoty
340 log_info_cont(*)
'# of transion region (x and y; one side) = ', itrans, jtrans
344 end subroutine copytopo_transgrid
349 subroutine copytopo_setalpha( &
360 real(rp),
intent(in) :: ctrx (
ia)
361 real(rp),
intent(in) :: ctry (
ja)
362 real(rp),
intent(out) :: alpha(
ia,
ja)
364 real(rp) :: coef_x, alpha_x1
365 real(rp) :: coef_y, alpha_y1
372 copytopo_fracx = max( min( copytopo_fracx, 1.0_rp ), eps )
373 copytopo_fracy = max( min( copytopo_fracy, 1.0_rp ), eps )
375 if ( copytopo_taux <= 0.0_rp )
then
378 coef_x = 1.0_rp / copytopo_taux
381 if ( copytopo_tauy <= 0.0_rp )
then
384 coef_y = 1.0_rp / copytopo_tauy
393 if ( ee1 <= 1.0_rp - copytopo_fracx )
then
396 ee1 = ( ee1 - 1.0_rp + copytopo_fracx ) / copytopo_fracx
399 if ( copytopo_linear_h )
then
400 alpha_x1 = coef_x * ee1
402 alpha_x1 = coef_x * ee1 * exp( -(1.0_rp-ee1) * copytopo_exp_h )
407 if ( ee1 <= 1.0_rp - copytopo_fracy )
then
410 ee1 = ( ee1 - 1.0_rp + copytopo_fracy ) / copytopo_fracy
413 if ( copytopo_linear_h )
then
414 alpha_y1 = coef_y * ee1
416 alpha_y1 = coef_y * ee1 * exp( -(1.0_rp-ee1) * copytopo_exp_h )
419 alpha(i,j) = max( alpha_x1, alpha_y1 )
423 call comm_vars8( alpha(:,:), 1 )
424 call comm_wait ( alpha(:,:), 1 )
427 end subroutine copytopo_setalpha
432 subroutine copytopo_input_data( &
455 mapprojection_lonlat2xy
464 real(rp),
intent(out) :: topo_parent(:,:)
466 integer :: ia_org, ja_org
468 real(rp),
allocatable :: lon_org (:,:), lat_org(:,:)
469 real(rp),
allocatable :: x_org (:,:), y_org (:,:)
470 real(rp),
allocatable :: topo_org(:,:)
475 real(rp) :: dummy(1,1)
476 integer :: idx_i(
ia,
ja,nest_interp_level)
477 integer :: idx_j(
ia,
ja,nest_interp_level)
478 real(rp) :: hfact(
ia,
ja,nest_interp_level)
479 integer :: interp_level
480 logical :: zonal, pole
482 real(rp) :: topo_sign(
ia,
ja)
483 real(rp) :: ocean_flag
488 select case ( copytopo_in_filetype )
490 call copytopo_get_size_scale( ia_org, ja_org )
499 log_error(
"COPYTOPO_input_data",*)
'COPYTOPO_IN_FILETYPE must be "SCALE" or "GrADS"'
504 call comm_bcast( ia_org )
505 call comm_bcast( ja_org )
508 allocate( lon_org(ia_org,ja_org) )
509 allocate( lat_org(ia_org,ja_org) )
510 allocate( topo_org(ia_org,ja_org) )
512 select case ( copytopo_in_filetype )
515 lon_org(:,:), lat_org(:,:), topo_org(:,:) )
519 lon_org(:,:), lat_org(:,:), topo_org(:,:) )
523 lon_org(:,:), lat_org(:,:), topo_org(:,:) )
527 call comm_bcast( ia_org, ja_org, lon_org(:,:) )
528 call comm_bcast( ia_org, ja_org, lat_org(:,:) )
529 call comm_bcast( ia_org, ja_org, topo_org(:,:) )
533 lon(:,:), lat(:,:), dummy(:,:), dummy(:,:), &
536 select case ( copytopo_intrp_type )
539 allocate( x_org(ia_org,ja_org) )
540 allocate( y_org(ia_org,ja_org) )
541 call mapprojection_lonlat2xy( ia_org, 1, ia_org, ja_org, 1, ja_org, &
542 lon_org(:,:), lat_org(:,:), &
543 x_org(:,:), y_org(:,:) )
545 zonal = ( maxval(lon_org) - minval(lon_org) ) > 2.0_rp * pi * 0.9_rp
546 pole = ( maxval(lat_org) > pi * 0.5_rp * 0.9_rp ) .or. ( minval(lat_org) < - pi * 0.5_rp * 0.9_rp )
547 call interp_factor2d( ia_org, ja_org,
ia,
ja, &
548 x_org(:,:), y_org(:,:), &
550 idx_i(:,:,:), idx_j(:,:,:), &
552 zonal = zonal, pole = pole )
553 deallocate( x_org, y_org )
556 case (
"DIST-WEIGHT" )
558 call interp_factor2d( nest_interp_level, &
559 ia_org, ja_org,
ia,
ja, &
560 lon_org(:,:), lat_org(:,:), &
561 lon(:,:), lat(:,:), &
562 idx_i(:,:,:), idx_j(:,:,:), &
564 interp_level = nest_interp_level
567 log_error(
"COPYTOPO_input_data",*)
"Unsupported type of COPYTOPO_INTRP_TYPE : ", trim(copytopo_intrp_type)
568 log_error_cont(*)
'It must be "LINEAR" or "DIST-WEIGHT"'
581 deallocate( topo_org, lon_org, lat_org )
583 if ( copytopo_hypdiff_niter > 1 )
then
589 * ( 0.5_rp + sign( 0.5_rp, eps - abs(topo_parent(i,j)) ) )
590 topo_sign(i,j) = sign( 1.0_rp, topo_parent(i,j) ) * ( 1.0_rp - ocean_flag )
595 topo_parent(:,:), copytopo_hypdiff_order, copytopo_hypdiff_niter, &
596 limiter_sign = topo_sign(:,:) )
600 call comm_vars8( topo_parent(:,:), 1 )
601 call comm_wait ( topo_parent(:,:), 1 )
605 end subroutine copytopo_input_data
609 subroutine copytopo_mix_data( &
615 real(rp),
intent(in) :: topo_parent(:,:)
616 real(rp),
intent(in) :: alpha (:,:)
617 real(rp),
intent(inout) :: topo_child (:,:)
622 if ( copytopo_entire_region )
then
626 topo_child(i,j) = topo_parent(i,j)
633 topo_child(i,j) = ( 1.0_rp-alpha(i,j) ) * topo_child(i,j) &
634 + ( alpha(i,j) ) * topo_parent(i,j)
640 end subroutine copytopo_mix_data
643 subroutine copytopo_get_size_scale( &
649 integer,
intent(out) :: ia_org
650 integer,
intent(out) :: ja_org
654 copytopo_in_basename, &
655 copytopo_in_prc_num_x, &
656 copytopo_in_prc_num_y, &
657 copytopo_latlon_catalogue )
664 end subroutine copytopo_get_size_scale
681 integer,
intent(in) :: IA_org, JA_org
682 real(RP),
intent(out) :: LON_org (IA_org,JA_org)
683 real(RP),
intent(out) :: LAT_org (IA_org,JA_org)
684 real(RP),
intent(out) :: TOPO_org(IA_org,JA_org)
687 integer,
allocatable :: tile_id(:)
689 real(RP),
allocatable :: read2D(:,:)
691 integer :: tilei, tilej
692 integer :: cxs, cxe, cys, cye
693 integer :: pxs, pxe, pys, pye
700 num_tile = num_tile )
701 allocate( tile_id(num_tile) )
717 allocate( read2d(tilei,tilej) )
721 aggregate=.false., rankid=rank )
723 call file_read( fid,
"lon", read2d(:,:) )
724 lon_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
725 call file_read( fid,
"lat", read2d(:,:) )
726 lat_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
727 call file_read( fid, copytopo_in_varname, read2d(:,:) )
728 topo_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye)
745 integer,
intent(out) :: IA_org, JA_org
753 call file_grads_get_shape( file_id, copytopo_in_varname, &
776 integer,
intent(in) :: IA_org, JA_org
777 real(RP),
intent(out) :: LON_org (IA_org,JA_org)
778 real(RP),
intent(out) :: LAT_org (IA_org,JA_org)
779 real(RP),
intent(out) :: TOPO_org(IA_org,JA_org)
781 integer :: file_id, var_id
782 real(RP) :: lon1D(IA_org), lat1D(JA_org)
793 call file_grads_read( file_id, var_id, &
798 lon_org(i,j) = lon1d(i) * d2r
802 call file_grads_read( file_id, var_id, &
804 postfix = copytopo_in_postfix )
808 lon_org(i,j) = lon_org(i,j) * d2r
817 call file_grads_read( file_id, var_id, &
822 lat_org(i,j) = lat1d(j) * d2r
826 call file_grads_read( file_id, var_id, &
828 postfix = copytopo_in_postfix )
832 lat_org(i,j) = lat_org(i,j) * d2r
840 call file_grads_read( file_id, var_id, &
842 postfix = copytopo_in_postfix )
860 integer,
intent(out) :: IA_org, JA_org
866 postfix = copytopo_in_postfix, &
886 integer,
intent(in) :: IA_org, JA_org
887 real(RP),
intent(out) :: LON_org (IA_org,JA_org)
888 real(RP),
intent(out) :: LAT_org (IA_org,JA_org)
889 real(RP),
intent(out) :: TOPO_org(IA_org,JA_org)
896 postfix = copytopo_in_postfix, &
899 call file_read( fid,
"XLAT", lat_org(:,:) )
901 call file_read( fid,
"XLONG", lon_org(:,:) )
906 lat_org(i,j) = lat_org(i,j) * d2r
907 lon_org(i,j) = lon_org(i,j) * d2r
911 call file_read( fid,
"HGT", topo_org(:,:) )