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_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 74 of file mod_copytopo.F90.

74  use scale_prc, only: &
75  prc_abort
76  implicit none
77 
78  real(RP), intent(inout) :: TOPO_child(:,:)
79 
80  real(RP) :: CTRX (IA)
81  real(RP) :: CTRY (JA)
82  real(RP) :: alpha (IA,JA)
83  real(RP) :: TOPO_parent(IA,JA)
84 
85  namelist / param_copytopo / &
86  copytopo_in_filetype, &
87  copytopo_in_basename, &
88  copytopo_in_postfix, &
89  copytopo_in_varname, &
90  copytopo_in_lonname, &
91  copytopo_in_latname, &
92  copytopo_transition_dx, &
93  copytopo_transition_dy, &
94  copytopo_fracx, &
95  copytopo_fracy, &
96  copytopo_taux, &
97  copytopo_tauy, &
98  copytopo_entire_region, &
99  copytopo_linear_h, &
100  copytopo_exp_h, &
101  copytopo_hypdiff_order, &
102  copytopo_hypdiff_niter
103 
104  integer :: ierr
105  !---------------------------------------------------------------------------
106 
107  log_newline
108  log_info("COPYTOPO",*) 'Setup'
109 
110  !--- read namelist
111  rewind(io_fid_conf)
112  read(io_fid_conf,nml=param_copytopo,iostat=ierr)
113  if( ierr < 0 ) then !--- missing
114  log_info("COPYTOPO",*) 'Not found namelist. Default used.'
115  elseif( ierr > 0 ) then !--- fatal error
116  log_error("COPYTOPO",*) 'Not appropriate names in namelist PARAM_COPYTOPO. Check!'
117  call prc_abort
118  endif
119  log_nml(param_copytopo)
120 
121  ! copy topography from parent domain to transition region
122 
123  call copytopo_transgrid( ctrx(:), ctry(:) ) ! [OUT]
124 
125  call copytopo_setalpha ( ctrx(:), ctry(:), & ! [IN]
126  alpha(:,:) ) ! [OUT]
127 
128  call copytopo_input_data( topo_parent(:,:) ) ! [OUT]
129 
130  call copytopo_mix_data( topo_parent(:,:), & ! [IN]
131  alpha(:,:), & ! [IN]
132  topo_child(:,:) ) ! [INOUT]
133 
134  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_interp_level, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_y, 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_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::jsb, scale_landuse::landuse_fact_ocean, scale_comm_cartesc_nest::parent_imax, scale_comm_cartesc_nest::parent_jmax, 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 657 of file mod_copytopo.F90.

657  use scale_const, only: &
658  d2r => const_d2r
659  use scale_comm_cartesc_nest, only: &
660  nest_domain_shape => comm_cartesc_nest_domain_shape, &
661  nest_tile_id => comm_cartesc_nest_tile_id
662  use scale_file, only: &
663  file_open, &
664  file_read, &
665  file_close
666  implicit none
667  integer, intent(in) :: IA_org, JA_org
668  real(RP), intent(out) :: LON_org (IA_org,JA_org)
669  real(RP), intent(out) :: LAT_org (IA_org,JA_org)
670  real(RP), intent(out) :: TOPO_org(IA_org,JA_org)
671 
672  real(RP), allocatable :: read2D(:,:)
673 
674  integer :: tilei, tilej
675  integer :: cxs, cxe, cys, cye ! for child domain
676  integer :: pxs, pxe, pys, pye ! for parent domain
677  integer :: fid
678 
679  integer :: rank
680  integer :: n
681 
682  do n = 1, size( nest_tile_id(:) )
683  ! read data from split files
684  rank = nest_tile_id(n)
685 
686  call nest_domain_shape( tilei, tilej, & ! [OUT]
687  cxs, cxe, & ! [OUT]
688  cys, cye, & ! [OUT]
689  pxs, pxe, & ! [OUT]
690  pys, pye, & ! [OUT]
691  n ) ! [IN]
692 
693  allocate( read2d(tilei,tilej) )
694 
695  call file_open( copytopo_in_basename, & ! [IN]
696  fid, & ! [OUT]
697  aggregate=.false., rankid=rank ) ! [IN]
698 
699  call file_read( fid, "lon", read2d(:,:) )
700  lon_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
701  call file_read( fid, "lat", read2d(:,:) )
702  lat_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
703  call file_read( fid, copytopo_in_varname, read2d(:,:) )
704  topo_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye)
705  deallocate( read2d )
706 
707  call file_close( fid )
708 
709  enddo
710 
711  return

References scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape(), scale_comm_cartesc_nest::comm_cartesc_nest_tile_id, 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 717 of file mod_copytopo.F90.

717  use scale_file_grads, only: &
718  file_grads_open, &
719  file_grads_get_shape
720  implicit none
721  integer, intent(out) :: IA_org, JA_org
722 
723  integer :: file_id
724  integer :: shape(2)
725 
726  call file_grads_open( copytopo_in_basename, & ! (in)
727  file_id ) ! (out)
728 
729  call file_grads_get_shape( file_id, copytopo_in_varname, & ! (in)
730  shape(:) ) ! (out)
731 
732  ia_org = shape(1)
733  ja_org = shape(2)
734 
735  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 743 of file mod_copytopo.F90.

