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  logical, private :: landuse_in_check_coordinates = .true.
61  character(len=H_LONG), private :: landuse_out_basename = ''
62  character(len=H_MID), private :: landuse_out_title = 'SCALE-RM LANDUSE'
63  character(len=H_SHORT), private :: landuse_out_dtype = 'DEFAULT'
64  logical, private :: landuse_allocean = .false.
65  logical, private :: landuse_allland = .false.
66  logical, private :: landuse_allurban = .false.
67  logical, private :: landuse_mosaicworld = .false.
68 
69  !-----------------------------------------------------------------------------
70 contains
71  !-----------------------------------------------------------------------------
73  subroutine landuse_setup
74  use scale_process, only: &
76  implicit none
77 
78  namelist / param_landuse / &
79  landuse_in_basename, &
80  landuse_in_check_coordinates, &
81  landuse_out_basename, &
82  landuse_out_dtype, &
85  landuse_allocean, &
86  landuse_allland, &
87  landuse_allurban, &
88  landuse_mosaicworld
89 
90  integer :: ierr
91  !---------------------------------------------------------------------------
92 
93  if( io_l ) write(io_fid_log,*)
94  if( io_l ) write(io_fid_log,*) '++++++ Module[LANDUSE] / Categ[COUPLER] / Origin[SCALElib]'
95 
96  !--- read namelist
97  rewind(io_fid_conf)
98  read(io_fid_conf,nml=param_landuse,iostat=ierr)
99  if( ierr < 0 ) then !--- missing
100  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
101  elseif( ierr > 0 ) then !--- fatal error
102  write(*,*) 'xxx Not appropriate names in namelist PARAM_LANDUSE. Check!'
103  call prc_mpistop
104  endif
105  if( io_nml ) write(io_fid_nml,nml=param_landuse)
106 
107  allocate( landuse_frac_land(ia,ja) )
108  allocate( landuse_frac_lake(ia,ja) )
109  allocate( landuse_frac_urban(ia,ja) )
110  landuse_frac_land(:,:) = 0.0_rp
111  landuse_frac_lake(:,:) = 0.0_rp
112  landuse_frac_urban(:,:) = 0.0_rp
113 
116  landuse_frac_pft(:,:,:) = 0.0_rp
117  landuse_frac_pft(:,:,1) = 1.0_rp ! tentative, mosaic is off
118  landuse_index_pft(:,:,:) = 1 ! default
119 
120  allocate( landuse_fact_ocean(ia,ja) )
121  allocate( landuse_fact_land(ia,ja) )
122  allocate( landuse_fact_urban(ia,ja) )
123  landuse_fact_ocean(:,:) = 0.0_rp
124  landuse_fact_land(:,:) = 0.0_rp
125  landuse_fact_urban(:,:) = 0.0_rp
126 
127 
128  if ( landuse_allocean ) then
129  if( io_l ) write(io_fid_log,*) '*** Assume all grids are ocean'
130  call landuse_calc_fact
131  elseif( landuse_allland ) then
132  if( io_l ) write(io_fid_log,*) '*** Assume all grids are land'
133  landuse_frac_land(:,:) = 1.0_rp
134  call landuse_calc_fact
135  elseif( landuse_allurban ) then
136  if( io_l ) write(io_fid_log,*) '*** Assume all grids are land'
137  landuse_frac_land(:,:) = 1.0_rp
138  if( io_l ) write(io_fid_log,*) '*** Assume all lands are urban'
139  landuse_frac_urban(:,:) = 1.0_rp
140  call landuse_calc_fact
141  elseif( landuse_mosaicworld ) then
142  if( io_l ) write(io_fid_log,*) '*** Assume all grids have ocean, land, and urban'
143  landuse_frac_land(:,:) = 0.5_rp
144  landuse_frac_urban(:,:) = 0.5_rp
145  call landuse_calc_fact
146  else
147  ! read from file
148  call landuse_read
149  endif
150 
151  return
152  end subroutine landuse_setup
153 
154  !-----------------------------------------------------------------------------
155  subroutine landuse_calc_fact
156  implicit none
157  !---------------------------------------------------------------------------
158 
159  if( io_l ) write(io_fid_log,*)
160  if( io_l ) write(io_fid_log,*) '+++ calculate landuse factor'
161 
162  ! tentative treatment: The area of the lake is treated as the ocean
163  landuse_frac_land(:,:) = landuse_frac_land(:,:) * ( 1.0_rp - landuse_frac_lake(:,:) )
164 
165  ! make factors
166  landuse_fact_ocean(:,:) = ( 1.0_rp - landuse_frac_land(:,:) )
167  landuse_fact_land(:,:) = ( landuse_frac_land(:,:) ) * ( 1.0_rp - landuse_frac_urban(:,:) )
168  landuse_fact_urban(:,:) = ( landuse_frac_land(:,:) ) * ( landuse_frac_urban(:,:) )
169 
170  return
171  end subroutine landuse_calc_fact
172 
173  !-----------------------------------------------------------------------------
175  subroutine landuse_read
176  use scale_fileio, only: &
177  fileio_open, &
178  fileio_read, &
179  fileio_flush, &
180  fileio_check_coordinates, &
182  use scale_comm, only: &
183  comm_vars8, &
184  comm_wait
185  implicit none
186 
187  real(RP) :: temp(ia,ja)
188 
189  character(len=H_SHORT) :: varname
190  integer :: fid
191  integer :: p, tag
192  !---------------------------------------------------------------------------
193 
194  if( io_l ) write(io_fid_log,*)
195  if( io_l ) write(io_fid_log,*) '*** Input landuse file ***'
196 
197  if ( landuse_in_basename /= '' ) then
198 
199  call fileio_open( fid, landuse_in_basename )
200 
201  call fileio_read( landuse_frac_land(:,:), & ! [OUT]
202  fid, 'FRAC_LAND', 'XY', step=1 ) ! [IN]
203  call fileio_read( landuse_frac_lake(:,:), & ! [OUT]
204  fid, 'FRAC_LAKE', 'XY', step=1 ) ! [IN]
205  call fileio_read( landuse_frac_urban(:,:), & ! [OUT]
206  fid, 'FRAC_URBAN', 'XY', step=1 ) ! [IN]
207 
208  call fileio_read( landuse_fact_ocean(:,:), & ! [OUT]
209  fid, 'FRAC_OCEAN_abs', 'XY', step=1 ) ! [IN]
210  call fileio_read( landuse_fact_land(:,:), & ! [OUT]
211  fid, 'FRAC_LAND_abs', 'XY', step=1 ) ! [IN]
212  call fileio_read( landuse_fact_urban(:,:), & ! [OUT]
213  fid, 'FRAC_URBAN_abs', 'XY', step=1 ) ! [IN]
214 
215  call fileio_flush( fid )
216 
217  if ( landuse_in_check_coordinates ) then
218  call fileio_check_coordinates( fid )
219  end if
220 
221  call comm_vars8( landuse_frac_land(:,:), 1 )
222  call comm_vars8( landuse_frac_lake(:,:), 2 )
223  call comm_vars8( landuse_frac_urban(:,:), 3 )
224  call comm_vars8( landuse_fact_ocean(:,:), 4 )
225  call comm_vars8( landuse_fact_land(:,:), 5 )
226  call comm_vars8( landuse_fact_urban(:,:), 6 )
227  call comm_wait ( landuse_frac_land(:,:), 1 )
228  call comm_wait ( landuse_frac_lake(:,:), 2 )
229  call comm_wait ( landuse_frac_urban(:,:), 3 )
230  call comm_wait ( landuse_fact_ocean(:,:), 4 )
231  call comm_wait ( landuse_fact_land(:,:), 5 )
232  call comm_wait ( landuse_fact_urban(:,:), 6 )
233 
234  do p = 1, landuse_pft_mosaic
235  write(varname,'(A8,I1.1)') 'FRAC_PFT', p
236 
237  call fileio_read( landuse_frac_pft(:,:,p), & ! [OUT]
238  fid, varname, 'XY', step=1 ) ! [IN]
239  call fileio_flush( fid )
240 
241  tag = 7 + (p-1)*landuse_pft_mosaic
242  call comm_vars8( landuse_frac_pft(:,:,p), tag )
243  call comm_wait ( landuse_frac_pft(:,:,p), tag )
244 
245  write(varname,'(A9,I1.1)') 'INDEX_PFT', p
246 
247  call fileio_read( temp(:,:), & ! [OUT]
248  fid, varname, 'XY', step=1 ) ! [IN]
249  call fileio_flush( fid )
250 
251  tag = 8 + (p-1)*landuse_pft_mosaic
252  call comm_vars8( temp(:,:), tag )
253  call comm_wait ( temp(:,:), tag )
254 
255  landuse_index_pft(:,:,p) = int(temp(:,:))
256  enddo
257 
258  call fileio_close( fid )
259 
260  else
261  if( io_l ) write(io_fid_log,*) '*** landuse file is not specified.'
262  if( io_l ) write(io_fid_log,*) '*** Assume all grids are ocean'
263  endif
264 
265  return
266  end subroutine landuse_read
267 
268  !-----------------------------------------------------------------------------
270  subroutine landuse_write
271  use scale_fileio, only: &
272  fileio_create, &
273  fileio_def_var, &
274  fileio_enddef, &
275  fileio_write_var, &
277  implicit none
278 
279  real(RP) :: temp(ia,ja)
280 
281  integer :: fid, vid(6+landuse_pft_mosaic*2)
282  character(len=H_SHORT) :: varname
283  integer :: p
284  !---------------------------------------------------------------------------
285 
286  if ( landuse_out_basename /= '' ) then
287 
288  if( io_l ) write(io_fid_log,*)
289  if( io_l ) write(io_fid_log,*) '*** Output landuse file ***'
290 
291  call fileio_create( fid, landuse_out_basename, landuse_out_title, &
292  landuse_out_dtype, nozcoord=.true. )
293 
294  call fileio_def_var( fid, vid(1), 'FRAC_LAND', 'LAND fraction', '1', 'XY', landuse_out_dtype )
295  call fileio_def_var( fid, vid(2), 'FRAC_LAKE', 'LAKE fraction', '1', 'XY', landuse_out_dtype )
296  call fileio_def_var( fid, vid(3), 'FRAC_URBAN', 'URBAN fraction', '1', 'XY', landuse_out_dtype )
297  call fileio_def_var( fid, vid(4), 'FRAC_OCEAN_abs', 'absolute OCEAN fraction', '1', 'XY', landuse_out_dtype )
298  call fileio_def_var( fid, vid(5), 'FRAC_LAND_abs', 'absolute LAND fraction', '1', 'XY', landuse_out_dtype )
299  call fileio_def_var( fid, vid(6), 'FRAC_URBAN_abs', 'absolute URBAN fraction', '1', 'XY', landuse_out_dtype )
300 
301  do p = 1, landuse_pft_mosaic
302  write(varname,'(A8,I1.1)') 'FRAC_PFT', p
303  call fileio_def_var( fid, vid(7+2*(p-1)), varname, 'PFT fraction', '1', 'XY', landuse_out_dtype )
304  write(varname,'(A9,I1.1)') 'INDEX_PFT', p
305  call fileio_def_var( fid, vid(8+2*(p-1)), varname, 'PFT index', '1', 'XY', landuse_out_dtype )
306  end do
307 
308  call fileio_enddef( fid )
309 
310  call fileio_write_var( fid, vid(1), landuse_frac_land(:,:), 'FRAC_LAND', 'XY' )
311  call fileio_write_var( fid, vid(2), landuse_frac_lake(:,:), 'FRAC_LAKE', 'XY' )
312  call fileio_write_var( fid, vid(3), landuse_frac_urban(:,:), 'FRAC_URBAN', 'XY' )
313  call fileio_write_var( fid, vid(4), landuse_fact_ocean(:,:), 'FRAC_OCEAN_abs', 'XY' )
314  call fileio_write_var( fid, vid(5), landuse_fact_land(:,:), 'FRAC_LAND_abs', 'XY' )
315  call fileio_write_var( fid, vid(6), landuse_fact_urban(:,:), 'FRAC_URBAN_abs', 'XY' )
316 
317  do p = 1, landuse_pft_mosaic
318  write(varname,'(A8,I1.1)') 'FRAC_PFT', p
319  call fileio_write_var( fid, vid(7+2*(p-1)), landuse_frac_pft(:,:,p), varname, 'XY' )
320  write(varname,'(A9,I1.1)') 'INDEX_PFT', p
321  temp(:,:) = real(LANDUSE_index_PFT(:,:,p),kind=rp)
322  call fileio_write_var( fid, vid(8+2*(p-1)), temp(:,:), varname, 'XY' )
323  end do
324 
325  call fileio_close( fid )
326 
327  end if
328 
329  return
330  end subroutine landuse_write
331 
332 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:61
subroutine, public landuse_setup
Setup.
subroutine, public fileio_flush(fid)
Flush all pending requests to a netCDF file (PnetCDF only)
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
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:62
integer, public ia
of whole cells: x, local, with HALO
module LANDUSE
subroutine, public fileio_create(fid, basename, title, datatype, date, subsec, append, nozcoord)
Create/open a netCDF file.
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
subroutine, public fileio_enddef(fid)
Exit netCDF file define mode.
module profiler
Definition: scale_prof.F90:10
real(rp), dimension(:,:), allocatable, public landuse_frac_lake
lake fraction
subroutine, public fileio_open(fid, basename)
open a netCDF file for read
module PRECISION
subroutine, public fileio_def_var(fid, vid, varname, desc, unit, axistype, datatype, timeintv, nsteps)
Define a variable to file.
integer, public landuse_pft_nmax
number of plant functional type(PFT)
subroutine, public fileio_close(fid)
Close a netCDF file.
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 io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
integer, public ja
of whole cells: y, local, with HALO