SCALE-RM
Functions/Subroutines | Variables
mod_cnvlanduse Module Reference

module Convert LandUseIndex More...

Functions/Subroutines

subroutine, public cnvlanduse_setup
 Setup. More...
 
subroutine, public cnvlanduse
 Driver. More...
 

Variables

logical, public cnvlanduse_donothing
 
logical, public cnvlanduse_useglccv2 = .false.
 
logical, public cnvlanduse_uselu100m = .false.
 
logical, public cnvlanduse_usejibis = .false.
 

Detailed Description

module Convert LandUseIndex

Description
subroutines for preparing landuse index data (convert from external file)
Author
Team SCALE
NAMELIST
  • PARAM_CNVLANDUSE
    nametypedefault valuecomment
    CNVLANDUSE_NAME character(len=H_SHORT) 'NONE' keep backward compatibility
    CNVLANDUSE_USEGLCCV2 logical .false.
    CNVLANDUSE_USELU100M logical .false.
    CNVLANDUSE_USEJIBIS logical .false.
    CNVLANDUSE_UNITTILE_DDEG real(RP) 0.0_RP dx for unit tile [deg]
    CNVLANDUSE_OVERSAMPLING_FACTOR real(RP) 2.0_RP factor of min. dx against the unit tile

  • PARAM_CNVLANDUSE_GLCCv2
    nametypedefault valuecomment
    GLCCV2_IN_CATALOGUE character(len=H_LONG) '' metadata files for GLCCv2
    GLCCV2_IN_DIR character(len=H_LONG) '' directory contains GLCCv2 files (GrADS format)
    LIMIT_URBAN_FRACTION real(RP) 1.0_RP fraction limiter for urban area

  • PARAM_CNVLANDUSE_LU100M
    nametypedefault valuecomment
    LU100M_IN_CATALOGUE character(len=H_LONG) '' metadata files for LU100M
    LU100M_IN_DIR character(len=H_LONG) '' directory contains LU100M files (GrADS format)
    LIMIT_URBAN_FRACTION real(RP) 1.0_RP fraction limiter for urban area

History Output
No history output

Function/Subroutine Documentation

◆ cnvlanduse_setup()

subroutine, public mod_cnvlanduse::cnvlanduse_setup ( )

Setup.

Definition at line 59 of file mod_cnvlanduse.f90.

References cnvlanduse_donothing, cnvlanduse_useglccv2, cnvlanduse_usejibis, cnvlanduse_uselu100m, scale_const::const_d2r, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_process::prc_mpistop(), scale_grid_real::real_dlat, and scale_grid_real::real_dlon.

Referenced by mod_convert::convert_setup().

59  use scale_process, only: &
61  use scale_const, only: &
62  d2r => const_d2r
63  use scale_comm, only: &
64  comm_horizontal_min
65  use scale_grid_real, only: &
66  real_dlat, &
67  real_dlon
68  implicit none
69 
70  character(len=H_SHORT) :: CNVLANDUSE_name = 'NONE' ! keep backward compatibility
71 
72  namelist / param_cnvlanduse / &
73  cnvlanduse_name, &
74  cnvlanduse_useglccv2, &
75  cnvlanduse_uselu100m, &
76  cnvlanduse_usejibis, &
77  cnvlanduse_unittile_ddeg, &
78  cnvlanduse_oversampling_factor
79 
80  real(RP) :: drad(IA,JA)
81  real(RP) :: drad_min
82 
83  integer :: ierr
84  !---------------------------------------------------------------------------
85 
86  if( io_l ) write(io_fid_log,*)
87  if( io_l ) write(io_fid_log,*) '++++++ Module[convert landuseindex] / Categ[preprocess] / Origin[SCALE-RM]'
88  !--- read namelist
89  rewind(io_fid_conf)
90  read(io_fid_conf,nml=param_cnvlanduse,iostat=ierr)
91  if( ierr < 0 ) then !--- missing
92  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
93  elseif( ierr > 0 ) then !--- fatal error
94  write(*,*) 'xxx Not appropriate names in namelist PARAM_CNVLANDUSE. Check!'
95  call prc_mpistop
96  endif
97  if( io_nml ) write(io_fid_nml,nml=param_cnvlanduse)
98 
99  select case(cnvlanduse_name)
100  case('NONE')
101  ! do nothing
102  case('GLCCv2')
103  cnvlanduse_useglccv2 = .true.
104  cnvlanduse_uselu100m = .false.
105  cnvlanduse_usejibis = .false.
106  case('LU100M')
107  cnvlanduse_useglccv2 = .false.
108  cnvlanduse_uselu100m = .true.
109  cnvlanduse_usejibis = .false.
110  case('COMBINE')
111  cnvlanduse_useglccv2 = .true.
112  cnvlanduse_uselu100m = .true.
113  cnvlanduse_usejibis = .true.
114  case default
115  write(*,*) 'xxx Unsupported TYPE:', trim(cnvlanduse_name)
116  call prc_mpistop
117  endselect
118 
119  cnvlanduse_donothing = .true.
120 
121  if ( cnvlanduse_useglccv2 ) then
122  cnvlanduse_donothing = .false.
123  if( io_l ) write(io_fid_log,*) '*** Use GLCC ver.2, global 30 arcsec. data'
124  if ( cnvlanduse_uselu100m ) then
125  if( io_l ) write(io_fid_log,*) '*** Use KSJ landuse 100m data for Japan region'
126  if( io_l ) write(io_fid_log,*) '*** Overwrite Japan region'
127  if ( cnvlanduse_usejibis ) then
128  if( io_l ) write(io_fid_log,*) '*** Use J-IBIS map 100m data for Japan region'
129  if( io_l ) write(io_fid_log,*) '*** Overwrite Japan region (PFT only)'
130  endif
131  endif
132  elseif ( cnvlanduse_uselu100m ) then
133  cnvlanduse_donothing = .false.
134  if( io_l ) write(io_fid_log,*) '*** Use KSJ landuse 100m data, Japan region only'
135  if ( cnvlanduse_usejibis ) then
136  if( io_l ) write(io_fid_log,*) '*** Use J-IBIS map 100m data for Japan region'
137  if( io_l ) write(io_fid_log,*) '*** Overwrite Japan region (PFT only)'
138  endif
139  endif
140 
141  if ( cnvlanduse_donothing ) then
142  if( io_l ) write(io_fid_log,*) '*** Do nothing for landuse index'
143  else
144  drad(:,:) = min( real_dlat(:,:), real_dlon(:,:) )
145  call comm_horizontal_min( drad_min, drad(:,:) )
146 
147  if ( cnvlanduse_unittile_ddeg > 0.0_rp ) then
148  cnvlanduse_oversampling_factor = ( drad_min / d2r ) / cnvlanduse_unittile_ddeg
149  endif
150  cnvlanduse_oversampling_factor = max( 1.0_rp, cnvlanduse_oversampling_factor )
151  cnvlanduse_unittile_ddeg = ( drad_min / d2r ) / cnvlanduse_oversampling_factor
152 
153  if( io_l ) write(io_fid_log,*) '*** The size of tile [deg] = ', cnvlanduse_unittile_ddeg
154  if( io_l ) write(io_fid_log,*) '*** oversampling factor = ', cnvlanduse_oversampling_factor
155  endif
156 
157  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
module GRID (real space)
real(rp), dimension(:,:), allocatable, public real_dlon
delta longitude
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
real(rp), dimension(:,:), allocatable, public real_dlat
delta latitude
module CONSTANT
Definition: scale_const.F90:14
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cnvlanduse()

