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

Last change on this file since 2226 was 2214, checked in by kanani, 7 years ago

last commit documented

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