| 1 | ! radsurf_properties.f90 - Derived type for surface properties |
|---|
| 2 | ! |
|---|
| 3 | ! (C) Copyright 2017- ECMWF. |
|---|
| 4 | ! |
|---|
| 5 | ! This software is licensed under the terms of the Apache Licence Version 2.0 |
|---|
| 6 | ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0. |
|---|
| 7 | ! |
|---|
| 8 | ! In applying this licence, ECMWF does not waive the privileges and immunities |
|---|
| 9 | ! granted to it by virtue of its status as an intergovernmental organisation |
|---|
| 10 | ! nor does it submit to any jurisdiction. |
|---|
| 11 | ! |
|---|
| 12 | ! Author: Robin Hogan |
|---|
| 13 | ! Email: r.j.hogan@ecmwf.int |
|---|
| 14 | ! |
|---|
| 15 | |
|---|
| 16 | module radsurf_properties |
|---|
| 17 | |
|---|
| 18 | use parkind1, only : jpim, jprb |
|---|
| 19 | |
|---|
| 20 | implicit none |
|---|
| 21 | |
|---|
| 22 | public |
|---|
| 23 | |
|---|
| 24 | ! Number of tile types |
|---|
| 25 | integer(kind=jpim), parameter :: NTileTypes = 3 |
|---|
| 26 | |
|---|
| 27 | ! Codes for the different type of tile |
|---|
| 28 | enum, bind(c) |
|---|
| 29 | enumerator :: ITileFlat = 1, & |
|---|
| 30 | & ITileVegetation, & |
|---|
| 31 | & ITileUrban3D |
|---|
| 32 | end enum |
|---|
| 33 | |
|---|
| 34 | character(len=*), parameter :: TileRepresentationName(NTileTypes) & |
|---|
| 35 | & = (/ 'Flat ', & |
|---|
| 36 | & 'HomogeneousVegetation', & |
|---|
| 37 | & 'Urban3D ' /) |
|---|
| 38 | |
|---|
| 39 | ! Number of surface facets and regions |
|---|
| 40 | integer(kind=jpim), parameter :: NTileFacets (NTileTypes) = (/ 1, 1, 3 /) |
|---|
| 41 | integer(kind=jpim), parameter :: NTileRegions(NTileTypes) = (/ 0, 1, 1 /) |
|---|
| 42 | |
|---|
| 43 | |
|---|
| 44 | !--------------------------------------------------------------------- |
|---|
| 45 | ! Derived type storing a physical description of the properties of |
|---|
| 46 | ! the surface tiles needed to compute the boundary conditions for |
|---|
| 47 | ! the radiation scheme. |
|---|
| 48 | type surface_type |
|---|
| 49 | ! Skin temperature (ncol,nfacet) |
|---|
| 50 | real(kind=jprb), allocatable :: skin_temperature(:,:) |
|---|
| 51 | |
|---|
| 52 | ! Air temperature in canopy (ncol,ntile) |
|---|
| 53 | real(kind=jprb), allocatable :: canopy_temperature(:,:) |
|---|
| 54 | |
|---|
| 55 | ! Shortwave albedo: if sw_albedo_direct is not allocated then |
|---|
| 56 | ! sw_albedo will be used for both direct and diffuse solar |
|---|
| 57 | ! radiation; otherwise, sw_albedo will be used for diffuse |
|---|
| 58 | ! radiation and sw_albedo_direct for direct |
|---|
| 59 | ! radiation. Dimensioning is (ncol,nalbedobands,nfacet). |
|---|
| 60 | real(kind=jprb), allocatable, dimension(:,:,:) :: & |
|---|
| 61 | & sw_albedo, sw_albedo_direct |
|---|
| 62 | |
|---|
| 63 | ! Longwave emissivity: dimensioning is (ncol,nemissbands,nfacet) |
|---|
| 64 | real(kind=jprb), allocatable, dimension(:,:,:) :: & |
|---|
| 65 | & lw_emissivity |
|---|
| 66 | |
|---|
| 67 | ! Representation codes (ITileFlat etc) for each tile: dimensioning |
|---|
| 68 | ! is (ntile) |
|---|
| 69 | integer(kind=jpim), allocatable, dimension(:) :: & |
|---|
| 70 | & i_representation |
|---|
| 71 | |
|---|
| 72 | ! Indices to the facets: dimensioning is (ncol,ntile) |
|---|
| 73 | integer(kind=jpim), allocatable, dimension(:,:) :: & |
|---|
| 74 | & i_ground_facet, & |
|---|
| 75 | & i_roof_facet, & |
|---|
| 76 | & i_wall_facet |
|---|
| 77 | |
|---|
| 78 | ! Indices to the regions: dimensioning is (ncol,ntile) |
|---|
| 79 | integer(kind=jpim), allocatable, dimension(:,:) :: & |
|---|
| 80 | & i_region_1 !, i_region_2 |
|---|
| 81 | |
|---|
| 82 | ! Physical properties dimensioned (ncol,ntile) |
|---|
| 83 | real(kind=jprb), allocatable, dimension(:,:) :: & |
|---|
| 84 | & tile_fraction, & ! Fraction of column occupied by each tile |
|---|
| 85 | & canopy_depth ! Depth of canopy (m) |
|---|
| 86 | |
|---|
| 87 | ! Building properties in urban tile dimensioned (ncol,ntile) |
|---|
| 88 | real(kind=jprb), allocatable, dimension(:,:) :: & |
|---|
| 89 | & building_fraction, & ! Fraction of canopy containing buildings |
|---|
| 90 | & building_normalized_perimeter ! Perimeter length of buildings per unit area (m) |
|---|
| 91 | |
|---|
| 92 | ! Vegetation physical properties dimensioned (ncol,ntile) |
|---|
| 93 | real(kind=jprb), allocatable, dimension(:,:) :: & |
|---|
| 94 | & vegetation_optical_depth, & ! Same as LAI? |
|---|
| 95 | & vegetation_fractional_std ! Fractional standard deviation |
|---|
| 96 | |
|---|
| 97 | ! Vegetation optical properties |
|---|
| 98 | real(kind=jprb), allocatable, dimension(:,:,:) :: & |
|---|
| 99 | & vegetation_sw_albedo, & ! (ncol,nalbedobands,ntile) |
|---|
| 100 | & vegetation_lw_emissivity! (ncol,nemissbands,ntile) |
|---|
| 101 | |
|---|
| 102 | ! When 3D effects in vegetation canopies are represented, these will |
|---|
| 103 | ! be available |
|---|
| 104 | ! & vegetation_normalized_perimeter |
|---|
| 105 | ! & vegetation_fraction, & |
|---|
| 106 | |
|---|
| 107 | ! Number of tiles, columns, facets, regions |
|---|
| 108 | integer(kind=jpim) :: ntile, ncol, nfacet, nregion |
|---|
| 109 | |
|---|
| 110 | ! Number of albedo and emissivity bands |
|---|
| 111 | integer(kind=jpim) :: nalbedobands, nemissbands |
|---|
| 112 | |
|---|
| 113 | ! Is this a simple surface containing one "flat" tile? |
|---|
| 114 | logical :: is_simple = .true. |
|---|
| 115 | |
|---|
| 116 | contains |
|---|
| 117 | procedure :: allocate => allocate_surface |
|---|
| 118 | procedure :: deallocate => deallocate_surface |
|---|
| 119 | procedure :: read => read_from_netcdf |
|---|
| 120 | procedure :: set_facet_indices |
|---|
| 121 | |
|---|
| 122 | end type surface_type |
|---|
| 123 | |
|---|
| 124 | |
|---|
| 125 | contains |
|---|
| 126 | |
|---|
| 127 | !--------------------------------------------------------------------- |
|---|
| 128 | subroutine allocate_surface(this, ncol, nalbedobands, nemissbands, & |
|---|
| 129 | & i_representation, & |
|---|
| 130 | & use_sw_albedo_direct) |
|---|
| 131 | |
|---|
| 132 | use yomhook, only : lhook, dr_hook |
|---|
| 133 | |
|---|
| 134 | class(surface_type), intent(inout) :: this |
|---|
| 135 | integer(kind=jpim), intent(in) :: ncol, nalbedobands, nemissbands |
|---|
| 136 | integer(kind=jpim), intent(in), optional :: i_representation(:) |
|---|
| 137 | logical(kind=jpim), intent(in), optional :: use_sw_albedo_direct |
|---|
| 138 | |
|---|
| 139 | ! Number of facets and regions |
|---|
| 140 | integer(kind=jpim) :: nfacet, nregion, ntile |
|---|
| 141 | |
|---|
| 142 | real(jprb) :: hook_handle |
|---|
| 143 | |
|---|
| 144 | if (lhook) call dr_hook('radiation_surface:allocate',0,hook_handle) |
|---|
| 145 | |
|---|
| 146 | call this%deallocate() |
|---|
| 147 | |
|---|
| 148 | this%ncol = ncol |
|---|
| 149 | if (present(i_representation)) then |
|---|
| 150 | nfacet = sum(NTileFacets (i_representation)) |
|---|
| 151 | nregion = sum(NTileRegions(i_representation)) |
|---|
| 152 | ntile = size(i_representation) |
|---|
| 153 | this%is_simple = .false. |
|---|
| 154 | else |
|---|
| 155 | nfacet = 1 |
|---|
| 156 | nregion = 0 |
|---|
| 157 | ntile = 1 |
|---|
| 158 | this%is_simple = .true. |
|---|
| 159 | end if |
|---|
| 160 | this%nfacet = nfacet |
|---|
| 161 | this%nregion= nregion |
|---|
| 162 | this%ntile = ntile |
|---|
| 163 | this%nalbedobands = nalbedobands |
|---|
| 164 | this%nemissbands = nemissbands |
|---|
| 165 | |
|---|
| 166 | allocate(this%skin_temperature(ncol,nfacet)) |
|---|
| 167 | allocate(this%sw_albedo(ncol, nalbedobands, nfacet)) |
|---|
| 168 | |
|---|
| 169 | if (present(use_sw_albedo_direct)) then |
|---|
| 170 | if (use_sw_albedo_direct) then |
|---|
| 171 | allocate(this%sw_albedo_direct(ncol, nalbedobands, nfacet)) |
|---|
| 172 | end if |
|---|
| 173 | end if |
|---|
| 174 | |
|---|
| 175 | allocate(this%lw_emissivity(ncol, nemissbands, nfacet)) |
|---|
| 176 | |
|---|
| 177 | if (this%ntile > 1) then |
|---|
| 178 | allocate(this%tile_fraction(ncol,ntile)) |
|---|
| 179 | end if |
|---|
| 180 | |
|---|
| 181 | if (nregion > 0) then |
|---|
| 182 | allocate(this%canopy_depth(ncol,ntile)) |
|---|
| 183 | allocate(this%canopy_temperature(ncol,ntile)) |
|---|
| 184 | this%canopy_depth = 0.0_jprb |
|---|
| 185 | this%canopy_temperature = -1.0_jprb |
|---|
| 186 | end if |
|---|
| 187 | |
|---|
| 188 | if (.not. this%is_simple) then |
|---|
| 189 | allocate(this%i_representation(ntile)) |
|---|
| 190 | this%i_representation = i_representation |
|---|
| 191 | |
|---|
| 192 | if (any(i_representation == ITileUrban3D)) then |
|---|
| 193 | allocate(this%building_fraction(ncol,ntile)) |
|---|
| 194 | allocate(this%building_normalized_perimeter(ncol,ntile)) |
|---|
| 195 | ! Initialize to default values |
|---|
| 196 | this%building_fraction = 0.0_jprb |
|---|
| 197 | this%building_normalized_perimeter = 0.0_jprb |
|---|
| 198 | end if |
|---|
| 199 | |
|---|
| 200 | if (any(i_representation == ITileVegetation)) then |
|---|
| 201 | allocate(this%vegetation_optical_depth(ncol,ntile)) |
|---|
| 202 | allocate(this%vegetation_fractional_std(ncol,ntile)) |
|---|
| 203 | allocate(this%vegetation_sw_albedo(ncol,nalbedobands,ntile)) |
|---|
| 204 | allocate(this%vegetation_lw_emissivity(ncol,nemissbands,ntile)) |
|---|
| 205 | ! Initialize to default values |
|---|
| 206 | this%vegetation_optical_depth = 0.0_jprb |
|---|
| 207 | this%vegetation_fractional_std = 0.0_jprb |
|---|
| 208 | this%vegetation_sw_albedo = 0.0_jprb |
|---|
| 209 | this%vegetation_lw_emissivity = 1.0_jprb |
|---|
| 210 | end if |
|---|
| 211 | |
|---|
| 212 | call this%set_facet_indices |
|---|
| 213 | |
|---|
| 214 | end if |
|---|
| 215 | |
|---|
| 216 | if (lhook) call dr_hook('radiation_surface:allocate',1,hook_handle) |
|---|
| 217 | |
|---|
| 218 | end subroutine allocate_surface |
|---|
| 219 | |
|---|
| 220 | |
|---|
| 221 | !--------------------------------------------------------------------- |
|---|
| 222 | ! Set the indices to facets |
|---|
| 223 | subroutine set_facet_indices(this) |
|---|
| 224 | |
|---|
| 225 | use radiation_io, only : nulerr, radiation_abort |
|---|
| 226 | |
|---|
| 227 | class(surface_type), intent(inout) :: this |
|---|
| 228 | |
|---|
| 229 | ! Indices to tiles and facets |
|---|
| 230 | integer(kind=jpim) :: ifacet, iregion, jtile |
|---|
| 231 | if (.not. this%is_simple) then |
|---|
| 232 | |
|---|
| 233 | allocate(this%i_ground_facet(this%ncol,this%ntile)) |
|---|
| 234 | this%i_ground_facet = 0 |
|---|
| 235 | |
|---|
| 236 | if (any(this%i_representation == ITileUrban3D)) then |
|---|
| 237 | allocate(this%i_roof_facet(this%ncol,this%ntile)) |
|---|
| 238 | allocate(this%i_wall_facet(this%ncol,this%ntile)) |
|---|
| 239 | ! Initialize to default values |
|---|
| 240 | this%i_roof_facet = 0 |
|---|
| 241 | this%i_wall_facet = 0 |
|---|
| 242 | end if |
|---|
| 243 | |
|---|
| 244 | if (this%nregion > 0) then |
|---|
| 245 | allocate(this%i_region_1(this%ncol,this%ntile)) |
|---|
| 246 | this%i_region_1 = 0 |
|---|
| 247 | end if |
|---|
| 248 | |
|---|
| 249 | ! Set the indices to the various facets |
|---|
| 250 | ifacet = 1 |
|---|
| 251 | iregion = 1 |
|---|
| 252 | do jtile = 1,this%ntile |
|---|
| 253 | this%i_ground_facet(:,jtile) = ifacet |
|---|
| 254 | ifacet = ifacet + 1 |
|---|
| 255 | if (this%i_representation(jtile) == ITileVegetation) then |
|---|
| 256 | this%i_region_1(:,jtile) = iregion |
|---|
| 257 | iregion = iregion+1 |
|---|
| 258 | else if (this%i_representation(jtile) == ITileUrban3D) then |
|---|
| 259 | this%i_roof_facet(:,jtile) = ifacet |
|---|
| 260 | this%i_wall_facet(:,jtile) = ifacet+1 |
|---|
| 261 | ifacet = ifacet+2 |
|---|
| 262 | this%i_region_1(:,jtile) = iregion |
|---|
| 263 | iregion = iregion+1 |
|---|
| 264 | else if (this%i_representation(jtile) /= ITileFlat) then |
|---|
| 265 | write(nulerr,'(a,i0,a)') '*** Error: tile representation ', & |
|---|
| 266 | & this%i_representation(jtile), ' not understood' |
|---|
| 267 | call radiation_abort() |
|---|
| 268 | end if |
|---|
| 269 | end do |
|---|
| 270 | end if |
|---|
| 271 | end subroutine set_facet_indices |
|---|
| 272 | |
|---|
| 273 | !--------------------------------------------------------------------- |
|---|
| 274 | subroutine deallocate_surface(this) |
|---|
| 275 | |
|---|
| 276 | use yomhook, only : lhook, dr_hook |
|---|
| 277 | |
|---|
| 278 | class(surface_type), intent(inout) :: this |
|---|
| 279 | |
|---|
| 280 | real(jprb) :: hook_handle |
|---|
| 281 | |
|---|
| 282 | if (lhook) call dr_hook('radiation_surface:deallocate',0,hook_handle) |
|---|
| 283 | |
|---|
| 284 | if (allocated(this%skin_temperature)) then |
|---|
| 285 | deallocate(this%skin_temperature) |
|---|
| 286 | end if |
|---|
| 287 | if (allocated(this%sw_albedo)) then |
|---|
| 288 | deallocate(this%sw_albedo) |
|---|
| 289 | end if |
|---|
| 290 | if (allocated(this%sw_albedo_direct)) then |
|---|
| 291 | deallocate(this%sw_albedo_direct) |
|---|
| 292 | end if |
|---|
| 293 | if (allocated(this%lw_emissivity)) then |
|---|
| 294 | deallocate(this%lw_emissivity) |
|---|
| 295 | end if |
|---|
| 296 | if (allocated(this%i_representation)) then |
|---|
| 297 | deallocate(this%i_representation) |
|---|
| 298 | end if |
|---|
| 299 | if (allocated(this%i_ground_facet)) then |
|---|
| 300 | deallocate(this%i_ground_facet) |
|---|
| 301 | end if |
|---|
| 302 | if (allocated(this%i_roof_facet)) then |
|---|
| 303 | deallocate(this%i_roof_facet) |
|---|
| 304 | end if |
|---|
| 305 | if (allocated(this%i_wall_facet)) then |
|---|
| 306 | deallocate(this%i_wall_facet) |
|---|
| 307 | end if |
|---|
| 308 | if (allocated(this%tile_fraction)) then |
|---|
| 309 | deallocate(this%tile_fraction) |
|---|
| 310 | end if |
|---|
| 311 | if (allocated(this%canopy_depth)) then |
|---|
| 312 | deallocate(this%canopy_depth) |
|---|
| 313 | end if |
|---|
| 314 | if (allocated(this%canopy_temperature)) then |
|---|
| 315 | deallocate(this%canopy_temperature) |
|---|
| 316 | end if |
|---|
| 317 | if (allocated(this%building_fraction)) then |
|---|
| 318 | deallocate(this%building_fraction) |
|---|
| 319 | end if |
|---|
| 320 | if (allocated(this%building_normalized_perimeter)) then |
|---|
| 321 | deallocate(this%building_normalized_perimeter) |
|---|
| 322 | end if |
|---|
| 323 | if (allocated(this%vegetation_optical_depth)) then |
|---|
| 324 | deallocate(this%vegetation_optical_depth) |
|---|
| 325 | end if |
|---|
| 326 | if (allocated(this%vegetation_fractional_std)) then |
|---|
| 327 | deallocate(this%vegetation_fractional_std) |
|---|
| 328 | end if |
|---|
| 329 | if (allocated(this%vegetation_sw_albedo)) then |
|---|
| 330 | deallocate(this%vegetation_sw_albedo) |
|---|
| 331 | end if |
|---|
| 332 | if (allocated(this%vegetation_lw_emissivity)) then |
|---|
| 333 | deallocate(this%vegetation_lw_emissivity) |
|---|
| 334 | end if |
|---|
| 335 | |
|---|
| 336 | this%ntile = 0 |
|---|
| 337 | this%ncol = 0 |
|---|
| 338 | this%nfacet = 0 |
|---|
| 339 | |
|---|
| 340 | if (lhook) call dr_hook('radiation_surface:deallocate',1,hook_handle) |
|---|
| 341 | |
|---|
| 342 | end subroutine deallocate_surface |
|---|
| 343 | |
|---|
| 344 | |
|---|
| 345 | !--------------------------------------------------------------------- |
|---|
| 346 | ! Print a description of the surface tile types |
|---|
| 347 | subroutine print_surface_representation(i_representation) |
|---|
| 348 | |
|---|
| 349 | use radiation_io, only : nulout |
|---|
| 350 | |
|---|
| 351 | integer(kind=jpim), dimension(:), allocatable, intent(in) :: i_representation |
|---|
| 352 | |
|---|
| 353 | integer :: ntile, jtile, ifacet, iregion |
|---|
| 354 | |
|---|
| 355 | write(nulout,'(a)') 'Surface tile representation:' |
|---|
| 356 | if (.not. allocated(i_representation)) then |
|---|
| 357 | write(nulout,'(a)') ' Simple (one flat tile)' |
|---|
| 358 | else |
|---|
| 359 | ntile = size(i_representation,1) |
|---|
| 360 | ifacet = 1 |
|---|
| 361 | iregion = 1; |
|---|
| 362 | do jtile = 1,ntile |
|---|
| 363 | write(nulout,'(a,i0,a,a)') ' Tile ', jtile, ': ', trim(TileRepresentationName(i_representation(jtile))) |
|---|
| 364 | write(nulout,'(a,i0,a)') ' Facet ', ifacet, ': ground' |
|---|
| 365 | ifacet = ifacet+1 |
|---|
| 366 | |
|---|
| 367 | select case (i_representation(jtile)) |
|---|
| 368 | |
|---|
| 369 | case (ITileVegetation) |
|---|
| 370 | write(nulout,'(a,i0,a)') ' Region ', iregion, ': vegetation canopy' |
|---|
| 371 | iregion = iregion+1 |
|---|
| 372 | |
|---|
| 373 | case (ITileUrban3D) |
|---|
| 374 | write(nulout,'(a,i0,a)') ' Facet ', ifacet, ': roof' |
|---|
| 375 | write(nulout,'(a,i0,a)') ' Facet ', ifacet+1, ': wall' |
|---|
| 376 | ifacet = ifacet+1 |
|---|
| 377 | write(nulout,'(a,i0,a)') ' Region ', iregion, ': street canyon' |
|---|
| 378 | iregion = iregion+1 |
|---|
| 379 | |
|---|
| 380 | end select |
|---|
| 381 | |
|---|
| 382 | end do |
|---|
| 383 | end if |
|---|
| 384 | |
|---|
| 385 | end subroutine print_surface_representation |
|---|
| 386 | |
|---|
| 387 | |
|---|
| 388 | !--------------------------------------------------------------------- |
|---|
| 389 | subroutine read_from_netcdf(this, file) |
|---|
| 390 | |
|---|
| 391 | use parkind1, only : jprb, jpim |
|---|
| 392 | use easy_netcdf, only : netcdf_file |
|---|
| 393 | |
|---|
| 394 | implicit none |
|---|
| 395 | |
|---|
| 396 | type(netcdf_file), intent(in) :: file |
|---|
| 397 | class(surface_type), intent(inout) :: this |
|---|
| 398 | |
|---|
| 399 | real(kind=jprb), allocatable :: data_1d(:) |
|---|
| 400 | |
|---|
| 401 | integer(kind=jpim), parameter :: ipermute(3) = [2,3,1] |
|---|
| 402 | |
|---|
| 403 | call this%deallocate |
|---|
| 404 | |
|---|
| 405 | this%is_simple = .false. |
|---|
| 406 | |
|---|
| 407 | call file%get('skin_temperature', this%skin_temperature, do_transp=.true.) |
|---|
| 408 | call file%get('canopy_temperature', this%canopy_temperature, do_transp=.true.) |
|---|
| 409 | call file%get('sw_albedo', this%sw_albedo, ipermute=ipermute) |
|---|
| 410 | if (file%exists('sw_albedo_direct')) then |
|---|
| 411 | call file%get('sw_albedo', this%sw_albedo_direct, ipermute=ipermute) |
|---|
| 412 | end if |
|---|
| 413 | call file%get('lw_emissivity', this%lw_emissivity, ipermute=ipermute) |
|---|
| 414 | call file%get('tile_representation', data_1d) |
|---|
| 415 | this%ntile = size(data_1d) |
|---|
| 416 | allocate(this%i_representation(this%ntile)) |
|---|
| 417 | this%i_representation = int(data_1d) |
|---|
| 418 | call file%get('tile_fraction', this%tile_fraction, do_transp=.true.) |
|---|
| 419 | call file%get('canopy_depth', this%canopy_depth, do_transp=.true.) |
|---|
| 420 | call file%get('building_fraction', this%building_fraction, do_transp=.true.) |
|---|
| 421 | if (file%exists('building_normalized_perimeter')) then |
|---|
| 422 | call file%get('building_normalized_perimeter', & |
|---|
| 423 | & this%building_normalized_perimeter, do_transp=.true.) |
|---|
| 424 | else |
|---|
| 425 | ! Convert building scale to normalized perimeter |
|---|
| 426 | call file%get('building_scale', this%building_normalized_perimeter, do_transp=.true.) |
|---|
| 427 | this%building_normalized_perimeter = 4.0_jprb * this%building_fraction & |
|---|
| 428 | & * (1.0_jprb-this%building_fraction) / max(1.0e-8_jprb,this%building_normalized_perimeter) |
|---|
| 429 | end if |
|---|
| 430 | call file%get('vegetation_optical_depth', this%vegetation_optical_depth, do_transp=.true.) |
|---|
| 431 | call file%get('vegetation_sw_albedo', this%vegetation_sw_albedo, ipermute=ipermute) |
|---|
| 432 | call file%get('vegetation_lw_emissivity', this%vegetation_lw_emissivity, ipermute=ipermute) |
|---|
| 433 | |
|---|
| 434 | this%nregion = sum(NTileRegions(this%i_representation)) |
|---|
| 435 | this%ncol = size(this%skin_temperature,1) |
|---|
| 436 | this%nfacet = size(this%skin_temperature,2) |
|---|
| 437 | |
|---|
| 438 | this%nalbedobands = size(this%sw_albedo,2) |
|---|
| 439 | this%nemissbands = size(this%lw_emissivity,2) |
|---|
| 440 | |
|---|
| 441 | call this%set_facet_indices |
|---|
| 442 | |
|---|
| 443 | end subroutine read_from_netcdf |
|---|
| 444 | |
|---|
| 445 | |
|---|
| 446 | end module radsurf_properties |
|---|