SCALE-RM
Functions/Subroutines
mod_copytopo Module Reference

module Copy topography More...

Functions/Subroutines

subroutine, public copytopo (TOPO_child)
 Setup and Main. More...
 
subroutine copytopo_get_data_scale (IA_org, JA_org, LON_org, LAT_org, TOPO_org)
 
subroutine copytopo_get_size_grads (IA_org, JA_org)
 
subroutine copytopo_get_data_grads (IA_org, JA_org, LON_org, LAT_org, TOPO_org)
 
subroutine copytopo_get_size_wrfarw (IA_org, JA_org)
 
subroutine copytopo_get_data_wrfarw (IA_org, JA_org, LON_org, LAT_org, TOPO_org)
 

Detailed Description

module Copy topography

Description
subroutines for preparing topography data (copy from parent domain)
Author
Team SCALE
NAMELIST
  • PARAM_COPYTOPO
    nametypedefault valuecomment
    COPYTOPO_IN_FILETYPE character(len=H_SHORT) '' !< 'SCALE', 'GrADS', or 'WRF-ARW'
    COPYTOPO_IN_BASENAME character(len=H_LONG) ''
    COPYTOPO_IN_PRC_NUM_X integer 0 for old files
    COPYTOPO_IN_PRC_NUM_Y integer 0 for old files
    COPYTOPO_LATLON_CATALOGUE character(len=H_LONG) '' for old files
    COPYTOPO_IN_POSTFIX character(len=H_LONG) ''
    COPYTOPO_IN_VARNAME character(len=H_SHORT) 'topo'
    COPYTOPO_IN_LONNAME character(len=H_SHORT) 'lon'
    COPYTOPO_IN_LATNAME character(len=H_SHORT) 'lat'
    COPYTOPO_TRANSITION_DX real(RP) -1.0_RP thickness of transition region [m]: x
    COPYTOPO_TRANSITION_DY real(RP) -1.0_RP thickness of transition region [m]: y
    COPYTOPO_FRACX real(RP) 1.0_RP fraction of transition region (x) (0-1)
    COPYTOPO_FRACY real(RP) 1.0_RP fraction of transition region (y) (0-1)
    COPYTOPO_TAUX real(RP) 1.0_RP maximum value for mixing tau (x) [s]
    COPYTOPO_TAUY real(RP) 1.0_RP maximum value for mixing tau (y) [s]
    COPYTOPO_ENTIRE_REGION logical .false. copy parent topo over an entire region
    COPYTOPO_LINEAR_H logical .true. linear or non-linear profile of relax region
    COPYTOPO_EXP_H real(RP) 2.0_RP factor of non-linear profile of relax region
    COPYTOPO_HYPDIFF_ORDER integer 4
    COPYTOPO_HYPDIFF_NITER integer 20

History Output
No history output

Function/Subroutine Documentation

◆ copytopo()

subroutine, public mod_copytopo::copytopo ( real(rp), dimension(:,:), intent(inout)  TOPO_child)

Setup and Main.

Parameters
[in,out]topo_childtopography of child domain

Definition at line 79 of file mod_copytopo.F90.

79  use scale_prc, only: &
80  prc_abort
81  implicit none
82 
83  real(RP), intent(inout) :: TOPO_child(:,:)
84 
85  real(RP) :: CTRX (IA)
86  real(RP) :: CTRY (JA)
87  real(RP) :: alpha (IA,JA)
88  real(RP) :: TOPO_parent(IA,JA)
89 
90  namelist / param_copytopo / &
91  copytopo_in_filetype, &
92  copytopo_in_basename, &
93  copytopo_in_prc_num_x, & ! for old files
94  copytopo_in_prc_num_y, & ! for old files
95  copytopo_latlon_catalogue, & ! for old files
96  copytopo_in_postfix, &
97  copytopo_in_varname, &
98  copytopo_in_lonname, &
99  copytopo_in_latname, &
100  copytopo_transition_dx, &
101  copytopo_transition_dy, &
102  copytopo_fracx, &
103  copytopo_fracy, &
104  copytopo_taux, &
105  copytopo_tauy, &
106  copytopo_entire_region, &
107  copytopo_linear_h, &
108  copytopo_exp_h, &
109  copytopo_hypdiff_order, &
110  copytopo_hypdiff_niter
111 
112  integer :: ierr
113  !---------------------------------------------------------------------------
114 
115  log_newline
116  log_info("COPYTOPO",*) 'Setup'
117 
118  !--- read namelist
119  rewind(io_fid_conf)
120  read(io_fid_conf,nml=param_copytopo,iostat=ierr)
121  if( ierr < 0 ) then !--- missing
122  log_info("COPYTOPO",*) 'Not found namelist. Default used.'
123  elseif( ierr > 0 ) then !--- fatal error
124  log_error("COPYTOPO",*) 'Not appropriate names in namelist PARAM_COPYTOPO. Check!'
125  call prc_abort
126  endif
127  log_nml(param_copytopo)
128 
129  ! copy topography from parent domain to transition region
130 
131  call copytopo_transgrid( ctrx(:), ctry(:) ) ! [OUT]
132 
133  call copytopo_setalpha ( ctrx(:), ctry(:), & ! [IN]
134  alpha(:,:) ) ! [OUT]
135 
136  call copytopo_input_data( topo_parent(:,:) ) ! [OUT]
137 
138  call copytopo_mix_data( topo_parent(:,:), & ! [IN]
139  alpha(:,:), & ! [IN]
140  topo_child(:,:) ) ! [INOUT]
141 
142  return

References scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfxg, scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfyg, scale_atmos_grid_cartesc::atmos_grid_cartesc_cx, scale_atmos_grid_cartesc::atmos_grid_cartesc_cxg, scale_atmos_grid_cartesc::atmos_grid_cartesc_cy, scale_atmos_grid_cartesc::atmos_grid_cartesc_cyg, scale_atmos_grid_cartesc::atmos_grid_cartesc_fxg, scale_atmos_grid_cartesc::atmos_grid_cartesc_fyg, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon, scale_comm_cartesc_nest::comm_cartesc_nest_domain_regist_file(), scale_comm_cartesc_nest::comm_cartesc_nest_interp_level, scale_comm_cartesc_nest::comm_cartesc_nest_parent_info(), scale_const::const_eps, scale_const::const_pi, copytopo_get_data_grads(), copytopo_get_data_scale(), copytopo_get_data_wrfarw(), copytopo_get_size_grads(), copytopo_get_size_wrfarw(), scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::iag, scale_atmos_grid_cartesc_index::ieb, scale_atmos_grid_cartesc_index::ihalo, scale_atmos_grid_cartesc_index::imax, scale_atmos_grid_cartesc_index::imaxg, scale_interp::interp_domain_compatibility(), scale_interp::interp_interp2d(), scale_io::io_fid_conf, scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::jag, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::jhalo, scale_atmos_grid_cartesc_index::jmax, scale_atmos_grid_cartesc_index::jmaxg, scale_atmos_grid_cartesc_index::jsb, scale_landuse::landuse_fact_ocean, scale_prc_cartesc::prc_2drank, scale_prc::prc_abort(), scale_prc::prc_ismaster, and scale_prc::prc_myrank.

Referenced by mod_cnvtopo::cnvtopo().

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

◆ copytopo_get_data_scale()

subroutine mod_copytopo::copytopo_get_data_scale ( integer, intent(in)  IA_org,
integer, intent(in)  JA_org,
real(rp), dimension (ia_org,ja_org), intent(out)  LON_org,
real(rp), dimension (ia_org,ja_org), intent(out)  LAT_org,
real(rp), dimension(ia_org,ja_org), intent(out)  TOPO_org 
)

Definition at line 671 of file mod_copytopo.F90.

671  use scale_const, only: &
672  d2r => const_d2r
673  use scale_comm_cartesc_nest, only: &
676  use scale_file, only: &
677  file_open, &
678  file_read, &
679  file_close
680  implicit none
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)
685 
686  integer :: num_tile
687  integer, allocatable :: tile_id(:)
688 
689  real(RP), allocatable :: read2D(:,:)
690 
691  integer :: tilei, tilej
692  integer :: cxs, cxe, cys, cye ! for child domain
693  integer :: pxs, pxe, pys, pye ! for parent domain
694  integer :: fid
695 
696  integer :: rank
697  integer :: n
698 
699  call comm_cartesc_nest_parent_info( copytopo_domid, & ! [IN]
700  num_tile = num_tile )
701  allocate( tile_id(num_tile) )
702  call comm_cartesc_nest_parent_info( copytopo_domid, & ! [IN]
703  tile_id = tile_id )
704 
705  do n = 1, num_tile
706  ! read data from split files
707  rank = tile_id(n)
708 
709  call comm_cartesc_nest_domain_shape( tilei, tilej, & ! [OUT]
710  cxs, cxe, & ! [OUT]
711  cys, cye, & ! [OUT]
712  pxs, pxe, & ! [OUT]
713  pys, pye, & ! [OUT]
714  copytopo_domid, & ! [IN]
715  n ) ! [IN]
716 
717  allocate( read2d(tilei,tilej) )
718 
719  call file_open( copytopo_in_basename, & ! [IN]
720  fid, & ! [OUT]
721  aggregate=.false., rankid=rank ) ! [IN]
722 
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)
729  deallocate( read2d )
730 
731  call file_close( fid )
732 
733  enddo
734 
735  return

References scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape(), scale_comm_cartesc_nest::comm_cartesc_nest_parent_info(), scale_const::const_d2r, scale_file::file_close(), and scale_file::file_open().

Referenced by copytopo().

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

◆ copytopo_get_size_grads()

subroutine mod_copytopo::copytopo_get_size_grads ( integer, intent(out)  IA_org,
integer, intent(out)  JA_org 
)

Definition at line 741 of file mod_copytopo.F90.

741  use scale_file_grads, only: &
742  file_grads_open, &
743  file_grads_get_shape
744  implicit none
745  integer, intent(out) :: IA_org, JA_org
746 
747  integer :: file_id
748  integer :: shape(2)
749 
750  call file_grads_open( copytopo_in_basename, & ! (in)
751  file_id ) ! (out)
752 
753  call file_grads_get_shape( file_id, copytopo_in_varname, & ! (in)
754  shape(:) ) ! (out)
755 
756  ia_org = shape(1)
757  ja_org = shape(2)
758 
759  return

References scale_file_grads::file_grads_open().

Referenced by copytopo().

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

◆ copytopo_get_data_grads()

subroutine mod_copytopo::copytopo_get_data_grads ( integer, intent(in)  IA_org,
integer, intent(in)  JA_org,
real(rp), dimension (ia_org,ja_org), intent(out)  LON_org,
real(rp), dimension (ia_org,ja_org), intent(out)  LAT_org,
real(rp), dimension(ia_org,ja_org), intent(out)  TOPO_org 
)

Definition at line 767 of file mod_copytopo.F90.

767  use scale_const, only: &
768  d2r => const_d2r
769  use scale_file_grads, only: &
770  file_grads_open, &
773  file_grads_read, &
775  implicit none
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)
780 
781  integer :: file_id, var_id
782  real(RP) :: lon1D(IA_org), lat1D(JA_org)
783 
784  integer :: i, j
785 
786  call file_grads_open( copytopo_in_basename, & ! (in)
787  file_id ) ! (out)
788 
789  ! lon
790  call file_grads_varid( file_id, copytopo_in_lonname, & ! (in)
791  var_id ) ! (out)
792  if ( file_grads_isoned( file_id, var_id ) ) then
793  call file_grads_read( file_id, var_id, & ! (in)
794  lon1d(:) ) ! (out)
795  !$omp parallel do collapse(2)
796  do j = 1, ja_org
797  do i = 1, ia_org
798  lon_org(i,j) = lon1d(i) * d2r
799  end do
800  end do
801  else
802  call file_grads_read( file_id, var_id, & ! (in)
803  lon_org(:,:), & ! (out)
804  postfix = copytopo_in_postfix ) ! (in)
805  !$omp parallel do collapse(2)
806  do j = 1, ja_org
807  do i = 1, ia_org
808  lon_org(i,j) = lon_org(i,j) * d2r
809  end do
810  end do
811  end if
812 
813  ! lat
814  call file_grads_varid( file_id, copytopo_in_latname, & ! (in)
815  var_id ) ! (out)
816  if ( file_grads_isoned( file_id, var_id ) ) then
817  call file_grads_read( file_id, var_id, & ! (in)
818  lat1d(:) ) ! (out)
819  !$omp parallel do collapse(2)
820  do j = 1, ja_org
821  do i = 1, ia_org
822  lat_org(i,j) = lat1d(j) * d2r
823  end do
824  end do
825  else
826  call file_grads_read( file_id, var_id, & ! (in)
827  lat_org(:,:), & ! (out)
828  postfix = copytopo_in_postfix ) ! (in)
829  !$omp parallel do collapse(2)
830  do j = 1, ja_org
831  do i = 1, ia_org
832  lat_org(i,j) = lat_org(i,j) * d2r
833  end do
834  end do
835  end if
836 
837  ! topo
838  call file_grads_varid( file_id, copytopo_in_varname, & ! (in)
839  var_id ) ! (out)
840  call file_grads_read( file_id, var_id, & ! (in)
841  topo_org(:,:), & ! (out)
842  postfix = copytopo_in_postfix ) ! (in)
843 
844  ! close
845  call file_grads_close( file_id )
846 
847  return

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

