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_allocean = .false.
64  logical, private :: landuse_allland = .false.
65  logical, private :: landuse_allurban = .false.
66  logical, private :: landuse_mosaicworld = .false.
67 
68  !-----------------------------------------------------------------------------
69 contains
70  !-----------------------------------------------------------------------------
72  subroutine landuse_setup
73  use scale_process, only: &
75  implicit none
76 
77  namelist / param_landuse / &
78  landuse_in_basename, &
79  landuse_out_basename, &
80  landuse_out_dtype, &
83  landuse_allocean, &
84  landuse_allland, &
85  landuse_allurban, &
86  landuse_mosaicworld
87 
88  integer :: ierr
89  !---------------------------------------------------------------------------
90 
91  if( io_l ) write(io_fid_log,*)
92  if( io_l ) write(io_fid_log,*) '++++++ Module[LANDUSE] / Categ[COUPLER] / Origin[SCALElib]'
93 
94  !--- read namelist
95  rewind(io_fid_conf)
96  read(io_fid_conf,nml=param_landuse,iostat=ierr)
97  if( ierr < 0 ) then !--- missing
98  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
99  elseif( ierr > 0 ) then !--- fatal error
100  write(*,*) 'xxx Not appropriate names in namelist PARAM_LANDUSE. Check!'
101  call prc_mpistop
102  endif
103  if( io_lnml ) write(io_fid_log,nml=param_landuse)
104 
105  allocate( landuse_frac_land(ia,ja) )
106  allocate( landuse_frac_lake(ia,ja) )
107  allocate( landuse_frac_urban(ia,ja) )
108  landuse_frac_land(:,:) = 0.0_rp
109  landuse_frac_lake(:,:) = 0.0_rp
110  landuse_frac_urban(:,:) = 0.0_rp
111 
114  landuse_frac_pft(:,:,:) = 0.0_rp
115  landuse_frac_pft(:,:,1) = 1.0_rp ! tentative, mosaic is off
116  landuse_index_pft(:,:,:) = 1 ! default
117 
118  allocate( landuse_fact_ocean(ia,ja) )
119  allocate( landuse_fact_land(ia,ja) )
120  allocate( landuse_fact_urban(ia,ja) )
121  landuse_fact_ocean(:,:) = 0.0_rp
122  landuse_fact_land(:,:) = 0.0_rp
123  landuse_fact_urban(:,:) = 0.0_rp
124 
125 
126  if ( landuse_allocean ) then
127  if( io_l ) write(io_fid_log,*) '*** Assume all grids are ocean'
128  call landuse_calc_fact
129  elseif( landuse_allland ) then
130  if( io_l ) write(io_fid_log,*) '*** Assume all grids are land'
131  landuse_frac_land(:,:) = 1.0_rp
132  call landuse_calc_fact
133  elseif( landuse_allurban ) then
134  if( io_l ) write(io_fid_log,*) '*** Assume all grids are land'
135  landuse_frac_land(:,:) = 1.0_rp
136  if( io_l ) write(io_fid_log,*) '*** Assume all lands are urban'
137  landuse_frac_urban(:,:) = 1.0_rp
138  call landuse_calc_fact
139  elseif( landuse_mosaicworld ) then
140  if( io_l ) write(io_fid_log,*) '*** Assume all grids have ocean, land, and urban'
141  landuse_frac_land(:,:) = 0.5_rp
142  landuse_frac_urban(:,:) = 0.5_rp
143  call landuse_calc_fact
144  else
145  ! read from file
146  call landuse_read
147  endif
148 
149  return
150  end subroutine landuse_setup
151 
152  !-----------------------------------------------------------------------------
153  subroutine landuse_calc_fact
154  implicit none
155  !---------------------------------------------------------------------------
156 
157  if( io_l ) write(io_fid_log,*)
158  if( io_l ) write(io_fid_log,*) '+++ calculate landuse factor'
159 
160  ! tentative treatment for lake fraction
161  landuse_frac_land(:,:) = landuse_frac_land(:,:) * ( 1.0_rp - landuse_frac_lake(:,:) )
162 
163  ! make factors
164  landuse_fact_ocean(:,:) = ( 1.0_rp - landuse_frac_land(:,:) )
165  landuse_fact_land(:,:) = ( landuse_frac_land(:,:) ) * ( 1.0_rp - landuse_frac_urban(:,:) )
166  landuse_fact_urban(:,:) = ( landuse_frac_land(:,:) ) * ( landuse_frac_urban(:,:) )
167 
168  return
169  end subroutine landuse_calc_fact
170 
171  !-----------------------------------------------------------------------------
173  subroutine landuse_read
174  use scale_fileio, only: &
175  fileio_read
176  use scale_comm, only: &
177  comm_vars8, &
178  comm_wait
179  implicit none
180 
181  real(RP) :: temp(ia,ja)
182 
183  character(len=H_SHORT) :: varname
184  integer :: p, tag
185  !---------------------------------------------------------------------------
186 
187  if( io_l ) write(io_fid_log,*)
188  if( io_l ) write(io_fid_log,*) '*** Input landuse file ***'
189 
190  if ( landuse_in_basename /= '' ) then
191 
192  call fileio_read( landuse_frac_land(:,:), & ! [OUT]
193  landuse_in_basename, 'FRAC_LAND', 'XY', step=1 ) ! [IN]
194  call fileio_read( landuse_frac_lake(:,:), & ! [OUT]
195  landuse_in_basename, 'FRAC_LAKE', 'XY', step=1 ) ! [IN]
196  call fileio_read( landuse_frac_urban(:,:), & ! [OUT]
197  landuse_in_basename, 'FRAC_URBAN', 'XY', step=1 ) ! [IN]
198 
199  call comm_vars8( landuse_frac_land(:,:), 1 )
200  call comm_vars8( landuse_frac_lake(:,:), 2 )
201  call comm_vars8( landuse_frac_urban(:,:), 3 )
202  call comm_wait ( landuse_frac_land(:,:), 1 )
203  call comm_wait ( landuse_frac_lake(:,:), 2 )
204  call comm_wait ( landuse_frac_urban(:,:), 3 )
205 
206  call fileio_read( landuse_fact_ocean(:,:), & ! [OUT]
207  landuse_in_basename, 'FRAC_OCEAN_abs', 'XY', step=1 ) ! [IN]
208  call fileio_read( landuse_fact_land(:,:), & ! [OUT]
209  landuse_in_basename, 'FRAC_LAND_abs', 'XY', step=1 ) ! [IN]
210  call fileio_read( landuse_fact_urban(:,:), & ! [OUT]
211  landuse_in_basename, 'FRAC_URBAN_abs', 'XY', step=1 ) ! [IN]
212 
213  call comm_vars8( landuse_fact_ocean(:,:), 4 )
214  call comm_vars8( landuse_fact_land(:,:), 5 )
215  call comm_vars8( landuse_fact_urban(:,:), 6 )
216  call comm_wait ( landuse_fact_ocean(:,:), 4 )
217  call comm_wait ( landuse_fact_land(:,:), 5 )
218  call comm_wait ( landuse_fact_urban(:,:), 6 )
219 
220  do p = 1, landuse_pft_mosaic
221  write(varname,'(A8,I1.1)') 'FRAC_PFT', p
222 
223  call fileio_read( landuse_frac_pft(:,:,p), & ! [OUT]
224  landuse_in_basename, varname, 'XY', step=1 ) ! [IN]
225 
226  tag = 7 + (p-1)*landuse_pft_mosaic
227  call comm_vars8( landuse_frac_pft(:,:,p), tag )
228  call comm_wait ( landuse_frac_pft(:,:,p), tag )
229 
230  write(varname,'(A9,I1.1)') 'INDEX_PFT', p
231 
232  call fileio_read( temp(:,:), & ! [OUT]
233  landuse_in_basename, varname, 'XY', step=1 ) ! [IN]
234 
235  tag = 8 + (p-1)*landuse_pft_mosaic
236  call comm_vars8( temp(:,:), tag )
237  call comm_wait ( temp(:,:), tag )
238 
239  landuse_index_pft(:,:,p) = int(temp(:,:))
240  enddo
241 
242  else
243  if( io_l ) write(io_fid_log,*) '*** landuse file is not specified.'
244  if( io_l ) write(io_fid_log,*) '*** Assume all grids are ocean'
245  endif
246 
247  return
248  end subroutine landuse_read
249 
250  !-----------------------------------------------------------------------------
252  subroutine landuse_write
253  use scale_fileio, only: &
254  fileio_write
255  implicit none
256 
257  real(RP) :: temp(ia,ja)
258 
259  character(len=H_SHORT) :: varname
260  integer :: p
261  !---------------------------------------------------------------------------
262 
263  if ( landuse_out_basename /= '' ) then
264 
265  if( io_l ) write(io_fid_log,*)
266  if( io_l ) write(io_fid_log,*) '*** Output landuse file ***'
267 
268  call fileio_write( landuse_frac_land(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
269  'FRAC_LAND', 'LAND fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
270  nozcoord=.true. )
271 
272  call fileio_write( landuse_frac_lake(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
273  'FRAC_LAKE', 'LAKE fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
274  nozcoord=.true. )
275  call fileio_write( landuse_frac_urban(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
276  'FRAC_URBAN', 'URBAN fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
277  nozcoord=.true. )
278 
279  do p = 1, landuse_pft_mosaic
280  write(varname,'(A8,I1.1)') 'FRAC_PFT', p
281 
282  call fileio_write( landuse_frac_pft(:,:,p), landuse_out_basename, landuse_out_title, & ! [IN]
283  varname, 'PFT fraction', '1', 'XY', landuse_out_dtype, & ! [IN]
284  nozcoord=.true. )
285 
286  write(varname,'(A9,I1.1)') 'INDEX_PFT', p
287  temp(:,:) = real(LANDUSE_index_PFT(:,:,p),kind=rp)
288 
289  call fileio_write( temp(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
290  varname, 'PFT index', '1', 'XY', landuse_out_dtype, & ! [IN]
291  nozcoord=.true. )
292  enddo
293 
294  call fileio_write( landuse_fact_ocean(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
295  'FRAC_OCEAN_abs', 'absolute OCEAN fraction', '0-1', 'XY', landuse_out_dtype, & ! [IN]
296  nozcoord=.true. )
297 
298  call fileio_write( landuse_fact_land(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
299  'FRAC_LAND_abs ', 'absolute LAND fraction', '0-1', 'XY', landuse_out_dtype, & ! [IN]
300  nozcoord=.true. )
301 
302  call fileio_write( landuse_fact_urban(:,:), landuse_out_basename, landuse_out_title, & ! [IN]
303  'FRAC_URBAN_abs', 'absolute URBAN fraction', '0-1', 'XY', landuse_out_dtype, & ! [IN]
304  nozcoord=.true. )
305  endif
306 
307  return
308  end subroutine landuse_write
309 
310 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)