743  use scale_const, only: &
744  d2r => const_d2r
745  use scale_file_grads, only: &
746  file_grads_open, &
749  file_grads_read, &
751  implicit none
752  integer, intent(in) :: IA_org, JA_org
753  real(RP), intent(out) :: LON_org (IA_org,JA_org)
754  real(RP), intent(out) :: LAT_org (IA_org,JA_org)
755  real(RP), intent(out) :: TOPO_org(IA_org,JA_org)
756 
757  integer :: file_id, var_id
758  real(RP) :: lon1D(IA_org), lat1D(JA_org)
759 
760  integer :: i, j
761 
762  call file_grads_open( copytopo_in_basename, & ! (in)
763  file_id ) ! (out)
764 
765  ! lon
766  call file_grads_varid( file_id, copytopo_in_lonname, & ! (in)
767  var_id ) ! (out)
768  if ( file_grads_isoned( file_id, var_id ) ) then
769  call file_grads_read( file_id, var_id, & ! (in)
770  lon1d(:) ) ! (out)
771  !$omp parallel do collapse(2)
772  do j = 1, ja_org
773  do i = 1, ia_org
774  lon_org(i,j) = lon1d(i) * d2r
775  end do
776  end do
777  else
778  call file_grads_read( file_id, var_id, & ! (in)
779  lon_org(:,:), & ! (out)
780  postfix = copytopo_in_postfix ) ! (in)
781  !$omp parallel do collapse(2)
782  do j = 1, ja_org
783  do i = 1, ia_org
784  lon_org(i,j) = lon_org(i,j) * d2r
785  end do
786  end do
787  end if
788 
789  ! lat
790  call file_grads_varid( file_id, copytopo_in_latname, & ! (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  lat1d(:) ) ! (out)
795  !$omp parallel do collapse(2)
796  do j = 1, ja_org
797  do i = 1, ia_org
798  lat_org(i,j) = lat1d(j) * d2r
799  end do
800  end do
801  else
802  call file_grads_read( file_id, var_id, & ! (in)
803  lat_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  lat_org(i,j) = lat_org(i,j) * d2r
809  end do
810  end do
811  end if
812 
813  ! topo
814  call file_grads_varid( file_id, copytopo_in_varname, & ! (in)
815  var_id ) ! (out)
816  call file_grads_read( file_id, var_id, & ! (in)
817  topo_org(:,:), & ! (out)
818  postfix = copytopo_in_postfix ) ! (in)
819 
820  ! close
821  call file_grads_close( file_id )
822 
823  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 829 of file mod_copytopo.F90.

829  use scale_file_grads, only: &
830  file_grads_open, &
831  file_grads_get_shape
832  use scale_file, only: &
833  file_open, &
835  implicit none
836  integer, intent(out) :: IA_org, JA_org
837 
838  integer :: fid
839 
840  call file_open( copytopo_in_basename, & ! (in)
841  fid, & ! (out)
842  postfix = copytopo_in_postfix, & ! (in)
843  single = .true. ) ! (in)
844 
845  call file_get_dimlength( fid, "west_east", ia_org )
846  call file_get_dimlength( fid, "south_north", ja_org )
847 
848  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 855 of file mod_copytopo.F90.

855  use scale_const, only: &
856  d2r => const_d2r
857  use scale_file, only: &
858  file_open, &
859  file_read, &
860  file_close
861  implicit none
862  integer, intent(in) :: IA_org, JA_org
863  real(RP), intent(out) :: LON_org (IA_org,JA_org)
864  real(RP), intent(out) :: LAT_org (IA_org,JA_org)
865  real(RP), intent(out) :: TOPO_org(IA_org,JA_org)
866 
867  integer :: fid
868  integer :: i, j
869 
870  call file_open( copytopo_in_basename, & ! (in)
871  fid, & ! (out)
872  postfix = copytopo_in_postfix, & ! (in)
873  single = .true. ) ! (in)
874 
875  call file_read( fid, "XLAT", lat_org(:,:) )
876 
877  call file_read( fid, "XLONG", lon_org(:,:) )
878 
879  !$omp parallel do collapse(2)
880  do j = 1, ja_org
881  do i = 1, ia_org
882  lat_org(i,j) = lat_org(i,j) * d2r
883  lon_org(i,j) = lon_org(i,j) * d2r
884  end do
885  end do
886 
887  call file_read( fid, "HGT", topo_org(:,:) )
888 
889  call file_close( fid )
890 

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_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_comm_cartesc_nest::comm_cartesc_nest_tile_id
integer, dimension(:), allocatable, public comm_cartesc_nest_tile_id
parent tile real id
Definition: scale_comm_cartesC_nest.F90:52
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:300
scale_file::file_get_dimlength
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:565
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:4738
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_file::file_open
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:487
scale_file_grads
module file_grads
Definition: scale_file_grads.F90:11
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, iloc)
Return shape of ParentDomain at the specified rank (for offline)
Definition: scale_comm_cartesC_nest.F90:1061
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
scale_file_grads::file_grads_close
subroutine, public file_grads_close(file_id)
Definition: scale_file_grads.F90:624