Referenced by copytopo().

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

◆ copytopo_get_size_wrfarw()

subroutine mod_copytopo::copytopo_get_size_wrfarw ( integer, intent(out)  IA_org,
integer, intent(out)  JA_org 
)

Definition at line 853 of file mod_copytopo.F90.

853  use scale_file_grads, only: &
854  file_grads_open, &
855  file_grads_get_shape
856  use scale_file, only: &
857  file_open, &
859  implicit none
860  integer, intent(out) :: IA_org, JA_org
861 
862  integer :: fid
863 
864  call file_open( copytopo_in_basename, & ! (in)
865  fid, & ! (out)
866  postfix = copytopo_in_postfix, & ! (in)
867  single = .true. ) ! (in)
868 
869  call file_get_dimlength( fid, "west_east", ia_org )
870  call file_get_dimlength( fid, "south_north", ja_org )
871 
872  return

References scale_file::file_get_dimlength(), scale_file_grads::file_grads_open(), and scale_file::file_open().

Referenced by copytopo().

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

◆ copytopo_get_data_wrfarw()

subroutine mod_copytopo::copytopo_get_data_wrfarw ( integer, intent(in)  IA_org,
integer, intent(in)  JA_org,
real(rp), dimension (ia_org,ja_org), intent(out)  LON_org,
real(rp), dimension (ia_org,ja_org), intent(out)  LAT_org,
real(rp), dimension(ia_org,ja_org), intent(out)  TOPO_org 
)

Definition at line 879 of file mod_copytopo.F90.

879  use scale_const, only: &
880  d2r => const_d2r
881  use scale_file, only: &
882  file_open, &
883  file_read, &
884  file_close
885  implicit none
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)
890 
891  integer :: fid
892  integer :: i, j
893 
894  call file_open( copytopo_in_basename, & ! (in)
895  fid, & ! (out)
896  postfix = copytopo_in_postfix, & ! (in)
897  single = .true. ) ! (in)
898 
899  call file_read( fid, "XLAT", lat_org(:,:) )
900 
901  call file_read( fid, "XLONG", lon_org(:,:) )
902 
903  !$omp parallel do collapse(2)
904  do j = 1, ja_org
905  do i = 1, ia_org
906  lat_org(i,j) = lat_org(i,j) * d2r
907  lon_org(i,j) = lon_org(i,j) * d2r
908  end do
909  end do
910 
911  call file_read( fid, "HGT", topo_org(:,:) )
912 
913  call file_close( fid )
914 

References scale_const::const_d2r, scale_file::file_close(), and scale_file::file_open().

Referenced by copytopo().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape
subroutine, public comm_cartesc_nest_domain_shape(tilei, tilej, cxs, cxe, cys, cye, pxs, pxe, pys, pye, dom_id, iloc, xstg, ystg)
Return shape of ParentDomain at the specified rank (for offline)
Definition: scale_comm_cartesC_nest.F90:1187
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_file::file_open
subroutine, public file_open(basename, fid, mode, single, allnodes, aggregate, rankid, postfix)
Definition: scale_file.F90:536
scale_comm_cartesc_nest
module Communication CartesianC nesting
Definition: scale_comm_cartesC_nest.F90:12
scale_file_grads::file_grads_varid
subroutine, public file_grads_varid(file_id, var_name, var_id)
Definition: scale_file_grads.F90:312
scale_file::file_get_dimlength
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:633
scale_comm_cartesc_nest::comm_cartesc_nest_parent_info
subroutine, public comm_cartesc_nest_parent_info(dom_id, KMAX, LKMAX, IMAXG, JMAXG, num_tile, tile_id)
Return infomation of parent domain (for offline)
Definition: scale_comm_cartesC_nest.F90:1130
scale_file
module file
Definition: scale_file.F90:15
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_file::file_close
subroutine, public file_close(fid, abort)
Definition: scale_file.F90:6234
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_file_grads
module file_grads
Definition: scale_file_grads.F90:11
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_file_grads::file_grads_close
subroutine, public file_grads_close(file_id)
Definition: scale_file_grads.F90:713