SCALE-RM
scale_landuse.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
18  use scale_prof
20  !-----------------------------------------------------------------------------
21  implicit none
22  private
23  !-----------------------------------------------------------------------------
24  !
25  !++ Public procedure
26  !
27  public :: landuse_setup
28  public :: landuse_calc_fact
29  public :: landuse_write
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35  real(RP), public, allocatable :: landuse_fact_ocean(:,:)
36  real(RP), public, allocatable :: landuse_fact_land (:,:)
37  real(RP), public, allocatable :: landuse_fact_urban(:,:)
38 
39  real(RP), public, allocatable :: landuse_frac_land (:,:)
40  real(RP), public, allocatable :: landuse_frac_lake (:,:)
41  real(RP), public, allocatable :: landuse_frac_urban(:,:)
42 
43  integer, public :: landuse_pft_mosaic = 2
44  integer, public :: landuse_pft_nmax = 15
45 
46  real(RP), public, allocatable :: landuse_frac_pft (:,:,:)
47  integer, public, allocatable :: landuse_index_pft(:,:,:)
48 
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private procedure
52  !
53  private :: landuse_read
54 
55  !-----------------------------------------------------------------------------
56  !
57  !++ Private parameters & variables
58  !
59  character(len=H_LONG), private :: landuse_in_basename = ''
60  character(len=H_LONG), private :: landuse_out_basename = ''
61  character(len=H_MID), private :: landuse_out_title = 'SCALE-RM LANDUSE'
62  character(len=H_MID), private :: landuse_out_dtype = 'DEFAULT'
63  logical, private :: landuse_allland = .false.
64  logical, private :: landuse_allurban = .false.
65  logical, private :: landuse_mosaicworld = .false.
66 
67  !-----------------------------------------------------------------------------
68 contains
69  !-----------------------------------------------------------------------------
71  subroutine landuse_setup
72  use scale_process, only: &
74  implicit none
75 
76  namelist / param_landuse / &
77  landuse_in_basename, &
78  landuse_out_basename, &
79  landuse_out_dtype, &
82  landuse_allland, &
83  landuse_allurban, &
84  landuse_mosaicworld
85 
86  integer :: ierr
87  !---------------------------------------------------------------------------
88 
89  if( io_l ) write(io_fid_log,*)
90  if( io_l ) write(io_fid_log,*) '++++++ Module[LANDUSE] / Categ[COUPLER] / Origin[SCALElib]'
91 
92  !--- read namelist
93  rewind(io_fid_conf)
94  read(io_fid_conf,nml=param_landuse,iostat=ierr)
95  if( ierr < 0 ) then !--- missing
96  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
97  elseif( ierr > 0 ) then !--- fatal error
98  write(*,*) 'xxx Not appropriate names in namelist PARAM_LANDUSE. Check!'
99  call prc_mpistop
100  endif
101  if( io_lnml ) write(io_fid_log,nml=param_landuse)
102 
103  allocate( landuse_frac_land(ia,ja) )
104  allocate( landuse_frac_lake(ia,ja) )
105  allocate( landuse_frac_urban(ia,ja) )
106  landuse_frac_land(:,:) = 0.0_rp
107  landuse_frac_lake(:,:) = 0.0_rp
108  landuse_frac_urban(:,:) = 0.0_rp
109 
112  landuse_frac_pft(:,:,:) = 0.0_rp
113  landuse_frac_pft(:,:,1) = 1.0_rp ! tentative, mosaic is off
114  landuse_index_pft(:,:,:) = 1 ! default
115 
116  allocate( landuse_fact_ocean(ia,ja) )
117  allocate( landuse_fact_land(ia,ja) )
118  allocate( landuse_fact_urban(ia,ja) )
119  landuse_fact_ocean(:,:) = 1.0_rp
120  landuse_fact_land(:,:) = 0.0_rp
121  landuse_fact_urban(:,:) = 0.0_rp
122 
123 
124  if ( landuse_allland ) then
125  if( io_l ) write(io_fid_log,*) '*** Assume all grids are land'
126  landuse_frac_land(:,:) = 1.0_rp
127  call landuse_calc_fact
128  elseif( landuse_allurban ) then
129  if( io_l ) write(io_fid_log,*) '*** Assume all grids are land'
130  landuse_frac_land(:,:) = 1.0_rp
131  if( io_l ) write(io_fid_log,*) '*** Assume all lands are urban'
132  landuse_frac_urban(:,:) = 1.0_rp
133  call landuse_calc_fact
134  elseif( landuse_mosaicworld ) then
135  if( io_l ) write(io_fid_log,*) '*** Assume all grids have ocean, land, and urban'
136  landuse_frac_land(:,:) = 0.5_rp
137  landuse_frac_urban(:,:) = 0.5_rp
138  call landuse_calc_fact
139  else
140  ! read from file
141  call landuse_read
142  endif
143 
144  return
145  end subroutine landuse_setup
146 
147  !-----------------------------------------------------------------------------
148  subroutine landuse_calc_fact
149  implicit none
150  !---------------------------------------------------------------------------
151 
152  if( io_l ) write(io_fid_log,*)
153  if( io_l ) write(io_fid_log,*) '+++ calculate landuse factor'
154 
155  ! tentative treatment for lake fraction
156  landuse_frac_land(:,:) = landuse_frac_land(:,:) * ( 1.0_rp - landuse_frac_lake(:,:) )
157 
158  ! make factors
159  landuse_fact_ocean(:,:) = ( 1.0_rp - landuse_frac_land(:,:) )
160  landuse_fact_land(:,:) = ( landuse_frac_land(:,:) ) * ( 1.0_rp - landuse_frac_urban(:,:) )
161  landuse_fact_urban(:,:) = ( landuse_frac_land(:,:) ) * ( landuse_frac_urban(:,:) )
162 
163  return
164  end subroutine landuse_calc_fact
165 
166  !-----------------------------------------------------------------------------
168  subroutine landuse_read
169  use scale_fileio, only: &
170  fileio_read
171  use scale_comm, only: &
172  comm_vars8, &
173  comm_wait
174  implicit none
175 
176  real(RP) :: temp(ia,ja)
177 
178  character(len=H_SHORT) :: varname
179  integer :: p, tag
180  !---------------------------------------------------------------------------
181 
182  if( io_l ) write(io_fid_log,*)
183  if( io_l ) write(io_fid_log,*) '*** Input landuse file ***'
184 
185  if ( landuse_in_basename /= '' ) then
186 
187  call fileio_read( landuse_frac_land(:,:), & ! [OUT]
188  landuse_in_basename, 'FRAC_LAND', 'XY', step=1 ) ! [IN]
189  call fileio_read( landuse_frac_lake(:,:), & ! [OUT]
190  landuse_in_basename, 'FRAC_LAKE', 'XY', step=1 ) ! [IN]
191  call fileio_read( landuse_frac_urban(:,:), & ! [OUT]
192  landuse_in_basename, 'FRAC_URBAN', 'XY', step=1 ) ! [IN]
193 
194  call comm_vars8( landuse_frac_land(:,:), 1 )
195  call comm_vars8( landuse_frac_lake(:,:), 2 )
196  call comm_vars8( landuse_frac_urban(:,:), 3 )
197  call comm_wait ( landuse_frac_land(:,:), 1 )
198  call comm_wait ( landuse_frac_lake(:,:), 2 )
199  call comm_wait ( landuse_frac_urban(:,:), 3 )
200 
201  call fileio_read( landuse_fact_ocean(:,:), & ! [OUT]
202  landuse_in_basename, 'FRAC_OCEAN_abs', 'XY', step=1 ) ! [IN]
203  call fileio_read( landuse_fact_land(:,:), & ! [OUT]
204  landuse_in_basename, 'FRAC_LAND_abs', 'XY', step=1 ) ! [IN]
205  call fileio_read( landuse_fact_urban(:,:), & ! [OUT]
206  landuse_in_basename, 'FRAC_URBAN_abs', 'XY', step=1 ) ! [IN]
207 
208  call comm_vars8( landuse_fact_ocean(:,:), 4 )
209  call comm_vars8( landuse_fact_land(:,:), 5 )
210  call comm_vars8( landuse_fact_urban(:,:), 6 )
211  call comm_wait ( landuse_fact_ocean(:,:), 4 )
212  call comm_wait ( landuse_fact_land(:,:), 5 )
213  call comm_wait ( landuse_fact_urban(:,:), 6 )
214 
215  do p = 1, landuse_pft_mosaic
216  write(varname,'(A8,I1.1)') 'FRAC_PFT', p
217 
218  call fileio_read( landuse_frac_pft(:,:,p), & ! [OUT]
219  landuse_in_basename, varname, 'XY', step=1 ) ! [IN]
220 
221  tag = 7 + (p-1)*landuse_pft_mosaic
222  call comm_vars8( landuse_frac_pft(:,:,p), tag )
223  call comm_wait ( landuse_frac_pft(:,:,p), tag )
224 
225  write(varname,'(A9,I1.1)') 'INDEX_PFT', p
226 
227  call fileio_read( temp(:,:), & ! [OUT]
228  landuse_in_basename, varname, 'XY', step=1 ) ! [IN]
229 
230  tag = 8 + (p-1)*landuse_pft_mosaic
231  call comm_vars8( temp(:,:), tag )
232  call comm_wait ( temp(:,:), tag )
233 
234  landuse_index_pft(:,:,p) = int(temp(:,:))
235  enddo
236 
237  else
238  if( io_l ) write(io_fid_log,*) '*** landuse file is not specified.'
239  if( io_l ) write(io_fid_log,*) '*** Assume all grids are ocean'
240  endif
241 
242  return
243  end subroutine landuse_read
244 
245  !-----------------------------------------------------------------------------
247  subroutine landuse_write
248  use scale_fileio, only: &
249  fileio_write
250  implicit none
251 
252  real(RP) :: temp(ia,ja)
253 
254  character(len=H_SHORT) :: varname
255  integer :: p
256  !---------------------------------------------------------------------------
257 
258  if ( landuse_out_basename /= '' ) then
259 
260  if( io_l ) write(io_fid_log,*)
261  if( io_l ) write(io_fid_log,*) '*** Output landuse file ***'
262 
263  call fileio_write( landuse_frac_land(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
264  'FRAC_LAND', 'LAND fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
265  nozcoord=.true. )
266 
267  call fileio_write( landuse_frac_lake(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
268  'FRAC_LAKE', 'LAKE fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
269  nozcoord=.true. )
270  call fileio_write( landuse_frac_urban(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
271  'FRAC_URBAN', 'URBAN fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
272  nozcoord=.true. )
273 
274  do p = 1, landuse_pft_mosaic
275  write(varname,'(A8,I1.1)') 'FRAC_PFT', p
276 
277  call fileio_write( landuse_frac_pft(:,:,p), landuse_out_basename, landuse_out_title, & ! [IN]
278  varname, 'PFT fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
279  nozcoord=.true. )
280 
281  write(varname,'(A9,I1.1)') 'INDEX_PFT', p
282  temp(:,:) = real(LANDUSE_index_PFT(:,:,p),kind=rp)
283 
284  call fileio_write( temp(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
285  varname, 'PFT index', '1', 'XY', landuse_out_dtype, & ! [IN]
286  nozcoord=.true. )
287  enddo
288 
289  call fileio_write( landuse_fact_ocean(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
290  'FRAC_OCEAN_abs', 'absolute OCEAN fraction', '0-1', 'XY', landuse_out_dtype, & ! [IN]
291  nozcoord=.true. )
292 
293  call fileio_write( landuse_fact_land(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
294  'FRAC_LAND_abs ', 'absolute LAND fraction', '0-1', 'XY', landuse_out_dtype, & ! [IN]
295  nozcoord=.true. )
296 
297  call fileio_write( landuse_fact_urban(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
298  'FRAC_URBAN_abs', 'absolute URBAN fraction', '0-1', 'XY', landuse_out_dtype, & ! [IN]
299  nozcoord=.true. )
300  endif
301 
302  return
303  end subroutine landuse_write
304 
305 end module scale_landuse
real(rp), dimension(:,:,:), allocatable, public landuse_frac_pft
fraction of PFT for each mosaic
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public landuse_setup
Setup.
module STDIO
Definition: scale_stdio.F90:12
subroutine, public landuse_calc_fact
real(rp), dimension(:,:), allocatable, public landuse_frac_urban
urban fraction
module FILE I/O (netcdf)
module grid index
integer, public ia
of x whole cells (local, with HALO)
module LANDUSE
subroutine, public landuse_write
Write landuse data.
module COMMUNICATION
Definition: scale_comm.F90:23
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
module PROCESS
integer, public landuse_pft_mosaic
number of PFT mosaic
integer, dimension(:,:,:), allocatable, public landuse_index_pft
index of PFT for each mosaic
module profiler
Definition: scale_prof.F90:10
real(rp), dimension(:,:), allocatable, public landuse_frac_lake
lake fraction
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
module PRECISION
integer, public landuse_pft_nmax
number of plant functional type(PFT)
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public rp
real(rp), dimension(:,:), allocatable, public landuse_frac_land
land fraction
integer, public ja
of y whole cells (local, with HALO)