source: palm/trunk/SOURCE/diffusivities.f90 @ 1356

Last change on this file since 1356 was 1354, checked in by heinze, 11 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 7.6 KB
RevLine 
[97]1 SUBROUTINE diffusivities( var, var_reference )
[1]2
[1036]3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later 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!
[1310]17! Copyright 1997-2014 Leibniz Universitaet Hannover
[1036]18!--------------------------------------------------------------------------------!
19!
[484]20! Current revisions:
[1]21! -----------------
[1341]22!
[1354]23!
[1321]24! Former revisions:
25! -----------------
26! $Id: diffusivities.f90 1354 2014-04-08 15:22:57Z heinze $
27!
[1354]28! 1353 2014-04-08 15:21:23Z heinze
29! REAL constants provided with KIND-attribute
30!
[1341]31! 1340 2014-03-25 19:45:13Z kanani
32! REAL constants defined as wp-kind
33!
[1323]34! 1322 2014-03-20 16:38:49Z raasch
35! REAL constants defined as wp-kind
36!
[1321]37! 1320 2014-03-20 08:40:49Z raasch
[1320]38! ONLY-attribute added to USE-statements,
39! kind-parameters added to all INTEGER and REAL declaration statements,
40! kinds are defined in new module kinds,
41! revision history before 2012 removed,
42! comment fields (!:) to be used for variable explanations added to
43! all variable declaration statements
[1017]44!
[1182]45! 1179 2013-06-14 05:57:58Z raasch
46! use_reference renamed use_single_reference_value
47!
[1037]48! 1036 2012-10-22 13:43:42Z raasch
49! code put under GPL (PALM 3.9)
50!
[1017]51! 1015 2012-09-27 09:23:24Z raasch
52! OpenACC statements added + code changes required for GPU optimization,
53! adjustment of mixing length to the Prandtl mixing length at first grid point
54! above ground removed
55!
[1]56! Revision 1.1  1997/09/19 07:41:10  raasch
57! Initial revision
58!
59!
60! Description:
61! ------------
62! Computation of the turbulent diffusion coefficients for momentum and heat
63! according to Prandtl-Kolmogorov
64!------------------------------------------------------------------------------!
65
[1320]66    USE arrays_3d,                                                             &
67        ONLY:  dd2zu, e, kh, km, l_grid, l_wall
68       
69    USE control_parameters,                                                    &
70        ONLY:  atmos_ocean_sign, e_min, g, outflow_l, outflow_n, outflow_r,    &
71                outflow_s, use_single_reference_value, wall_adjustment
72               
73    USE indices,                                                               &
74        ONLY:  nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb_s_inner, nzb, nzt
75    USE kinds
76   
[1]77    USE pegrid
[1320]78   
79    USE statistics,                                                            &
80        ONLY :  rmask, statistic_regions, sums_l_l
[1]81
82    IMPLICIT NONE
83
[1320]84    INTEGER(iwp) ::  i                   !:
85    INTEGER(iwp) ::  j                   !:
86    INTEGER(iwp) ::  k                   !:
87    INTEGER(iwp) ::  omp_get_thread_num  !:
88    INTEGER(iwp) ::  sr                  !:
89    INTEGER(iwp) ::  tn                  !:
[1]90
[1320]91    REAL(wp)     ::  dvar_dz             !:
92    REAL(wp)     ::  l                   !:
93    REAL(wp)     ::  ll                  !:
94    REAL(wp)     ::  l_stable            !:
95    REAL(wp)     ::  sqrt_e              !:
96    REAL(wp)     ::  var_reference       !:
[1]97
[1320]98    REAL(wp)     ::  var(nzb:nzt+1,nysg:nyng,nxlg:nxrg)  !:
[97]99
[1]100
101!
102!-- Default thread number in case of one thread
103    tn = 0
104
105!
106!-- Initialization for calculation of the mixing length profile
[1340]107    sums_l_l = 0.0_wp
[1]108
109!
110!-- Compute the turbulent diffusion coefficient for momentum
[1015]111    !$OMP PARALLEL PRIVATE (dvar_dz,i,j,k,l,ll,l_stable,sqrt_e,sr,tn)
[1]112!$  tn = omp_get_thread_num()
113
[1015]114!
115!-- Data declerations for accelerators
116    !$acc data present( dd2zu, e, km, kh, l_grid, l_wall, nzb_s_inner, rif, var )
117    !$acc kernels
118
119!
120!-- Introduce an optional minimum tke
[1340]121    IF ( e_min > 0.0_wp )  THEN
[1015]122       !$OMP DO
123       !$acc loop
124       DO  i = nxlg, nxrg
125          DO  j = nysg, nyng
126             !$acc loop vector( 32 )
127             DO  k = 1, nzt
128                IF ( k > nzb_s_inner(j,i) )  THEN
129                   e(k,j,i) = MAX( e(k,j,i), e_min )
130                ENDIF
131             ENDDO
132          ENDDO
133       ENDDO
134    ENDIF
135
[1]136    !$OMP DO
[1015]137    !$acc loop
[667]138    DO  i = nxlg, nxrg
139       DO  j = nysg, nyng
[1015]140          !$acc loop vector( 32 )
141          DO  k = 1, nzt
[1]142
[1015]143             IF ( k > nzb_s_inner(j,i) )  THEN
[1]144
[1015]145                sqrt_e = SQRT( e(k,j,i) )
[1]146!
[1015]147!--             Determine the mixing length
148                dvar_dz = atmos_ocean_sign * &  ! inverse effect of pt/rho gradient
149                          ( var(k+1,j,i) - var(k-1,j,i) ) * dd2zu(k)
[1340]150                IF ( dvar_dz > 0.0_wp ) THEN
[1179]151                   IF ( use_single_reference_value )  THEN
[1322]152                      l_stable = 0.76_wp * sqrt_e /                            &
[1340]153                                    SQRT( g / var_reference * dvar_dz ) + 1E-5_wp
[1015]154                   ELSE
[1322]155                      l_stable = 0.76_wp * sqrt_e /                            &
[1340]156                                    SQRT( g / var(k,j,i) * dvar_dz ) + 1E-5_wp
[1015]157                   ENDIF
158                ELSE
159                   l_stable = l_grid(k)
160                ENDIF
[1]161!
[1015]162!--             Adjustment of the mixing length
163                IF ( wall_adjustment )  THEN
164                   l  = MIN( l_wall(k,j,i), l_grid(k), l_stable )
165                   ll = MIN( l_wall(k,j,i), l_grid(k) )
[57]166                ELSE
[1015]167                   l  = MIN( l_grid(k), l_stable )
168                   ll = l_grid(k)
[57]169                ENDIF
[1015]170
171      !
172      !--       Compute diffusion coefficients for momentum and heat
[1322]173                km(k,j,i) = 0.1_wp * l * sqrt_e
174                kh(k,j,i) = ( 1.0_wp + 2.0_wp * l / ll ) * km(k,j,i)
[1015]175
[1]176             ENDIF
177
[1015]178#if ! defined( __openacc )
[1]179!
[1015]180!++          Statistics still have to be realized for accelerators
181!--          Summation for averaged profile (cf. flow_statistics)
182             DO  sr = 0, statistic_regions
183                sums_l_l(k,sr,tn) = sums_l_l(k,sr,tn) + l * rmask(j,i,sr)
184             ENDDO
185#endif
[1]186
187          ENDDO
188       ENDDO
189    ENDDO
190
[1015]191#if ! defined( __openacc )
192!
193!++ Statistics still have to be realized for accelerators
[1]194    sums_l_l(nzt+1,:,tn) = sums_l_l(nzt,:,tn)   ! quasi boundary-condition for
[667]195                                                  ! data output
[1015]196#endif
[1]197    !$OMP END PARALLEL
198
199!
200!-- Set vertical boundary values (Neumann conditions both at bottom and top).
201!-- Horizontal boundary conditions at vertical walls are not set because
202!-- so far vertical walls require usage of a Prandtl-layer where the boundary
203!-- values of the diffusivities are not needed
204    !$OMP PARALLEL DO
[1015]205    !$acc loop
[667]206    DO  i = nxlg, nxrg
207       DO  j = nysg, nyng
[1]208          km(nzb_s_inner(j,i),j,i) = km(nzb_s_inner(j,i)+1,j,i)
209          km(nzt+1,j,i)            = km(nzt,j,i)
210          kh(nzb_s_inner(j,i),j,i) = kh(nzb_s_inner(j,i)+1,j,i)
211          kh(nzt+1,j,i)            = kh(nzt,j,i)
212       ENDDO
213    ENDDO
214
215!
216!-- Set Neumann boundary conditions at the outflow boundaries in case of
217!-- non-cyclic lateral boundaries
218    IF ( outflow_l )  THEN
219       km(:,:,nxl-1) = km(:,:,nxl)
220       kh(:,:,nxl-1) = kh(:,:,nxl)
221    ENDIF
222    IF ( outflow_r )  THEN
223       km(:,:,nxr+1) = km(:,:,nxr)
224       kh(:,:,nxr+1) = kh(:,:,nxr)
225    ENDIF
226    IF ( outflow_s )  THEN
227       km(:,nys-1,:) = km(:,nys,:)
228       kh(:,nys-1,:) = kh(:,nys,:)
229    ENDIF
230    IF ( outflow_n )  THEN
231       km(:,nyn+1,:) = km(:,nyn,:)
232       kh(:,nyn+1,:) = kh(:,nyn,:)
233    ENDIF
234
[1015]235    !$acc end kernels
236    !$acc end data
[1]237
238 END SUBROUTINE diffusivities
Note: See TracBrowser for help on using the repository browser.