source: palm/trunk/SOURCE/urban_surface_mod.f90 @ 2012

Last change on this file since 2012 was 2012, checked in by kanani, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 199.3 KB
Line 
1!> @file urban_surface_mod.f90
2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 2015-2016 Czech Technical University in Prague
18! Copyright 1997-2016 Leibniz Universitaet Hannover
19!--------------------------------------------------------------------------------!
20!
21! Current revisions:
22! ------------------
23!
24!
25! Former revisions:
26! -----------------
27! $Id: urban_surface_mod.f90 2012 2016-09-19 17:31:38Z kanani $
28!
29! 2011 2016-09-19 17:29:57Z kanani
30! Major reformatting according to PALM coding standard (comments, blanks,
31! alphabetical ordering, etc.),
32! removed debug_prints,
33! removed auxiliary SUBROUTINE get_usm_info, instead, USM flag urban_surface is
34! defined in MODULE control_parameters (modules.f90) to avoid circular
35! dependencies,
36! renamed canopy_heat_flux to pc_heating_rate, as meaning of quantity changed.
37!
38! 2007 2016-08-24 15:47:17Z kanani
39! Initial revision
40!
41!
42! Description:
43! ------------
44! 2016/6/9 - Initial version of the USM (Urban Surface Model)
45!            authors: Jaroslav Resler, Pavel Krc (CTU in Prague, ICS AS in Prague)
46!            with contributions: Michal Belda, Nina Benesova, Ondrej Vlcek
47!            partly inspired by PALM LSM (B. Maronga)
48!            parameterizations of Ra checked with TUF3D (E. S. Krayenhoff)
49!> Module for Urban Surface Model (USM)
50!> The module includes:
51!>    1. radiation model with direct/diffuse radiation, shading, reflections
52!>       and integration with plant canopy
53!>    2. wall and wall surface model
54!>    3. surface layer energy balance
55!>    4. anthropogenic heat (only from transportation so far)
56!>    5. necessary auxiliary subroutines (reading inputs, writing outputs,
57!>       restart simulations, ...)
58!> It also make use of standard radiation and integrates it into
59!> urban surface model.
60!>
61!> Further work:
62!> -------------
63!> 1. Reduce number of shape view factors by merging factors for distant surfaces
64!>    under shallow angles. Idea: Iteratively select the smallest shape view
65!>    factor by value (among all sources and targets) which has a similarly
66!>    oriented source neighbor (or near enough) SVF and merge them by adding
67!>    value of the smaller SVF to the larger one and deleting the smaller one.
68!>    This will allow for better scaling at higher resolutions.
69!>
70!> 2. Remove global arrays surfouts, surfoutl and only keep track of radiosity
71!>    from surfaces that are visible from local surfaces (i.e. there is a SVF
72!>    where target is local). To do that, radiosity will be exchanged after each
73!>    reflection step using MPI_Alltoall instead of current MPI_Allgather.
74!>
75!------------------------------------------------------------------------------!
76 MODULE urban_surface_mod
77
78    USE arrays_3d,                                                             &
79        ONLY:  zu, pt, pt_1, pt_2, p, ol, shf, ts, us, u, v, w, hyp, tend
80
81    USE cloud_parameters,                                                      &
82        ONLY:  cp, r_d
83
84    USE constants,                                                             &
85        ONLY:  pi
86   
87    USE control_parameters,                                                    &
88        ONLY:  dz, topography, dt_3d, intermediate_timestep_count,             &
89               initializing_actions, intermediate_timestep_count_max,          &
90               simulated_time, end_time, timestep_scheme, tsc,                 &
91               coupling_char, io_blocks, io_group, message_string,             &
92               time_since_reference_point, surface_pressure,                   &
93               g, pt_surface, large_scale_forcing, lsf_surf,                   &
94               time_do3d, dt_do3d, average_count_3d, urban_surface
95
96    USE cpulog,                                                                &
97        ONLY:  cpu_log, log_point, log_point_s
98     
99    USE grid_variables,                                                        &
100        ONLY:  dx, dy, ddx, ddy, ddx2, ddy2
101   
102    USE indices,                                                               &
103        ONLY:  nx, ny, nnx, nny, nnz, nxl, nxlg, nxr, nxrg, nyn, nyng, nys,    &
104               nysg, nzb_s_inner, nzb_s_outer, nzb, nzt, nbgp
105
106    USE, INTRINSIC :: iso_c_binding 
107
108    USE kinds
109             
110    USE pegrid
111   
112    USE plant_canopy_model_mod,                                                &
113        ONLY:  plant_canopy, pch_index,                                        &
114               pc_heating_rate, lad_s
115   
116    USE radiation_model_mod,                                                   &
117        ONLY:  radiation, calc_zenith, zenith, day_init, time_utc_init,        &
118               rad_net, rad_sw_in, rad_lw_in, rad_sw_out, rad_lw_out,          &
119               sigma_sb, sun_direction, sun_dir_lat, sun_dir_lon,              &
120               force_radiation_call
121
122    USE statistics,                                                            &
123        ONLY:  hom, statistic_regions
124
125               
126
127    IMPLICIT NONE
128
129!-- configuration parameters (they can be setup in PALM config)
130    LOGICAL                                        ::  split_diffusion_radiation = .TRUE. !< split direct and diffusion dw radiation
131                                                                                          !< (.F. in case the radiation model already does it)   
132    LOGICAL                                        ::  usm_energy_balance_land = .TRUE.   !< flag parameter indicating wheather the energy balance is calculated for land and roofs
133    LOGICAL                                        ::  usm_energy_balance_wall = .TRUE.   !< flag parameter indicating wheather the energy balance is calculated for land and roofs
134    LOGICAL                                        ::  usm_material_model = .TRUE.        !< flag parameter indicating wheather the  model of heat in materials is used
135    LOGICAL                                        ::  usm_anthropogenic_heat = .FALSE.   !< flag parameter indicating wheather the anthropogenic heat sources (e.g.transportation) are used
136    LOGICAL                                        ::  force_radiation_call_l = .FALSE.   !< flag parameter for unscheduled radiation model calls
137    LOGICAL                                        ::  mrt_factors = .FALSE.              !< whether to generate MRT factor files during init
138    LOGICAL                                        ::  write_svf_on_init = .FALSE.
139    LOGICAL                                        ::  read_svf_on_init = .FALSE.
140    LOGICAL                                        ::  usm_lad_rma = .TRUE.               !< use MPI RMA to access LAD for raytracing (instead of global array)
141   
142    INTEGER(iwp)                                   ::  nrefsteps = 0                      !< number of reflection steps to perform
143   
144    INTEGER(iwp)                                   ::  land_category = 2                  !< default category for land surface
145    INTEGER(iwp)                                   ::  wall_category = 2                  !< default category for wall surface over pedestrian zone
146    INTEGER(iwp)                                   ::  pedestrant_category = 2            !< default category for wall surface in pedestrian zone
147    INTEGER(iwp)                                   ::  roof_category = 2                  !< default category for root surface
148    REAL(wp)                                       ::  roof_height_limit = 4._wp          !< height for distinguish between land surfaces and roofs
149
150    REAL(wp), PARAMETER                            ::  ext_coef = 0.6_wp                  !< extinction coefficient (a.k.a. alpha)
151    REAL(wp)                                       ::  ra_horiz_coef = 5.0_wp             !< mysterious coefficient for correction of overestimation
152                                                                                          !< of r_a for horizontal surfaces -> TODO
153   
154!-- parameters of urban surface model
155    INTEGER(iwp), PARAMETER                        ::  usm_version_len = 10               !< length of identification string of usm version
156    CHARACTER(usm_version_len), PARAMETER          ::  usm_version = 'USM v. 1.0'         !< identification of version of binary svf and restart files
157    INTEGER(iwp), PARAMETER                        ::  svf_code_len = 15                  !< length of code for verification of the end of svf file
158    CHARACTER(svf_code_len), PARAMETER             ::  svf_code = '*** end svf ***'       !< code for verification of the end of svf file
159    INTEGER(iwp)                                   ::  nzu                                !< number of layers of urban surface (will be calculated)
160    INTEGER(iwp)                                   ::  nzub,nzut                          !< bottom and top layer of urban surface (will be calculated)
161    INTEGER(iwp), PARAMETER                        ::  nzut_free = 3                      !< number of free layers in urban surface layer above top of buildings
162    INTEGER(iwp), PARAMETER                        ::  ndsvf = 2                          !< number of dimensions of real values in SVF
163    INTEGER(iwp), PARAMETER                        ::  ndcsf = 2                          !< number of dimensions of real values in CSF
164    INTEGER(iwp), PARAMETER                        ::  kdcsf = 4                          !< number of dimensions of integer values in CSF
165    INTEGER(iwp), PARAMETER                        ::  id = 1                             !< position of d-index in surfl and surf
166    INTEGER(iwp), PARAMETER                        ::  iz = 2                             !< position of k-index in surfl and surf
167    INTEGER(iwp), PARAMETER                        ::  iy = 3                             !< position of j-index in surfl and surf
168    INTEGER(iwp), PARAMETER                        ::  ix = 4                             !< position of i-index in surfl and surf
169    INTEGER(iwp), PARAMETER                        ::  iroof = 0                          !< 0 - index of ground or roof
170    INTEGER(iwp), PARAMETER                        ::  isouth = 1                         !< 1 - index of south facing wall
171    INTEGER(iwp), PARAMETER                        ::  inorth = 2                         !< 2 - index of north facing wall
172    INTEGER(iwp), PARAMETER                        ::  iwest  = 3                         !< 3 - index of west facing wall
173    INTEGER(iwp), PARAMETER                        ::  ieast  = 4                         !< 4 - index of east facing wall
174    INTEGER(iwp), PARAMETER                        ::  isky = 5                           !< 5 - index of top border of the urban surface layer ("urban sky")
175    INTEGER(iwp), PARAMETER                        ::  inorthb = 6                        !< 6 - index of free north border of the domain (south facing)
176    INTEGER(iwp), PARAMETER                        ::  isouthb = 7                        !< 7 - index of north south border of the domain (north facing)
177    INTEGER(iwp), PARAMETER                        ::  ieastb  = 8                        !< 8 - index of east border of the domain (west facing)
178    INTEGER(iwp), PARAMETER                        ::  iwestb  = 9                        !< 9 - index of wast border of the domain (east facing)
179    INTEGER(iwp), DIMENSION(0:9), PARAMETER        ::  idir = (/0,0,0,-1,1,0,0,0,-1,1/)   !< surface normal direction x indices
180    INTEGER(iwp), DIMENSION(0:9), PARAMETER        ::  jdir = (/0,-1,1,0,0,0,-1,1,0,0/)   !< surface normal direction y indices
181    INTEGER(iwp), DIMENSION(0:9), PARAMETER        ::  kdir = (/1,0,0,0,0,-1,0,0,0,0/)    !< surface normal direction z indices
182    REAL(wp), DIMENSION(1:4)                       ::  ddxy2                              !< 1/dx^2 or 1/dy^2 (in surface normal direction)
183    INTEGER(iwp), DIMENSION(1:4,6:9)               ::  ijdb                               !< start and end of the local domain border coordinates (set in code)
184    LOGICAL, DIMENSION(6:9)                        ::  isborder                           !< is PE on the border of the domain in four corresponding directions
185                                                                                          !< parameter but set in the code
186
187!-- indices and sizes of urban surface model
188    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE      ::  surfl            !< coordinates of i-th local surface in local grid - surfl[:,k] = [d, z, y, x]
189    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE      ::  surf             !< coordinates of i-th surface in grid - surf[:,k] = [d, z, y, x]
190    INTEGER(iwp)                                   ::  nsurfl           !< number of all surfaces in local processor
191    INTEGER(iwp), DIMENSION(:), ALLOCATABLE        ::  nsurfs           !< array of number of all surfaces in individual processors
192    INTEGER(iwp)                                   ::  startsky         !< start index of block of sky
193    INTEGER(iwp)                                   ::  endsky           !< end index of block of sky
194    INTEGER(iwp)                                   ::  nskys            !< number of sky surfaces in local processor
195    INTEGER(iwp)                                   ::  startland        !< start index of block of land and roof surfaces
196    INTEGER(iwp)                                   ::  endland          !< end index of block of land and roof surfaces
197    INTEGER(iwp)                                   ::  nlands           !< number of land and roof surfaces in local processor
198    INTEGER(iwp)                                   ::  startwall        !< start index of block of wall surfaces
199    INTEGER(iwp)                                   ::  endwall          !< end index of block of wall surfaces
200    INTEGER(iwp)                                   ::  nwalls           !< number of wall surfaces in local processor
201    INTEGER(iwp)                                   ::  startenergy      !< start index of block of real surfaces (land, walls and roofs)
202    INTEGER(iwp)                                   ::  endenergy        !< end index of block of real surfaces (land, walls and roofs)
203    INTEGER(iwp)                                   ::  nenergy          !< number of real surfaces in local processor
204    INTEGER(iwp)                                   ::  nsurf            !< global number of surfaces in index array of surfaces (nsurf = Σproc nsurfs)
205    INTEGER(iwp), DIMENSION(:), ALLOCATABLE        ::  surfstart        !< starts of blocks of surfaces for individual processors in array surf
206                                                                        !< respective block for particular processor is surfstart[iproc]+1 : surfstart[iproc+1]
207    INTEGER(iwp)                                   ::  nsvfl            !< number of svf (excluding csf) for local processor
208    INTEGER(iwp)                                   ::  ncsfl            !< no. of csf in local processor
209                                                                        !< needed only during calc_svf but must be here because it is
210                                                                        !< shared between subroutines usm_calc_svf and usm_raytrace                                                                         
211    INTEGER(iwp)                                   ::  nsvfcsfl         !< sum of svf+csf for local processor
212
213
214!-- type for calculation of svf
215    TYPE t_svf
216        INTEGER(iwp)                               :: isurflt           !<
217        INTEGER(iwp)                               :: isurfs            !<
218        REAL(wp)                                   :: rsvf              !<
219        REAL(wp)                                   :: rtransp           !<
220    END TYPE
221
222!-- type for calculation of csf
223    TYPE t_csf
224        INTEGER(iwp)                               :: ip                !<
225        INTEGER(iwp)                               :: itx               !<
226        INTEGER(iwp)                               :: ity               !<
227        INTEGER(iwp)                               :: itz               !<
228        INTEGER(iwp)                               :: isurfs            !<
229        REAL(wp)                                   :: rsvf              !<
230        REAL(wp)                                   :: rtransp           !<
231    END TYPE
232
233!-- arrays for calculation of svf and csf
234    TYPE(t_svf), DIMENSION(:), POINTER             ::  asvf             !< pointer to growing svc array
235    TYPE(t_csf), DIMENSION(:), POINTER             ::  acsf             !< pointer to growing csf array
236    TYPE(t_svf), DIMENSION(:), ALLOCATABLE, TARGET ::  asvf1, asvf2     !< realizations of svf array
237    TYPE(t_csf), DIMENSION(:), ALLOCATABLE, TARGET ::  acsf1, acsf2     !< realizations of csf array
238    INTEGER(iwp)                                   ::  nsvfla           !< dimmension of array allocated for storage of svf in local processor
239    INTEGER(iwp)                                   ::  ncsfla           !< dimmension of array allocated for storage of csf in local processor
240    INTEGER(iwp)                                   ::  msvf, mcsf       !< mod for swapping the growing array
241    INTEGER(iwp), PARAMETER                        ::  gasize = 10000   !< initial size of growing arrays
242
243!-- arrays storing the values of USM
244    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE      ::  svfsurf          !< svfsurf[:,isvf] = index of source and target surface for svf[isvf]
245    REAL(wp), DIMENSION(:,:), ALLOCATABLE          ::  svf              !< array of shape view factors+direct irradiation factors
246                                                                        !< for individual local surfaces and plant canopy sinks
247    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfins          !< array of sw radiation falling to local surface after i-th reflection
248    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinl          !< array of lw radiation for local surface after i-th reflection
249   
250                                                                        !< Inward radiation is also valid for virtual surfaces (radiation leaving domain)
251    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinsw         !< array of sw radiation falling to local surface including radiation from reflections
252    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinlw         !< array of lw radiation falling to local surface including radiation from reflections
253    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinswdir      !< array of direct sw radiation falling to local surface
254    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinswdif      !< array of diffuse sw radiation from sky and model boundary falling to local surface
255    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinlwdif      !< array of diffuse lw radiation from sky and model boundary falling to local surface
256   
257                                                                        !< Outward radiation is only valid for nonvirtual surfaces
258    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfoutsl        !< array of reflected sw radiation for local surface in i-th reflection
259    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfoutll        !< array of reflected + emitted lw radiation for local surface in i-th reflection
260    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfouts         !< array of reflected sw radiation for all surfaces in i-th reflection
261    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfoutl         !< array of reflected + emitted lw radiation for all surfaces in i-th reflection
262    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfoutsw        !< array of total sw radiation outgoing from nonvirtual surfaces surfaces after all reflection
263    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfoutlw        !< array of total lw radiation outgoing from nonvirtual surfaces surfaces after all reflection
264    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfhf           !< array of total radiation flux incoming to minus outgoing from local surface
265    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  rad_net_l        !< local copy of rad_net (net radiation at surface)
266
267!-- arrays for time averages
268    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  rad_net_av       !< average of rad_net_l
269    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinsw_av      !< average of sw radiation falling to local surface including radiation from reflections
270    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinlw_av      !< average of lw radiation falling to local surface including radiation from reflections
271    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinswdir_av   !< average of direct sw radiation falling to local surface
272    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinswdif_av   !< average of diffuse sw radiation from sky and model boundary falling to local surface
273    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinlwdif_av   !< average of diffuse lw radiation from sky and model boundary falling to local surface
274    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinswref_av   !< average of sw radiation falling to surface from reflections
275    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfinlwref_av   !< average of lw radiation falling to surface from reflections
276    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfoutsw_av     !< average of total sw radiation outgoing from nonvirtual surfaces surfaces after all reflection
277    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfoutlw_av     !< average of total lw radiation outgoing from nonvirtual surfaces surfaces after all reflection
278    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  surfhf_av        !< average of total radiation flux incoming to minus outgoing from local surface   
279   
280!-- block variables needed for calculation of the plant canopy model inside the urban surface model
281    REAL(wp), DIMENSION(:,:), ALLOCATABLE          ::  pcbl             !< z,y,x coordinates of i-th local plant canopy box pcbl[:,i] = [z, y, x]
282    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE    ::  gridpcbl         !< index of local pcb[z,y,x]
283    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  pcbinsw          !< array of absorbed sw radiation for local plant canopy box
284    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  pcbinlw          !< array of absorbed lw radiation for local plant canopy box
285    INTEGER(iwp)                                   ::  npcbl            !< number of the plant canopy gridboxes in local processor
286    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE      ::  pch              !< heights of the plant canopy
287    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE      ::  pct              !< top layer of the plant canopy
288    REAL(wp), DIMENSION(:,:,:), POINTER            ::  usm_lad          !< subset of lad_s within urban surface, transformed to plain Z coordinate
289    REAL(wp), DIMENSION(:), POINTER                ::  usm_lad_g        !< usm_lad globalized (used to avoid MPI RMA calls in raytracing)
290    REAL(wp)                                       ::  prototype_lad    !< prototype leaf area density for computing effective optical depth
291    INTEGER(iwp), DIMENSION(:), ALLOCATABLE        ::  nzterr, plantt   !< temporary global arrays for raytracing
292   
293!-- radiation related arrays (it should be better in interface of radiation module of PALM
294    REAL(wp), DIMENSION(:,:), ALLOCATABLE          ::  rad_sw_in_dir    !< direct sw radiation
295    REAL(wp), DIMENSION(:,:), ALLOCATABLE          ::  rad_sw_in_diff   !< diffusion sw radiation
296    REAL(wp), DIMENSION(:,:), ALLOCATABLE          ::  rad_lw_in_diff   !< diffusion lw radiation
297
298!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
299!-- anthropogenic heat sources
300!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
301    REAL(wp), DIMENSION(:,:), ALLOCATABLE          ::  aheat             !< daily average of anthropogenic heat (W/m2)
302    REAL(wp), DIMENSION(:), ALLOCATABLE            ::  aheatprof         !< diurnal profile of anthropogenic heat
303
304!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
305!-- wall surface model
306!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
307!-- wall surface model constants
308    INTEGER(iwp), PARAMETER                        :: nzb_wall = 0       !< inner side of the wall model (to be switched)
309    INTEGER(iwp), PARAMETER                        :: nzt_wall = 3       !< outer side of the wall model (to be switched)
310    INTEGER(iwp), PARAMETER                        :: nzw = 4            !< number of wall layers (fixed for now)
311
312    REAL(wp), DIMENSION(nzb_wall:nzt_wall)         :: zwn_default = (/0.0242_wp, 0.0969_wp, 0.346_wp, 1.0_wp /)
313                                                                         !< normalized soil, wall and roof layer depths (m/m)
314                                                                       
315    REAL(wp)                                       ::   wall_inner_temperature = 296.0_wp    !< temperature of the inner wall surface (~23 degrees C) (K)
316    REAL(wp)                                       ::   roof_inner_temperature = 296.0_wp    !< temperature of the inner roof surface (~23 degrees C) (K)
317    REAL(wp)                                       ::   soil_inner_temperature = 283.0_wp    !< temperature of the deep soil (~10 degrees C) (K)
318
319!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
320!-- surface and material model variables for walls, ground, roofs
321!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
322    INTEGER(iwp), DIMENSION(:), ALLOCATABLE        :: surface_types      !< array of types of wall parameters
323    REAL(wp), DIMENSION(:), ALLOCATABLE            :: zwn                !< normalized wall layer depths (m)
324    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: ddz_wall           !< 1/dz_wall
325    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: ddz_wall_stag      !< 1/dz_wall_stag
326    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: dz_wall            !< wall grid spacing (center-center)
327    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: dz_wall_stag       !< wall grid spacing (edge-edge)
328    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: zw                 !< wall layer depths (m)
329
330#if defined( __nopointer )
331    REAL(wp), DIMENSION(:), ALLOCATABLE, TARGET    :: t_surf             !< wall surface temperature (K)
332    REAL(wp), DIMENSION(:), ALLOCATABLE, TARGET    :: t_surf_p           !< progn. wall surface temperature (K)
333#else
334    REAL(wp), DIMENSION(:), POINTER                :: t_surf
335    REAL(wp), DIMENSION(:), POINTER                :: t_surf_p 
336
337    REAL(wp), DIMENSION(:), ALLOCATABLE, TARGET    :: t_surf_1
338    REAL(wp), DIMENSION(:), ALLOCATABLE, TARGET    :: t_surf_2
339#endif
340    REAL(wp), DIMENSION(:), ALLOCATABLE, TARGET    :: t_surf_av          !< average of wall surface temperature (K)
341
342!-- Temporal tendencies for time stepping           
343    REAL(wp), DIMENSION(:), ALLOCATABLE            :: tt_surface_m       !< surface temperature tendency (K)
344
345!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
346!-- Energy balance variables
347!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
348!-- parameters of the land, roof and wall surfaces
349    LOGICAL,  DIMENSION(:), ALLOCATABLE            :: isroof_surf        !< is the surface the part of a roof
350    REAL(wp), DIMENSION(:), ALLOCATABLE            :: albedo_surf        !< albedo of the surface
351!-- parameters of the wall surfaces
352    REAL(wp), DIMENSION(:), ALLOCATABLE            :: c_surface          !< heat capacity of the wall surface skin ( J m−2 K−1 )
353    REAL(wp), DIMENSION(:), ALLOCATABLE            :: emiss_surf         !< emissivity of the wall surface
354    REAL(wp), DIMENSION(:), ALLOCATABLE            :: lambda_surf        !< heat conductivity λS between air and surface ( W m−2 K−1 )
355   
356!-- parameters of the walls material
357    REAL(wp), DIMENSION(:), ALLOCATABLE            :: thickness_wall     !< thickness of the wall, roof and soil layers
358    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: rho_c_wall         !< volumetric heat capacity of the material ( J m-3 K-1 ) (= 2.19E6)
359    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: lambda_h           !< heat conductivity λT of the material ( W m-1 K-1 )
360    REAL(wp), DIMENSION(:), ALLOCATABLE            :: roughness_wall     !< roughness relative to concrete
361   
362!-- output wall heat flux arrays
363    REAL(wp), DIMENSION(:), ALLOCATABLE            :: wshf               !< kinematic wall heat flux of sensible heat (needed for diffusion_s!<)
364    REAL(wp), DIMENSION(:), ALLOCATABLE            :: wshf_eb            !< wall heat flux of sensible heat in wall normal direction
365    REAL(wp), DIMENSION(:), ALLOCATABLE            :: wshf_eb_av         !< average of wshf_eb
366    REAL(wp), DIMENSION(:), ALLOCATABLE            :: wghf_eb            !< wall ground heat flux
367    REAL(wp), DIMENSION(:), ALLOCATABLE            :: wghf_eb_av         !< average of wghf_eb
368
369#if defined( __nopointer )
370    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET  :: t_wall             !< Wall temperature (K)
371    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET  :: t_wall_av          !< Average of t_wall
372    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET  :: t_wall_p           !< Prog. wall temperature (K)
373#else
374    REAL(wp), DIMENSION(:,:), POINTER              :: t_wall, t_wall_p
375    REAL(wp), DIMENSION(:,:), ALLOCATABLE, TARGET  :: t_wall_av, t_wall_1, t_wall_2
376#endif
377
378!-- Wall temporal tendencies for time stepping
379    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: tt_wall_m          !< t_wall prognostic array
380
381!-- Surface and material parameters classes (surface_type)
382!-- albedo, emissivity, lambda_surf, roughness, thickness, volumetric heat capacity, thermal conductivity
383    INTEGER(iwp)                                   :: n_surface_types      !< number of the wall type categories
384    INTEGER(iwp), PARAMETER                        :: n_surface_params = 8 !< number of parameters for each type of the wall
385    INTEGER(iwp), PARAMETER                        :: ialbedo  = 1         !< albedo of the surface
386    INTEGER(iwp), PARAMETER                        :: iemiss   = 2         !< emissivity of the surface
387    INTEGER(iwp), PARAMETER                        :: ilambdas = 3         !< heat conductivity λS between air and surface ( W m−2 K−1 )
388    INTEGER(iwp), PARAMETER                        :: irough   = 4         !< roughness relative to concrete
389    INTEGER(iwp), PARAMETER                        :: icsurf   = 5         !< Surface skin layer heat capacity (J m−2 K−1 )
390    INTEGER(iwp), PARAMETER                        :: ithick   = 6         !< thickness of the surface (wall, roof, land)  ( m )
391    INTEGER(iwp), PARAMETER                        :: irhoC    = 7         !< volumetric heat capacity rho*C of the material ( J m−3 K−1 )
392    INTEGER(iwp), PARAMETER                        :: ilambdah = 8         !< thermal conductivity λH of the wall (W m−1 K−1 )
393    CHARACTER(12), DIMENSION(:), ALLOCATABLE       :: surface_type_names   !< names of wall types (used only for reports)
394    INTEGER(iwp), DIMENSION(:), ALLOCATABLE        :: surface_type_codes   !< codes of wall types
395    REAL(wp), DIMENSION(:,:), ALLOCATABLE          :: surface_params       !< parameters of wall types
396   
397    CHARACTER(len=*), PARAMETER                    :: svf_file_name='usm_svf'
398   
399!-- interfaces of subroutines accessed from outside of this module
400    INTERFACE usm_check_data_output
401       MODULE PROCEDURE usm_check_data_output
402    END INTERFACE usm_check_data_output
403   
404    INTERFACE usm_check_parameters
405       MODULE PROCEDURE usm_check_parameters
406    END INTERFACE usm_check_parameters
407   
408    INTERFACE usm_data_output_3d
409       MODULE PROCEDURE usm_data_output_3d
410    END INTERFACE usm_data_output_3d
411   
412    INTERFACE usm_define_netcdf_grid
413       MODULE PROCEDURE usm_define_netcdf_grid
414    END INTERFACE usm_define_netcdf_grid
415
416    INTERFACE usm_init_urban_surface
417       MODULE PROCEDURE usm_init_urban_surface
418    END INTERFACE usm_init_urban_surface
419
420    INTERFACE usm_material_heat_model
421       MODULE PROCEDURE usm_material_heat_model
422    END INTERFACE usm_material_heat_model
423   
424    INTERFACE usm_parin
425       MODULE PROCEDURE usm_parin
426    END INTERFACE usm_parin
427
428    INTERFACE usm_radiation
429       MODULE PROCEDURE usm_radiation
430    END INTERFACE usm_radiation
431   
432    INTERFACE usm_read_restart_data 
433       MODULE PROCEDURE usm_read_restart_data
434    END INTERFACE usm_read_restart_data
435
436    INTERFACE usm_surface_energy_balance
437       MODULE PROCEDURE usm_surface_energy_balance
438    END INTERFACE usm_surface_energy_balance
439   
440    INTERFACE usm_swap_timelevel
441       MODULE PROCEDURE usm_swap_timelevel
442    END INTERFACE usm_swap_timelevel
443   
444    INTERFACE usm_wall_heat_flux
445       MODULE PROCEDURE usm_wall_heat_flux
446       MODULE PROCEDURE usm_wall_heat_flux_ij
447    END INTERFACE usm_wall_heat_flux
448   
449    INTERFACE usm_write_restart_data
450       MODULE PROCEDURE usm_write_restart_data
451    END INTERFACE usm_write_restart_data
452   
453    SAVE
454
455    PRIVATE
456   
457!-- Public parameters, constants and initial values
458    PUBLIC split_diffusion_radiation,                                          &
459           usm_anthropogenic_heat, usm_material_model, mrt_factors,            &
460           usm_check_parameters,                                               &
461           usm_energy_balance_land, usm_energy_balance_wall, nrefsteps,        &
462           usm_init_urban_surface, usm_radiation, usm_read_restart_data,       &
463           usm_wall_heat_flux,                                                 &
464           usm_surface_energy_balance, usm_material_heat_model,                &
465           usm_swap_timelevel, usm_check_data_output, usm_average_3d_data,     &
466           usm_data_output_3d, usm_define_netcdf_grid, usm_parin,              &
467           usm_write_restart_data,                                             &
468           nzub, nzut, ra_horiz_coef, usm_lad_rma,                             &
469           land_category, pedestrant_category, wall_category, roof_category,   &
470           write_svf_on_init, read_svf_on_init
471
472
473 CONTAINS
474
475 
476!------------------------------------------------------------------------------!
477! Description:
478! ------------
479!> This subroutine creates the necessary indices of the urban surfaces
480!> and plant canopy and it allocates the needed arrays for USM
481!------------------------------------------------------------------------------!
482    SUBROUTINE usm_allocate_urban_surface
483   
484        IMPLICIT NONE
485       
486        INTEGER(iwp)                            :: i, j, k, d, l, ir, jr, ids
487        INTEGER(iwp)                            :: nzubl, nzutl, isurf, ipcgb
488        INTEGER                                 :: procid
489
490       
491
492       
493!--     auxiliary vars
494        ddxy2 = (/ddy2,ddy2,ddx2,ddx2/)      !< 1/dx^2 or 1/dy^2 (in surface normal direction)
495       
496        CALL location_message( '', .TRUE. )
497        CALL location_message( '    allocation of needed arrays', .TRUE. )
498!--     find nzub, nzut, nzu
499        nzubl = minval(nzb_s_inner(nys:nyn,nxl:nxr))
500        nzutl = maxval(nzb_s_inner(nys:nyn,nxl:nxr))
501        nzubl = max(nzubl,nzb)
502       
503        IF ( plant_canopy )  THEN
504!--         allocate needed arrays
505            ALLOCATE( pct(nys:nyn,nxl:nxr) )
506            ALLOCATE( pch(nys:nyn,nxl:nxr) )
507
508!--         calculate plant canopy height
509            npcbl = 0
510            pct = 0.0_wp
511            pch = 0.0_wp
512            DO i = nxl, nxr
513                DO j = nys, nyn
514                    DO k = nzt+1, 0, -1
515                        IF ( lad_s(k,j,i) /= 0.0_wp )  THEN
516!--                         we are at the top of the pcs
517                            pct(j,i) = k + nzb_s_inner(j,i)
518                            pch(j,i) = k
519                            npcbl = npcbl + pch(j,i)
520                            EXIT
521                        ENDIF
522                    ENDDO
523                ENDDO
524            ENDDO
525           
526            nzutl = max(nzutl, maxval(pct))
527!--         code of plant canopy model uses parameter pch_index
528!--         we need to setup it here to right value
529!--         (pch_index, lad_s and other arrays in PCM are defined flat)
530            pch_index = maxval(pch)
531
532            prototype_lad = maxval(lad_s) * .9_wp  !< better be *1.0 if lad is either 0 or maxval(lad) everywhere
533            IF ( prototype_lad <= 0._wp ) prototype_lad = .3_wp
534            !WRITE(message_string, '(a,f6.3)') 'Precomputing effective box optical ' &
535            !    // 'depth using prototype leaf area density = ', prototype_lad
536            !CALL message('usm_init_urban_surface', 'PA0520', 0, 0, -1, 6, 0)
537        ENDIF
538       
539        nzutl = min(nzutl+nzut_free, nzt)
540                 
541#if defined( __parallel )
542        CALL MPI_AllReduce(nzubl,nzub,1,MPI_INTEGER,MPI_MIN,comm2d,ierr);
543        CALL MPI_AllReduce(nzutl,nzut,1,MPI_INTEGER,MPI_MAX,comm2d,ierr);
544#else
545        nzub = nzubl
546        nzut = nzutl
547#endif
548
549!--     global number of urban layers
550        nzu = nzut - nzub + 1
551       
552!--     allocate urban surfaces grid
553!--     calc number of surfaces in local proc
554        CALL location_message( '    calculation of indices for surfaces', .TRUE. )
555        nsurfl = 0
556!--     calculate land surface and roof
557        startland = nsurfl+1
558        nsurfl = nsurfl+(nxr-nxl+1)*(nyn-nys+1)
559        endland = nsurfl
560        nlands = endland-startland+1
561
562!--     calculation of the walls
563        startwall = nsurfl+1
564        DO i = nxl, nxr
565            DO j = nys, nyn
566!--             test for walls
567!--             (we don't use array flags because it isn't calculated in case of masking_method=.T.)
568                DO ids = 1, 4  !-- four wall directions
569                    jr = min(max(j-jdir(ids),0),ny)
570                    ir = min(max(i-idir(ids),0),nx)
571                    nsurfl = nsurfl + max(0, nzb_s_inner(jr,ir)-nzb_s_inner(j,i))
572                ENDDO
573            ENDDO
574        ENDDO
575        endwall = nsurfl
576        nwalls = endwall-startwall+1
577       
578!--     range of energy balance surfaces
579        nenergy = 0
580        IF ( usm_energy_balance_land )  THEN
581            startenergy = startland
582            nenergy = nenergy + nlands
583        ELSE
584            startenergy = startwall
585        ENDIF
586        IF ( usm_energy_balance_wall )  THEN
587            endenergy = endwall
588            nenergy = nenergy + nwalls
589        ELSE
590            endenergy = endland
591        ENDIF
592
593!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
594!--     block of virtual surfaces
595!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
596!--     calculate sky surfaces
597        startsky = nsurfl+1
598        nsurfl = nsurfl+(nxr-nxl+1)*(nyn-nys+1)
599        endsky = nsurfl
600        nskys = endsky-startsky+1
601       
602!--     border flags
603#if defined( __parallel )
604        isborder = (/ north_border_pe, south_border_pe, right_border_pe, left_border_pe /)
605#else
606        isborder = (/.TRUE.,.TRUE.,.TRUE.,.TRUE./)
607#endif
608!--     fill array of the limits of the local domain borders
609        ijdb = RESHAPE( (/ nxl,nxr,nyn,nyn,nxl,nxr,nys,nys,nxr,nxr,nys,nyn,nxl,nxl,nys,nyn /), (/4, 4/) )
610!--     calulation of the free borders of the domain
611        DO ids = 6,9
612            IF ( isborder(ids) )  THEN
613!--             free border of the domain in direction ids
614                DO i = ijdb(1,ids), ijdb(2,ids)
615                    DO j = ijdb(3,ids), ijdb(4,ids)
616                        k = nzut - max(nzb_s_inner(j,i), nzb_s_inner(j-jdir(ids),i-idir(ids)))
617                        nsurfl = nsurfl + k
618                    ENDDO
619                ENDDO
620            ENDIF
621        ENDDO
622       
623!--     fill gridpcbl and pcbl
624        IF ( plant_canopy )  THEN
625            ALLOCATE( pcbl(iz:ix, 1:npcbl) )
626            ALLOCATE( gridpcbl(nzub:nzut,nys:nyn,nxl:nxr) )
627            gridpcbl(:,:,:) = 0
628            ipcgb = 0
629            DO i = nxl, nxr
630                DO j = nys, nyn
631                    DO k = nzb_s_inner(j,i)+1, pct(j,i)
632                        ipcgb = ipcgb + 1
633                        gridpcbl(k,j,i) = ipcgb
634                        pcbl(:,ipcgb) = (/ k, j, i /)
635                    ENDDO
636                ENDDO
637            ENDDO
638
639            ALLOCATE( pcbinsw( 1:npcbl ) )
640            ALLOCATE( pcbinlw( 1:npcbl ) )
641        ENDIF
642
643!--     fill surfl
644        ALLOCATE(surfl(4,nsurfl))
645        isurf = 0
646       
647!--     add land surfaces or roofs
648        DO i = nxl, nxr
649            DO j = nys, nyn
650                isurf = isurf + 1
651                k = nzb_s_inner(j,i)+1
652                surfl(:,isurf) = (/iroof,k,j,i/)
653            ENDDO
654        ENDDO
655
656!--     add walls
657        DO i = nxl, nxr
658            DO j = nys, nyn
659                DO ids = 1, 4  !> four wall directions
660                    jr = min(max(j-jdir(ids),0),ny)
661                    ir = min(max(i-idir(ids),0),nx)
662                    DO k = nzb_s_inner(j,i)+1, nzb_s_inner(jr,ir)
663                        isurf = isurf + 1
664                        surfl(:,isurf) = (/ids,k,j,i/)
665                    ENDDO
666                ENDDO
667            ENDDO
668        ENDDO
669
670!--     add sky
671        DO i = nxl, nxr
672            DO j = nys, nyn
673                isurf = isurf + 1
674                k = nzut
675                surfl(:,isurf) = (/isky,k,j,i/)
676            ENDDO
677        ENDDO
678       
679!--     calulation of the free borders of the domain
680        DO ids = 6,9
681            IF ( isborder(ids) )  THEN
682!--             free border of the domain in direction ids
683                DO i = ijdb(1,ids), ijdb(2,ids)
684                    DO j = ijdb(3,ids), ijdb(4,ids)
685                        DO k = max(nzb_s_inner(j,i),nzb_s_inner(j-jdir(ids),i-idir(ids)))+1, nzut
686                            isurf = isurf + 1
687                            surfl(:,isurf) = (/ids,k,j,i/)
688                        ENDDO
689                    ENDDO
690                ENDDO
691            ENDIF
692        ENDDO
693       
694!--     global array surf of indices of surfaces and displacement index array surfstart
695        ALLOCATE(nsurfs(0:numprocs-1))
696       
697#if defined( __parallel )
698        CALL MPI_Allgather(nsurfl,1,MPI_INTEGER,nsurfs,1,MPI_INTEGER,comm2d,ierr)
699#else
700        nsurfs(0) = nsurfl
701#endif
702        ALLOCATE(surfstart(0:numprocs))
703        k = 0
704        DO i=0,numprocs-1
705            surfstart(i) = k
706            k = k+nsurfs(i)
707        ENDDO
708        surfstart(numprocs) = k
709        nsurf = k
710        ALLOCATE(surf(4,nsurf))
711       
712#if defined( __parallel )
713        CALL MPI_AllGatherv(surfl, nsurfl*4, MPI_INTEGER, surf, nsurfs*4, surfstart*4, MPI_INTEGER, comm2d, ierr)
714#else
715        surf = surfl
716#endif
717       
718!--
719!--     allocation of the arrays for direct and diffusion radiation
720        CALL location_message( '    allocation of radiation arrays', .TRUE. )
721!--     rad_sw_in, rad_lw_in are computed in radiation model,
722!--     splitting of direct and diffusion part is done
723!--     in usm_calc_diffusion_radiation for now
724        ALLOCATE( rad_sw_in_dir(nysg:nyng,nxlg:nxrg) )
725        ALLOCATE( rad_sw_in_diff(nysg:nyng,nxlg:nxrg) )
726        ALLOCATE( rad_lw_in_diff(nysg:nyng,nxlg:nxrg) )
727       
728!--     allocate radiation arrays
729        ALLOCATE( surfins(nsurfl) )
730        ALLOCATE( surfinl(nsurfl) )
731        ALLOCATE( surfinsw(nsurfl) )
732        ALLOCATE( surfinlw(nsurfl) )
733        ALLOCATE( surfinswdir(nsurfl) )
734        ALLOCATE( surfinswdif(nsurfl) )
735        ALLOCATE( surfinlwdif(nsurfl) )
736        ALLOCATE( surfoutsl(startenergy:endenergy) )
737        ALLOCATE( surfoutll(startenergy:endenergy) )
738        ALLOCATE( surfoutsw(startenergy:endenergy) )
739        ALLOCATE( surfoutlw(startenergy:endenergy) )
740        ALLOCATE( surfouts(nsurf) ) !TODO: global surfaces without virtual
741        ALLOCATE( surfoutl(nsurf) ) !TODO: global surfaces without virtual
742        ALLOCATE( surfhf(startenergy:endenergy) )
743        ALLOCATE( rad_net_l(startenergy:endenergy) )
744
745!--     Wall surface model
746!--     allocate arrays for wall surface model and define pointers
747       
748!--     allocate array of wall types and wall parameters
749        ALLOCATE ( surface_types(startenergy:endenergy) )
750       
751!--     broadband albedo of the land, roof and wall surface
752!--     for domain border and sky set artifically to 1.0
753!--     what allows us to calculate heat flux leaving over
754!--     side and top borders of the domain
755        ALLOCATE ( albedo_surf(nsurfl) )
756        albedo_surf = 1.0_wp
757       
758!--     wall and roof surface parameters
759        ALLOCATE ( isroof_surf(startenergy:endenergy) )
760        ALLOCATE ( emiss_surf(startenergy:endenergy) )
761        ALLOCATE ( lambda_surf(startenergy:endenergy) )
762        ALLOCATE ( c_surface(startenergy:endenergy) )
763        ALLOCATE ( roughness_wall(startenergy:endenergy) )
764       
765!--     allocate wall and roof material parameters
766        ALLOCATE ( thickness_wall(startenergy:endenergy) )
767        ALLOCATE ( lambda_h(nzb_wall:nzt_wall,startenergy:endenergy) )
768        ALLOCATE ( rho_c_wall(nzb_wall:nzt_wall,startenergy:endenergy) )
769
770!--     allocate wall and roof layers sizes
771        ALLOCATE ( zwn(nzb_wall:nzt_wall) )
772        ALLOCATE ( dz_wall(nzb_wall:nzt_wall+1, startenergy:endenergy) )
773        ALLOCATE ( ddz_wall(nzb_wall:nzt_wall+1, startenergy:endenergy) )
774        ALLOCATE ( dz_wall_stag(nzb_wall:nzt_wall, startenergy:endenergy) )
775        ALLOCATE ( ddz_wall_stag(nzb_wall:nzt_wall, startenergy:endenergy) )
776        ALLOCATE ( zw(nzb_wall:nzt_wall, startenergy:endenergy) )
777
778!--     allocate wall and roof temperature arrays
779#if defined( __nopointer )
780        ALLOCATE ( t_surf(startenergy:endenergy) )
781        ALLOCATE ( t_surf_p(startenergy:endenergy) )
782        ALLOCATE ( t_wall(nzb_wall:nzt_wall+1,startenergy:endenergy) )
783        ALLOCATE ( t_wall_p(nzb_wall:nzt_wall+1,startenergy:endenergy) )
784#else
785        ALLOCATE ( t_surf_1(startenergy:endenergy) )
786        ALLOCATE ( t_surf_2(startenergy:endenergy) )
787        ALLOCATE ( t_wall_1(nzb_wall:nzt_wall+1,startenergy:endenergy) )
788        ALLOCATE ( t_wall_2(nzb_wall:nzt_wall+1,startenergy:endenergy) )
789
790!--     initial assignment of the pointers
791        t_wall    => t_wall_1;    t_wall_p    => t_wall_2
792        t_surf => t_surf_1; t_surf_p => t_surf_2
793#endif
794
795!--     allocate intermediate timestep arrays
796        ALLOCATE ( tt_surface_m(startenergy:endenergy) )
797        ALLOCATE ( tt_wall_m(nzb_wall:nzt_wall+1,startenergy:endenergy) )
798
799!--     allocate wall heat flux output array
800        ALLOCATE ( wshf(startwall:endwall) )
801        ALLOCATE ( wshf_eb(startenergy:endenergy) )
802        ALLOCATE ( wghf_eb(startenergy:endenergy) )
803
804!--     set inital values for prognostic quantities
805        tt_surface_m = 0.0_wp
806        tt_wall_m    = 0.0_wp
807
808        wshf = 0.0_wp
809        wshf_eb = 0.0_wp
810        wghf_eb = 0.0_wp
811       
812    END SUBROUTINE usm_allocate_urban_surface
813
814
815
816!------------------------------------------------------------------------------!
817! Description:
818! ------------
819!> Sum up and time-average urban surface output quantities as well as allocate
820!> the array necessary for storing the average.
821!------------------------------------------------------------------------------!
822    SUBROUTINE usm_average_3d_data( mode, variable )
823
824        IMPLICIT NONE
825
826        CHARACTER (len=*), INTENT(IN) ::  mode
827        CHARACTER (len=*), INTENT(IN) :: variable
828 
829        INTEGER(iwp)                                       :: i, j, k, l, ids, iwl,istat
830        CHARACTER (len=20)                                 :: var, surfid
831        INTEGER(iwp), PARAMETER                            :: nd = 5
832        CHARACTER(len=6), DIMENSION(0:nd-1), PARAMETER     :: dirname = (/ '_roof ', '_south', '_north', '_west ', '_east ' /)
833
834!--     find the real name of the variable
835        var = TRIM(variable)
836        DO i = 0, nd-1
837            k = len(TRIM(var))
838            j = len(TRIM(dirname(i)))
839            IF ( var(k-j+1:k) == dirname(i) )  THEN
840                ids = i
841                var = var(:k-j)
842                EXIT
843            ENDIF
844        ENDDO
845        IF ( ids == -1 )  THEN
846            var = TRIM(variable)
847        ENDIF
848        IF ( var(1:11) == 'usm_t_wall_'  .AND.  len(TRIM(var)) >= 12 )  THEN
849!--          wall layers
850            READ(var(12:12), '(I1)', iostat=istat ) iwl
851            IF ( istat == 0  .AND.  iwl >= nzb_wall  .AND.  iwl <= nzt_wall )  THEN
852                var = var(1:10)
853            ELSE
854!--             wrong wall layer index
855                RETURN
856            ENDIF
857        ENDIF
858
859        IF ( mode == 'allocate' )  THEN
860           
861           SELECT CASE ( TRIM( var ) )
862               
863                CASE ( 'usm_radnet' )
864!--                 array of complete radiation balance
865                    IF ( .NOT.  ALLOCATED(rad_net_av) )  THEN
866                        ALLOCATE( rad_net_av(startenergy:endenergy) )
867                        rad_net_av = 0.0_wp
868                    ENDIF
869                   
870                CASE ( 'usm_rad_insw' )
871!--                 array of sw radiation falling to surface after i-th reflection
872                    IF ( .NOT.  ALLOCATED(surfinsw_av) )  THEN
873                        ALLOCATE( surfinsw_av(startenergy:endenergy) )
874                        surfinsw_av = 0.0_wp
875                    ENDIF
876                                   
877                CASE ( 'usm_rad_inlw' )
878!--                 array of lw radiation falling to surface after i-th reflection
879                    IF ( .NOT.  ALLOCATED(surfinlw_av) )  THEN
880                        ALLOCATE( surfinlw_av(startenergy:endenergy) )
881                        surfinlw_av = 0.0_wp
882                    ENDIF
883
884                CASE ( 'usm_rad_inswdir' )
885!--                 array of direct sw radiation falling to surface from sun
886                    IF ( .NOT.  ALLOCATED(surfinswdir_av) )  THEN
887                        ALLOCATE( surfinswdir_av(startenergy:endenergy) )
888                        surfinswdir_av = 0.0_wp
889                    ENDIF
890
891                CASE ( 'usm_rad_inswdif' )
892!--                 array of difusion sw radiation falling to surface from sky and borders of the domain
893                    IF ( .NOT.  ALLOCATED(surfinswdif_av) )  THEN
894                        ALLOCATE( surfinswdif_av(startenergy:endenergy) )
895                        surfinswdif_av = 0.0_wp
896                    ENDIF
897
898                CASE ( 'usm_rad_inswref' )
899!--                 array of sw radiation falling to surface from reflections
900                    IF ( .NOT.  ALLOCATED(surfinswref_av) )  THEN
901                        ALLOCATE( surfinswref_av(startenergy:endenergy) )
902                        surfinswref_av = 0.0_wp
903                    ENDIF
904
905                CASE ( 'usm_rad_inlwdif' )
906!--                 array of sw radiation falling to surface after i-th reflection
907                    IF ( .NOT.  ALLOCATED(surfinlwdif_av) )  THEN
908                        ALLOCATE( surfinlwdif_av(startenergy:endenergy) )
909                        surfinlwdif_av = 0.0_wp
910                    ENDIF
911
912                CASE ( 'usm_rad_inlwref' )
913!--                 array of lw radiation falling to surface from reflections
914                    IF ( .NOT.  ALLOCATED(surfinlwref_av) )  THEN
915                        ALLOCATE( surfinlwref_av(startenergy:endenergy) )
916                        surfinlwref_av = 0.0_wp
917                    ENDIF
918
919                CASE ( 'usm_rad_outsw' )
920!--                 array of sw radiation emitted from surface after i-th reflection
921                    IF ( .NOT.  ALLOCATED(surfoutsw_av) )  THEN
922                        ALLOCATE( surfoutsw_av(startenergy:endenergy) )
923                        surfoutsw_av = 0.0_wp
924                    ENDIF
925
926                CASE ( 'usm_rad_outlw' )
927!--                 array of lw radiation emitted from surface after i-th reflection
928                    IF ( .NOT.  ALLOCATED(surfoutlw_av) )  THEN
929                        ALLOCATE( surfoutlw_av(startenergy:endenergy) )
930                        surfoutlw_av = 0.0_wp
931                    ENDIF
932
933                CASE ( 'usm_rad_hf' )
934!--                 array of heat flux from radiation for surfaces after i-th reflection
935                    IF ( .NOT.  ALLOCATED(surfhf_av) )  THEN
936                        ALLOCATE( surfhf_av(startenergy:endenergy) )
937                        surfhf_av = 0.0_wp
938                    ENDIF
939
940                CASE ( 'usm_wshf' )
941!--                 array of sensible heat flux from surfaces
942!--                 land surfaces
943                    IF ( .NOT.  ALLOCATED(wshf_eb_av) )  THEN
944                        ALLOCATE( wshf_eb_av(startenergy:endenergy) )
945                        wshf_eb_av = 0.0_wp
946                    ENDIF
947
948                CASE ( 'usm_wghf' )
949!--                 array of heat flux from ground (wall, roof, land)
950                    IF ( .NOT.  ALLOCATED(wghf_eb_av) )  THEN
951                        ALLOCATE( wghf_eb_av(startenergy:endenergy) )
952                        wghf_eb_av = 0.0_wp
953                    ENDIF
954
955                CASE ( 'usm_t_surf' )
956!--                 surface temperature for surfaces
957                    IF ( .NOT.  ALLOCATED(t_surf_av) )  THEN
958                        ALLOCATE( t_surf_av(startenergy:endenergy) )
959                        t_surf_av = 0.0_wp
960                    ENDIF
961                   
962                CASE ( 'usm_t_wall' )
963!--                 wall temperature for iwl layer of walls and land
964                    IF ( .NOT.  ALLOCATED(t_wall_av) )  THEN
965                        ALLOCATE( t_wall_av(nzb_wall:nzt_wall,startenergy:endenergy) )
966                        t_wall_av = 0.0_wp
967                    ENDIF
968
969               CASE DEFAULT
970                   CONTINUE
971
972           END SELECT
973
974        ELSEIF ( mode == 'sum' )  THEN
975           
976           SELECT CASE ( TRIM( var ) )
977               
978                CASE ( 'usm_radnet' )
979!--                 array of complete radiation balance
980                    DO l = startenergy, endenergy
981                        IF ( surfl(id,l) == ids )  THEN
982                            rad_net_av(l) = rad_net_av(l) + rad_net_l(l)
983                        ENDIF
984                    ENDDO
985                   
986                CASE ( 'usm_rad_insw' )
987!--                 array of sw radiation falling to surface after i-th reflection
988                    DO l = startenergy, endenergy
989                        IF ( surfl(id,l) == ids )  THEN
990                            surfinsw_av(l) = surfinsw_av(l) + surfinsw(l)
991                        ENDIF
992                    ENDDO
993                             
994                CASE ( 'usm_rad_inlw' )
995!--                 array of lw radiation falling to surface after i-th reflection
996                    DO l = startenergy, endenergy
997                        IF ( surfl(id,l) == ids )  THEN
998                            surfinlw_av(l) = surfinlw_av(l) + surfinlw(l)
999                        ENDIF
1000                    ENDDO
1001                   
1002                CASE ( 'usm_rad_inswdir' )
1003!--                 array of direct sw radiation falling to surface from sun
1004                    DO l = startenergy, endenergy
1005                        IF ( surfl(id,l) == ids )  THEN
1006                            surfinswdir_av(l) = surfinswdir_av(l) + surfinswdir(l)
1007                        ENDIF
1008                    ENDDO
1009                   
1010                CASE ( 'usm_rad_inswdif' )
1011!--                 array of difusion sw radiation falling to surface from sky and borders of the domain
1012                    DO l = startenergy, endenergy
1013                        IF ( surfl(id,l) == ids )  THEN
1014                            surfinswdif_av(l) = surfinswdif_av(l) + surfinswdif(l)
1015                        ENDIF
1016                    ENDDO
1017                   
1018                CASE ( 'usm_rad_inswref' )
1019!--                 array of sw radiation falling to surface from reflections
1020                    DO l = startenergy, endenergy
1021                        IF ( surfl(id,l) == ids )  THEN
1022                            surfinswref_av(l) = surfinswref_av(l) + surfinsw(l) - &
1023                                                surfinswdir(l) - surfinswdif(l)
1024                        ENDIF
1025                    ENDDO
1026                   
1027                CASE ( 'usm_rad_inlwdif' )
1028!--                 array of sw radiation falling to surface after i-th reflection
1029                    DO l = startenergy, endenergy
1030                        IF ( surfl(id,l) == ids )  THEN
1031                            surfinlwdif_av(l) = surfinlwdif_av(l) + surfinlwdif(l)
1032                        ENDIF
1033                    ENDDO
1034                   
1035                CASE ( 'usm_rad_inlwref' )
1036!--                 array of lw radiation falling to surface from reflections
1037                    DO l = startenergy, endenergy
1038                        IF ( surfl(id,l) == ids )  THEN
1039                            surfinlwref_av(l) = surfinlwref_av(l) + &
1040                                                surfinlw(l) - surfinlwdif(l)
1041                        ENDIF
1042                    ENDDO
1043                   
1044                CASE ( 'usm_rad_outsw' )
1045!--                 array of sw radiation emitted from surface after i-th reflection
1046                    DO l = startenergy, endenergy
1047                        IF ( surfl(id,l) == ids )  THEN
1048                            surfoutsw_av(l) = surfoutsw_av(l) + surfoutsw(l)
1049                        ENDIF
1050                    ENDDO
1051                   
1052                CASE ( 'usm_rad_outlw' )
1053!--                 array of lw radiation emitted from surface after i-th reflection
1054                    DO l = startenergy, endenergy
1055                        IF ( surfl(id,l) == ids )  THEN
1056                            surfoutlw_av(l) = surfoutlw_av(l) + surfoutlw(l)
1057                        ENDIF
1058                    ENDDO
1059                   
1060                CASE ( 'usm_rad_hf' )
1061!--                 array of heat flux from radiation for surfaces after i-th reflection
1062                    DO l = startenergy, endenergy
1063                        IF ( surfl(id,l) == ids )  THEN
1064                            surfhf_av(l) = surfhf_av(l) + surfhf(l)
1065                        ENDIF
1066                    ENDDO
1067                   
1068                CASE ( 'usm_wshf' )
1069!--                 array of sensible heat flux from surfaces (land, roof, wall)
1070                    DO l = startenergy, endenergy
1071                        IF ( surfl(id,l) == ids )  THEN
1072                            wshf_eb_av(l) = wshf_eb_av(l) + wshf_eb(l)
1073                        ENDIF
1074                    ENDDO
1075                   
1076                CASE ( 'usm_wghf' )
1077!--                 array of heat flux from ground (wall, roof, land)
1078                    DO l = startenergy, endenergy
1079                        IF ( surfl(id,l) == ids )  THEN
1080                            wghf_eb_av(l) = wghf_eb_av(l) + wghf_eb(l)
1081                        ENDIF
1082                    ENDDO
1083                   
1084                CASE ( 'usm_t_surf' )
1085!--                 surface temperature for surfaces
1086                    DO l = startenergy, endenergy
1087                        IF ( surfl(id,l) == ids )  THEN
1088                            t_surf_av(l) = t_surf_av(l) + t_surf(l)
1089                        ENDIF
1090                    ENDDO
1091                   
1092                CASE ( 'usm_t_wall' )
1093!--                 wall temperature for  iwl layer of walls and land
1094                    DO l = startenergy, endenergy
1095                        IF ( surfl(id,l) == ids )  THEN
1096                            t_wall_av(iwl, l) = t_wall_av(iwl,l) + t_wall(iwl, l)
1097                        ENDIF
1098                    ENDDO
1099                   
1100                CASE DEFAULT
1101                    CONTINUE
1102
1103           END SELECT
1104
1105        ELSEIF ( mode == 'average' )  THEN
1106           
1107           SELECT CASE ( TRIM( var ) )
1108               
1109                CASE ( 'usm_radnet' )
1110!--                 array of complete radiation balance
1111                    DO l = startenergy, endenergy
1112                        IF ( surfl(id,l) == ids )  THEN
1113                            rad_net_av(l) = rad_net_av(l) / REAL( average_count_3d, kind=wp )
1114                        ENDIF
1115                    ENDDO
1116                   
1117                CASE ( 'usm_rad_insw' )
1118!--                 array of sw radiation falling to surface after i-th reflection
1119                    DO l = startenergy, endenergy
1120                        IF ( surfl(id,l) == ids )  THEN
1121                            surfinsw_av(l) = surfinsw_av(l) / REAL( average_count_3d, kind=wp )
1122                        ENDIF
1123                    ENDDO
1124                                   
1125                CASE ( 'usm_rad_inlw' )
1126!--                 array of lw radiation falling to surface after i-th reflection
1127                    DO l = startenergy, endenergy
1128                        IF ( surfl(id,l) == ids )  THEN
1129                            surfinlw_av(l) = surfinlw_av(l) / REAL( average_count_3d, kind=wp )
1130                        ENDIF
1131                    ENDDO
1132
1133                CASE ( 'usm_rad_inswdir' )
1134!--                 array of direct sw radiation falling to surface from sun
1135                    DO l = startenergy, endenergy
1136                        IF ( surfl(id,l) == ids )  THEN
1137                            surfinswdir_av(l) = surfinswdir_av(l) / REAL( average_count_3d, kind=wp )
1138                        ENDIF
1139                    ENDDO
1140
1141                CASE ( 'usm_rad_inswdif' )
1142!--                 array of difusion sw radiation falling to surface from sky and borders of the domain
1143                    DO l = startenergy, endenergy
1144                        IF ( surfl(id,l) == ids )  THEN
1145                            surfinswdif_av(l) = surfinswdif_av(l) / REAL( average_count_3d, kind=wp )
1146                        ENDIF
1147                    ENDDO
1148
1149                CASE ( 'usm_rad_inswref' )
1150!--                 array of sw radiation falling to surface from reflections
1151                    DO l = startenergy, endenergy
1152                        IF ( surfl(id,l) == ids )  THEN
1153                            surfinswref_av(l) = surfinswref_av(l) / REAL( average_count_3d, kind=wp )
1154                        ENDIF
1155                    ENDDO
1156
1157                CASE ( 'usm_rad_inlwdif' )
1158!--                 array of sw radiation falling to surface after i-th reflection
1159                    DO l = startenergy, endenergy
1160                        IF ( surfl(id,l) == ids )  THEN
1161                            surfinlwdif_av(l) = surfinlwdif_av(l) / REAL( average_count_3d, kind=wp )
1162                        ENDIF
1163                    ENDDO
1164
1165                CASE ( 'usm_rad_inlwref' )
1166!--                 array of lw radiation falling to surface from reflections
1167                    DO l = startenergy, endenergy
1168                        IF ( surfl(id,l) == ids )  THEN
1169                            surfinlwref_av(l) = surfinlwref_av(l) / REAL( average_count_3d, kind=wp )
1170                        ENDIF
1171                    ENDDO
1172
1173                CASE ( 'usm_rad_outsw' )
1174!--                 array of sw radiation emitted from surface after i-th reflection
1175                    DO l = startenergy, endenergy
1176                        IF ( surfl(id,l) == ids )  THEN
1177                            surfoutsw_av(l) = surfoutsw_av(l) / REAL( average_count_3d, kind=wp )
1178                        ENDIF
1179                    ENDDO
1180
1181                CASE ( 'usm_rad_outlw' )
1182!--                 array of lw radiation emitted from surface after i-th reflection
1183                    DO l = startenergy, endenergy
1184                        IF ( surfl(id,l) == ids )  THEN
1185                            surfoutlw_av(l) = surfoutlw_av(l) / REAL( average_count_3d, kind=wp )
1186                        ENDIF
1187                    ENDDO
1188
1189                CASE ( 'usm_rad_hf' )
1190!--                 array of heat flux from radiation for surfaces after i-th reflection
1191                    DO l = startenergy, endenergy
1192                        IF ( surfl(id,l) == ids )  THEN
1193                            surfhf_av(l) = surfhf_av(l) / REAL( average_count_3d, kind=wp )
1194                        ENDIF
1195                    ENDDO
1196
1197                CASE ( 'usm_wshf' )
1198!--                 array of sensible heat flux from surfaces (land, roof, wall)
1199                    DO l = startenergy, endenergy
1200                        IF ( surfl(id,l) == ids )  THEN
1201                            wshf_eb_av(l) = wshf_eb_av(l) / REAL( average_count_3d, kind=wp )
1202                        ENDIF
1203                    ENDDO
1204
1205                CASE ( 'usm_wghf' )
1206!--                 array of heat flux from ground (wall, roof, land)
1207                    DO l = startenergy, endenergy
1208                        IF ( surfl(id,l) == ids )  THEN
1209                            wghf_eb_av(l) = wghf_eb_av(l) / REAL( average_count_3d, kind=wp )
1210                        ENDIF
1211                    ENDDO
1212
1213                CASE ( 'usm_t_surf' )
1214!--                 surface temperature for surfaces
1215                    DO l = startenergy, endenergy
1216                        IF ( surfl(id,l) == ids )  THEN
1217                            t_surf_av(l) = t_surf_av(l) / REAL( average_count_3d, kind=wp )
1218                        ENDIF
1219                    ENDDO
1220                   
1221                CASE ( 'usm_t_wall' )
1222!--                 wall temperature for  iwl layer of walls and land
1223                    DO l = startenergy, endenergy
1224                        IF ( surfl(id,l) == ids )  THEN
1225                            t_wall_av(iwl, l) = t_wall_av(iwl,l) / REAL( average_count_3d, kind=wp )
1226                        ENDIF
1227                    ENDDO
1228               
1229           END SELECT
1230
1231        ENDIF
1232
1233    END SUBROUTINE usm_average_3d_data
1234
1235
1236!------------------------------------------------------------------------------!
1237!> Calculates radiation absorbed by box with given size and LAD.
1238!>
1239!> Simulates resol**2 rays (by equally spacing a bounding horizontal square
1240!> conatining all possible rays that would cross the box) and calculates
1241!> average transparency per ray. Returns fraction of absorbed radiation flux
1242!> and area for which this fraction is effective.
1243!------------------------------------------------------------------------------!
1244    PURE SUBROUTINE usm_box_absorb(boxsize, resol, dens, uvec, area, absorb)
1245        IMPLICIT NONE
1246
1247        REAL(wp), DIMENSION(3), INTENT(in) :: &
1248            boxsize, &      !< z, y, x size of box in m
1249            uvec            !< z, y, x unit vector of incoming flux
1250        INTEGER(iwp), INTENT(in) :: &
1251            resol           !< No. of rays in x and y dimensions
1252        REAL(wp), INTENT(in) :: &
1253            dens            !< box density (e.g. Leaf Area Density)
1254        REAL(wp), INTENT(out) :: &
1255            area, &         !< horizontal area for flux absorbtion
1256            absorb          !< fraction of absorbed flux
1257        REAL(wp) :: &
1258            xshift, yshift, &
1259            xmin, xmax, ymin, ymax, &
1260            xorig, yorig, &
1261            dx1, dy1, dz1, dx2, dy2, dz2, &
1262            crdist, &
1263            transp
1264        INTEGER(iwp) :: &
1265            i, j
1266
1267        xshift = uvec(3) / uvec(1) * boxsize(1)
1268        xmin = min(0._wp, -xshift)
1269        xmax = boxsize(3) + max(0._wp, -xshift)
1270        yshift = uvec(2) / uvec(1) * boxsize(1)
1271        ymin = min(0._wp, -yshift)
1272        ymax = boxsize(2) + max(0._wp, -yshift)
1273
1274        transp = 0._wp
1275        DO i = 1, resol
1276            xorig = xmin + (xmax-xmin) * (i-.5_wp) / resol
1277            DO j = 1, resol
1278                yorig = ymin + (ymax-ymin) * (j-.5_wp) / resol
1279
1280                dz1 = 0._wp
1281                dz2 = boxsize(1)/uvec(1)
1282
1283                IF ( uvec(2) > 0._wp )  THEN
1284                    dy1 = -yorig             / uvec(2) !< crossing with y=0
1285                    dy2 = (boxsize(2)-yorig) / uvec(2) !< crossing with y=boxsize(2)
1286                ELSE IF ( uvec(2) < 0._wp )  THEN
1287                    dy1 = (boxsize(2)-yorig) / uvec(2) !< crossing with y=boxsize(2)
1288                    dy2 = -yorig             / uvec(2) !< crossing with y=0
1289                ELSE !uvec(2)==0
1290                    dy1 = -huge(1._wp)
1291                    dy2 = huge(1._wp)
1292                ENDIF
1293
1294                IF ( uvec(3) > 0._wp )  THEN
1295                    dx1 = -xorig             / uvec(3) !< crossing with x=0
1296                    dx2 = (boxsize(3)-xorig) / uvec(3) !< crossing with x=boxsize(3)
1297                ELSE IF ( uvec(3) < 0._wp )  THEN
1298                    dx1 = (boxsize(3)-xorig) / uvec(3) !< crossing with x=boxsize(3)
1299                    dx2 = -xorig             / uvec(3) !< crossing with x=0
1300                ELSE !uvec(1)==0
1301                    dx1 = -huge(1._wp)
1302                    dx2 = huge(1._wp)
1303                ENDIF
1304
1305                crdist = max(0._wp, (min(dz2, dy2, dx2) - max(dz1, dy1, dx1)))
1306                transp = transp + exp(-ext_coef * dens * crdist)
1307            ENDDO
1308        ENDDO
1309        transp = transp / resol**2
1310        area = (boxsize(3)+xshift)*(boxsize(2)+yshift)
1311        absorb = 1._wp - transp
1312       
1313    END SUBROUTINE usm_box_absorb
1314   
1315   
1316!------------------------------------------------------------------------------!
1317! Description:
1318! ------------
1319!> This subroutine splits direct and diffusion dw radiation
1320!> It sould not be called in case the radiation model already does it
1321!> It follows <CITATION>
1322!------------------------------------------------------------------------------!
1323    SUBROUTINE usm_calc_diffusion_radiation 
1324   
1325        REAL(wp), PARAMETER                          ::  sol_const = 1367.0_wp   !< solar conbstant
1326        REAL(wp), PARAMETER                          :: lowest_solarUp = 0.1_wp  !< limit the sun elevation to protect stability of the calculation
1327        INTEGER(iwp)                                 :: i, j
1328        REAL(wp), PARAMETER                          ::  year_seconds = 86400._wp * 365._wp
1329        REAL(wp)                                     ::  year_angle              !< angle
1330        REAL(wp)                                     ::  etr                     !< extraterestrial radiation
1331        REAL(wp)                                     ::  corrected_solarUp       !< corrected solar up radiation
1332        REAL(wp)                                     ::  horizontalETR           !< horizontal extraterestrial radiation
1333        REAL(wp)                                     ::  clearnessIndex          !< clearness index
1334        REAL(wp)                                     ::  diff_frac               !< diffusion fraction of the radiation
1335
1336       
1337!--     Calculate current day and time based on the initial values and simulation time
1338        year_angle = ((day_init*86400) + time_utc_init+time_since_reference_point) &
1339                       / year_seconds * 2.0_wp * pi
1340       
1341        etr = sol_const * (1.00011_wp +                                            &
1342                          0.034221_wp * cos(year_angle) +                          &
1343                          0.001280_wp * sin(year_angle) +                          &
1344                          0.000719_wp * cos(2.0_wp * year_angle) +                 &
1345                          0.000077_wp * sin(2.0_wp * year_angle))
1346       
1347!--   
1348!--     Under a very low angle, we keep extraterestrial radiation at
1349!--     the last small value, therefore the clearness index will be pushed
1350!--     towards 0 while keeping full continuity.
1351!--   
1352        IF ( zenith(0) <= lowest_solarUp )  THEN
1353            corrected_solarUp = lowest_solarUp
1354        ELSE
1355            corrected_solarUp = zenith(0)
1356        ENDIF
1357       
1358        horizontalETR = etr * corrected_solarUp
1359       
1360        DO i = nxlg, nxrg
1361            DO j = nysg, nyng
1362                clearnessIndex = rad_sw_in(0,j,i) / horizontalETR
1363                diff_frac = 1.0_wp / (1.0_wp + exp(-5.0033_wp + 8.6025_wp * clearnessIndex))
1364                rad_sw_in_diff(j,i) = rad_sw_in(0,j,i) * diff_frac
1365                rad_sw_in_dir(j,i)  = rad_sw_in(0,j,i) * (1.0_wp - diff_frac)
1366                rad_lw_in_diff(j,i) = rad_lw_in(0,j,i)
1367            ENDDO
1368        ENDDO
1369       
1370    END SUBROUTINE usm_calc_diffusion_radiation
1371   
1372
1373!------------------------------------------------------------------------------!
1374! Description:
1375! ------------
1376!> Calculates shape view factors SVF and plant sink canopy factors PSCF
1377!> !!!!!DESCRIPTION!!!!!!!!!!
1378!------------------------------------------------------------------------------!
1379    SUBROUTINE usm_calc_svf
1380   
1381        IMPLICIT NONE
1382       
1383        INTEGER(iwp)                                :: i, j, k, l, d, ip, jp
1384        INTEGER(iwp)                                :: isvf, ksvf, icsf, kcsf, npcsfl, isvf_surflt, imrtt, imrtf
1385        INTEGER(iwp)                                :: sd, td, ioln, iproc
1386        REAL(wp),     DIMENSION(0:9)                :: facearea
1387        INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE   :: nzterrl, planthl
1388        REAL(wp),     DIMENSION(:,:), ALLOCATABLE   :: csflt, pcsflt
1389        INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE   :: kcsflt,kpcsflt
1390        INTEGER(iwp), DIMENSION(:), ALLOCATABLE     :: icsflt,dcsflt,ipcsflt,dpcsflt
1391        REAL(wp), DIMENSION(3)                      :: uv
1392        LOGICAL                                     :: visible
1393        REAL(wp)                                    :: sz, sy, sx, tz, ty, tx, transparency, rirrf, sqdist, svfsum
1394        !REAL(wp)                                    :: rsvf
1395        INTEGER(iwp)                                :: isurflt, isurfs, isurflt_prev
1396        INTEGER(iwp)                                :: itx, ity, itz
1397        CHARACTER(len=7)                            :: pid_char = ''
1398        INTEGER(iwp)                                :: win_lad, minfo
1399        REAL(wp), DIMENSION(:,:,:), POINTER         :: lad_s_rma       !< fortran pointer, but lower bounds are 1
1400        TYPE(c_ptr)                                 :: lad_s_rma_p     !< allocated c pointer
1401        INTEGER(kind=MPI_ADDRESS_KIND)              :: size_lad_rma
1402   
1403!--     calculation of the SVF
1404        CALL location_message( '    calculation of SVF and CSF', .TRUE. )
1405
1406!--     precalculate face areas for different face directions using normal vector
1407        DO d = 0, 9
1408            facearea(d) = 1._wp
1409            IF ( idir(d) == 0 ) facearea(d) = facearea(d) * dx
1410            IF ( jdir(d) == 0 ) facearea(d) = facearea(d) * dy
1411            IF ( kdir(d) == 0 ) facearea(d) = facearea(d) * dz
1412        ENDDO
1413
1414!--     initialize variables and temporary arrays for calculation of svf and csf
1415        nsvfl  = 0
1416        ncsfl  = 0
1417        nsvfla = gasize
1418        msvf   = 1
1419        ALLOCATE( asvf1(nsvfla) )
1420        asvf => asvf1
1421        IF ( plant_canopy )  THEN
1422            ncsfla = gasize
1423            mcsf   = 1
1424            ALLOCATE( acsf1(ncsfla) )
1425            acsf => acsf1
1426        ENDIF
1427       
1428!--     initialize temporary terrain and plant canopy height arrays (global 2D array!)
1429        ALLOCATE( nzterr(0:(nx+1)*(ny+1)-1) )
1430#if defined( __parallel )
1431        ALLOCATE( nzterrl(nys:nyn,nxl:nxr) )
1432        nzterrl = nzb_s_inner(nys:nyn,nxl:nxr)
1433        CALL MPI_AllGather( nzterrl, nnx*nny, MPI_INTEGER, &
1434                            nzterr, nnx*nny, MPI_INTEGER, comm2d, ierr )
1435        DEALLOCATE(nzterrl)
1436#else
1437        nzterr = RESHAPE( nzb_s_inner(nys:nyn,nxl:nxr), (/(nx+1)*(ny+1)/) )
1438#endif
1439        IF ( plant_canopy )  THEN
1440            ALLOCATE( plantt(0:(nx+1)*(ny+1)-1) )
1441#if defined( __parallel )
1442            ALLOCATE( planthl(nys:nyn,nxl:nxr) )
1443            planthl = pch(nys:nyn,nxl:nxr)
1444       
1445            CALL MPI_AllGather( planthl, nnx*nny, MPI_INTEGER, &
1446                                plantt, nnx*nny, MPI_INTEGER, comm2d, ierr )
1447            DEALLOCATE(planthl)
1448
1449            IF ( usm_lad_rma )  THEN
1450                CALL MPI_Info_create(minfo, ierr)
1451                CALL MPI_Info_set(minfo, 'accumulate_ordering', '', ierr)
1452                CALL MPI_Info_set(minfo, 'accumulate_ops', 'same_op', ierr)
1453                CALL MPI_Info_set(minfo, 'same_size', 'true', ierr)
1454                CALL MPI_Info_set(minfo, 'same_disp_unit', 'true', ierr)
1455
1456!--             Allocate and initialize the MPI RMA window
1457!--             must be in accordance with allocation of lad_s in plant_canopy_model
1458!--             optimization of memory should be done
1459!--             Argument X of function c_sizeof(X) needs arbitrary REAL(wp) value, set to 1.0_wp for now
1460                size_lad_rma = c_sizeof(1.0_wp)*nnx*nny*nzu
1461                CALL MPI_Win_allocate(size_lad_rma, c_sizeof(1.0_wp), minfo, comm2d, &
1462                                        lad_s_rma_p, win_lad, ierr)
1463                CALL c_f_pointer(lad_s_rma_p, lad_s_rma, (/ nzu, nny, nnx /))
1464                usm_lad(nzub:, nys:, nxl:) => lad_s_rma(:,:,:)
1465            ELSE
1466                ALLOCATE(usm_lad(nzub:nzut, nys:nyn, nxl:nxr))
1467            ENDIF
1468#else
1469            plantt = RESHAPE( pct(nys:nyn,nxl:nxr), (/(nx+1)*(ny+1)/) )
1470            ALLOCATE(usm_lad(nzub:nzut, nys:nyn, nxl:nxr))
1471#endif
1472            usm_lad(:,:,:) = 0._wp
1473            DO i = nxl, nxr
1474                DO j = nys, nyn
1475                    k = nzb_s_inner(j, i)
1476                    usm_lad(k:nzut, j, i) = lad_s(0:nzut-k, j, i)
1477                ENDDO
1478            ENDDO
1479
1480#if defined( __parallel )
1481            IF ( usm_lad_rma )  THEN
1482                CALL MPI_Info_free(minfo, ierr)
1483                CALL MPI_Win_lock_all(0, win_lad, ierr)
1484            ELSE
1485                ALLOCATE( usm_lad_g(0:(nx+1)*(ny+1)*nzu-1) )
1486                CALL MPI_AllGather( usm_lad, nnx*nny*nzu, MPI_REAL, &
1487                                    usm_lad_g, nnx*nny*nzu, MPI_REAL, comm2d, ierr )
1488            ENDIF
1489#endif
1490        ENDIF
1491
1492        IF ( mrt_factors )  THEN
1493            OPEN(153, file='MRT_TARGETS', access='SEQUENTIAL', &
1494                    action='READ', status='OLD', form='FORMATTED', err=524)
1495            OPEN(154, file='MRT_FACTORS'//myid_char, access='DIRECT', recl=(5*4+2*8), &
1496                    action='WRITE', status='REPLACE', form='UNFORMATTED', err=525)
1497            imrtf = 1
1498            DO
1499                READ(153, *, end=526, err=524) imrtt, i, j, k
1500                IF ( i < nxl  .OR.  i > nxr &
1501                     .OR.  j < nys  .OR.  j > nyn ) CYCLE
1502                tx = REAL(i)
1503                ty = REAL(j)
1504                tz = REAL(k)
1505
1506                DO isurfs = 1, nsurf
1507                    IF ( .NOT.  usm_facing(i, j, k, -1, &
1508                        surf(ix, isurfs), surf(iy, isurfs), &
1509                        surf(iz, isurfs), surf(id, isurfs)) )  THEN
1510                        CYCLE
1511                    ENDIF
1512                     
1513                    sd = surf(id, isurfs)
1514                    sz = REAL(surf(iz, isurfs), wp) - 0.5_wp * kdir(sd)
1515                    sy = REAL(surf(iy, isurfs), wp) - 0.5_wp * jdir(sd)
1516                    sx = REAL(surf(ix, isurfs), wp) - 0.5_wp * idir(sd)
1517
1518!--                 unit vector source -> target
1519                    uv = (/ (tz-sz)*dz, (ty-sy)*dy, (tx-sx)*dx /)
1520                    sqdist = SUM(uv(:)**2)
1521                    uv = uv / SQRT(sqdist)
1522
1523!--                 irradiance factor - see svf. Here we consider that target face is always normal,
1524!--                 i.e. the second dot product equals 1
1525                    rirrf = dot_product((/ kdir(sd), jdir(sd), idir(sd) /), uv) &
1526                        / (pi * sqdist) * facearea(sd)
1527
1528!--                 raytrace while not creating any canopy sink factors
1529                    CALL usm_raytrace((/sz,sy,sx/), (/tz,ty,tx/), isurfs, rirrf, 1._wp, .FALSE., &
1530                            visible, transparency, win_lad)
1531                    IF ( .NOT.  visible ) CYCLE
1532
1533                    !rsvf = rirrf * transparency
1534                    WRITE(154, rec=imrtf, err=525) INT(imrtt, kind=4), &
1535                        INT(surf(id, isurfs), kind=4), &
1536                        INT(surf(iz, isurfs), kind=4), &
1537                        INT(surf(iy, isurfs), kind=4), &
1538                        INT(surf(ix, isurfs), kind=4), &
1539                        REAL(rirrf, kind=8), REAL(transparency, kind=8)
1540                    imrtf = imrtf + 1
1541
1542                ENDDO !< isurfs
1543            ENDDO !< MRT_TARGETS record
1544
1545524         message_string = 'error reading file MRT_TARGETS'
1546            CALL message( 'usm_calc_svf', 'PA0524', 1, 2, 0, 6, 0 )
1547
1548525         message_string = 'error writing file MRT_FACTORS'//myid_char
1549            CALL message( 'usm_calc_svf', 'PA0525', 1, 2, 0, 6, 0 )
1550
1551526         CLOSE(153)
1552            CLOSE(154)
1553        ENDIF  !< mrt_factors
1554
1555       
1556        DO isurflt = 1, nsurfl
1557!--         determine face centers
1558            td = surfl(id, isurflt)
1559            IF ( td >= isky  .AND.  .NOT.  plant_canopy ) CYCLE
1560            tz = REAL(surfl(iz, isurflt), wp) - 0.5_wp * kdir(td)
1561            ty = REAL(surfl(iy, isurflt), wp) - 0.5_wp * jdir(td)
1562            tx = REAL(surfl(ix, isurflt), wp) - 0.5_wp * idir(td)
1563            DO isurfs = 1, nsurf
1564                IF ( .NOT.  usm_facing(surfl(ix, isurflt), surfl(iy, isurflt), &
1565                    surfl(iz, isurflt), surfl(id, isurflt), &
1566                    surf(ix, isurfs), surf(iy, isurfs), &
1567                    surf(iz, isurfs), surf(id, isurfs)) )  THEN
1568                    CYCLE
1569                ENDIF
1570                 
1571                sd = surf(id, isurfs)
1572                sz = REAL(surf(iz, isurfs), wp) - 0.5_wp * kdir(sd)
1573                sy = REAL(surf(iy, isurfs), wp) - 0.5_wp * jdir(sd)
1574                sx = REAL(surf(ix, isurfs), wp) - 0.5_wp * idir(sd)
1575
1576!--             unit vector source -> target
1577                uv = (/ (tz-sz)*dz, (ty-sy)*dy, (tx-sx)*dx /)
1578                sqdist = SUM(uv(:)**2)
1579                uv = uv / SQRT(sqdist)
1580               
1581!--             irradiance factor (our unshaded shape view factor) = view factor per differential target area * source area
1582                rirrf = dot_product((/ kdir(sd), jdir(sd), idir(sd) /), uv) & ! cosine of source normal and direction
1583                    * dot_product((/ kdir(td), jdir(td), idir(td) /), -uv) &  ! cosine of target normal and reverse direction
1584                    / (pi * sqdist) & ! square of distance between centers
1585                    * facearea(sd)
1586
1587!--             raytrace + process plant canopy sinks within
1588                CALL usm_raytrace((/sz,sy,sx/), (/tz,ty,tx/), isurfs, rirrf, facearea(td), .TRUE., &
1589                        visible, transparency, win_lad)
1590               
1591                IF ( .NOT.  visible ) CYCLE
1592                IF ( td >= isky ) CYCLE !< we calculated these only for raytracing
1593                                        !< to find plant canopy sinks, we don't need svf for them
1594                ! rsvf = rirrf * transparency
1595
1596!--             write to the svf array
1597                nsvfl = nsvfl + 1
1598!--             check dimmension of asvf array and enlarge it if needed
1599                IF ( nsvfla < nsvfl )  THEN
1600                    k = nsvfla * 2
1601                    IF ( msvf == 0 )  THEN
1602                        msvf = 1
1603                        ALLOCATE( asvf1(k) )
1604                        asvf => asvf1
1605                        asvf1(1:nsvfla) = asvf2
1606                        DEALLOCATE( asvf2 )
1607                    ELSE
1608                        msvf = 0
1609                        ALLOCATE( asvf2(k) )
1610                        asvf => asvf2
1611                        asvf2(1:nsvfla) = asvf1
1612                        DEALLOCATE( asvf1 )
1613                    ENDIF
1614                    nsvfla = k
1615                ENDIF
1616!--             write svf values into the array
1617                asvf(nsvfl)%isurflt = isurflt
1618                asvf(nsvfl)%isurfs = isurfs
1619                asvf(nsvfl)%rsvf = rirrf !we postopne multiplication by transparency
1620                asvf(nsvfl)%rtransp = transparency !a.k.a. Direct Irradiance Factor
1621            ENDDO
1622        ENDDO
1623
1624        CALL location_message( '    waiting for completion of SVF and CSF calculation in all processes', .TRUE. )
1625
1626#if defined( __parallel )
1627        IF ( plant_canopy )  THEN
1628            IF ( usm_lad_rma )  THEN
1629                CALL MPI_Win_flush_all(win_lad, ierr)
1630!--             unlock MPI window
1631                CALL MPI_Win_unlock_all(win_lad, ierr)
1632!--             free MPI window
1633                CALL MPI_Win_free(win_lad, ierr)
1634            ELSE
1635                DEALLOCATE(usm_lad)
1636                DEALLOCATE(usm_lad_g)
1637            ENDIF
1638        ENDIF
1639#else
1640        DEALLOCATE(usm_lad)
1641#endif
1642
1643!--     deallocate temporary global arrays
1644        IF ( ALLOCATED(nzterr) ) DEALLOCATE(nzterr)
1645        IF ( ALLOCATED(plantt) ) DEALLOCATE(plantt)
1646       
1647!--     sort svf ( a version of quicksort )
1648        CALL quicksort_svf(asvf,1,nsvfl)
1649       
1650        npcsfl = 0
1651        IF ( plant_canopy )  THEN
1652!--         sort and merge csf for the last time, keeping the array size to minimum
1653            CALL usm_merge_and_grow_csf(-1)
1654           
1655!--         aggregate csb among processors
1656!--         allocate necessary arrays
1657            ALLOCATE( csflt(ndcsf,max(ncsfl,ndcsf)) )
1658            ALLOCATE( kcsflt(kdcsf,max(ncsfl,kdcsf)) )
1659            ALLOCATE( icsflt(0:numprocs-1) )
1660            ALLOCATE( dcsflt(0:numprocs-1) )
1661            ALLOCATE( ipcsflt(0:numprocs-1) )
1662            ALLOCATE( dpcsflt(0:numprocs-1) )
1663           
1664!--         fill out arrays of csf values and
1665!--         arrays of number of elements and displacements
1666!--         for particular precessors
1667            icsflt = 0
1668            dcsflt = 0
1669            ip = -1
1670            j = -1
1671            d = 0
1672            DO kcsf = 1, ncsfl
1673                j = j+1
1674                IF ( acsf(kcsf)%ip /= ip )  THEN
1675!--                 new block of the processor
1676!--                 number of elements of previous block
1677                    IF ( ip>=0) icsflt(ip) = j
1678                    d = d+j
1679!--                 blank blocks
1680                    DO jp = ip+1, acsf(kcsf)%ip-1
1681!--                     number of elements is zero, displacement is equal to previous
1682                        icsflt(jp) = 0
1683                        dcsflt(jp) = d
1684                    ENDDO
1685!--                 the actual block
1686                    ip = acsf(kcsf)%ip
1687                    dcsflt(ip) = d
1688                    j = 0
1689                ENDIF
1690!--             fill out real values of rsvf, rtransp
1691                csflt(1,kcsf) = acsf(kcsf)%rsvf
1692                csflt(2,kcsf) = acsf(kcsf)%rtransp
1693!--             fill out integer values of itz,ity,itx,isurfs
1694                kcsflt(1,kcsf) = acsf(kcsf)%itz
1695                kcsflt(2,kcsf) = acsf(kcsf)%ity
1696                kcsflt(3,kcsf) = acsf(kcsf)%itx
1697                kcsflt(4,kcsf) = acsf(kcsf)%isurfs
1698            ENDDO
1699!--         last blank blocks at the end of array
1700            j = j+1
1701            IF ( ip>=0 ) icsflt(ip) = j
1702            d = d+j
1703            DO jp = ip+1, numprocs-1
1704!--             number of elements is zero, displacement is equal to previous
1705                icsflt(jp) = 0
1706                dcsflt(jp) = d
1707            ENDDO
1708           
1709!--         deallocate temporary acsf array
1710!--         DEALLOCATE(acsf) - ifort has a problem with deallocation of allocatable target
1711!--         via pointing pointer - we need to test original targets
1712            IF ( ALLOCATED(acsf1) )  THEN
1713                DEALLOCATE(acsf1)
1714            ENDIF
1715            IF ( ALLOCATED(acsf2) )  THEN
1716                DEALLOCATE(acsf2)
1717            ENDIF
1718                   
1719#if defined( __parallel )
1720!--         scatter and gather the number of elements to and from all processor
1721!--         and calculate displacements
1722            CALL mpi_alltoall(icsflt,1,MPI_INTEGER,ipcsflt,1,MPI_INTEGER,comm2d, ierr)
1723           
1724            npcsfl = SUM(ipcsflt)
1725            d = 0
1726            DO i = 0, numprocs-1
1727                dpcsflt(i) = d
1728                d = d + ipcsflt(i)
1729            ENDDO
1730       
1731!--         exchange csf fields between processors
1732            ALLOCATE( pcsflt(ndcsf,max(npcsfl,ndcsf)) )
1733            ALLOCATE( kpcsflt(kdcsf,max(npcsfl,kdcsf)) )
1734            CALL MPI_AlltoAllv(csflt, ndcsf*icsflt, ndcsf*dcsflt, MPI_REAL, &
1735                pcsflt, ndcsf*ipcsflt, ndcsf*dpcsflt, MPI_REAL, comm2d, ierr)
1736            CALL MPI_AlltoAllv(kcsflt, kdcsf*icsflt, kdcsf*dcsflt, MPI_INTEGER, &
1737                kpcsflt, kdcsf*ipcsflt, kdcsf*dpcsflt, MPI_INTEGER, comm2d, ierr)
1738           
1739#else
1740            npcsfl = ncsfl
1741            ALLOCATE( pcsflt(ndcsf,max(npcsfl,ndcsf)) )
1742            ALLOCATE( kpcsflt(kdcsf,max(npcsfl,kdcsf)) )
1743            pcsflt = csflt
1744            kpcsflt = kcsflt
1745#endif
1746
1747!--         deallocate temporary arrays
1748            DEALLOCATE( csflt )
1749            DEALLOCATE( kcsflt )
1750            DEALLOCATE( icsflt )
1751            DEALLOCATE( dcsflt )
1752            DEALLOCATE( ipcsflt )
1753            DEALLOCATE( dpcsflt )
1754
1755!--         sort csf ( a version of quicksort )
1756            CALL quicksort_csf2(kpcsflt, pcsflt, 1, npcsfl)
1757
1758!--         aggregate canopy sink factor records with identical box & source
1759!--         againg across all values from all processors
1760            IF ( npcsfl > 0 )  THEN
1761                icsf = 1 !< reading index
1762                kcsf = 1 !< writing index
1763                DO while (icsf < npcsfl)
1764!--                 here kpcsf(kcsf) already has values from kpcsf(icsf)
1765                    IF ( kpcsflt(3,icsf) == kpcsflt(3,icsf+1)  .AND.  &
1766                         kpcsflt(2,icsf) == kpcsflt(2,icsf+1)  .AND.  &
1767                         kpcsflt(1,icsf) == kpcsflt(1,icsf+1)  .AND.  &
1768                         kpcsflt(4,icsf) == kpcsflt(4,icsf+1) )  THEN
1769!--                     We could simply take either first or second rtransp, both are valid. As a very simple heuristic about which ray
1770!--                     probably passes nearer the center of the target box, we choose DIF from the entry with greater CSF, since that
1771!--                     might mean that the traced beam passes longer through the canopy box.
1772                        IF ( pcsflt(1,kcsf) < pcsflt(1,icsf+1) )  THEN
1773                            pcsflt(2,kcsf) = pcsflt(2,icsf+1)
1774                        ENDIF
1775                        pcsflt(1,kcsf) = pcsflt(1,kcsf) + pcsflt(1,icsf+1)
1776
1777!--                     advance reading index, keep writing index
1778                        icsf = icsf + 1
1779                    ELSE
1780!--                     not identical, just advance and copy
1781                        icsf = icsf + 1
1782                        kcsf = kcsf + 1
1783                        kpcsflt(:,kcsf) = kpcsflt(:,icsf)
1784                        pcsflt(:,kcsf) = pcsflt(:,icsf)
1785                    ENDIF
1786                ENDDO
1787!--             last written item is now also the last item in valid part of array
1788                npcsfl = kcsf
1789            ENDIF
1790           
1791        ENDIF  !< plant_canopy
1792       
1793        CALL location_message( '    calculation of the complete SVF array', .TRUE. )
1794        nsvfcsfl = nsvfl + npcsfl
1795       
1796        ALLOCATE( svf(ndsvf,nsvfcsfl) )
1797        ALLOCATE( svfsurf(ndsvf,nsvfcsfl) )
1798
1799        !< load svf from the structure array to plain arrays
1800        isurflt_prev = -1
1801        ksvf = 1
1802        svfsum = 0._wp
1803        DO isvf = 1, nsvfl
1804!--         normalize svf per target face
1805            IF ( asvf(ksvf)%isurflt /= isurflt_prev )  THEN
1806                IF ( isurflt_prev /= -1  .AND.  svfsum /= 0._wp )  THEN
1807!--                 TODO detect and log when normalization differs too much from 1
1808                    svf(1, isvf_surflt:isvf-1) = svf(1, isvf_surflt:isvf-1) / svfsum
1809                ENDIF
1810                isurflt_prev = asvf(ksvf)%isurflt
1811                isvf_surflt = isvf
1812                svfsum = asvf(ksvf)%rsvf !?? / asvf(ksvf)%rtransp
1813            ELSE
1814                svfsum = svfsum + asvf(ksvf)%rsvf !?? / asvf(ksvf)%rtransp
1815            ENDIF
1816           
1817            svf(:, isvf) = (/ asvf(ksvf)%rsvf, asvf(ksvf)%rtransp /)
1818            svfsurf(:, isvf) = (/ asvf(ksvf)%isurflt, asvf(ksvf)%isurfs /)
1819                   
1820!--         next element
1821            ksvf = ksvf + 1   
1822        ENDDO
1823       
1824        IF ( isurflt_prev /= -1  .AND.  svfsum /= 0._wp )  THEN
1825!--         TODO detect and log when normalization differs too much from 1
1826            svf(1, isvf_surflt:nsvfl) = svf(1, isvf_surflt:nsvfl) / svfsum
1827        ENDIF
1828       
1829!--     deallocate temporary asvf array
1830!--     DEALLOCATE(asvf) - ifort has a problem with deallocation of allocatable target
1831!--     via pointing pointer - we need to test original targets
1832        IF ( ALLOCATED(asvf1) )  THEN
1833            DEALLOCATE(asvf1)
1834        ENDIF
1835        IF ( ALLOCATED(asvf2) )  THEN
1836            DEALLOCATE(asvf2)
1837        ENDIF
1838
1839        IF ( plant_canopy )  THEN
1840            CALL location_message( '    calculation of the complete CSF part of the array', .TRUE. )
1841            IF ( npcsfl > 0 )  THEN
1842                DO isvf = 1, npcsfl
1843                    svf(:,nsvfl+isvf) = pcsflt(:,isvf)
1844                    svfsurf(1,nsvfl+isvf) =  gridpcbl(kpcsflt(1,isvf),kpcsflt(2,isvf),kpcsflt(3,isvf))
1845                    svfsurf(2,nsvfl+isvf) =  kpcsflt(4,isvf)
1846                ENDDO
1847            ENDIF
1848           
1849!--         deallocation of temporary arrays
1850            DEALLOCATE( pcsflt )
1851            DEALLOCATE( kpcsflt )
1852           
1853        ENDIF
1854       
1855        RETURN
1856       
1857301     WRITE( message_string, * )  &
1858            'I/O error when processing shape view factors / ',  &
1859            'plant canopy sink factors / direct irradiance factors.'
1860        CALL message( 'init_urban_surface', 'PA0502', 2, 2, 0, 6, 0 )
1861       
1862    END SUBROUTINE usm_calc_svf
1863
1864
1865!------------------------------------------------------------------------------!
1866!
1867! Description:
1868! ------------
1869!> Subroutine checks variables and assigns units.
1870!> It is caaled out from subroutine check_parameters.
1871!------------------------------------------------------------------------------!
1872    SUBROUTINE usm_check_data_output( variable, unit )
1873       
1874        IMPLICIT NONE
1875 
1876        CHARACTER (len=*),INTENT(IN)    ::  variable !:
1877        CHARACTER (len=*),INTENT(OUT)   ::  unit     !:
1878       
1879        CHARACTER (len=20)              :: var
1880
1881        var = TRIM(variable)
1882        IF ( var(1:11)=='usm_radnet_'  .OR.  var(1:13)=='usm_rad_insw_'  .OR.           &
1883             var(1:13)=='usm_rad_inlw_'  .OR.  var(1:16)=='usm_rad_inswdir_'  .OR.      &
1884             var(1:16)=='usm_rad_inswdif_'  .OR.  var(1:16)=='usm_rad_inswref_'  .OR.   &
1885             var(1:16)=='usm_rad_inlwdif_'  .OR.  var(1:16)=='usm_rad_inlwref_'  .OR.   &
1886             var(1:14)=='usm_rad_outsw_'  .OR.  var(1:14)=='usm_rad_outlw_'  .OR.       &
1887             var(1:11)=='usm_rad_hf_'  .OR.                                             &
1888             var(1:9) =='usm_wshf_'  .OR.  var(1:9)=='usm_wghf_' )  THEN
1889            unit = 'W/m2'
1890        ELSE IF ( var(1:10) =='usm_t_surf'  .OR.  var(1:10) =='usm_t_wall' )  THEN
1891            unit = 'K'
1892        ELSE IF ( var(1:9) =='usm_surfz'  .OR.  var(1:7) =='usm_svf'  .OR.              & 
1893                  var(1:7) =='usm_dif'  .OR.  var(1:11) =='usm_surfcat'  .OR.           &
1894                  var(1:11) =='usm_surfalb'  .OR.  var(1:12) =='usm_surfemis')  THEN
1895            unit = '1'
1896        ELSE IF ( plant_canopy  .AND.  var(1:7) =='usm_lad' )  THEN
1897            unit = 'm2/m3'
1898        ELSE IF ( plant_canopy  .AND.  var(1:14) == 'usm_canopy_khf' )  THEN
1899            unit = 'K/s'
1900        ELSE
1901            unit = 'illegal'
1902        ENDIF
1903
1904    END SUBROUTINE usm_check_data_output
1905
1906
1907!------------------------------------------------------------------------------!
1908! Description:
1909! ------------
1910!> Check parameters routine for urban surface model
1911!------------------------------------------------------------------------------!
1912    SUBROUTINE usm_check_parameters
1913   
1914       USE control_parameters,                                                 &
1915           ONLY:  bc_pt_b, bc_q_b, constant_flux_layer, large_scale_forcing,   &
1916                  lsf_surf, topography
1917
1918!
1919!--    Dirichlet boundary conditions are required as the surface fluxes are
1920!--    calculated from the temperature/humidity gradients in the urban surface
1921!--    model
1922       IF ( bc_pt_b == 'neumann'   .OR.   bc_q_b == 'neumann' )  THEN
1923          message_string = 'urban surface model requires setting of '//        &
1924                           'bc_pt_b = "dirichlet" and '//                      &
1925                           'bc_q_b  = "dirichlet"'
1926          CALL message( 'check_parameters', 'PA0590', 1, 2, 0, 6, 0 )
1927       ENDIF
1928
1929       IF ( .NOT.  constant_flux_layer )  THEN
1930          message_string = 'urban surface model requires '//                   &
1931                           'constant_flux_layer = .T.'
1932          CALL message( 'check_parameters', 'PA0591', 1, 2, 0, 6, 0 )
1933       ENDIF
1934!       
1935!--    Surface forcing has to be disabled for LSF in case of enabled
1936!--    urban surface module
1937       IF ( large_scale_forcing )  THEN
1938          lsf_surf = .FALSE.
1939       ENDIF
1940!
1941!--    Topography
1942       IF ( topography == 'flat' )  THEN
1943          message_string = 'topography /= "flat" is required '//               &
1944                           'when using the urban surface model'
1945          CALL message( 'check_parameters', 'PA0592', 1, 2, 0, 6, 0 )
1946       ENDIF
1947
1948
1949    END SUBROUTINE usm_check_parameters
1950
1951
1952!------------------------------------------------------------------------------!
1953!
1954! Description:
1955! ------------
1956!> Output of the 3D-arrays in netCDF and/or AVS format
1957!> for variables of urban_surface model.
1958!> It resorts the urban surface module output quantities from surf style
1959!> indexing into temporary 3D array with indices (i,j,k).
1960!> It is called from subroutine data_output_3d.
1961!------------------------------------------------------------------------------!
1962    SUBROUTINE usm_data_output_3d( av, variable, found, local_pf, nzb_do, nzt_do )
1963       
1964        IMPLICIT NONE
1965
1966        INTEGER(iwp), INTENT(IN)       ::  av        !<
1967        CHARACTER (len=*), INTENT(IN)  ::  variable  !<
1968        INTEGER(iwp), INTENT(IN)       ::  nzb_do    !< lower limit of the data output (usually 0)
1969        INTEGER(iwp), INTENT(IN)       ::  nzt_do    !< vertical upper limit of the data output (usually nz_do3d)
1970        LOGICAL, INTENT(OUT)           ::  found     !<
1971        REAL(sp), DIMENSION(nxlg:nxrg,nysg:nyng,nzb_do:nzt_do) ::  local_pf   !< sp - it has to correspond to module data_output_3d
1972        REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg)     ::  temp_pf    !< temp array for urban surface output procedure
1973       
1974        CHARACTER (len=20)                                     :: var, surfid
1975        INTEGER(iwp), PARAMETER                                :: nd = 5
1976        CHARACTER(len=6), DIMENSION(0:nd-1), PARAMETER         :: dirname = (/ '_roof ', '_south', '_north', '_west ', '_east ' /)
1977        INTEGER(iwp), DIMENSION(0:nd-1), PARAMETER             :: dirint = (/ iroof, isouth, inorth, iwest, ieast /)
1978        INTEGER(iwp), DIMENSION(0:nd-1)                        :: dirstart
1979        INTEGER(iwp), DIMENSION(0:nd-1)                        :: dirend
1980        INTEGER(iwp)                                           :: ids,isurf,isvf,isurfs,isurflt
1981        INTEGER(iwp)                                           :: is,js,ks,i,j,k,iwl,istat
1982
1983        dirstart = (/ startland, startwall, startwall, startwall, startwall /)
1984        dirend = (/ endland, endwall, endwall, endwall, endwall /)
1985
1986        found = .TRUE.
1987        temp_pf = -1._wp
1988       
1989        ids = -1
1990        var = TRIM(variable)
1991        DO i = 0, nd-1
1992            k = len(TRIM(var))
1993            j = len(TRIM(dirname(i)))
1994            IF ( var(k-j+1:k) == dirname(i) )  THEN
1995                ids = i
1996                var = var(:k-j)
1997                EXIT
1998            ENDIF
1999        ENDDO
2000        IF ( ids == -1 )  THEN
2001            var = TRIM(variable)
2002        ENDIF
2003        IF ( var(1:11) == 'usm_t_wall_'  .AND.  len(TRIM(var)) >= 12 )  THEN
2004!--         wall layers
2005            READ(var(12:12), '(I1)', iostat=istat ) iwl
2006            IF ( istat == 0  .AND.  iwl >= nzb_wall  .AND.  iwl <= nzt_wall )  THEN
2007                var = var(1:10)
2008            ENDIF
2009        ENDIF
2010        IF ( (var(1:8) == 'usm_svf_'  .OR.  var(1:8) == 'usm_dif_')  .AND.  len(TRIM(var)) >= 13 )  THEN
2011!--         svf values to particular surface
2012            surfid = var(9:)
2013            i = index(surfid,'_')
2014            j = index(surfid(i+1:),'_')
2015            READ(surfid(1:i-1),*, iostat=istat ) is
2016            IF ( istat == 0 )  THEN
2017                READ(surfid(i+1:i+j-1),*, iostat=istat ) js
2018            ENDIF
2019            IF ( istat == 0 )  THEN
2020                READ(surfid(i+j+1:),*, iostat=istat ) ks
2021            ENDIF
2022            IF ( istat == 0 )  THEN
2023                var = var(1:7)
2024            ENDIF
2025        ENDIF
2026       
2027        SELECT CASE ( TRIM(var) )
2028
2029          CASE ( 'usm_surfz' )
2030!--           array of lw radiation falling to local surface after i-th reflection
2031              DO isurf = dirstart(ids), dirend(ids)
2032                  IF ( surfl(id,isurf) == ids )  THEN
2033                      IF ( surfl(id,isurf) == iroof )  THEN
2034                          temp_pf(0,surfl(iy,isurf),surfl(ix,isurf)) =             &
2035                                  max(temp_pf(0,surfl(iy,isurf),surfl(ix,isurf)),  &
2036                                      REAL(surfl(iz,isurf),wp))
2037                      ELSE
2038                          temp_pf(0,surfl(iy,isurf),surfl(ix,isurf)) =             &
2039                                  max(temp_pf(0,surfl(iy,isurf),surfl(ix,isurf)),  &
2040                                      REAL(surfl(iz,isurf),wp)+1.0_wp)
2041                      ENDIF
2042                  ENDIF
2043              ENDDO
2044
2045          CASE ( 'usm_surfcat' )
2046!--           surface category
2047              DO isurf = dirstart(ids), dirend(ids)
2048                 IF ( surfl(id,isurf) == ids )  THEN
2049                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surface_types(isurf)
2050                 ENDIF
2051              ENDDO
2052             
2053          CASE ( 'usm_surfalb' )
2054!--           surface albedo
2055              DO isurf = dirstart(ids), dirend(ids)
2056                 IF ( surfl(id,isurf) == ids )  THEN
2057                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = albedo_surf(isurf)
2058                 ENDIF
2059              ENDDO
2060             
2061          CASE ( 'usm_surfemis' )
2062!--           surface albedo
2063              DO isurf = dirstart(ids), dirend(ids)
2064                 IF ( surfl(id,isurf) == ids )  THEN
2065                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = emiss_surf(isurf)
2066                 ENDIF
2067              ENDDO
2068             
2069          CASE ( 'usm_svf', 'usm_dif' )
2070!--           shape view factors or iradiance factors to selected surface
2071              IF ( TRIM(var)=='usm_svf' )  THEN
2072                  k = 1
2073              ELSE
2074                  k = 2
2075              ENDIF
2076              DO isvf = 1, nsvfl
2077                  isurflt = svfsurf(1, isvf)
2078                  isurfs = svfsurf(2, isvf)
2079                             
2080                  IF ( surf(ix,isurfs) == is  .AND.  surf(iy,isurfs) == js  .AND.       &
2081                       surf(iz,isurfs) == ks  .AND.  surf(id,isurfs) == ids )  THEN
2082  !--                 correct source surface
2083                      temp_pf(surfl(iz,isurflt),surfl(iy,isurflt),surfl(ix,isurflt)) = svf(k,isvf)
2084                  ENDIF
2085              ENDDO
2086
2087          CASE ( 'usm_radnet' )
2088!--           array of complete radiation balance
2089              DO isurf = dirstart(ids), dirend(ids)
2090                 IF ( surfl(id,isurf) == ids )  THEN
2091                   IF ( av == 0 )  THEN
2092                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = rad_net_l(isurf)
2093                   ELSE
2094                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = rad_net_av(isurf)
2095                   ENDIF
2096                 ENDIF
2097              ENDDO
2098
2099          CASE ( 'usm_rad_insw' )
2100!--           array of sw radiation falling to surface after i-th reflection
2101              DO isurf = dirstart(ids), dirend(ids)
2102                 IF ( surfl(id,isurf) == ids )  THEN
2103                   IF ( av == 0 )  THEN
2104                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinsw(isurf)
2105                   ELSE
2106                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinsw_av(isurf)
2107                   ENDIF
2108                 ENDIF
2109              ENDDO
2110
2111          CASE ( 'usm_rad_inlw' )
2112!--           array of lw radiation falling to surface after i-th reflection
2113              DO isurf = dirstart(ids), dirend(ids)
2114                 IF ( surfl(id,isurf) == ids )  THEN
2115                   IF ( av == 0 )  THEN
2116                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinlw(isurf)
2117                   ELSE
2118                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinlw_av(isurf)
2119                   ENDIF
2120                 ENDIF
2121              ENDDO
2122
2123          CASE ( 'usm_rad_inswdir' )
2124!--           array of direct sw radiation falling to surface from sun
2125              DO isurf = dirstart(ids), dirend(ids)
2126                 IF ( surfl(id,isurf) == ids )  THEN
2127                   IF ( av == 0 )  THEN
2128                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinswdir(isurf)
2129                   ELSE
2130                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinswdir_av(isurf)
2131                   ENDIF
2132                 ENDIF
2133              ENDDO
2134
2135          CASE ( 'usm_rad_inswdif' )
2136!--           array of difusion sw radiation falling to surface from sky and borders of the domain
2137              DO isurf = dirstart(ids), dirend(ids)
2138                 IF ( surfl(id,isurf) == ids )  THEN
2139                   IF ( av == 0 )  THEN
2140                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinswdif(isurf)
2141                   ELSE
2142                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinswdif_av(isurf)
2143                   ENDIF
2144                 ENDIF
2145              ENDDO
2146
2147          CASE ( 'usm_rad_inswref' )
2148!--           array of sw radiation falling to surface from reflections
2149              DO isurf = dirstart(ids), dirend(ids)
2150                 IF ( surfl(id,isurf) == ids )  THEN
2151                   IF ( av == 0 )  THEN
2152                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = &
2153                       surfinsw(isurf) - surfinswdir(isurf) - surfinswdif(isurf)
2154                   ELSE
2155                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinswref_av(isurf)
2156                   ENDIF
2157                 ENDIF
2158              ENDDO
2159
2160          CASE ( 'usm_rad_inlwdif' )
2161!--           array of sw radiation falling to surface after i-th reflection
2162              DO isurf = dirstart(ids), dirend(ids)
2163                 IF ( surfl(id,isurf) == ids )  THEN
2164                   IF ( av == 0 )  THEN
2165                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinlwdif(isurf)
2166                   ELSE
2167                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinlwdif_av(isurf)
2168                   ENDIF
2169                 ENDIF
2170              ENDDO
2171
2172          CASE ( 'usm_rad_inlwref' )
2173!--           array of lw radiation falling to surface from reflections
2174              DO isurf = dirstart(ids), dirend(ids)
2175                 IF ( surfl(id,isurf) == ids )  THEN
2176                   IF ( av == 0 )  THEN
2177                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinlw(isurf) - surfinlwdif(isurf)
2178                   ELSE
2179                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfinlwref_av(isurf)
2180                   ENDIF
2181                 ENDIF
2182              ENDDO
2183
2184          CASE ( 'usm_rad_outsw' )
2185!--           array of sw radiation emitted from surface after i-th reflection
2186              DO isurf = dirstart(ids), dirend(ids)
2187                 IF ( surfl(id,isurf) == ids )  THEN
2188                   IF ( av == 0 )  THEN
2189                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfoutsw(isurf)
2190                   ELSE
2191                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfoutsw_av(isurf)
2192                   ENDIF
2193                 ENDIF
2194              ENDDO
2195
2196          CASE ( 'usm_rad_outlw' )
2197!--           array of lw radiation emitted from surface after i-th reflection
2198              DO isurf = dirstart(ids), dirend(ids)
2199                 IF ( surfl(id,isurf) == ids )  THEN
2200                   IF ( av == 0 )  THEN
2201                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfoutlw(isurf)
2202                   ELSE
2203                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfoutlw_av(isurf)
2204                   ENDIF
2205                 ENDIF
2206              ENDDO
2207
2208          CASE ( 'usm_rad_hf' )
2209!--           array of heat flux from radiation for surfaces after all reflections
2210              DO isurf = dirstart(ids), dirend(ids)
2211                 IF ( surfl(id,isurf) == ids )  THEN
2212                   IF ( av == 0 )  THEN
2213                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfhf(isurf)
2214                   ELSE
2215                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = surfhf_av(isurf)
2216                   ENDIF
2217                 ENDIF
2218              ENDDO
2219
2220          CASE ( 'usm_wshf' )
2221!--           array of sensible heat flux from surfaces
2222!--           horizontal surfaces
2223              DO isurf = dirstart(ids), dirend(ids)
2224                 IF ( surfl(id,isurf) == ids )  THEN
2225                   IF ( av == 0 )  THEN
2226                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = wshf_eb(isurf)
2227                   ELSE
2228                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = wshf_eb_av(isurf)
2229                   ENDIF
2230                 ENDIF
2231              ENDDO
2232
2233          CASE ( 'usm_wghf' )
2234!--           array of heat flux from ground (land, wall, roof)
2235              DO isurf = dirstart(ids), dirend(ids)
2236                 IF ( surfl(id,isurf) == ids )  THEN
2237                   IF ( av == 0 )  THEN
2238                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = wghf_eb(isurf)
2239                   ELSE
2240                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = wghf_eb_av(isurf)
2241                   ENDIF
2242                 ENDIF
2243              ENDDO
2244
2245          CASE ( 'usm_t_surf' )
2246!--           surface temperature for surfaces
2247              DO isurf = max(startenergy,dirstart(ids)), min(endenergy,dirend(ids))
2248                 IF ( surfl(id,isurf) == ids )  THEN
2249                   IF ( av == 0 )  THEN
2250                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = t_surf(isurf)
2251                   ELSE
2252                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = t_surf_av(isurf)
2253                   ENDIF
2254                 ENDIF
2255              ENDDO
2256             
2257          CASE ( 'usm_t_wall' )
2258!--           wall temperature for  iwl layer of walls and land
2259              DO isurf = dirstart(ids), dirend(ids)
2260                 IF ( surfl(id,isurf) == ids )  THEN
2261                   IF ( av == 0 )  THEN
2262                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = t_wall(iwl,isurf)
2263                   ELSE
2264                     temp_pf(surfl(iz,isurf),surfl(iy,isurf),surfl(ix,isurf)) = t_wall_av(iwl,isurf)
2265                   ENDIF
2266                 ENDIF
2267              ENDDO
2268
2269          CASE ( 'usm_lad' )
2270!--           leaf area density
2271              DO i = nxl, nxr
2272                 DO j = nys, nyn
2273                     DO k = nzb_s_inner(j,i), nzut
2274                         temp_pf(k,j,i) = lad_s(k-nzb_s_inner(j,i),j,i)
2275                     ENDDO
2276                 ENDDO
2277              ENDDO
2278             
2279          CASE ( 'usm_canopy_khf' )
2280!--           canopy kinematic heat flux
2281              DO i = nxl, nxr
2282                 DO j = nys, nyn
2283                     DO k = nzb_s_inner(j,i), nzut
2284                         temp_pf(k,j,i) = pc_heating_rate(k-nzb_s_inner(j,i),j,i)
2285                     ENDDO
2286                 ENDDO
2287              ENDDO
2288             
2289          CASE DEFAULT
2290              found = .FALSE.
2291             
2292        END SELECT
2293       
2294!--     fill out array local_pf which is subsequently treated by data_output_3d
2295        CALL exchange_horiz( temp_pf, nbgp )
2296        DO j = nysg,nyng
2297            DO i = nxlg,nxrg
2298                DO k = nzb_do, nzt_do
2299                    local_pf(i,j,k) = temp_pf(k,j,i)
2300                ENDDO
2301            ENDDO
2302        ENDDO
2303       
2304    END SUBROUTINE usm_data_output_3d
2305   
2306
2307!------------------------------------------------------------------------------!
2308!
2309! Description:
2310! ------------
2311!> Soubroutine defines appropriate grid for netcdf variables.
2312!> It is called out from subroutine netcdf.
2313!------------------------------------------------------------------------------!
2314    SUBROUTINE usm_define_netcdf_grid( variable, found, grid_x, grid_y, grid_z )
2315   
2316        IMPLICIT NONE
2317
2318        CHARACTER (len=*), INTENT(IN)  ::  variable    !<
2319        LOGICAL, INTENT(OUT)           ::  found       !<
2320        CHARACTER (len=*), INTENT(OUT) ::  grid_x      !<
2321        CHARACTER (len=*), INTENT(OUT) ::  grid_y      !<
2322        CHARACTER (len=*), INTENT(OUT) ::  grid_z      !<
2323
2324        CHARACTER (len=20)            :: var
2325
2326        var = TRIM(variable)
2327        IF ( var(1:11)=='usm_radnet_'  .OR.  var(1:13) =='usm_rad_insw_'  .OR.            &
2328             var(1:13) =='usm_rad_inlw_'  .OR.  var(1:16) =='usm_rad_inswdir_'  .OR.      &
2329             var(1:16) =='usm_rad_inswdif_'  .OR.  var(1:16) =='usm_rad_inswref_'  .OR.   &
2330             var(1:16) =='usm_rad_inlwdif_'  .OR.  var(1:16) =='usm_rad_inlwref_'  .OR.   &
2331             var(1:14) =='usm_rad_outsw_'  .OR.  var(1:14) =='usm_rad_outlw_'  .OR.       &
2332             var(1:11) =='usm_rad_hf_'  .OR.                                              &
2333             var(1:9) == 'usm_wshf_'  .OR.  var(1:9)== 'usm_wghf_'  .OR.                  &
2334             var(1:10) == 'usm_t_surf'  .OR.  var(1:10) == 'usm_t_wall'  .OR.             &
2335             var(1:9) == 'usm_surfz'  .OR.  var(1:7) == 'usm_svf'  .OR.                   & 
2336             var(1:7) =='usm_dif'  .OR.  var(1:11) =='usm_surfcat'  .OR.                  &
2337             var(1:11) =='usm_surfalb'  .OR.  var(1:12) =='usm_surfemis'  .OR.            &
2338             var(1:7) == 'usm_lad'  .OR.  var(1:14) == 'usm_canopy_khf' )  THEN
2339
2340            found = .TRUE.
2341            grid_x = 'x'
2342            grid_y = 'y'
2343            grid_z = 'zu'
2344        ELSE
2345            found  = .FALSE.
2346            grid_x = 'none'
2347            grid_y = 'none'
2348            grid_z = 'none'
2349        ENDIF
2350
2351    END SUBROUTINE usm_define_netcdf_grid
2352   
2353   
2354!------------------------------------------------------------------------------!
2355!> Finds first model boundary crossed by a ray
2356!------------------------------------------------------------------------------!
2357    PURE SUBROUTINE usm_find_boundary_face(origin, uvect, bdycross)
2358        IMPLICIT NONE
2359        REAL(wp), DIMENSION(3), INTENT(in)      :: origin    !< ray origin
2360        REAL(wp), DIMENSION(3), INTENT(in)      :: uvect     !< ray unit vector
2361        INTEGER(iwp), DIMENSION(4), INTENT(out) :: bdycross  !< found boundary crossing (d, z, y, x)
2362        REAL(wp), DIMENSION(3)                  :: crossdist !< crossing distance
2363        INTEGER(iwp), DIMENSION(3)              :: bdyd      !< boundary direction
2364        REAL(wp)                                :: bdydim    !<
2365        REAL(wp)                                :: dist      !<
2366        INTEGER(iwp)                            :: seldim    !< found fist crossing index
2367        INTEGER(iwp)                            :: d         !<
2368
2369        bdydim = nzut + .5_wp !< top boundary
2370        bdyd(1) = isky
2371        crossdist(1) = (bdydim - origin(1)) / uvect(1)
2372
2373        IF ( uvect(2) >= 0._wp )  THEN
2374            bdydim = ny + .5_wp !< north global boundary
2375            bdyd(2) = inorthb
2376        ELSE
2377            bdydim = -.5_wp !< south global boundary
2378            bdyd(2) = isouthb
2379        ENDIF
2380        crossdist(2) = (bdydim - origin(2)) / uvect(2)
2381
2382        IF ( uvect(3) >= 0._wp )  THEN
2383            bdydim = nx + .5_wp !< east global boundary
2384            bdyd(3) = ieastb
2385        ELSE
2386            bdydim = -.5_wp !< west global boundary
2387            bdyd(3) = iwestb
2388        ENDIF
2389        crossdist(3) = (bdydim - origin(3)) / uvect(3)
2390
2391        seldim = minloc(crossdist, 1)
2392        dist = crossdist(seldim)
2393        d = bdyd(seldim)
2394
2395        bdycross(1) = d
2396        bdycross(2:4) = NINT( origin(:) + uvect(:)*dist &
2397                        + .5_wp * (/ kdir(d), jdir(d), idir(d) /) )
2398    END SUBROUTINE
2399
2400
2401!------------------------------------------------------------------------------!
2402!> Determines whether two faces are oriented towards each other
2403!------------------------------------------------------------------------------!
2404    PURE LOGICAL FUNCTION usm_facing(x, y, z, d, x2, y2, z2, d2)
2405        IMPLICIT NONE
2406        INTEGER(iwp),   INTENT(in)  :: x, y, z, d, x2, y2, z2, d2
2407     
2408        usm_facing = .FALSE.
2409        IF ( d==iroof  .AND.  d2==iroof ) RETURN
2410        IF ( d==isky  .AND.  d2==isky ) RETURN
2411        IF ( (d==isouth  .OR.  d==inorthb)  .AND.  (d2==isouth.OR.d2==inorthb) ) RETURN
2412        IF ( (d==inorth  .OR.  d==isouthb)  .AND.  (d2==inorth.OR.d2==isouthb) ) RETURN
2413        IF ( (d==iwest  .OR.  d==ieastb)  .AND.  (d2==iwest.OR.d2==ieastb) ) RETURN
2414        IF ( (d==ieast  .OR.  d==iwestb)  .AND.  (d2==ieast.OR.d2==iwestb) ) RETURN
2415
2416        SELECT CASE (d)
2417            CASE (iroof)                   !< ground, roof
2418                IF ( z2 < z ) RETURN
2419            CASE (isky)                    !< sky
2420                IF ( z2 > z ) RETURN
2421            CASE (isouth, inorthb)         !< south facing
2422                IF ( y2 > y ) RETURN
2423            CASE (inorth, isouthb)         !< north facing
2424                IF ( y2 < y ) RETURN
2425            CASE (iwest, ieastb)           !< west facing
2426                IF ( x2 > x ) RETURN
2427            CASE (ieast, iwestb)           !< east facing
2428                IF ( x2 < x ) RETURN
2429        END SELECT
2430
2431        SELECT CASE (d2)
2432            CASE (iroof)                   !< ground, roof
2433                IF ( z < z2 ) RETURN
2434            CASE (isky)                    !< sky
2435                IF ( z > z2 ) RETURN
2436            CASE (isouth, inorthb)         !< south facing
2437                IF ( y > y2 ) RETURN
2438            CASE (inorth, isouthb)         !< north facing
2439                IF ( y < y2 ) RETURN
2440            CASE (iwest, ieastb)           !< west facing
2441                IF ( x > x2 ) RETURN
2442            CASE (ieast, iwestb)           !< east facing
2443                IF ( x < x2 ) RETURN
2444            CASE (-1)
2445                CONTINUE
2446        END SELECT
2447
2448        usm_facing = .TRUE.
2449       
2450    END FUNCTION usm_facing
2451   
2452
2453!------------------------------------------------------------------------------!
2454! Description:
2455! ------------
2456!> Initialization of the wall surface model
2457!------------------------------------------------------------------------------!
2458    SUBROUTINE usm_init_material_model
2459
2460        IMPLICIT NONE
2461
2462        INTEGER(iwp) ::  k, l            !< running indices
2463       
2464        CALL location_message( '    initialization of wall surface model', .TRUE. )
2465       
2466!--     Calculate wall grid spacings.
2467!--     Temperature is defined at the center of the wall layers,
2468!--     whereas gradients/fluxes are defined at the edges (_stag)
2469        DO l = nzb_wall, nzt_wall
2470           zwn(l) = zwn_default(l)
2471        ENDDO
2472       
2473!--     apply for all particular wall grids
2474        DO l = startenergy, endenergy
2475           zw(:,l) = zwn(:) * thickness_wall(l)
2476           dz_wall(nzb_wall,l) = zw(nzb_wall,l)
2477           DO k = nzb_wall+1, nzt_wall
2478               dz_wall(k,l) = zw(k,l) - zw(k-1,l)
2479           ENDDO
2480           
2481           dz_wall(nzt_wall+1,l) = dz_wall(nzt_wall,l)
2482
2483           DO k = nzb_wall, nzt_wall-1
2484               dz_wall_stag(k,l) = 0.5 * (dz_wall(k+1,l) + dz_wall(k,l))
2485           ENDDO
2486           dz_wall_stag(nzt_wall,l) = dz_wall(nzt_wall,l)
2487        ENDDO
2488       
2489        ddz_wall      = 1.0_wp / dz_wall
2490        ddz_wall_stag = 1.0_wp / dz_wall_stag
2491       
2492        CALL location_message( '    wall structures filed out', .TRUE. )
2493
2494        CALL location_message( '    initialization of wall surface model finished', .TRUE. )
2495
2496    END SUBROUTINE usm_init_material_model
2497
2498 
2499!------------------------------------------------------------------------------!
2500! Description:
2501! ------------
2502!> Initialization of the urban surface model
2503!------------------------------------------------------------------------------!
2504    SUBROUTINE usm_init_urban_surface
2505   
2506        IMPLICIT NONE
2507
2508        INTEGER(iwp) ::  i, j, k, l            !< running indices
2509        REAL(wp)     ::  c, d, tin, exn
2510       
2511        CALL cpu_log( log_point_s(78), 'usm_init', 'start' )
2512!--     surface forcing have to be disabled for LSF
2513!--     in case of enabled urban surface module
2514        IF ( large_scale_forcing )  THEN
2515            lsf_surf = .FALSE.
2516        ENDIF
2517       
2518!--     init anthropogenic sources of heat
2519        CALL usm_allocate_urban_surface()
2520       
2521!--     read the surface_types array somewhere
2522        CALL usm_read_urban_surface_types()
2523       
2524!--     init material heat model
2525        CALL usm_init_material_model()
2526       
2527        IF ( usm_anthropogenic_heat )  THEN
2528!--         init anthropogenic sources of heat (from transportation for now)
2529            CALL usm_read_anthropogenic_heat()
2530        ENDIF
2531       
2532        IF ( read_svf_on_init )  THEN
2533!--         read svf and svfsurf data from file
2534            CALL location_message( '    Start reading SVF from file', .TRUE. )
2535            CALL usm_read_svf_from_file()
2536            CALL location_message( '    Reading SVF from file has finished', .TRUE. )
2537        ELSE
2538!--         calculate SFV and CSF
2539            CALL location_message( '    Start calculation of SVF', .TRUE. )
2540            CALL cpu_log( log_point_s(79), 'usm_calc_svf', 'start' )
2541            CALL usm_calc_svf()
2542            CALL cpu_log( log_point_s(79), 'usm_calc_svf', 'stop' )
2543            CALL location_message( '    Calculation of SVF has finished', .TRUE. )
2544        ENDIF
2545
2546        IF ( write_svf_on_init )  THEN
2547!--         write svf and svfsurf data to file
2548            CALL usm_write_svf_to_file()
2549        ENDIF
2550       
2551        IF ( plant_canopy )  THEN
2552!--         gridpcbl was only necessary for initialization
2553            DEALLOCATE( gridpcbl )
2554            IF ( .NOT.  ALLOCATED(pc_heating_rate) )  THEN
2555!--             then pc_heating_rate is allocated in init_plant_canopy
2556!--             in case of cthf /= 0 => we need to allocate it for our use here
2557                ALLOCATE( pc_heating_rate(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2558            ENDIF
2559        ENDIF
2560
2561!--     Intitialization of the surface and wall/ground/roof temperature
2562
2563!--     Initialization for restart runs
2564        IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
2565
2566!--         restore data from restart file
2567            CALL usm_read_restart_data()
2568        ELSE
2569       
2570!--         Calculate initial surface temperature
2571            exn = ( surface_pressure / 1000.0_wp )**0.286_wp
2572
2573            DO l = startenergy, endenergy
2574                k = surfl(iz,l)
2575                j = surfl(iy,l)
2576                i = surfl(ix,l)
2577
2578!--              Initial surface temperature set from pt of adjacent gridbox
2579                t_surf(l) = pt(k,j,i) * exn
2580            ENDDO
2581     
2582!--         initial values for t_wall
2583!--         outer value is set to surface temperature
2584!--         inner value is set to wall_inner_temperature
2585!--         and profile is logaritmic (linear in nz)
2586            DO l = startenergy, endenergy
2587                IF ( isroof_surf(l) )  THEN
2588                    tin = roof_inner_temperature
2589                ELSE IF ( surf(id,l)==iroof )  THEN
2590                    tin = soil_inner_temperature
2591                ELSE
2592                    tin = wall_inner_temperature
2593                ENDIF
2594                DO k = nzb_wall, nzt_wall+1
2595                    c = REAL(k-nzb_wall,wp)/REAL(nzt_wall+1-nzb_wall,wp)
2596                    t_wall(k,:) = (1.0_wp-c)*t_surf(:) + c*tin
2597                ENDDO
2598            ENDDO
2599        ENDIF
2600       
2601!--   
2602!--        Possibly DO user-defined actions (e.g. define heterogeneous wall surface)
2603        CALL user_init_urban_surface
2604
2605!--     initialize prognostic values for the first timestep
2606        t_surf_p = t_surf
2607        t_wall_p = t_wall
2608       
2609!--     Adjust radiative fluxes for urban surface at model start
2610        CALL usm_radiation
2611       
2612        CALL cpu_log( log_point_s(78), 'usm_init', 'stop' )
2613       
2614       
2615    END SUBROUTINE usm_init_urban_surface
2616
2617
2618!------------------------------------------------------------------------------!
2619! Description:
2620! ------------
2621!
2622!> Wall model as part of the urban surface model. The model predicts wall
2623!> temperature.
2624!------------------------------------------------------------------------------!
2625    SUBROUTINE usm_material_heat_model
2626
2627
2628        IMPLICIT NONE
2629
2630        INTEGER(iwp) ::  i,j,k,l,kw                      !< running indices
2631
2632        REAL(wp), DIMENSION(nzb_wall:nzt_wall) :: wtend  !< tendency
2633
2634                                               
2635        DO l = startenergy, endenergy
2636!--         calculate frequently used parameters
2637            k = surfl(iz,l)
2638            j = surfl(iy,l)
2639            i = surfl(ix,l)
2640
2641            !
2642!--         prognostic equation for ground/wall/roof temperature t_wall
2643            wtend(:) = 0.0_wp
2644            wtend(nzb_wall) = (1.0_wp/rho_c_wall(nzb_wall,l)) *                     &
2645                       ( lambda_h(nzb_wall,l) * ( t_wall(nzb_wall+1,l)              &
2646                         - t_wall(nzb_wall,l) ) * ddz_wall(nzb_wall+1,l)            &
2647                         + wghf_eb(l) ) * ddz_wall_stag(nzb_wall,l)
2648           
2649            DO  kw = nzb_wall+1, nzt_wall
2650                wtend(kw) = (1.0_wp/rho_c_wall(kw,l))                               &
2651                              * (   lambda_h(kw,l)                                  &
2652                                 * ( t_wall(kw+1,l) - t_wall(kw,l) )                &
2653                                 * ddz_wall(kw+1,l)                                 &
2654                              - lambda_h(kw-1,l)                                    &
2655                                 * ( t_wall(kw,l) - t_wall(kw-1,l) )                &
2656                                 * ddz_wall(kw,l)                                   &
2657                              ) * ddz_wall_stag(kw,l)
2658            ENDDO
2659
2660            t_wall_p(nzb_wall:nzt_wall,l) = t_wall(nzb_wall:nzt_wall,l)             &
2661                                             + dt_3d * ( tsc(2)                     &
2662                                             * wtend(nzb_wall:nzt_wall) + tsc(3)    &
2663                                             * tt_wall_m(nzb_wall:nzt_wall,l) )   
2664           
2665            !
2666!--         calculate t_wall tendencies for the next Runge-Kutta step
2667            IF ( timestep_scheme(1:5) == 'runge' )  THEN
2668                IF ( intermediate_timestep_count == 1 )  THEN
2669                   DO  kw = nzb_wall, nzt_wall
2670                      tt_wall_m(kw,l) = wtend(kw)
2671                   ENDDO
2672                ELSEIF ( intermediate_timestep_count <                              &
2673                         intermediate_timestep_count_max )  THEN
2674                    DO  kw = nzb_wall, nzt_wall
2675                        tt_wall_m(kw,l) = -9.5625_wp * wtend(kw) + 5.3125_wp        &
2676                                         * tt_wall_m(kw,l)
2677                    ENDDO
2678                ENDIF
2679            ENDIF
2680        ENDDO
2681
2682    END SUBROUTINE usm_material_heat_model
2683
2684
2685!------------------------------------------------------------------------------!
2686! Description:
2687! ------------
2688!> This subroutine calculates interaction of the solar radiation
2689!> with urban surface and updates surface, roofs and walls heatfluxes.
2690!> It also updates rad_sw_out and rad_lw_out.
2691!------------------------------------------------------------------------------!
2692    SUBROUTINE usm_radiation
2693   
2694        IMPLICIT NONE
2695       
2696        INTEGER(iwp)               :: i, j, k, kk, is, js, d, ku, refstep
2697        INTEGER(iwp)               :: nzubl, nzutl, isurf, isurfsrc, isurf1, isvf, ipcgb
2698        INTEGER(iwp), DIMENSION(4) :: bdycross
2699        REAL(wp), DIMENSION(3,3)   :: mrot            !< grid rotation matrix (xyz)
2700        REAL(wp), DIMENSION(3,0:9) :: vnorm           !< face direction normal vectors (xyz)
2701        REAL(wp), DIMENSION(3)     :: sunorig         !< grid rotated solar direction unit vector (xyz)
2702        REAL(wp), DIMENSION(3)     :: sunorig_grid    !< grid squashed solar direction unit vector (zyx)
2703        REAL(wp), DIMENSION(0:9)   :: costheta        !< direct irradiance factor of solar angle
2704        REAL(wp), DIMENSION(nzub:nzut) :: pchf_prep   !< precalculated factor for canopy temp tendency
2705        REAL(wp), PARAMETER        :: alpha = 0._wp   !< grid rotation (TODO: add to namelist or remove)
2706        REAL(wp)                   :: rx, ry, rz
2707        REAL(wp)                   :: pc_box_area, pc_abs_frac, pc_abs_eff
2708        INTEGER(iwp)               :: pc_box_dimshift !< transform for best accuracy
2709       
2710       
2711        IF ( plant_canopy )  THEN
2712            pchf_prep(:) = r_d * (hyp(nzub:nzut) / 100000.0_wp)**0.286_wp &
2713                        / (cp * hyp(nzub:nzut) * dx*dy*dz) !< equals to 1 / (rho * c_p * Vbox * T)
2714        ENDIF
2715
2716        sun_direction = .TRUE.
2717        CALL calc_zenith  !< required also for diffusion radiation
2718
2719!--     prepare rotated normal vectors and irradiance factor
2720        vnorm(1,:) = idir(:)
2721        vnorm(2,:) = jdir(:)
2722        vnorm(3,:) = kdir(:)
2723        mrot(1, :) = (/ cos(alpha), -sin(alpha), 0._wp /)
2724        mrot(2, :) = (/ sin(alpha),  cos(alpha), 0._wp /)
2725        mrot(3, :) = (/ 0._wp,       0._wp,      1._wp /)
2726        sunorig = (/ sun_dir_lon, sun_dir_lat, zenith(0) /)
2727        sunorig = matmul(mrot, sunorig)
2728        DO d = 0, 9
2729            costheta(d) = dot_product(sunorig, vnorm(:,d))
2730        ENDDO
2731       
2732        IF ( zenith(0) > 0 )  THEN
2733!--         now we will "squash" the sunorig vector by grid box size in
2734!--         each dimension, so that this new direction vector will allow us
2735!--         to traverse the ray path within grid coordinates directly
2736            sunorig_grid = (/ sunorig(3)/dz, sunorig(2)/dy, sunorig(1)/dx /)
2737!--         sunorig_grid = sunorig_grid / norm2(sunorig_grid)
2738            sunorig_grid = sunorig_grid / SQRT(SUM(sunorig_grid**2))
2739
2740            IF ( plant_canopy )  THEN
2741!--            precompute effective box depth with prototype Leaf Area Density
2742               pc_box_dimshift = maxloc(sunorig, 1) - 1
2743               CALL usm_box_absorb(cshift((/dx,dy,dz/), pc_box_dimshift),      &
2744                                   60, prototype_lad,                          &
2745                                   cshift(sunorig, pc_box_dimshift),           &
2746                                   pc_box_area, pc_abs_frac)
2747               pc_box_area = pc_box_area * sunorig(pc_box_dimshift+1) / sunorig(3)
2748               pc_abs_eff = log(1._wp - pc_abs_frac) / prototype_lad
2749            ENDIF
2750        ENDIF
2751       
2752!--     split diffusion and direct part of the solar downward radiation
2753!--     comming from radiation model and store it in 2D arrays
2754!--     rad_sw_in_diff, rad_sw_in_dir and rad_lw_in_diff
2755        IF ( split_diffusion_radiation )  THEN
2756            CALL usm_calc_diffusion_radiation
2757        ELSE
2758            rad_sw_in_diff = 0.0_wp
2759            rad_sw_in_dir(:,:)  = rad_sw_in(0,:,:)
2760            rad_lw_in_diff(:,:) = rad_lw_in(0,:,:)
2761        ENDIF
2762
2763!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2764!--     First pass: direct + diffuse irradiance
2765!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2766        surfinswdir   = 0._wp
2767        surfinswdif   = 0._wp
2768        surfinlwdif   = 0._wp
2769        surfins   = 0._wp
2770        surfinl   = 0._wp
2771        surfoutsl    = 0._wp
2772        surfoutll    = 0._wp
2773       
2774!--     Set up thermal radiation from surfaces
2775!--     emiss_surf is defined only for surfaces for which energy balance is calculated
2776        surfoutll(startenergy:endenergy) = emiss_surf(startenergy:endenergy) * sigma_sb   &
2777                                           * t_surf(startenergy:endenergy)**4
2778       
2779#if defined( __parallel )
2780!--     might be optimized and gather only values relevant for current processor
2781        CALL MPI_AllGatherv(surfoutll, nenergy, MPI_REAL, &
2782                            surfoutl, nsurfs, surfstart, MPI_REAL, comm2d, ierr)
2783#else
2784        surfoutl(:) = surfoutll(:)
2785#endif
2786       
2787        isurf1 = -1   !< previous processed surface
2788        DO isvf = 1, nsvfl
2789            isurf = svfsurf(1, isvf)
2790            k = surfl(iz, isurf)
2791            j = surfl(iy, isurf)
2792            i = surfl(ix, isurf)
2793            isurfsrc = svfsurf(2, isvf)
2794            IF ( zenith(0) > 0  .AND.  isurf /= isurf1 )  THEN
2795!--             locate the virtual surface where the direct solar ray crosses domain boundary
2796!--             (once per target surface)
2797                d = surfl(id, isurf)
2798                rz = REAL(k, wp) - 0.5_wp * kdir(d)
2799                ry = REAL(j, wp) - 0.5_wp * jdir(d)
2800                rx = REAL(i, wp) - 0.5_wp * idir(d)
2801               
2802                CALL usm_find_boundary_face( (/ rz, ry, rx /), sunorig_grid, bdycross)
2803               
2804                isurf1 = isurf
2805            ENDIF
2806
2807            IF ( surf(id, isurfsrc) >= isky )  THEN
2808!--             diffuse rad from boundary surfaces. Since it is a simply
2809!--             calculated value, it is not assigned to surfref(s/l),
2810!--             instead it is used directly here
2811!--             we consider the radiation from the radiation model falling on surface
2812!--             as the radiation falling on the top of urban layer into the place of the source surface
2813!--             we consider it as a very reasonable simplification which allow as avoid
2814!--             necessity of other global range arrays and some all to all mpi communication
2815                surfinswdif(isurf) = surfinswdif(isurf) + rad_sw_in_diff(j,i) * svf(1,isvf) * svf(2,isvf)
2816                                                                !< canopy shading is applied only to shortwave
2817                surfinlwdif(isurf) = surfinlwdif(isurf) + rad_lw_in_diff(j,i) * svf(1,isvf)
2818            ELSE
2819!--             for surface-to-surface factors we calculate thermal radiation in 1st pass
2820                surfinl(isurf) = surfinl(isurf) + svf(1,isvf) * surfoutl(isurfsrc)
2821            ENDIF
2822           
2823            IF ( zenith(0) > 0  .AND.  all( surf(:, isurfsrc) == bdycross ) )  THEN
2824!--             found svf between model boundary and the face => face isn't shaded
2825                surfinswdir(isurf) = rad_sw_in_dir(j, i) &
2826                    * costheta(surfl(id, isurf)) * svf(2,isvf) / zenith(0)
2827
2828            ENDIF
2829        ENDDO
2830
2831        IF ( plant_canopy )  THEN
2832       
2833            pcbinsw(:) = 0._wp
2834            pcbinlw(:) = 0._wp  !< will stay always 0 since we don't absorb lw anymore
2835            !
2836!--         pcsf first pass
2837            isurf1 = -1  !< previous processed pcgb
2838            DO isvf = nsvfl+1, nsvfcsfl
2839                ipcgb = svfsurf(1, isvf)
2840                i = pcbl(ix,ipcgb)
2841                j = pcbl(iy,ipcgb)
2842                k = pcbl(iz,ipcgb)
2843                isurfsrc = svfsurf(2, isvf)
2844
2845                IF ( zenith(0) > 0  .AND.  ipcgb /= isurf1 )  THEN
2846!--                 locate the virtual surface where the direct solar ray crosses domain boundary
2847!--                 (once per target PC gridbox)
2848                    rz = REAL(k, wp)
2849                    ry = REAL(j, wp)
2850                    rx = REAL(i, wp)
2851                    CALL usm_find_boundary_face( (/ rz, ry, rx /), &
2852                        sunorig_grid, bdycross)
2853
2854                    isurf1 = ipcgb
2855                ENDIF
2856
2857                IF ( surf(id, isurfsrc) >= isky )  THEN
2858!--                 Diffuse rad from boundary surfaces. See comments for svf above.
2859                    pcbinsw(ipcgb) = pcbinsw(ipcgb) + svf(1,isvf) * svf(2,isvf) * rad_sw_in_diff(j,i)
2860!--                 canopy shading is applied only to shortwave, therefore no absorbtion for lw
2861!--                 pcbinlw(ipcgb) = pcbinlw(ipcgb) + svf(1,isvf) * rad_lw_in_diff(j,i)
2862                !ELSE
2863!--                 Thermal radiation in 1st pass
2864!--                 pcbinlw(ipcgb) = pcbinlw(ipcgb) + svf(1,isvf) * surfoutl(isurfsrc)
2865                ENDIF
2866
2867                IF ( zenith(0) > 0  .AND.  all( surf(:, isurfsrc) == bdycross ) )  THEN
2868!--                 found svf between model boundary and the pcgb => pcgb isn't shaded
2869                    pc_abs_frac = 1._wp - exp(pc_abs_eff * lad_s(k,j,i))
2870                    pcbinsw(ipcgb) = pcbinsw(ipcgb) &
2871                        + rad_sw_in_dir(j, i) * pc_box_area * svf(2,isvf) * pc_abs_frac
2872                ENDIF
2873            ENDDO
2874        ENDIF
2875        surfins(startenergy:endenergy) = surfinswdir(startenergy:endenergy) + surfinswdif(startenergy:endenergy)
2876        surfinl(startenergy:endenergy) = surfinl(startenergy:endenergy) + surfinlwdif(startenergy:endenergy)
2877        surfinsw(:) = surfins(:)
2878        surfinlw(:) = surfinl(:)
2879        surfoutsw(:) = 0.0_wp
2880        surfoutlw(:) = surfoutll(:)
2881        surfhf(startenergy:endenergy) = surfinsw(startenergy:endenergy) + surfinlw(startenergy:endenergy) &
2882                                      - surfoutsw(startenergy:endenergy) - surfoutlw(startenergy:endenergy)
2883       
2884!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2885!--     Next passes - reflections
2886!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2887        DO refstep = 1, nrefsteps
2888       
2889            surfoutsl(startenergy:endenergy) = albedo_surf(startenergy:endenergy) * surfins(startenergy:endenergy)
2890!--         for non-transparent surfaces, longwave albedo is 1 - emissivity
2891            surfoutll(startenergy:endenergy) = (1._wp - emiss_surf(startenergy:endenergy)) * surfinl(startenergy:endenergy)
2892
2893#if defined( __parallel )
2894            CALL MPI_AllGatherv(surfoutsl, nsurfl, MPI_REAL, &
2895                surfouts, nsurfs, surfstart, MPI_REAL, comm2d, ierr)
2896            CALL MPI_AllGatherv(surfoutll, nsurfl, MPI_REAL, &
2897                surfoutl, nsurfs, surfstart, MPI_REAL, comm2d, ierr)
2898#else
2899            surfouts(:) = surfoutsl(:)
2900            surfoutl(:) = surfoutll(:)
2901#endif
2902
2903!--         reset for next pass input
2904            surfins(:) = 0._wp
2905            surfinl(:) = 0._wp
2906           
2907!--         reflected radiation
2908            DO isvf = 1, nsvfl
2909                isurf = svfsurf(1, isvf)
2910                isurfsrc = svfsurf(2, isvf)
2911
2912!--             TODO: to remove if, use start+end for isvf
2913                IF ( surf(id, isurfsrc) < isky )  THEN
2914                    surfins(isurf) = surfins(isurf) + svf(1,isvf) * svf(2,isvf) * surfouts(isurfsrc)
2915                    surfinl(isurf) = surfinl(isurf) + svf(1,isvf) * surfoutl(isurfsrc)
2916                ENDIF
2917            ENDDO
2918
2919!--         radiation absorbed by plant canopy
2920            DO isvf = nsvfl+1, nsvfcsfl
2921                ipcgb = svfsurf(1, isvf)
2922                isurfsrc = svfsurf(2, isvf)
2923
2924                IF ( surf(id, isurfsrc) < isky )  THEN
2925                    pcbinsw(ipcgb) = pcbinsw(ipcgb) + svf(1,isvf) * svf(2,isvf) * surfouts(isurfsrc)
2926!--                 pcbinlw(ipcgb) = pcbinlw(ipcgb) + svf(1,isvf) * surfoutl(isurfsrc)
2927                ENDIF
2928            ENDDO
2929           
2930            surfinsw(:) = surfinsw(:)  + surfins(:)
2931            surfinlw(:) = surfinlw(:)  + surfinl(:)
2932            surfoutsw(startenergy:endenergy) = surfoutsw(startenergy:endenergy) + surfoutsl(startenergy:endenergy)
2933            surfoutlw(startenergy:endenergy) = surfoutlw(startenergy:endenergy) + surfoutll(startenergy:endenergy)
2934            surfhf(startenergy:endenergy) = surfinsw(startenergy:endenergy) + surfinlw(startenergy:endenergy) &
2935                                          - surfoutsw(startenergy:endenergy) - surfoutlw(startenergy:endenergy)
2936       
2937        ENDDO
2938
2939!--     push heat flux absorbed by plant canopy to respective 3D arrays
2940        IF ( plant_canopy )  THEN
2941            pc_heating_rate(:,:,:) = 0._wp
2942            DO ipcgb = 1, npcbl
2943                j = pcbl(iy, ipcgb)
2944                i = pcbl(ix, ipcgb)
2945                k = pcbl(iz, ipcgb)
2946                kk = k - nzb_s_inner(j,i)  !- lad arrays are defined flat
2947                pc_heating_rate(kk, j, i) = (pcbinsw(ipcgb) + pcbinlw(ipcgb)) &
2948                    * pchf_prep(k) * pt(k, j, i) !-- = dT/dt
2949            ENDDO
2950        ENDIF
2951
2952!--     return surface radiation to horizontal surfaces
2953!--     to rad_sw_in, rad_lw_in and rad_net for outputs
2954        !!!!!!!!!!
2955!--     we need the original radiation on urban top layer
2956!--     for calculation of MRT so we can't do adjustment here for now
2957        !!!!!!!!!!
2958        !!!DO isurf = 1, nsurfl
2959        !!!    i = surfl(ix,isurf)
2960        !!!    j = surfl(iy,isurf)
2961        !!!    k = surfl(iz,isurf)
2962        !!!    d = surfl(id,isurf)
2963        !!!    IF ( d==iroof )  THEN
2964        !!!        rad_sw_in(:,j,i) = surfinsw(isurf)
2965        !!!        rad_lw_in(:,j,i) = surfinlw(isurf)
2966        !!!        rad_net(j,i) = rad_sw_in(k,j,i) - rad_sw_out(k,j,i) + rad_lw_in(k,j,i) - rad_lw_out(k,j,i)
2967        !!!    ENDIF
2968        !!!ENDDO
2969
2970    END SUBROUTINE usm_radiation
2971
2972   
2973!------------------------------------------------------------------------------!
2974! Description:
2975! ------------
2976!> Raytracing for detecting obstacles and calculating compound canopy sink
2977!> factors. (A simple obstacle detection would only need to process faces in
2978!> 3 dimensions without any ordering.)
2979!> Assumtions:
2980!> -----------
2981!> 1. The ray always originates from a face midpoint (only one coordinate equals
2982!>    *.5, i.e. wall) and doesn't travel parallel to the surface (that would mean
2983!>    shape factor=0). Therefore, the ray may never travel exactly along a face
2984!>    or an edge.
2985!> 2. From grid bottom to urban surface top the grid has to be *equidistant*
2986!>    within each of the dimensions, including vertical (but the resolution
2987!>    doesn't need to be the same in all three dimensions).
2988!------------------------------------------------------------------------------!
2989    SUBROUTINE usm_raytrace(src, targ, isrc, rirrf, atarg, create_csf, visible, transparency, win_lad)
2990        IMPLICIT NONE
2991
2992        REAL(wp), DIMENSION(3), INTENT(in)     :: src, targ    !< real coordinates z,y,x
2993        INTEGER(iwp), INTENT(in)               :: isrc         !< index of source face for csf
2994        REAL(wp), INTENT(in)                   :: rirrf        !< irradiance factor for csf
2995        REAL(wp), INTENT(in)                   :: atarg        !< target surface area for csf
2996        LOGICAL, INTENT(in)                    :: create_csf   !< whether to generate new CSFs during raytracing
2997        LOGICAL, INTENT(out)                   :: visible
2998        REAL(wp), INTENT(out)                  :: transparency !< along whole path
2999        INTEGER(iwp), INTENT(in)               :: win_lad
3000        INTEGER(iwp)                           :: k, d
3001        INTEGER(iwp)                           :: seldim       !< dimension to be incremented
3002        INTEGER(iwp)                           :: ncsb         !< no of written plant canopy sinkboxes
3003        INTEGER(iwp)                           :: maxboxes     !< max no of gridboxes visited
3004        REAL(wp)                               :: distance     !< euclidean along path
3005        REAL(wp)                               :: crlen        !< length of gridbox crossing
3006        REAL(wp)                               :: lastdist     !< beginning of current crossing
3007        REAL(wp)                               :: nextdist     !< end of current crossing
3008        REAL(wp)                               :: realdist     !< distance in meters per unit distance
3009        REAL(wp)                               :: crmid        !< midpoint of crossing
3010        REAL(wp)                               :: cursink      !< sink factor for current canopy box
3011        REAL(wp), DIMENSION(3)                 :: delta        !< path vector
3012        REAL(wp), DIMENSION(3)                 :: uvect        !< unit vector
3013        REAL(wp), DIMENSION(3)                 :: dimnextdist  !< distance for each dimension increments
3014        INTEGER(iwp), DIMENSION(3)             :: box          !< gridbox being crossed
3015        INTEGER(iwp), DIMENSION(3)             :: dimnext      !< next dimension increments along path
3016        INTEGER(iwp), DIMENSION(3)             :: dimdelta     !< dimension direction = +- 1
3017        INTEGER(iwp)                           :: px, py       !< number of processors in x and y dir before
3018                                                               !< the processor in the question
3019        INTEGER(iwp)                           :: ig, ip
3020        REAL(wp)                               :: lad_s_target
3021        INTEGER(kind=MPI_ADDRESS_KIND)         :: lad_disp
3022        REAL(wp), PARAMETER                    :: grow_factor = 1.5_wp
3023
3024!--     Maximum number of gridboxes visited equals to maximum number of boundaries crossed in each dimension plus one. That's also
3025!--     the maximum number of plant canopy boxes written. We grow the acsf array accordingly using exponential factor.
3026        maxboxes = SUM(ABS(NINT(targ) - NINT(src))) + 1
3027        IF ( plant_canopy  .AND.  ncsfl + maxboxes > ncsfla )  THEN
3028!--         use this code for growing by fixed exponential increments (equivalent to case where ncsfl always increases by 1)
3029!--         k = CEILING(grow_factor ** real(CEILING(log(real(ncsfl + maxboxes, kind=wp)) &
3030!--                                                / log(grow_factor)), kind=wp))
3031!--         or use this code to simply always keep some extra space after growing
3032            k = CEILING(REAL(ncsfl + maxboxes, kind=wp) * grow_factor)
3033
3034            CALL usm_merge_and_grow_csf(k)
3035        ENDIF
3036       
3037        transparency = 1._wp
3038        ncsb = 0
3039
3040        delta(:) = targ(:) - src(:)
3041        distance = SQRT(SUM(delta(:)**2))
3042        IF ( distance == 0._wp )  THEN
3043            visible = .TRUE.
3044            RETURN
3045        ENDIF
3046        uvect(:) = delta(:) / distance
3047        realdist = SQRT(SUM( (uvect(:)*(/dz,dy,dx/))**2 ))
3048
3049        lastdist = 0._wp
3050
3051!--     Since all face coordinates have values *.5 and we'd like to use
3052!--     integers, all these have .5 added
3053        DO d = 1, 3
3054            IF ( uvect(d) == 0._wp )  THEN
3055                dimnext(d) = 999999999
3056                dimdelta(d) = 999999999
3057                dimnextdist(d) = 1.0E20_wp
3058            ELSE IF ( uvect(d) > 0._wp )  THEN
3059                dimnext(d) = CEILING(src(d) + .5_wp)
3060                dimdelta(d) = 1
3061                dimnextdist(d) = (dimnext(d) - .5_wp - src(d)) / uvect(d)
3062            ELSE
3063                dimnext(d) = FLOOR(src(d) + .5_wp)
3064                dimdelta(d) = -1
3065                dimnextdist(d) = (dimnext(d) - .5_wp - src(d)) / uvect(d)
3066            ENDIF
3067        ENDDO
3068
3069        DO
3070!--         along what dimension will the next wall crossing be?
3071            seldim = minloc(dimnextdist, 1)
3072            nextdist = dimnextdist(seldim)
3073            IF ( nextdist > distance ) nextdist = distance
3074
3075            crlen = nextdist - lastdist
3076            IF ( crlen > .001_wp )  THEN
3077                crmid = (lastdist + nextdist) * .5_wp
3078                box = NINT(src(:) + uvect(:) * crmid)
3079
3080!--             calculate index of the grid with global indices (box(2),box(3))
3081!--             in the array nzterr and plantt and id of the coresponding processor
3082                px = box(3)/nnx
3083                py = box(2)/nny
3084                ip = px*pdims(2)+py
3085                ig = ip*nnx*nny + (box(3)-px*nnx)*nny + box(2)-py*nny
3086                IF ( box(1) <= nzterr(ig) )  THEN
3087                    visible = .FALSE.
3088                    IF ( ncsb > 0 )  THEN
3089!--                     rewind written plant canopy sink factors - they are invalid
3090                        ncsfl = ncsfl - ncsb
3091                    ENDIF
3092                    RETURN
3093                ENDIF
3094
3095                IF ( plant_canopy )  THEN
3096                    IF ( box(1) <= plantt(ig) )  THEN
3097#if defined( __parallel )
3098                        lad_disp = (box(3)-px*nnx)*(nny*nzu) + (box(2)-py*nny)*nzu + box(1)-nzub
3099                        IF ( usm_lad_rma )  THEN
3100!--                         Read LAD using MPI RMA
3101                            CALL cpu_log( log_point_s(77), 'usm_init_rma', 'start' )
3102                            CALL MPI_Get(lad_s_target, 1, MPI_REAL, ip, lad_disp, 1, MPI_REAL, &
3103                                win_lad, ierr)
3104                            IF ( ierr /= 0 )  THEN
3105                                WRITE(message_string, *) 'MPI error ', ierr, ' at MPI_Get'
3106                                CALL message( 'usm_raytrace', 'PA0519', 1, 2, 0, 6, 0 )
3107                            ENDIF
3108                           
3109                            CALL MPI_Win_flush_local(ip, win_lad, ierr)
3110                            IF ( ierr /= 0 )  THEN
3111                                WRITE(message_string, *) 'MPI error ', ierr, ' at MPI_Win_flush_local'
3112                                CALL message( 'usm_raytrace', 'PA0519', 1, 2, 0, 6, 0 )
3113                            ENDIF
3114                            CALL cpu_log( log_point_s(77), 'usm_init_rma', 'stop' )
3115                        ELSE
3116                            lad_s_target = usm_lad_g(ip*nnx*nny*nzu + lad_disp)
3117                        ENDIF
3118#else
3119                        lad_s_target = usm_lad(box(1),box(2),box(3))
3120#endif
3121                        cursink = 1._wp - exp(-ext_coef * lad_s_target &
3122                            * crlen*realdist)
3123
3124                        IF ( create_csf )  THEN
3125!--                         write svf values into the array
3126                            ncsb = ncsb + 1
3127                            ncsfl = ncsfl + 1
3128                            acsf(ncsfl)%ip = ip
3129                            acsf(ncsfl)%itx = box(3)
3130                            acsf(ncsfl)%ity = box(2)
3131                            acsf(ncsfl)%itz = box(1)
3132                            acsf(ncsfl)%isurfs = isrc
3133                            acsf(ncsfl)%rsvf = REAL(cursink*rirrf*atarg, wp) !-- we postpone multiplication by transparency
3134                            acsf(ncsfl)%rtransp = REAL(transparency, wp)
3135                        ENDIF  !< create_csf
3136
3137                        transparency = transparency * (1._wp - cursink)
3138                       
3139                    ENDIF
3140                ENDIF
3141            ENDIF
3142
3143            IF ( nextdist >= distance ) EXIT
3144            lastdist = nextdist
3145            dimnext(seldim) = dimnext(seldim) + dimdelta(seldim)
3146            dimnextdist(seldim) = (dimnext(seldim) - .5_wp - src(seldim)) / uvect(seldim)
3147        ENDDO
3148       
3149        visible = .TRUE.
3150    END SUBROUTINE usm_raytrace
3151   
3152 
3153!------------------------------------------------------------------------------!
3154! Description:
3155! ------------
3156!
3157!> This subroutine is part of the urban surface model.
3158!> It reads daily heat produced by anthropogenic sources
3159!> and the diurnal cycle of the heat.
3160!------------------------------------------------------------------------------!
3161    SUBROUTINE usm_read_anthropogenic_heat
3162   
3163        INTEGER(iwp)                  :: i,j,ii
3164        REAL(wp)                      :: heat
3165       
3166!--     allocation of array of sources of anthropogenic heat and their diural profile
3167        ALLOCATE( aheat(nys:nyn,nxl:nxr) )
3168        ALLOCATE( aheatprof(0:24) )
3169
3170!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3171!--     read daily amount of heat and its daily cycle
3172!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3173        aheat = 0.0_wp
3174        DO  ii = 0, io_blocks-1
3175            IF ( ii == io_group )  THEN
3176
3177!--             open anthropogenic heat file
3178                OPEN( 151, file='ANTHROPOGENIC_HEAT'//TRIM(coupling_char), action='read', &
3179                           status='old', form='formatted', err=11 )
3180                i = 0
3181                j = 0
3182                DO
3183                    READ( 151, *, err=12, end=13 )  i, j, heat
3184                    IF ( i >= nxl  .AND.  i <= nxr  .AND.  j >= nys  .AND.  j <= nyn )  THEN
3185!--                     write heat into the array
3186                        aheat(j,i) = heat
3187                    ENDIF
3188                    CYCLE
3189 12                 WRITE(message_string,'(a,2i4)') 'error in file ANTHROPOGENIC_HEAT'//TRIM(coupling_char)//' after line ',i,j
3190                    CALL message( 'usm_read_anthropogenic_heat', 'PA0515', 0, 1, 0, 6, 0 )
3191                ENDDO
3192 13             CLOSE(151)
3193                CYCLE
3194 11             message_string = 'file ANTHROPOGENIC_HEAT'//TRIM(coupling_char)//' does not exist'
3195                CALL message( 'usm_read_anthropogenic_heat', 'PA0516', 1, 2, 0, 6, 0 )
3196            ENDIF
3197           
3198#if defined( __parallel ) && ! defined ( __check )
3199            CALL mpi_barrier( comm2d, ierr )
3200#endif
3201        ENDDO
3202       
3203!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3204!--     read diurnal profiles of heat sources
3205!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3206        aheatprof = 0.0_wp
3207        DO  ii = 0, io_blocks-1
3208            IF ( ii == io_group )  THEN
3209
3210!--             open anthropogenic heat profile file
3211                OPEN( 151, file='ANTHROPOGENIC_HEAT_PROFILE'//TRIM(coupling_char), action='read', &
3212                           status='old', form='formatted', err=21 )
3213                i = 0
3214                DO
3215                    READ( 151, *, err=22, end=23 )  i, heat
3216                    IF ( i >= 0  .AND.  i <= 24 )  THEN
3217!--                     write heat into the array
3218                        aheatprof(i) = heat
3219                    ENDIF
3220                    CYCLE
3221 22                 WRITE(message_string,'(a,i4)') 'error in file ANTHROPOGENIC_HEAT_PROFILE'// &
3222                                                     TRIM(coupling_char)//' after line ',i
3223                    CALL message( 'usm_read_anthropogenic_heat', 'PA0517', 0, 1, 0, 6, 0 )
3224                ENDDO
3225                aheatprof(24) = aheatprof(0)
3226 23             CLOSE(151)
3227                CYCLE
3228 21             message_string = 'file ANTHROPOGENIC_HEAT_PROFILE'//TRIM(coupling_char)//' does not exist'
3229                CALL message( 'usm_read_anthropogenic_heat', 'PA0518', 1, 2, 0, 6, 0 )
3230            ENDIF
3231           
3232#if defined( __parallel ) && ! defined ( __check )
3233            CALL mpi_barrier( comm2d, ierr )
3234#endif
3235        ENDDO
3236       
3237    END SUBROUTINE usm_read_anthropogenic_heat
3238   
3239
3240!------------------------------------------------------------------------------!
3241!
3242! Description:
3243! ------------
3244!> Soubroutine reads t_surf and t_wall data from restart files
3245!kanani: Renamed this routine according to corresponging routines in PALM
3246!kanani: Modified the routine to match read_var_list, from where usm_read_restart_data
3247!        shall be called in the future. This part has not been tested yet. (see virtual_flight_mod)
3248!        Also, I had some trouble with the allocation of t_surf, since this is a pointer.
3249!        So, I added some directives here.
3250!------------------------------------------------------------------------------!
3251    SUBROUTINE usm_read_restart_data
3252
3253
3254       IMPLICIT NONE
3255       
3256       CHARACTER (LEN=30) ::  variable_chr  !< dummy variable to read string
3257       
3258       INTEGER            ::  i             !< running index     
3259       
3260
3261       DO  i = 0, io_blocks-1
3262          IF ( i == io_group )  THEN
3263             READ ( 13 )  variable_chr
3264             DO   WHILE ( TRIM( variable_chr ) /= '*** end usm ***' )
3265
3266                SELECT CASE ( TRIM( variable_chr ) )
3267               
3268                   CASE ( 't_surf' )
3269#if defined( __nopointer )                   
3270                      IF ( .NOT.  ALLOCATED( t_surf ) )                         &
3271                         ALLOCATE( t_surf(startenergy:endenergy) )
3272                      READ ( 13 )  t_surf
3273#else                     
3274                      IF ( .NOT.  ALLOCATED( t_surf_1 ) )                         &
3275                         ALLOCATE( t_surf_1(startenergy:endenergy) )
3276                      READ ( 13 )  t_surf_1
3277#endif
3278
3279                   CASE ( 't_wall' )
3280#if defined( __nopointer )
3281                      IF ( .NOT.  ALLOCATED( t_wall ) )                         &
3282                         ALLOCATE( t_wall(nzb_wall:nzt_wall+1,startenergy:endenergy) )
3283                      READ ( 13 )  t_wall
3284#else
3285                      IF ( .NOT.  ALLOCATED( t_wall_1 ) )                         &
3286                         ALLOCATE( t_wall_1(nzb_wall:nzt_wall+1,startenergy:endenergy) )
3287                      READ ( 13 )  t_wall_1
3288#endif
3289
3290                   CASE DEFAULT
3291                      WRITE ( message_string, * )  'unknown variable named "', &
3292                                        TRIM( variable_chr ), '" found in',    &
3293                                        '&data from prior run on PE ', myid
3294                      CALL message( 'user_read_restart_data', 'UI0012', 1, 2, 0, 6, 0 )
3295
3296                END SELECT
3297
3298                READ ( 13 )  variable_chr
3299
3300             ENDDO
3301          ENDIF
3302#if defined( __parallel )
3303          CALL MPI_BARRIER( comm2d, ierr )
3304#endif
3305       ENDDO
3306
3307    END SUBROUTINE usm_read_restart_data
3308
3309
3310!------------------------------------------------------------------------------!
3311!
3312! Description:
3313! ------------
3314!> Soubroutine reads svf and svfsurf data from saved file
3315!------------------------------------------------------------------------------!
3316    SUBROUTINE usm_read_svf_from_file
3317
3318        IMPLICIT NONE
3319        INTEGER                      :: fsvf = 89
3320        INTEGER                      :: i
3321        CHARACTER(usm_version_len)   :: usm_version_field
3322        CHARACTER(svf_code_len)      :: svf_code_field
3323
3324        DO  i = 0, io_blocks-1
3325            IF ( i == io_group )  THEN
3326                OPEN ( fsvf, file=TRIM(svf_file_name)//TRIM(coupling_char)//myid_char,               &
3327                    form='unformatted', status='old' )
3328
3329!--             read and check version
3330                READ ( fsvf ) usm_version_field
3331                IF ( TRIM(usm_version_field) /= TRIM(usm_version) )  THEN
3332                    WRITE( message_string, * ) 'Version of binary SVF file "',           &
3333                                            TRIM(usm_version_field), '" does not match ',            &
3334                                            'the version of model "', TRIM(usm_version), '"'
3335                    CALL message( 'usm_read_svf_from_file', 'UI0012', 1, 2, 0, 6, 0 )
3336                ENDIF
3337               
3338!--             read nsvfcsfl, nsvfl
3339                READ ( fsvf ) nsvfcsfl, nsvfl
3340                IF ( nsvfcsfl <= 0 )  THEN
3341                    WRITE( message_string, * ) 'Wrong number of SVF or CSF'
3342                    CALL message( 'usm_read_svf_from_file', 'UI0012', 1, 2, 0, 6, 0 )
3343                ELSE
3344                    WRITE(message_string,*) '    Number of SVF and CSF to read', nsvfcsfl, nsvfl
3345                    CALL location_message( message_string, .TRUE. )
3346                ENDIF
3347               
3348                ALLOCATE(svf(ndsvf,nsvfcsfl))
3349                ALLOCATE(svfsurf(ndsvf,nsvfcsfl))
3350               
3351                READ(fsvf) svf
3352                READ(fsvf) svfsurf
3353                READ ( fsvf ) svf_code_field
3354               
3355                IF ( TRIM(svf_code_field) /= TRIM(svf_code) )  THEN
3356                    WRITE( message_string, * ) 'Wrong structure of binary svf file'
3357                    CALL message( 'usm_read_svf_from_file', 'UI0012', 1, 2, 0, 6, 0 )
3358                ENDIF
3359               
3360                CLOSE (fsvf)
3361               
3362            ENDIF
3363#if defined( __parallel )
3364            CALL mpi_barrier( comm2d, ierr )
3365#endif
3366        ENDDO
3367
3368    END SUBROUTINE usm_read_svf_from_file
3369
3370   
3371!------------------------------------------------------------------------------!
3372! Description:
3373! ------------
3374!
3375!> This subroutine reads walls, roofs and land categories and it parameters
3376!> from input files.
3377!------------------------------------------------------------------------------!
3378    SUBROUTINE usm_read_urban_surface_types
3379   
3380        CHARACTER(12)                                         :: wtn
3381        INTEGER                                               :: wtc
3382        REAL(wp), DIMENSION(n_surface_params)                 :: wtp
3383   
3384        INTEGER(iwp), DIMENSION(0:17, nysg:nyng, nxlg:nxrg)   :: usm_par
3385        REAL(wp), DIMENSION(1:14, nysg:nyng, nxlg:nxrg)       :: usm_val
3386        INTEGER(iwp)                                          :: k, l, d, iw, jw, kw, it, ip, ii, ij
3387        INTEGER(iwp)                                          :: i, j
3388        INTEGER(iwp)                                          :: nz, roof, dirwe, dirsn
3389        INTEGER(iwp)                                          :: category
3390        INTEGER(iwp)                                          :: weheight1, wecat1, snheight1, sncat1
3391        INTEGER(iwp)                                          :: weheight2, wecat2, snheight2, sncat2
3392        INTEGER(iwp)                                          :: weheight3, wecat3, snheight3, sncat3
3393        REAL(wp)                                              :: height, albedo, thick
3394        REAL(wp)                                              :: wealbedo1, wethick1, snalbedo1, snthick1
3395        REAL(wp)                                              :: wealbedo2, wethick2, snalbedo2, snthick2
3396        REAL(wp)                                              :: wealbedo3, wethick3, snalbedo3, snthick3
3397       
3398!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3399!--     read categories of walls and their parameters
3400!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3401        DO  ii = 0, io_blocks-1
3402            IF ( ii == io_group )  THEN
3403
3404!--             open urban surface file
3405                OPEN( 151, file='SURFACE_PARAMETERS'//coupling_char, action='read', &
3406                           status='old', form='formatted', err=15 ) 
3407!--             first test and get n_surface_types
3408                k = 0
3409                l = 0
3410                DO
3411                    l = l+1
3412                    READ( 151, *, err=11, end=12 )  wtc, wtp, wtn
3413                    k = k+1
3414                    CYCLE
3415 11                 CONTINUE
3416                ENDDO
3417 12             n_surface_types = k
3418                ALLOCATE( surface_type_names(n_surface_types) )
3419                ALLOCATE( surface_type_codes(n_surface_types) )
3420                ALLOCATE( surface_params(n_surface_params, n_surface_types) )
3421!--             real reading
3422                rewind( 151 )
3423                k = 0
3424                DO
3425                    READ( 151, *, err=13, end=14 )  wtc, wtp, wtn
3426                    k = k+1
3427                    surface_type_codes(k) = wtc
3428                    surface_params(:,k) = wtp
3429                    surface_type_names(k) = wtn
3430                    CYCLE
343113                  WRITE(6,'(i3,a,2i5)') myid, 'readparams2 error k=', k
3432                    FLUSH(6)
3433                    CONTINUE
3434                ENDDO
3435 14             CLOSE(151)
3436                CYCLE
3437 15             message_string = 'file SURFACE_PARAMETERS'//TRIM(coupling_char)//' does not exist'
3438                CALL message( 'usm_read_urban_surface_types', 'PA0513', 1, 2, 0, 6, 0 )
3439            ENDIF
3440        ENDDO
3441   
3442!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3443!--     read types of surfaces
3444!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3445        usm_par = 0
3446        DO  ii = 0, io_blocks-1
3447            IF ( ii == io_group )  THEN
3448
3449                !
3450!--             open csv urban surface file
3451                OPEN( 151, file='URBAN_SURFACE'//TRIM(coupling_char), action='read', &
3452                      status='old', form='formatted', err=23 )
3453               
3454                l = 0
3455                DO
3456                    l = l+1
3457!--                 i, j, height, nz, roof, dirwe, dirsn, category, soilcat,
3458!--                 weheight1, wecat1, snheight1, sncat1, weheight2, wecat2, snheight2, sncat2,
3459!--                 weheight3, wecat3, snheight3, sncat3
3460                    READ( 151, *, err=21, end=25 )  i, j, height, nz, roof, dirwe, dirsn,            &
3461                                            category, albedo, thick,                                 &
3462                                            weheight1, wecat1, wealbedo1, wethick1,                  &
3463                                            weheight2, wecat2, wealbedo2, wethick2,                  &
3464                                            weheight3, wecat3, wealbedo3, wethick3,                  &
3465                                            snheight1, sncat1, snalbedo1, snthick1,                  &
3466                                            snheight2, sncat2, snalbedo2, snthick2,                  &
3467                                            snheight3, sncat3, snalbedo3, snthick3
3468
3469                    IF ( i >= nxlg  .AND.  i <= nxrg  .AND.  j >= nysg  .AND.  j <= nyng )  THEN
3470!--                     write integer variables into array
3471                        usm_par(:,j,i) = (/1, nz, roof, dirwe, dirsn, category,                      &
3472                                          weheight1, wecat1, weheight2, wecat2, weheight3, wecat3,   &
3473                                          snheight1, sncat1, snheight2, sncat2, snheight3, sncat3 /)
3474!--                     write real values into array
3475                        usm_val(:,j,i) = (/ albedo, thick,                                           &
3476                                           wealbedo1, wethick1, wealbedo2, wethick2,                 &
3477                                           wealbedo3, wethick3, snalbedo1, snthick1,                 &
3478                                           snalbedo2, snthick2, snalbedo3, snthick3 /)
3479                    ENDIF
3480                    CYCLE
3481 21                 WRITE (message_string, "(A,I5)") 'errors in file URBAN_SURFACE'//TRIM(coupling_char)//' on line ', l
3482                    CALL message( 'usm_read_urban_surface_types', 'PA0512', 0, 1, 0, 6, 0 )
3483                ENDDO
3484         
3485 23             message_string = 'file URBAN_SURFACE'//TRIM(coupling_char)//' does not exist'
3486                CALL message( 'usm_read_urban_surface_types', 'PA0514', 1, 2, 0, 6, 0 )
3487
3488 25             CLOSE( 90 )
3489
3490            ENDIF
3491#if defined( __parallel ) && ! defined ( __check )
3492            CALL mpi_barrier( comm2d, ierr )
3493#endif
3494        ENDDO
3495       
3496        !
3497!--     check completeness and formal correctness of the data
3498        DO i = nxlg, nxrg
3499            DO j = nysg, nyng
3500                IF ( usm_par(0,j,i) /= 0  .AND.  (        &  !< incomplete data,supply default values later
3501                     usm_par(1,j,i) < nzb  .OR.           &
3502                     usm_par(1,j,i) > nzt  .OR.           &  !< incorrect height (nz < nzb  .OR.  nz > nzt)
3503                     usm_par(2,j,i) < 0  .OR.             &
3504                     usm_par(2,j,i) > 1  .OR.             &  !< incorrect roof sign
3505                     usm_par(3,j,i) < nzb-nzt  .OR.       & 
3506                     usm_par(3,j,i) > nzt-nzb  .OR.       &  !< incorrect west-east wall direction sign
3507                     usm_par(4,j,i) < nzb-nzt  .OR.       &
3508                     usm_par(4,j,i) > nzt-nzb  .OR.       &  !< incorrect south-north wall direction sign
3509                     usm_par(6,j,i) < nzb  .OR.           & 
3510                     usm_par(6,j,i) > nzt  .OR.           &  !< incorrect pedestrian level height for west-east wall
3511                     usm_par(8,j,i) > nzt  .OR.           &
3512                     usm_par(10,j,i) > nzt  .OR.          &  !< incorrect wall or roof level height for west-east wall
3513                     usm_par(12,j,i) < nzb  .OR.          & 
3514                     usm_par(12,j,i) > nzt  .OR.          &  !< incorrect pedestrian level height for south-north wall
3515                     usm_par(14,j,i) > nzt  .OR.          &
3516                     usm_par(16,j,i) > nzt                &  !< incorrect wall or roof level height for south-north wall
3517                    ) )  THEN
3518!--                 incorrect input data
3519                    WRITE (message_string, "(A,2I5)") 'missing or incorrect data in file URBAN_SURFACE'// &
3520                                                       TRIM(coupling_char)//' for i,j=', i,j
3521                    CALL message( 'usm_read_urban_surface', 'PA0504', 1, 2, 0, 6, 0 )
3522                ENDIF
3523               
3524            ENDDO
3525        ENDDO
3526       
3527!--     assign the surface types to local surface array
3528        DO  l = startenergy, endenergy
3529           
3530            d = surfl(id,l)
3531            kw = surfl(iz,l)
3532            j = surfl(iy,l)
3533            i = surfl(ix,l)
3534            IF ( d == iroof )  THEN
3535!--             horizontal surface - land or roof
3536                iw = i
3537                jw = j
3538                IF ( usm_par(5,jw,iw) == 0 )  THEN
3539                    IF ( zu(kw) >= roof_height_limit )  THEN
3540                        isroof_surf(l) = .TRUE.
3541                        surface_types(l) = roof_category         !< default category for root surface
3542                    ELSE
3543                        isroof_surf(l) = .FALSE.
3544                        surface_types(l) = land_category         !< default category for land surface
3545                    ENDIF
3546                    albedo_surf(l) = -1.0_wp
3547                    thickness_wall(l) = -1.0_wp
3548                ELSE
3549                    IF ( usm_par(2,jw,iw)==0 )  THEN
3550                        isroof_surf(l) = .FALSE.
3551                        thickness_wall(l) = -1.0_wp
3552                    ELSE
3553                        isroof_surf(l) = .TRUE.
3554                        thickness_wall(l) = usm_val(2,jw,iw)
3555                    ENDIF
3556                    surface_types(l) = usm_par(5,jw,iw)
3557                    albedo_surf(l) = usm_val(1,jw,iw)
3558                ENDIF
3559            ELSE
3560                SELECT CASE (d)
3561                    CASE (iwest)
3562                        iw = i
3563                        jw = j
3564                        ii = 6
3565                        ij = 3
3566                    CASE (ieast)
3567                        iw = i-1
3568                        jw = j
3569                        ii = 6
3570                        ij = 3
3571                    CASE (isouth)
3572                        iw = i
3573                        jw = j
3574                        ii = 12
3575                        ij = 9
3576                    CASE (inorth)
3577                        iw = i
3578                        jw = j-1
3579                        ii = 12
3580                        ij = 9
3581                END SELECT
3582               
3583                IF ( kw <= usm_par(ii,jw,iw) )  THEN
3584!--                 pedestrant zone
3585                    isroof_surf(l) = .FALSE.
3586                    IF ( usm_par(ii+1,jw,iw) == 0 )  THEN
3587                        surface_types(l) = pedestrant_category   !< default category for wall surface in pedestrant zone
3588                        albedo_surf(l) = -1.0_wp
3589                        thickness_wall(l) = -1.0_wp
3590                    ELSE
3591                        surface_types(l) = usm_par(ii+1,jw,iw)
3592                        albedo_surf(l) = usm_val(ij,jw,iw)
3593                        thickness_wall(l) = usm_val(ij+1,jw,iw)
3594                    ENDIF
3595                ELSE IF ( kw <= usm_par(ii+2,jw,iw) )  THEN
3596!--                 wall zone
3597                    isroof_surf(l) = .FALSE.
3598                    IF ( usm_par(ii+3,jw,iw) == 0 )  THEN
3599                        surface_types(l) = wall_category         !< default category for wall surface
3600                        albedo_surf(l) = -1.0_wp
3601                        thickness_wall(l) = -1.0_wp
3602                    ELSE
3603                        surface_types(l) = usm_par(ii+3,jw,iw)
3604                        albedo_surf(l) = usm_val(ij+2,jw,iw)
3605                        thickness_wall(l) = usm_val(ij+3,jw,iw)
3606                    ENDIF
3607                ELSE IF ( kw <= usm_par(ii+4,jw,iw) )  THEN
3608!--                 roof zone
3609                    isroof_surf(l) = .TRUE.
3610                    IF ( usm_par(ii+5,jw,iw) == 0 )  THEN
3611                        surface_types(l) = roof_category         !< default category for roof surface
3612                        albedo_surf(l) = -1.0_wp
3613                        thickness_wall(l) = -1.0_wp
3614                    ELSE
3615                        surface_types(l) = usm_par(ii+5,jw,iw)
3616                        albedo_surf(l) = usm_val(ij+4,jw,iw)
3617                        thickness_wall(l) = usm_val(ij+5,jw,iw)
3618                    ENDIF
3619                ELSE
3620!--                 something wrong
3621                    CALL message( 'usm_read_urban_surface', 'PA0505', 1, 2, 0, 6, 0 )
3622                ENDIF
3623            ENDIF
3624           
3625!--         find the type position
3626            it = surface_types(l)
3627            ip = -99999
3628            DO k = 1, n_surface_types
3629                IF ( surface_type_codes(k) == it )  THEN
3630                    ip = k
3631                    EXIT
3632                ENDIF
3633            ENDDO
3634            IF ( ip == -99999 )  THEN
3635!--             wall category not found
3636                WRITE (message_string, "(A,I5,A,3I5)") 'wall category ', it, ' not found  for i,j,k=', iw,jw,kw
3637                CALL message( 'usm_read_urban_surface', 'PA0506', 1, 2, 0, 6, 0 )
3638            ENDIF
3639           
3640!--         Fill out the parameters of the wall
3641!--         wall surface:
3642           
3643!--         albedo
3644            IF ( albedo_surf(l) < 0.0_wp )  THEN
3645                albedo_surf(l) = surface_params(ialbedo, ip)
3646            ENDIF
3647           
3648!--         emissivity of the wall
3649            emiss_surf(l) = surface_params(iemiss, ip)
3650           
3651!--         heat conductivity λS between air and wall ( W m−2 K−1 )
3652            lambda_surf(l) = surface_params(ilambdas, ip)
3653           
3654!--         roughness relative to concrete
3655            roughness_wall(l) = surface_params(irough, ip)
3656           
3657!--         Surface skin layer heat capacity (J m−2 K−1 )
3658            c_surface(l) = surface_params(icsurf, ip)
3659           
3660!--         wall material parameters:
3661           
3662!--         thickness of the wall (m)
3663!--         missing values are replaced by default value for category
3664            IF ( thickness_wall(l) <= 0.001_wp )  THEN
3665                thickness_wall(l) = surface_params(ithick, ip)
3666            ENDIF
3667           
3668!--         volumetric heat capacity rho*C of the wall ( J m−3 K−1 )
3669            rho_c_wall(:,l) = surface_params(irhoC, ip)
3670           
3671!--         thermal conductivity λH of the wall (W m−1 K−1 )
3672            lambda_h(:,l) = surface_params(ilambdah, ip)
3673           
3674        ENDDO
3675
3676        CALL location_message( '    types and parameters of urban surfaces read', .TRUE. )
3677   
3678    END SUBROUTINE usm_read_urban_surface_types
3679
3680
3681!------------------------------------------------------------------------------!
3682! Description:
3683! ------------
3684!> Solver for the energy balance at the ground/roof/wall surface.
3685!> It follows basic ideas and structure of lsm_energy_balance
3686!> with many simplifications and adjustments.
3687!> TODO better description
3688!------------------------------------------------------------------------------!
3689    SUBROUTINE usm_surface_energy_balance
3690
3691        IMPLICIT NONE
3692
3693        INTEGER(iwp)                          :: i, j, k, l, d      !< running indices
3694       
3695        REAL(wp)                              :: pt1                !< temperature at first grid box adjacent to surface
3696        REAL(wp)                              :: u1,v1,w1           !< near wall u,v,w
3697        REAL(wp)                              :: stend              !< surface tendency
3698        REAL(wp)                              :: coef_1             !< first coeficient for prognostic equation
3699        REAL(wp)                              :: coef_2             !< second  coeficient for prognostic equation
3700        REAL(wp)                              :: rho_cp             !< rho_wall_surface * cp
3701        REAL(wp)                              :: r_a                !< aerodynamic resistance for horizontal and vertical surfaces
3702        REAL(wp)                              :: f_shf              !< factor for shf_eb
3703        REAL(wp)                              :: lambda_surface     !< current value of lambda_surface (heat conductivity between air and wall)
3704        REAL(wp)                              :: Ueff               !< effective wind speed for calculation of heat transfer coefficients
3705        REAL(wp)                              :: httc               !< heat transfer coefficient
3706        REAL(wp), DIMENSION(nzub:nzut)        :: exn                !< value of the Exner function in layers
3707       
3708        REAL(wp), DIMENSION(0:4)              :: dxdir              !< surface normal direction gridbox length
3709        REAL(wp)                              :: dtime              !< simulated time of day (in UTC)
3710        INTEGER(iwp)                          :: dhour              !< simulated hour of day (in UTC)
3711        REAL(wp)                              :: acoef              !< actual coefficient of diurnal profile of anthropogenic heat
3712
3713        dxdir = (/dz,dy,dy,dx,dx/)
3714       
3715        exn(:) = (hyp(nzub:nzut) / 100000.0_wp )**0.286_wp          !< Exner function
3716           
3717!--   
3718        DO l = startenergy, endenergy
3719!--         Calculate frequently used parameters
3720            d = surfl(id,l)
3721            k = surfl(iz,l)
3722            j = surfl(iy,l)
3723            i = surfl(ix,l)
3724
3725!--         TODO - how to calculate lambda_surface for horizontal surfaces
3726!--         (lambda_surface is set according to stratification in land surface model)
3727            IF ( ol(j,i) >= 0.0_wp )  THEN
3728                lambda_surface = lambda_surf(l)
3729            ELSE
3730                lambda_surface = lambda_surf(l)
3731            ENDIF
3732           
3733            pt1  = pt(k,j,i)
3734
3735!--         calculate rho * cp coefficient at surface layer
3736            rho_cp  = cp * hyp(k) / ( r_d * pt1 * exn(k) )
3737
3738!--         calculate aerodyamic resistance.
3739            IF ( d == iroof )  THEN
3740!--             calculation for horizontal surfaces follows LSM formulation
3741!--             pt, us, ts are not available for the prognostic time step,
3742!--             data from the last time step is used here.
3743               
3744                r_a = (pt1 - t_surf(l)/exn(k)) / (ts(j,i) * us(j,i) + 1.0E-10_wp)
3745               
3746!--             make sure that the resistance does not drop to zero
3747                IF ( ABS(r_a) < 1.0E-10_wp )  r_a = 1.0E-10_wp
3748               
3749!--             the parameterization is developed originally for larger scales
3750!--             (compare with remark in TUF-3D)
3751!--             our first experiences show that the parameterization underestimates
3752!--             r_a in meter resolution.
3753!--             temporary solution - multiplication by magic constant :-(.
3754                r_a = r_a * ra_horiz_coef
3755               
3756!--             factor for shf_eb
3757                f_shf  = rho_cp / r_a
3758            ELSE
3759!--             calculation of r_a for vertical surfaces
3760!--
3761!--             heat transfer coefficient for forced convection along vertical walls
3762!--             follows formulation in TUF3d model (Krayenhoff & Voogt, 2006)
3763!--           
3764!--             H = httc (Tsfc - Tair)
3765!--             httc = rw * (11.8 + 4.2 * Ueff) - 4.0
3766!--           
3767!--                   rw: wall patch roughness relative to 1.0 for concrete
3768!--                   Ueff: effective wind speed
3769!--                   - 4.0 is a reduction of Rowley et al (1930) formulation based on
3770!--                   Cole and Sturrock (1977)
3771!--           
3772!--                   Ucan: Canyon wind speed
3773!--                   wstar: convective velocity
3774!--                   Qs: surface heat flux
3775!--                   zH: height of the convective layer
3776!--                   wstar = (g/Tcan*Qs*zH)**(1./3.)
3777               
3778!--             staggered grid needs to be taken into consideration
3779                IF ( d == inorth )  THEN
3780                    u1 = (u(k,j,i)+u(k,j,i+1))*0.5_wp
3781                    v1 = v(k,j+1,i)
3782                ELSE IF ( d == isouth )  THEN
3783                    u1 = (u(k,j,i)+u(k,j,i+1))*0.5_wp
3784                    v1 = v(k,j,i)
3785                ELSE IF ( d == ieast )  THEN
3786                    u1 = u(k,j,i+1)
3787                    v1 = (v(k,j,i)+v(k,j+1,i))*0.5_wp
3788                ELSE IF ( d == iwest )  THEN
3789                    u1 = u(k,j,i)
3790                    v1 = (v(k,j,i)+v(k,j+1,i))*0.5_wp
3791                ELSE
3792                    STOP
3793                ENDIF
3794                w1 = (w(k,j,i)+w(k-1,j,i))*0.5_wp
3795               
3796                Ueff = SQRT(u1**2 + v1**2 + w1**2)
3797                httc = roughness_wall(l) * (11.8 + 4.2 * Ueff) - 4.0
3798                f_shf  = httc
3799            ENDIF
3800       
3801!--         add LW up so that it can be removed in prognostic equation
3802            rad_net_l(l) = surfinsw(l) - surfoutsw(l) + surfinlw(l) - surfoutlw(l)
3803
3804!--         numerator of the prognostic equation
3805            coef_1 = rad_net_l(l) +    &    ! coef +1 corresponds to -lwout included in calculation of radnet_l
3806                     (3.0_wp+1.0_wp) * emiss_surf(l) * sigma_sb * t_surf(l) ** 4 +      & 
3807                     f_shf  * pt1 +                                                     &
3808                     lambda_surface * t_wall(nzb_wall,l)
3809
3810!--         denominator of the prognostic equation
3811            coef_2 = 4.0_wp * emiss_surf(l) * sigma_sb * t_surf(l) ** 3                 &
3812                         + lambda_surface + f_shf / exn(k)
3813
3814!--         implicit solution when the surface layer has no heat capacity,
3815!--         otherwise use RK3 scheme.
3816            t_surf_p(l) = ( coef_1 * dt_3d * tsc(2) + c_surface(l) * t_surf(l) ) /      & 
3817                              ( c_surface(l) + coef_2 * dt_3d * tsc(2) ) 
3818
3819!--         add RK3 term
3820            t_surf_p(l) = t_surf_p(l) + dt_3d * tsc(3) * tt_surface_m(l)
3821           
3822!--         calculate true tendency
3823            stend = (t_surf_p(l) - t_surf(l) - dt_3d * tsc(3) * tt_surface_m(l)) / (dt_3d  * tsc(2))
3824
3825!--         calculate t_surf tendencies for the next Runge-Kutta step
3826            IF ( timestep_scheme(1:5) == 'runge' )  THEN
3827                IF ( intermediate_timestep_count == 1 )  THEN
3828                    tt_surface_m(l) = stend
3829                ELSEIF ( intermediate_timestep_count <                                  &
3830                         intermediate_timestep_count_max )  THEN
3831                    tt_surface_m(l) = -9.5625_wp * stend + 5.3125_wp                    &
3832                                       * tt_surface_m(l)
3833                ENDIF
3834            ENDIF
3835
3836!--         in case of fast changes in the skin temperature, it is required to
3837!--         update the radiative fluxes in order to keep the solution stable
3838            IF ( ABS( t_surf_p(l) - t_surf(l) ) > 1.0_wp )  THEN
3839               force_radiation_call_l = .TRUE.
3840            ENDIF
3841           
3842!--         for horizontal surfaces is pt(nzb_s_inner(j,i),j,i) = pt_surf.
3843!--         there is no equivalent surface gridpoint for vertical surfaces.
3844!--         pt(k,j,i) is calculated for all directions in diffusion_s
3845!--         using surface and wall heat fluxes
3846            IF ( d == iroof )  THEN
3847               pt(nzb_s_inner(j,i),j,i) = t_surf_p(l) / exn(k)
3848            ENDIF
3849
3850!--         calculate fluxes
3851!--         rad_net_l is never used!           
3852            rad_net_l(l)     = rad_net_l(l) + 3.0_wp * sigma_sb                         &
3853                                * t_surf(l)**4 - 4.0_wp * sigma_sb                      &
3854                                * t_surf(l)**3 * t_surf_p(l)
3855            wghf_eb(l)       = lambda_surface * (t_surf_p(l) - t_wall(nzb_wall,l))
3856
3857!--         ground/wall/roof surface heat flux
3858            wshf_eb(l)  = - f_shf  * ( pt1 - t_surf_p(l) )
3859           
3860!--         store kinematic surface heat fluxes for utilization in other processes
3861!--         diffusion_s, surface_layer_fluxes,...
3862            IF ( d == iroof )  THEN
3863!--             shf is used in diffusion_s and also
3864!--             for calculation of surface layer fluxes
3865!--             update for horizontal surfaces
3866                shf(j,i) = wshf_eb(l) / rho_cp
3867            ELSE
3868!--             surface heat flux for vertical surfaces
3869!--             used in diffusion_s
3870                wshf(l) = wshf_eb(l) / rho_cp
3871            ENDIF
3872
3873        ENDDO
3874       
3875       
3876        IF ( usm_anthropogenic_heat  .AND.  &
3877             intermediate_timestep_count == intermediate_timestep_count_max )  THEN
3878!--         application of the additional anthropogenic heat sources
3879!--         we considere the traffic for now so all heat is absorbed
3880!--         to the first layer, generalization would be worth
3881           
3882!--         calculation of actual profile coefficient
3883!--         ??? check time_since_reference_point ???
3884            dtime = mod(simulated_time + time_utc_init, 24.0_wp*3600.0_wp)
3885            dhour = INT(dtime/3600.0_wp)
3886!--         linear interpolation of coeficient
3887            acoef = (REAL(dhour+1,wp)-dtime/3600.0_wp)*aheatprof(dhour) + (dtime/3600.0_wp-REAL(dhour,wp))*aheatprof(dhour+1)
3888            DO i = nxl, nxr
3889                DO j = nys, nyn
3890                    IF ( aheat(j,i) > 0.0_wp )  THEN
3891!--                     TODO the increase of pt in box i,j,nzb_s_inner(j,i)+1 in time dt_3d
3892!--                     given to anthropogenic heat aheat*acoef (W*m-2)
3893!--                     k = nzb_s_inner(j,i)+1
3894!--                     pt(k,j,i) = pt(k,j,i) + aheat(j,i)*acoef*dt_3d/(exn(k)*rho_cp*dz)
3895!--                     Instead of this, we can adjust shf in case AH only at surface
3896                        shf(j,i) = shf(j,i) + aheat(j,i)*acoef * ddx * ddy / rho_cp
3897                    ENDIF
3898                ENDDO
3899            ENDDO
3900        ENDIF
3901       
3902!--     pt and shf are defined on nxlg:nxrg,nysg:nyng
3903!--     get the borders from neighbours
3904        CALL exchange_horiz( pt, nbgp )
3905        CALL exchange_horiz_2d( shf )
3906
3907
3908!--    calculation of force_radiation_call:
3909!--    Make logical OR for all processes.
3910!--    Force radiation call if at least one processor forces it.
3911       IF ( intermediate_timestep_count == intermediate_timestep_count_max-1 )          &
3912       THEN
3913#if defined( __parallel )
3914          IF ( collective_wait )  CALL mpi_barrier( comm2d, ierr )
3915              CALL mpi_allreduce( force_radiation_call_l, force_radiation_call,         &
3916                                  1, MPI_LOGICAL, MPI_LOR, comm2d, ierr )
3917#else
3918          force_radiation_call = force_radiation_call_l
3919#endif
3920          force_radiation_call_l = .FALSE.
3921       ENDIF
3922
3923    END SUBROUTINE usm_surface_energy_balance
3924
3925
3926!------------------------------------------------------------------------------!
3927! Description:
3928! ------------
3929!> Swapping of timelevels for t_surf and t_wall
3930!> called out from subroutine swap_timelevel
3931!------------------------------------------------------------------------------!
3932    SUBROUTINE usm_swap_timelevel ( mod_count )
3933
3934       IMPLICIT NONE
3935
3936       INTEGER, INTENT(IN) :: mod_count
3937       INTEGER :: i
3938     
3939#if defined( __nopointer )
3940       t_surf    = t_surf_p
3941       t_wall    = t_wall_p
3942#else
3943       SELECT CASE ( mod_count )
3944          CASE ( 0 )
3945             t_surf  => t_surf_1; t_surf_p  => t_surf_2
3946             t_wall     => t_wall_1;    t_wall_p     => t_wall_2
3947          CASE ( 1 )
3948             t_surf  => t_surf_2; t_surf_p  => t_surf_1
3949             t_wall     => t_wall_2;    t_wall_p     => t_wall_1
3950       END SELECT
3951#endif
3952       
3953    END SUBROUTINE usm_swap_timelevel
3954
3955
3956!------------------------------------------------------------------------------!
3957! Description:
3958! ------------
3959!
3960!> This function applies the kinematic wall heat fluxes
3961!> for walls in four directions for all gridboxes in urban layer.
3962!> It is called out from subroutine prognostic_equations.
3963!> TODO Compare performance with cycle runnig l=startwall,endwall...
3964!------------------------------------------------------------------------------!
3965    SUBROUTINE usm_wall_heat_flux
3966   
3967        IMPLICIT NONE
3968
3969        INTEGER(iwp)              ::  i,j,k,d,l             !< running indices
3970       
3971        DO l = startenergy, endenergy
3972            j = surfl(iy,l)
3973            i = surfl(ix,l)
3974            k = surfl(iz,l)
3975            d = surfl(id,l)
3976            tend(k,j,i) = tend(k,j,i) + wshf(l) * ddxy2(d)
3977        ENDDO
3978
3979    END SUBROUTINE usm_wall_heat_flux
3980 
3981 
3982!------------------------------------------------------------------------------!
3983! Description:
3984! ------------
3985!
3986!> This function applies the kinematic wall heat fluxes
3987!> for walls in four directions around the gridbox i,j.
3988!> It is called out from subroutine prognostic_equations.
3989!------------------------------------------------------------------------------!
3990    SUBROUTINE usm_wall_heat_flux_ij(i,j) 
3991   
3992        IMPLICIT NONE
3993
3994        INTEGER(iwp), INTENT(in)  ::  i,j                   !< indices of grid box
3995        INTEGER(iwp)              ::  ii,jj,k,d,l
3996       
3997        DO l = startenergy, endenergy
3998            jj = surfl(iy,l)
3999            ii = surfl(ix,l)
4000            IF ( ii == i  .AND.  jj == j ) THEN
4001               k = surfl(iz,l)
4002               IF ( k >=  nzb_s_inner(j,i)+1  .AND.  k <=  nzb_s_outer(j,i) ) THEN
4003                  d = surfl(id,l)
4004                  IF ( d >= 1 .and. d <= 4 )   THEN
4005                     tend(k,j,i) = tend(k,j,i) + wshf(l) * ddxy2(d)
4006                  ENDIF
4007               ENDIF
4008            ENDIF
4009        ENDDO
4010
4011    END SUBROUTINE usm_wall_heat_flux_ij
4012 
4013
4014!------------------------------------------------------------------------------!
4015!
4016! Description:
4017! ------------
4018!> Subroutine writes t_surf and t_wall data into restart files
4019!kanani: Renamed this routine according to corresponging routines in PALM
4020!kanani: Modified the routine to match write_var_list, from where usm_write_restart_data
4021!        shall be called in the future. This part has not been tested yet. (see virtual_flight_mod)
4022!        Also, I had some trouble with the allocation of t_surf, since this is a pointer.
4023!        So, I added some directives here.
4024!------------------------------------------------------------------------------!
4025    SUBROUTINE usm_write_restart_data
4026   
4027       IMPLICIT NONE
4028       
4029       INTEGER ::  i
4030
4031       DO  i = 0, io_blocks-1
4032          IF ( i == io_group )  THEN
4033             WRITE ( 14 )  't_surf                        '
4034#if defined( __nopointer )             
4035             WRITE ( 14 )  t_surf
4036#else
4037             WRITE ( 14 )  t_surf_1
4038#endif
4039             WRITE ( 14 )  't_wall                        '
4040#if defined( __nopointer )             
4041             WRITE ( 14 )  t_wall
4042#else
4043             WRITE ( 14 )  t_wall_1
4044#endif
4045             WRITE ( 14 )  '*** end usm ***               '
4046          ENDIF
4047#if defined( __parallel )
4048          CALL MPI_BARRIER( comm2d, ierr )
4049#endif
4050       ENDDO
4051
4052       
4053    END SUBROUTINE usm_write_restart_data
4054
4055
4056!------------------------------------------------------------------------------!
4057!
4058! Description:
4059! ------------
4060!> Subroutine stores svf and svfsurf data to files
4061!------------------------------------------------------------------------------!
4062    SUBROUTINE usm_write_svf_to_file
4063   
4064        IMPLICIT NONE
4065        INTEGER             :: fsvf = 89
4066        INTEGER             :: i
4067       
4068        DO  i = 0, io_blocks-1
4069            IF ( i == io_group )  THEN
4070                OPEN ( fsvf, file=TRIM(svf_file_name)//TRIM(coupling_char)//myid_char,               &
4071                    form='unformatted', status='new' )
4072
4073                WRITE ( fsvf )  usm_version
4074                WRITE ( fsvf )  nsvfcsfl, nsvfl
4075                WRITE ( fsvf )  svf
4076                WRITE ( fsvf )  svfsurf
4077                WRITE ( fsvf )  TRIM(svf_code)
4078               
4079                CLOSE (fsvf)
4080#if defined( __parallel )
4081                CALL mpi_barrier( comm2d, ierr )
4082#endif
4083            ENDIF
4084        ENDDO
4085    END SUBROUTINE usm_write_svf_to_file
4086
4087
4088!------------------------------------------------------------------------------!
4089! Description:
4090! ------------
4091!> Parin for &usm_par for urban surface model
4092!------------------------------------------------------------------------------!
4093    SUBROUTINE usm_parin
4094
4095       IMPLICIT NONE
4096
4097       CHARACTER (LEN=80) ::  line  !< string containing current line of file PARIN
4098
4099       NAMELIST /urban_surface_par/                                            & 
4100                           land_category,                                      &
4101                           mrt_factors,                                        &
4102                           nrefsteps,                                          &
4103                           pedestrant_category,                                &
4104                           ra_horiz_coef,                                      &
4105                           read_svf_on_init,                                   &
4106                           roof_category,                                      &
4107                           split_diffusion_radiation,                          &
4108                           urban_surface,                                      &
4109                           usm_anthropogenic_heat,                             &
4110                           usm_energy_balance_land,                            &
4111                           usm_energy_balance_wall,                            &
4112                           usm_material_model,                                 &
4113                           usm_lad_rma,                                        &
4114                           wall_category,                                      &
4115                           write_svf_on_init
4116
4117       line = ' '
4118       
4119!
4120!--    Try to find urban surface model package
4121       REWIND ( 11 )
4122       line = ' '
4123       DO   WHILE ( INDEX( line, '&urban_surface_par' ) == 0 )
4124          READ ( 11, '(A)', END=10 )  line
4125       ENDDO
4126       BACKSPACE ( 11 )
4127
4128!
4129!--    Read user-defined namelist
4130       READ ( 11, urban_surface_par )
4131
4132!
4133!--    Set flag that indicates that the land surface model is switched on
4134       urban_surface = .TRUE.
4135       
4136
4137 10    CONTINUE
4138   
4139    END SUBROUTINE usm_parin   
4140
4141   
4142!------------------------------------------------------------------------------!
4143!
4144! Description:
4145! ------------
4146!> Block of auxiliary subroutines:
4147!> 1. quicksort and corresponding comparison
4148!> 2. usm_merge_and_grow_csf for implementation of "dynamical growing" array
4149!>    for svf and csf
4150!------------------------------------------------------------------------------!   
4151    PURE FUNCTION svf_lt(svf1,svf2) result (res)
4152      TYPE (t_svf), INTENT(in) :: svf1,svf2
4153      LOGICAL                  :: res
4154      IF ( svf1%isurflt < svf2%isurflt  .OR.    &
4155          (svf1%isurflt == svf2%isurflt  .AND.  svf1%isurfs < svf2%isurfs) )  THEN
4156          res = .TRUE.
4157      ELSE
4158          res = .FALSE.
4159      ENDIF
4160    END FUNCTION svf_lt
4161   
4162 
4163!-- quicksort.f -*-f90-*-
4164!-- Author: t-nissie, adaptation J.Resler
4165!-- License: GPLv3
4166!-- Gist: https://gist.github.com/t-nissie/479f0f16966925fa29ea
4167    RECURSIVE SUBROUTINE quicksort_svf(svfl, first, last)
4168        IMPLICIT NONE
4169        TYPE(t_svf), DIMENSION(:), INTENT(INOUT)  :: svfl
4170        INTEGER, INTENT(IN)                       :: first, last
4171        TYPE(t_svf)                               :: x, t
4172        INTEGER                                   :: i, j
4173
4174        IF ( first>=last ) RETURN
4175        x = svfl( (first+last) / 2 )
4176        i = first
4177        j = last
4178        DO
4179            DO while ( svf_lt(svfl(i),x) )
4180                i=i+1
4181            ENDDO
4182            DO while ( svf_lt(x,svfl(j)) )
4183                j=j-1
4184            ENDDO
4185            IF ( i >= j ) EXIT
4186            t = svfl(i);  svfl(i) = svfl(j);  svfl(j) = t
4187            i=i+1
4188            j=j-1
4189        ENDDO
4190        IF ( first < i-1 ) CALL quicksort_svf(svfl, first, i-1)
4191        IF ( j+1 < last )  CALL quicksort_svf(svfl, j+1, last)
4192    END SUBROUTINE quicksort_svf
4193
4194   
4195    PURE FUNCTION csf_lt(csf1,csf2) result (res)
4196      TYPE (t_csf), INTENT(in) :: csf1,csf2
4197      LOGICAL                  :: res
4198      IF ( csf1%ip < csf2%ip  .OR.    &
4199           (csf1%ip == csf2%ip  .AND.  csf1%itx < csf2%itx)  .OR.  &
4200           (csf1%ip == csf2%ip  .AND.  csf1%itx == csf2%itx  .AND.  csf1%ity < csf2%ity)  .OR.  &
4201           (csf1%ip == csf2%ip  .AND.  csf1%itx == csf2%itx  .AND.  csf1%ity == csf2%ity  .AND.   &
4202            csf1%itz < csf2%itz)  .OR.  &
4203           (csf1%ip == csf2%ip  .AND.  csf1%itx == csf2%itx  .AND.  csf1%ity == csf2%ity  .AND.   &
4204            csf1%itz == csf2%itz  .AND.  csf1%isurfs < csf2%isurfs) )  THEN
4205          res = .TRUE.
4206      ELSE
4207          res = .FALSE.
4208      ENDIF
4209    END FUNCTION csf_lt
4210
4211
4212!-- quicksort.f -*-f90-*-
4213!-- Author: t-nissie, adaptation J.Resler
4214!-- License: GPLv3
4215!-- Gist: https://gist.github.com/t-nissie/479f0f16966925fa29ea
4216    RECURSIVE SUBROUTINE quicksort_csf(csfl, first, last)
4217        IMPLICIT NONE
4218        TYPE(t_csf), DIMENSION(:), INTENT(INOUT)  :: csfl
4219        INTEGER, INTENT(IN)                       :: first, last
4220        TYPE(t_csf)                               :: x, t
4221        INTEGER                                   :: i, j
4222
4223        IF ( first>=last ) RETURN
4224        x = csfl( (first+last)/2 )
4225        i = first
4226        j = last
4227        DO
4228            DO while ( csf_lt(csfl(i),x) )
4229                i=i+1
4230            ENDDO
4231            DO while ( csf_lt(x,csfl(j)) )
4232                j=j-1
4233            ENDDO
4234            IF ( i >= j ) EXIT
4235            t = csfl(i);  csfl(i) = csfl(j);  csfl(j) = t
4236            i=i+1
4237            j=j-1
4238        ENDDO
4239        IF ( first < i-1 ) CALL quicksort_csf(csfl, first, i-1)
4240        IF ( j+1 < last )  CALL quicksort_csf(csfl, j+1, last)
4241    END SUBROUTINE quicksort_csf
4242
4243   
4244    SUBROUTINE usm_merge_and_grow_csf(newsize)
4245        INTEGER(iwp), INTENT(in)            :: newsize  !< new array size after grow, must be >= ncsfl
4246                                                        !< or -1 to shrink to minimum
4247        INTEGER(iwp)                        :: iread, iwrite
4248        TYPE(t_csf), DIMENSION(:), POINTER  :: acsfnew
4249
4250        IF ( newsize == -1 )  THEN
4251!--         merge in-place
4252            acsfnew => acsf
4253        ELSE
4254!--         allocate new array
4255            IF ( mcsf == 0 )  THEN
4256                ALLOCATE( acsf1(newsize) )
4257                acsfnew => acsf1
4258            ELSE
4259                ALLOCATE( acsf2(newsize) )
4260                acsfnew => acsf2
4261            ENDIF
4262        ENDIF
4263
4264        IF ( ncsfl >= 1 )  THEN
4265!--         sort csf in place (quicksort)
4266            CALL quicksort_csf(acsf,1,ncsfl)
4267
4268!--         while moving to a new array, aggregate canopy sink factor records with identical box & source
4269            acsfnew(1) = acsf(1)
4270            iwrite = 1
4271            DO iread = 2, ncsfl
4272!--             here acsf(kcsf) already has values from acsf(icsf)
4273                IF ( acsfnew(iwrite)%itx == acsf(iread)%itx &
4274                         .AND.  acsfnew(iwrite)%ity == acsf(iread)%ity &
4275                         .AND.  acsfnew(iwrite)%itz == acsf(iread)%itz &
4276                         .AND.  acsfnew(iwrite)%isurfs == acsf(iread)%isurfs )  THEN
4277!--                 We could simply take either first or second rtransp, both are valid. As a very simple heuristic about which ray
4278!--                 probably passes nearer the center of the target box, we choose DIF from the entry with greater CSF, since that
4279!--                 might mean that the traced beam passes longer through the canopy box.
4280                    IF ( acsfnew(iwrite)%rsvf < acsf(iread)%rsvf )  THEN
4281                        acsfnew(iwrite)%rtransp = acsf(iread)%rtransp
4282                    ENDIF
4283                    acsfnew(iwrite)%rsvf = acsfnew(iwrite)%rsvf + acsf(iread)%rsvf
4284!--                 advance reading index, keep writing index
4285                ELSE
4286!--                 not identical, just advance and copy
4287                    iwrite = iwrite + 1
4288                    acsfnew(iwrite) = acsf(iread)
4289                ENDIF
4290            ENDDO
4291            ncsfl = iwrite
4292        ENDIF
4293
4294        IF ( newsize == -1 )  THEN
4295!--         allocate new array and copy shrinked data
4296            IF ( mcsf == 0 )  THEN
4297                ALLOCATE( acsf1(ncsfl) )
4298                acsf1(1:ncsfl) = acsf2(1:ncsfl)
4299            ELSE
4300                ALLOCATE( acsf2(ncsfl) )
4301                acsf2(1:ncsfl) = acsf1(1:ncsfl)
4302            ENDIF
4303        ENDIF
4304
4305!--     deallocate old array
4306        IF ( mcsf == 0 )  THEN
4307            mcsf = 1
4308            acsf => acsf1
4309            DEALLOCATE( acsf2 )
4310        ELSE
4311            mcsf = 0
4312            acsf => acsf2
4313            DEALLOCATE( acsf1 )
4314        ENDIF
4315        ncsfla = newsize
4316    END SUBROUTINE usm_merge_and_grow_csf
4317
4318   
4319!-- quicksort.f -*-f90-*-
4320!-- Author: t-nissie, adaptation J.Resler
4321!-- License: GPLv3
4322!-- Gist: https://gist.github.com/t-nissie/479f0f16966925fa29ea
4323    RECURSIVE SUBROUTINE quicksort_csf2(kpcsflt, pcsflt, first, last)
4324        IMPLICIT NONE
4325        INTEGER(iwp), DIMENSION(:,:), INTENT(INOUT)  :: kpcsflt
4326        REAL(wp), DIMENSION(:,:), INTENT(INOUT)      :: pcsflt
4327        INTEGER, INTENT(IN)                          :: first, last
4328        REAL(wp), DIMENSION(ndcsf)                   :: t2
4329        INTEGER(iwp), DIMENSION(kdcsf)               :: x, t1
4330        INTEGER                                      :: i, j
4331
4332        IF ( first>=last ) RETURN
4333        x = kpcsflt(:, (first+last)/2 )
4334        i = first
4335        j = last
4336        DO
4337            DO while ( csf_lt2(kpcsflt(:,i),x) )
4338                i=i+1
4339            ENDDO
4340            DO while ( csf_lt2(x,kpcsflt(:,j)) )
4341                j=j-1
4342            ENDDO
4343            IF ( i >= j ) EXIT
4344            t1 = kpcsflt(:,i);  kpcsflt(:,i) = kpcsflt(:,j);  kpcsflt(:,j) = t1
4345            t2 = pcsflt(:,i);  pcsflt(:,i) = pcsflt(:,j);  pcsflt(:,j) = t2
4346            i=i+1
4347            j=j-1
4348        ENDDO
4349        IF ( first < i-1 ) CALL quicksort_csf2(kpcsflt, pcsflt, first, i-1)
4350        IF ( j+1 < last )  CALL quicksort_csf2(kpcsflt, pcsflt, j+1, last)
4351    END SUBROUTINE quicksort_csf2
4352   
4353
4354    PURE FUNCTION csf_lt2(item1, item2) result(res)
4355        INTEGER(iwp), DIMENSION(kdcsf), INTENT(in)  :: item1, item2
4356        LOGICAL                                     :: res
4357        res = ( (item1(3) < item2(3))                                                        &
4358             .OR.  (item1(3) == item2(3)  .AND.  item1(2) < item2(2))                            &
4359             .OR.  (item1(3) == item2(3)  .AND.  item1(2) == item2(2)  .AND.  item1(1) < item2(1)) &
4360             .OR.  (item1(3) == item2(3)  .AND.  item1(2) == item2(2)  .AND.  item1(1) == item2(1) &
4361                 .AND.  item1(4) < item2(4)) )
4362    END FUNCTION csf_lt2
4363
4364   
4365 END MODULE urban_surface_mod
Note: See TracBrowser for help on using the repository browser.