subroutine, public mod_cnvlanduse::cnvlanduse ( )

Driver.

Definition at line 163 of file mod_cnvlanduse.f90.

References cnvlanduse_donothing, cnvlanduse_useglccv2, cnvlanduse_usejibis, cnvlanduse_uselu100m, scale_const::const_d2r, scale_const::const_eps, scale_const::const_pi, scale_const::const_radius, scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_nml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_landuse::landuse_calc_fact(), scale_landuse::landuse_frac_lake, scale_landuse::landuse_frac_land, scale_landuse::landuse_frac_pft, scale_landuse::landuse_frac_urban, scale_landuse::landuse_index_pft, scale_landuse::landuse_pft_mosaic, scale_landuse::landuse_pft_nmax, scale_landuse::landuse_write(), scale_process::prc_mpistop(), scale_grid_real::real_laty, and scale_grid_real::real_lonx.

Referenced by mod_convert::convert().

163  use scale_process, only: &
165  use scale_landuse, only: &
168  implicit none
169  !---------------------------------------------------------------------------
170 
171  if ( cnvlanduse_donothing ) then
172  if( io_l ) write(io_fid_log,*)
173  if( io_l ) write(io_fid_log,*) '++++++ SKIP CONVERT LANDUSE DATA ++++++'
174  else
175  if( io_l ) write(io_fid_log,*)
176  if( io_l ) write(io_fid_log,*) '++++++ START CONVERT LANDUSE DATA ++++++'
177 
178  if ( cnvlanduse_useglccv2 ) then
179  call cnvlanduse_glccv2
180  endif
181 
182  if ( cnvlanduse_uselu100m ) then
183  call cnvlanduse_lu100m
184  endif
185 
186  if ( cnvlanduse_usejibis ) then
187  call cnvlanduse_jibis
188  endif
189 
190  ! calculate landuse factors
191  call landuse_calc_fact
192 
193  if( io_l ) write(io_fid_log,*) '++++++ END CONVERT LANDUSE DATA ++++++'
194 
195  ! output landuse file
196  call landuse_write
197  endif
198 
199  return
subroutine, public prc_mpistop
Abort MPI.
subroutine, public landuse_calc_fact
module LANDUSE
subroutine, public landuse_write
Write landuse data.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ cnvlanduse_donothing

logical, public mod_cnvlanduse::cnvlanduse_donothing

Definition at line 34 of file mod_cnvlanduse.f90.

Referenced by cnvlanduse(), and cnvlanduse_setup().

34  logical, public :: CNVLANDUSE_DoNothing

◆ cnvlanduse_useglccv2

logical, public mod_cnvlanduse::cnvlanduse_useglccv2 = .false.

Definition at line 35 of file mod_cnvlanduse.f90.

Referenced by cnvlanduse(), and cnvlanduse_setup().

35  logical, public :: CNVLANDUSE_UseGLCCv2 = .false.

◆ cnvlanduse_uselu100m

logical, public mod_cnvlanduse::cnvlanduse_uselu100m = .false.

Definition at line 36 of file mod_cnvlanduse.f90.

Referenced by cnvlanduse(), and cnvlanduse_setup().

36  logical, public :: CNVLANDUSE_UseLU100M = .false.

◆ cnvlanduse_usejibis

logical, public mod_cnvlanduse::cnvlanduse_usejibis = .false.

Definition at line 37 of file mod_cnvlanduse.f90.

Referenced by cnvlanduse(), and cnvlanduse_setup().

37  logical, public :: CNVLANDUSE_UseJIBIS = .false.