source: palm/trunk/SOURCE/surface_layer_fluxes.f90 @ 1748

Last change on this file since 1748 was 1748, checked in by raasch, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 42.3 KB
Line 
1!> @file surface_layer_fluxes.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 terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2015 Leibniz Universitaet Hannover
17!
18!--------------------------------------------------------------------------------!
19! Current revisions:
20! ------------------
21!
22!
23! Former revisions:
24! -----------------
25! $Id: surface_layer_fluxes.f90 1748 2016-02-08 12:38:17Z raasch $
26!
27! 1747 2016-02-08 12:25:53Z raasch
28! adjustments for OpenACC usage
29!
30! 1709 2015-11-04 14:47:01Z maronga
31! Bugfix: division by zero could occur when calculating rib at low wind speeds
32! Bugfix: calculation of uv_total for neutral = .T., initial value for ol for
33! neutral = .T.
34!
35! 1705 2015-11-02 14:28:56Z maronga
36! Typo removed
37!
38! 1697 2015-10-28 17:14:10Z raasch
39! FORTRAN and OpenMP errors removed
40!
41! 1696 2015-10-27 10:03:34Z maronga
42! Modularized and completely re-written version of prandtl_fluxes.f90. In the
43! course of the re-writing two additional methods have been implemented. See
44! updated description.
45!
46! 1551 2015-03-03 14:18:16Z maronga
47! Removed land surface model part. The surface fluxes are now always calculated
48! within prandtl_fluxes, based on the given surface temperature/humidity (which
49! is either provided by the land surface model, by large scale forcing data, or
50! directly prescribed by the user.
51!
52! 1496 2014-12-02 17:25:50Z maronga
53! Adapted for land surface model
54!
55! 1494 2014-11-21 17:14:03Z maronga
56! Bugfixes: qs is now calculated before calculation of Rif. calculation of
57! buoyancy flux in Rif corrected (added missing humidity term), allow use of
58! topography for coupled runs (not tested)
59!
60! 1361 2014-04-16 15:17:48Z hoffmann
61! Bugfix: calculation of turbulent fluxes of rain water content (qrsws) and rain
62! drop concentration (nrsws) added
63!
64! 1340 2014-03-25 19:45:13Z kanani
65! REAL constants defined as wp-kind
66!
67! 1320 2014-03-20 08:40:49Z raasch
68! ONLY-attribute added to USE-statements,
69! kind-parameters added to all INTEGER and REAL declaration statements,
70! kinds are defined in new module kinds,
71! old module precision_kind is removed,
72! revision history before 2012 removed,
73! comment fields (!:) to be used for variable explanations added to
74! all variable declaration statements
75!
76! 1276 2014-01-15 13:40:41Z heinze
77! Use LSF_DATA also in case of Dirichlet bottom boundary condition for scalars
78!
79! 1257 2013-11-08 15:18:40Z raasch
80! openACC "kernels do" replaced by "kernels loop", "loop independent" added
81!
82! 1036 2012-10-22 13:43:42Z raasch
83! code put under GPL (PALM 3.9)
84!
85! 1015 2012-09-27 09:23:24Z raasch
86! OpenACC statements added
87!
88! 978 2012-08-09 08:28:32Z fricke
89! roughness length for scalar quantities z0h added
90!
91! Revision 1.1  1998/01/23 10:06:06  raasch
92! Initial revision
93!
94!
95! Description:
96! ------------
97!> Diagnostic computation of vertical fluxes in the constant flux layer from the
98!> values of the variables at grid point k=1. Three different methods are
99!> available:
100!> 1) the "old" version (most_method = 'circular') which is fast, but inaccurate
101!> 2) a Newton iteration method (most_method = 'newton'), which is accurate, but
102!>    slower
103!> 3) a method using a lookup table which is fast and accurate. Note, however,
104!>    that this method cannot be used in case of roughness heterogeneity
105!>
106!> @todo (re)move large_scale_forcing actions
107!> @todo check/optimize OpenMP and OpenACC directives
108!------------------------------------------------------------------------------!
109 MODULE surface_layer_fluxes_mod
110
111    USE arrays_3d,                                                             &
112        ONLY:  e, kh, nr, nrs, nrsws, ol, pt, q, ql, qr, qrs, qrsws, qs, qsws, &
113               shf, ts, u, us, usws, v, vpt, vsws, zu, zw, z0, z0h
114
115    USE cloud_parameters,                                                      &
116        ONLY:  l_d_cp, pt_d_t
117
118    USE constants,                                                             &
119        ONLY:  pi
120
121    USE cpulog
122
123    USE control_parameters,                                                    &
124        ONLY:  cloud_physics, constant_heatflux, constant_waterflux,           &
125               coupling_mode, g, humidity, ibc_e_b, ibc_pt_b, icloud_scheme,   &
126               initializing_actions, kappa, intermediate_timestep_count,       &
127               intermediate_timestep_count_max, large_scale_forcing, lsf_surf, &
128               message_string, most_method, neutral, passive_scalar,           &
129               precipitation, pt_surface, q_surface, run_coupled,              &
130               surface_pressure, simulated_time, terminate_run, zeta_max,      &
131               zeta_min 
132
133    USE indices,                                                               &
134        ONLY:  nxl, nxlg, nxr, nxrg, nys, nysg, nyn, nyng, nzb_s_inner,        &
135               nzb_u_inner, nzb_v_inner
136
137    USE kinds
138
139    USE pegrid
140
141    USE land_surface_model_mod,                                                &
142        ONLY:  land_surface, skip_time_do_lsm
143
144    IMPLICIT NONE
145
146    INTEGER(iwp) ::  i            !< loop index x direction
147    INTEGER(iwp) ::  j            !< loop index y direction
148    INTEGER(iwp) ::  k            !< loop index z direction
149
150    INTEGER(iwp), PARAMETER     :: num_steps = 15000  !< number of steps in the lookup table
151
152    LOGICAL      ::  coupled_run  !< Flag for coupled atmosphere-ocean runs
153
154    REAL(wp), DIMENSION(:,:), ALLOCATABLE :: pt1,      & !< Potential temperature at first grid level (required for cloud_physics = .T.)
155                                             qv1,      & !< Specific humidity at first grid level (required for cloud_physics = .T.)
156                                             uv_total    !< Total velocity at first grid level
157
158    REAL(wp), DIMENSION(0:num_steps-1) :: rib_tab,  & !< Lookup table bulk Richardson number
159                                          ol_tab      !< Lookup table values of L
160
161    REAL(wp)     ::  e_s,               & !< Saturation water vapor pressure
162                     l_bnd  = 7500,     & !< Lookup table index of the last time step
163                     ol_max = 1.0E6_wp, & !< Maximum Obukhov length
164                     rib_max,           & !< Maximum Richardson number in lookup table
165                     rib_min,           & !< Minimum Richardson number in lookup table
166                     z_mo                 !< Height of the constant flux layer where MOST is assumed
167
168
169    SAVE
170
171    PRIVATE
172
173    PUBLIC init_surface_layer_fluxes, pt1, qv1, surface_layer_fluxes, uv_total
174
175    INTERFACE init_surface_layer_fluxes
176       MODULE PROCEDURE init_surface_layer_fluxes
177    END INTERFACE init_surface_layer_fluxes
178
179    INTERFACE surface_layer_fluxes
180       MODULE PROCEDURE surface_layer_fluxes
181    END INTERFACE surface_layer_fluxes
182
183
184 CONTAINS
185
186
187!------------------------------------------------------------------------------!
188! Description:
189! ------------
190!> Main routine to compute the surface fluxes
191!------------------------------------------------------------------------------!
192    SUBROUTINE surface_layer_fluxes
193
194       IMPLICIT NONE
195
196!
197!--    In case cloud physics is used, it is required to derive potential
198!--    temperature and specific humidity at first grid level from the fields pt
199!--    and q
200       IF ( cloud_physics )  THEN
201          CALL calc_pt_q
202       ENDIF
203
204!
205!--    First, calculate the new Obukhov length, then new friction velocity,
206!--    followed by the new scaling parameters (th*, q*, etc.), and the new
207!--    surface fluxes if required. The old routine ("circular") requires a
208!--    different order of calls as the scaling parameters from the previous time
209!--    steps are used to calculate the Obukhov length
210
211!
212!--    Depending on setting of most_method use the "old" routine
213       IF ( most_method == 'circular' )  THEN
214
215          CALL calc_scaling_parameters
216
217          CALL calc_uv_total
218
219          IF ( .NOT. neutral )  THEN
220             CALL calc_ol
221          ENDIF
222
223          CALL calc_us
224
225          CALL calc_surface_fluxes
226
227!
228!--    Use either Newton iteration or a lookup table for the bulk Richardson
229!--    number to calculate the Obukhov length
230       ELSEIF ( most_method == 'newton' .OR. most_method == 'lookup' )  THEN
231
232          CALL calc_uv_total
233
234          IF ( .NOT. neutral )  THEN
235             CALL calc_ol
236          ENDIF
237
238          CALL calc_us
239
240          CALL calc_scaling_parameters
241
242          CALL calc_surface_fluxes
243
244       ENDIF
245
246    END SUBROUTINE surface_layer_fluxes
247
248
249!------------------------------------------------------------------------------!
250! Description:
251! ------------
252!> Initializing actions for the surface layer routine. Basically, this involves
253!> the preparation of a lookup table for the the bulk Richardson number vs
254!> Obukhov length L when using the lookup table method.
255!------------------------------------------------------------------------------!
256    SUBROUTINE init_surface_layer_fluxes
257
258       IMPLICIT NONE
259
260       INTEGER(iwp) :: l,          & !< Index for loop to create lookup table
261                       num_steps_n   !< Number of non-stretched zeta steps
262
263       LOGICAL :: terminate_run_l = .FALSE.    !< Flag to terminate run (global)
264
265       REAL(wp), PARAMETER ::  zeta_stretch = -10.0_wp !< Start of stretching in the free convection limit
266                               
267       REAL(wp), DIMENSION(:), ALLOCATABLE :: zeta_tmp
268
269
270       REAL(wp) :: zeta_step,            & !< Increment of zeta
271                   regr      = 1.01_wp,  & !< Stretching factor of zeta_step in the free convection limit
272                   regr_old  = 1.0E9_wp, & !< Stretching factor of last iteration step
273                   z0h_min   = 0.0_wp,   & !< Minimum value of z0h to create table
274                   z0_min    = 0.0_wp      !< Minimum value of z0 to create table
275!
276!--    When cloud physics is used, arrays for storing potential temperature and
277!--    specific humidity at first grid level are required
278       IF ( cloud_physics )  THEN
279          ALLOCATE ( pt1(nysg:nyng,nxlg:nxrg) )
280          ALLOCATE ( qv1(nysg:nyng,nxlg:nxrg) )
281       ENDIF
282
283!
284!--    Allocate field for storing the horizontal velocity
285       ALLOCATE ( uv_total(nysg:nyng,nxlg:nxrg) )
286
287
288!
289!--    In case of runs with neutral statification, set Obukhov length to a
290!--    large value
291       IF ( neutral ) ol = 1.0E10_wp
292
293       IF ( most_method == 'lookup' )  THEN
294
295!
296!--       Check for roughness heterogeneity. In that case terminate run and
297!--       inform user
298          IF ( MINVAL( z0h ) /= MAXVAL( z0h ) .OR.                             &
299               MINVAL( z0  ) /= MAXVAL( z0  ) )  THEN
300             terminate_run_l = .TRUE.
301          ENDIF
302
303#if defined( __parallel )
304!
305!--       Make a logical OR for all processes. Force termiation of model if result
306!--       is TRUE
307          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
308          CALL MPI_ALLREDUCE( terminate_run_l, terminate_run, 1, MPI_LOGICAL,  &
309                              MPI_LOR, comm2d, ierr )
310#else
311          terminate_run = terminate_run_l
312#endif
313
314          IF ( terminate_run )  THEN
315             message_string = 'most_method = "lookup" cannot be used in ' //   &
316                              'combination with a prescribed roughness '  //   &
317                              'heterogeneity'
318             CALL message( 'surface_layer_fluxes', 'PA0417', 1, 2, 0, 6, 0 )
319          ENDIF
320
321          ALLOCATE(  zeta_tmp(0:num_steps-1) )
322
323!
324!--       Use the lowest possible value for z_mo
325          k    = MINVAL(nzb_s_inner)
326          z_mo = zu(k+1) - zw(k)
327
328!
329!--       Calculate z/L range from zeta_stretch to zeta_max using 90% of the
330!--       available steps (num_steps). The calculation is done with negative
331!--       values of zeta in order to simplify the stretching in the free
332!--       convection limit for the remaining 10% of steps.
333          zeta_tmp(0) = - zeta_max
334          num_steps_n = ( num_steps * 9 / 10 ) - 1
335          zeta_step   = (zeta_max - zeta_stretch) / REAL(num_steps_n)
336
337          DO l = 1, num_steps_n
338             zeta_tmp(l) = zeta_tmp(l-1) + zeta_step
339          ENDDO
340
341!
342!--       Calculate stretching factor for the free convection range
343          DO WHILE ( ABS( (regr-regr_old) / regr_old ) > 1.0E-10_wp )
344             regr_old = regr
345             regr = ( 1.0_wp - ( -zeta_min / zeta_step ) * ( 1.0_wp - regr )   &
346                    )**( 10.0_wp / REAL(num_steps) )
347          ENDDO
348
349!
350!--       Calculate z/L range from zeta_min to zeta_stretch
351          DO l = num_steps_n+1, num_steps-1
352             zeta_tmp(l) = zeta_tmp(l-1) + zeta_step
353             zeta_step = zeta_step * regr
354          ENDDO
355
356!
357!--       Invert array and switch sign, then calculate Obukhov length and bulk
358!--       Richardson number
359          z0h_min = MINVAL(z0h)
360          z0_min  = MINVAL(z0)
361         
362!
363!--       Calculate lookup table for the Richardson number versus Obukhov length
364!--       The Richardson number (rib) is defined depending on the choice of
365!--       boundary conditions for temperature
366          IF ( ibc_pt_b == 1 )  THEN
367             DO l = 0, num_steps-1
368                ol_tab(l)  = - z_mo / zeta_tmp(num_steps-1-l)
369                rib_tab(l) = z_mo / ol_tab(l)  / ( LOG( z_mo / z0_min )        &
370                                                - psi_m( z_mo / ol_tab(l) )    &
371                                                + psi_m( z0_min / ol_tab(l) )  &
372                                                  )**3
373             ENDDO 
374          ELSE
375             DO l = 0, num_steps-1
376                ol_tab(l)  = - z_mo / zeta_tmp(num_steps-1-l)
377                rib_tab(l) = z_mo / ol_tab(l)  * ( LOG( z_mo / z0h_min )       &
378                                              - psi_h( z_mo / ol_tab(l) )      &
379                                              + psi_h( z0h_min / ol_tab(l) )   &
380                                            )                                  &
381                                          / ( LOG( z_mo / z0_min )             &
382                                              - psi_m( z_mo / ol_tab(l) )      &
383                                              + psi_m( z0_min / ol_tab(l) )    &
384                                            )**2
385             ENDDO
386          ENDIF
387
388!
389!--       Determine minimum values of rib in the lookup table. Set upper limit
390!--       to critical Richardson number (0.25)
391          rib_min  = MINVAL(rib_tab)
392          rib_max  = 0.25 !MAXVAL(rib_tab)
393
394          DEALLOCATE( zeta_tmp )
395       ENDIF
396
397    END SUBROUTINE init_surface_layer_fluxes
398
399
400!------------------------------------------------------------------------------!
401! Description:
402! ------------
403!> Compute the absolute value of the horizontal velocity (relative to the
404!> surface). This is required by all methods
405!------------------------------------------------------------------------------!
406    SUBROUTINE calc_uv_total
407
408       IMPLICIT NONE
409
410
411       !$OMP PARALLEL DO PRIVATE( k )
412       !$acc kernels loop present( nzb_s_inner, u, uv_total, v ) private( j, k )
413       DO  i = nxl, nxr
414          DO  j = nys, nyn
415
416             k   = nzb_s_inner(j,i)
417             uv_total(j,i) = SQRT( ( 0.5_wp * ( u(k+1,j,i) + u(k+1,j,i+1)      &
418                                         - u(k,j,i)   - u(k,j,i+1) ) )**2 +    &
419                              ( 0.5_wp * ( v(k+1,j,i) + v(k+1,j+1,i)           &
420                                         - v(k,j,i)   - v(k,j+1,i) ) )**2 )
421
422!
423!--          For too small values of the local wind, MOST does not work. A
424!--          threshold value is thus set if required
425!            uv_total(j,i) = MAX(0.01_wp,uv_total(j,i))
426
427          ENDDO
428       ENDDO
429
430!
431!--    Values of uv_total need to be exchanged at the ghost boundaries
432       !$acc update host( uv_total )
433       CALL exchange_horiz_2d( uv_total )
434       !$acc update device( uv_total )
435
436    END SUBROUTINE calc_uv_total
437
438
439!------------------------------------------------------------------------------!
440! Description:
441! ------------
442!> Calculate the Obukhov length (L) and Richardson flux number (z/L)
443!------------------------------------------------------------------------------!
444    SUBROUTINE calc_ol
445
446       IMPLICIT NONE
447
448       INTEGER(iwp) :: iter,  & !< Newton iteration step
449                       l        !< look index
450
451       REAL(wp), DIMENSION(nysg:nyng,nxlg:nxrg) :: rib !< Bulk Richardson number
452
453       REAL(wp)     :: f,      & !< Function for Newton iteration: f = Ri - [...]/[...]^2 = 0
454                       f_d_ol, & !< Derivative of f
455                       ol_l,   & !< Lower bound of L for Newton iteration
456                       ol_m,   & !< Previous value of L for Newton iteration
457                       ol_old, & !< Previous time step value of L
458                       ol_u      !< Upper bound of L for Newton iteration
459
460       IF ( TRIM( most_method ) /= 'circular' )  THEN
461     
462          !$acc data present( nzb_s_inner, pt, q, qsws, rib, shf, uv_total, vpt, zu, zw )
463
464          !$OMP PARALLEL DO PRIVATE( k, z_mo )
465          !$acc kernels loop private( j, k, z_mo )
466          DO  i = nxl, nxr
467             DO  j = nys, nyn
468
469                k   = nzb_s_inner(j,i)
470                z_mo = zu(k+1) - zw(k)
471
472!
473!--             Evaluate bulk Richardson number (calculation depends on
474!--             definition based on setting of boundary conditions
475                IF ( ibc_pt_b /= 1 )  THEN
476                   IF ( humidity )  THEN
477                      rib(j,i) = g * z_mo * ( vpt(k+1,j,i) - vpt(k,j,i) )      &
478                           / ( uv_total(j,i)**2 * vpt(k+1,j,i) + 1.0E-20_wp )
479                   ELSE
480                      rib(j,i) = g * z_mo * ( pt(k+1,j,i) - pt(k,j,i) )        &
481                           / ( uv_total(j,i)**2 * pt(k+1,j,i)  + 1.0E-20_wp )
482                   ENDIF     
483                ELSE
484!
485!--                When using Neumann boundary conditions, the buoyancy flux
486!--                is required but cannot be calculated at the surface, as pt
487!--                and q are not known at the surface. Hence the values at
488!--                first grid level are used to estimate the buoyancy flux
489                   IF ( humidity )  THEN
490                      rib(j,i) = - g * z_mo * ( ( 1.0_wp + 0.61_wp             &
491                                 * q(k+1,j,i) ) * shf(j,i) + 0.61_wp           &
492                                 * pt(k+1,j,i) * qsws(j,i) )   &
493                                 / ( uv_total(j,i)**3 * vpt(k+1,j,i) * kappa**2&
494                                 + 1.0E-20_wp)
495                   ELSE
496                      rib(j,i) = - g * z_mo * shf(j,i)                         &
497                           / ( uv_total(j,i)**3 * pt(k+1,j,i) * kappa**2       &
498                           + 1.0E-20_wp )
499                   ENDIF
500                ENDIF 
501     
502             ENDDO
503          ENDDO 
504          !$acc end data
505
506       ENDIF
507
508!
509!--    Calculate the Obukhov length either using a Newton iteration
510!--    method, via a lookup table, or using the old circular way
511       IF ( TRIM( most_method ) == 'newton' )  THEN
512
513          !$OMP PARALLEL DO PRIVATE( k, z_mo )
514          !# WARNING: does not work on GPU so far
515          !!!!!!$acc kernels loop
516          DO  i = nxl, nxr
517             DO  j = nys, nyn
518
519                k   = nzb_s_inner(j,i)
520                z_mo = zu(k+1) - zw(k)
521
522!
523!--             Store current value in case the Newton iteration fails
524                ol_old = ol(j,i)
525
526!
527!--             Ensure that the bulk Richardson number and the Obukhov
528!--             lengtH have the same sign
529                IF ( rib(j,i) * ol(j,i) < 0.0_wp .OR.                        &
530                     ABS( ol(j,i) ) == ol_max )  THEN
531                   IF ( rib(j,i) > 0.0_wp ) ol(j,i) =  0.01_wp
532                   IF ( rib(j,i) < 0.0_wp ) ol(j,i) = -0.01_wp
533                ENDIF
534!
535!--             Iteration to find Obukhov length
536                iter = 0
537                DO
538                   iter = iter + 1
539!
540!--                In case of divergence, use the value of the previous time step
541                   IF ( iter > 1000 )  THEN
542                      ol(j,i) = ol_old
543                      EXIT
544                   ENDIF
545
546                   ol_m = ol(j,i)
547                   ol_l = ol_m - 0.001_wp * ol_m
548                   ol_u = ol_m + 0.001_wp * ol_m
549
550
551                   IF ( ibc_pt_b /= 1 )  THEN
552!
553!--                   Calculate f = Ri - [...]/[...]^2 = 0
554                      f = rib(j,i) - ( z_mo / ol_m ) * ( LOG( z_mo / z0h(j,i) )&
555                                                    - psi_h( z_mo / ol_m )     &
556                                                    + psi_h( z0h(j,i) / ol_m ) &
557                                                   )                           &
558                                                 / ( LOG( z_mo / z0(j,i) )     &
559                                                    - psi_m( z_mo / ol_m )     &
560                                                    + psi_m( z0(j,i) / ol_m )  &
561                                                    )**2
562
563!
564!--                    Calculate df/dL
565                       f_d_ol = ( - ( z_mo / ol_u ) * ( LOG( z_mo / z0h(j,i) ) &
566                                                   - psi_h( z_mo / ol_u )      &
567                                                   + psi_h( z0h(j,i) / ol_u )  &
568                                                 )                             &
569                                               / ( LOG( z_mo / z0(j,i) )       &
570                                                   - psi_m( z_mo / ol_u )      &
571                                                   + psi_m( z0(j,i) / ol_u )   &
572                                                 )**2                          &
573                              + ( z_mo / ol_l ) * ( LOG( z_mo / z0h(j,i) )     &
574                                                   - psi_h( z_mo / ol_l )      &
575                                                   + psi_h( z0h(j,i) / ol_l )  &
576                                                 )                             &
577                                               / ( LOG( z_mo / z0(j,i) )       &
578                                                   - psi_m( z_mo / ol_l )      &
579                                                   + psi_m( z0(j,i) / ol_l )   &
580                                                 )**2                          &
581                                ) / ( ol_u - ol_l )
582                   ELSE
583!
584!--                   Calculate f = Ri - 1 /[...]^3 = 0
585                      f = rib(j,i) - ( z_mo / ol_m ) / ( LOG( z_mo / z0(j,i) )&
586                                                    - psi_m( z_mo / ol_m )    &
587                                                    + psi_m( z0(j,i) / ol_m ) &
588                                                       )**3
589
590!
591!--                   Calculate df/dL
592                      f_d_ol = ( - ( z_mo / ol_u ) / ( LOG( z_mo / z0(j,i) )  &
593                                                   - psi_m( z_mo / ol_u )     &
594                                                   + psi_m( z0(j,i) / ol_u )  &
595                                                 )**3                         &
596                              + ( z_mo / ol_l ) / ( LOG( z_mo / z0(j,i) )     &
597                                                   - psi_m( z_mo / ol_l )     &
598                                                   + psi_m( z0(j,i) / ol_l )  &
599                                                 )**3                         &
600                                     ) / ( ol_u - ol_l )
601                   ENDIF
602!
603!--                Calculate new L
604                   ol(j,i) = ol_m - f / f_d_ol
605
606!
607!--                Ensure that the bulk Richardson number and the Obukhov
608!--                length have the same sign and ensure convergence.
609                   IF ( ol(j,i) * ol_m < 0.0_wp )  ol(j,i) = ol_m * 0.5_wp
610
611!
612!--                If unrealistic value occurs, set L to the maximum
613!--                value that is allowed
614                   IF ( ABS( ol(j,i) ) > ol_max )  THEN
615                      ol(j,i) = ol_max
616                      EXIT
617                   ENDIF
618!
619!--                Check for convergence
620                   IF ( ABS( ( ol(j,i) - ol_m ) / ol(j,i) ) < 1.0E-4_wp )  THEN
621                      EXIT
622                   ELSE
623                      CYCLE
624                   ENDIF
625
626                ENDDO
627                       
628             ENDDO
629          ENDDO
630
631       ELSEIF ( TRIM( most_method ) == 'lookup' )  THEN
632
633          !$OMP PARALLEL DO PRIVATE( k, z_mo )
634          !# WARNING: does not work on GPU so far
635          !!!!!!$acc kernels loop
636          DO  i = nxl, nxr
637             DO  j = nys, nyn
638
639                k   = nzb_s_inner(j,i)
640
641!
642!--             If the bulk Richardson number is outside the range of the lookup
643!--             table, set it to the exceeding threshold value
644                IF ( rib(j,i) < rib_min )  rib(j,i) = rib_min
645                IF ( rib(j,i) > rib_max )  rib(j,i) = rib_max
646
647!
648!--             Find the correct index bounds for linear interpolation. As the
649!--             Richardson number will not differ very much from time step to
650!--             time step , use the index from the last step and search in the
651!--             correct direction
652                l = l_bnd
653                IF ( rib_tab(l) - rib(j,i) > 0.0_wp )  THEN
654                   DO WHILE ( rib_tab(l-1) - rib(j,i) > 0.0_wp .AND. l > 0 )
655                      l = l-1
656                   ENDDO
657                ELSE
658                   DO WHILE ( rib_tab(l) - rib(j,i) < 0.0_wp                &
659                              .AND. l < num_steps-1 )
660                      l = l+1
661                   ENDDO
662                ENDIF
663                l_bnd = l
664
665!
666!--             Linear interpolation to find the correct value of z/L
667                ol(j,i) = ( ol_tab(l-1) + ( ol_tab(l) - ol_tab(l-1) )       &
668                            / (  rib_tab(l) - rib_tab(l-1) )                &
669                            * ( rib(j,i) - rib_tab(l-1) ) )
670
671             ENDDO
672          ENDDO
673
674       ELSEIF ( TRIM( most_method ) == 'circular' )  THEN
675
676          !$OMP PARALLEL DO PRIVATE( k, z_mo )
677          !$acc kernels loop present( nzb_s_inner, ol, pt, pt1, q, ql, qs, qv1, ts, us, vpt, zu, zw ) private( j, k, z_mo )
678          DO  i = nxl, nxr
679             DO  j = nys, nyn
680
681                k   = nzb_s_inner(j,i)
682                z_mo = zu(k+1) - zw(k)
683
684                IF ( .NOT. humidity )  THEN
685                   ol(j,i) =  ( pt(k+1,j,i) *  us(j,i)**2 ) / ( kappa * g      &
686                              * ts(j,i) + 1E-30_wp )
687                ELSEIF ( cloud_physics )  THEN
688
689                   ol(j,i) =  ( vpt(k+1,j,i) * us(j,i)**2 ) / ( kappa * g      &
690                              * ( ts(j,i) + 0.61_wp * pt1(j,i) * qs(j,i)       &
691                              + 0.61_wp * qv1(j,i) * ts(j,i) - ts(j,i)         &
692                              * ql(k+1,j,i) ) + 1E-30_wp )
693                ELSE
694                   ol(j,i) =  ( vpt(k+1,j,i) *  us(j,i)**2 ) / ( kappa * g     &
695                              * ( ts(j,i) + 0.61_wp * pt(k+1,j,i) * qs(j,i)    &
696                                  + 0.61_wp * q(k+1,j,i) * ts(j,i) ) + 1E-30_wp )
697                ENDIF
698!
699!--             Limit the value range of the Obukhov length.
700!--             This is necessary for very small velocities (u,v --> 0), because
701!--             the absolute value of ol can then become very small, which in
702!--             consequence would result in very large shear stresses and very
703!--             small momentum fluxes (both are generally unrealistic).
704                IF ( ( z_mo / ol(j,i) ) < zeta_min )  ol(j,i) = z_mo / zeta_min
705                IF ( ( z_mo / ol(j,i) ) > zeta_max )  ol(j,i) = z_mo / zeta_max
706
707             ENDDO
708          ENDDO
709
710       ENDIF
711
712!
713!--    Values of ol at ghost point locations are needed for the evaluation
714!--    of usws and vsws.
715       !$acc update host( ol )
716       CALL exchange_horiz_2d( ol )
717       !$acc update device( ol )
718
719    END SUBROUTINE calc_ol
720
721!
722!-- Calculate friction velocity u*
723    SUBROUTINE calc_us
724
725       IMPLICIT NONE
726
727       !$OMP PARALLEL DO PRIVATE( k, z_mo )
728       ! ### WARNING: does not work on GPU because of function calls
729       !!!!!$acc kernels loop
730       !$acc update host( ol, uv_total )
731       DO  i = nxlg, nxrg
732          DO  j = nysg, nyng
733
734             k   = nzb_s_inner(j,i)+1
735             z_mo = zu(k+1) - zw(k)
736
737!
738!--          Compute u* at the scalars' grid points
739             us(j,i) = kappa * uv_total(j,i) / ( LOG( z_mo / z0(j,i) )         &
740                                          - psi_m( z_mo / ol(j,i) )            &
741                                          + psi_m( z0(j,i) / ol(j,i) ) )
742          ENDDO
743       ENDDO
744       !$acc update device( us )
745
746    END SUBROUTINE calc_us
747
748!
749!-- Calculate potential temperature and specific humidity at first grid level
750    SUBROUTINE calc_pt_q
751
752       IMPLICIT NONE
753
754       !$acc kernels loop present( nzb_s_inner, pt, pt1, pt_d_t, q, ql, qv1 ) private( j, k )
755       DO  i = nxlg, nxrg
756          DO  j = nysg, nyng
757             k   = nzb_s_inner(j,i)+1
758             pt1(j,i) = pt(k,j,i) + l_d_cp * pt_d_t(k) * ql(k,j,i)
759             qv1(j,i) = q(k,j,i) - ql(k,j,i)
760          ENDDO
761       ENDDO
762
763    END SUBROUTINE calc_pt_q
764
765!
766!-- Calculate the other MOST scaling parameters theta*, q*, (qr*, nr*)
767    SUBROUTINE calc_scaling_parameters
768
769       IMPLICIT NONE
770
771!
772!--    Data information for accelerators
773       !$acc data present( e, nrsws, nzb_u_inner, nzb_v_inner, nzb_s_inner, pt )  &
774       !$acc      present( q, qs, qsws, qrsws, shf, ts, u, us, usws, v )     &
775       !$acc      present( vpt, vsws, zu, zw, z0, z0h )
776!
777!--    Compute theta*
778       IF ( constant_heatflux )  THEN
779
780!
781!--       For a given heat flux in the surface layer:
782          !$OMP PARALLEL DO
783          !$acc kernels loop private( j )
784          DO  i = nxlg, nxrg
785             DO  j = nysg, nyng
786                ts(j,i) = -shf(j,i) / ( us(j,i) + 1E-30_wp )
787!
788!--             ts must be limited, because otherwise overflow may occur in case
789!--             of us=0 when computing ol further below
790                IF ( ts(j,i) < -1.05E5_wp )  ts(j,i) = -1.0E5_wp
791                IF ( ts(j,i) >   1.0E5_wp )  ts(j,i) =  1.0E5_wp
792             ENDDO
793          ENDDO
794
795       ELSE
796!
797!--       For a given surface temperature:
798          IF ( large_scale_forcing .AND. lsf_surf )  THEN
799             !$OMP PARALLEL DO
800             !$acc kernels loop private( j, k )
801             DO  i = nxlg, nxrg
802                DO  j = nysg, nyng
803                   k = nzb_s_inner(j,i)
804                   pt(k,j,i) = pt_surface
805                ENDDO
806             ENDDO
807          ENDIF
808
809          !$OMP PARALLEL DO PRIVATE( k, z_mo )
810          ! ### WARNING: leads to severe compile error
811          !!!!$acc kernels loop
812          DO  i = nxlg, nxrg
813             DO  j = nysg, nyng
814
815                k   = nzb_s_inner(j,i)
816                z_mo = zu(k+1) - zw(k)
817
818                IF ( cloud_physics )  THEN
819                   ts(j,i) = kappa * ( pt1(j,i) - pt(k,j,i) )                  &
820                                     / ( LOG( z_mo / z0h(j,i) )                &
821                                         - psi_h( z_mo / ol(j,i) )             &
822                                         + psi_h( z0h(j,i) / ol(j,i) ) )
823                ELSE
824                   ts(j,i) = kappa * ( pt(k+1,j,i) - pt(k,j,i) )               &
825                                     / ( LOG( z_mo / z0h(j,i) )                &
826                                         - psi_h( z_mo / ol(j,i) )             &
827                                         + psi_h( z0h(j,i) / ol(j,i) ) )
828                ENDIF
829
830             ENDDO
831          ENDDO
832       ENDIF
833
834!
835!--    If required compute q*
836       IF ( humidity  .OR.  passive_scalar )  THEN
837          IF ( constant_waterflux )  THEN
838!
839!--          For a given water flux in the Prandtl layer:
840             !$OMP PARALLEL DO
841             !$acc kernels loop private( j )
842             DO  i = nxlg, nxrg
843                DO  j = nysg, nyng
844                   qs(j,i) = -qsws(j,i) / ( us(j,i) + 1E-30_wp )
845                ENDDO
846             ENDDO
847
848          ELSE
849             coupled_run = ( coupling_mode == 'atmosphere_to_ocean' .AND.      &
850                             run_coupled )
851
852             IF ( large_scale_forcing .AND. lsf_surf )  THEN
853                !$OMP PARALLEL DO
854                !$acc kernels loop private( j, k )
855                DO  i = nxlg, nxrg
856                   DO  j = nysg, nyng
857                      k = nzb_s_inner(j,i)
858                      q(k,j,i) = q_surface
859                   ENDDO
860                ENDDO
861             ENDIF
862
863             !$OMP PARALLEL DO PRIVATE( e_s, k, z_mo )
864             ! ### WARNING: leads to severe compile error
865             !!!!$acc kernels loop independent
866             DO  i = nxlg, nxrg
867                !!!!$acc loop independent
868                DO  j = nysg, nyng
869
870                   k   = nzb_s_inner(j,i)
871                   z_mo = zu(k+1) - zw(k)
872
873!
874!--                Assume saturation for atmosphere coupled to ocean (but not
875!--                in case of precursor runs)
876                   IF ( coupled_run )  THEN
877                      e_s = 6.1_wp * &
878                              EXP( 0.07_wp * ( MIN(pt(k,j,i),pt(k+1,j,i))      &
879                                               - 273.15_wp ) )
880                      q(k,j,i) = 0.622_wp * e_s / ( surface_pressure - e_s )
881                   ENDIF
882
883                   IF ( cloud_physics )  THEN
884                      qs(j,i) = kappa * ( qv1(j,i) - q(k,j,i) )                &
885                                        / ( LOG( z_mo / z0h(j,i) )             &
886                                            - psi_h( z_mo / ol(j,i) )          &
887                                            + psi_h( z0h(j,i) / ol(j,i) ) )
888
889                   ELSE
890                      qs(j,i) = kappa * ( q(k+1,j,i) - q(k,j,i) )              &
891                                        / ( LOG( z_mo / z0h(j,i) )             &
892                                            - psi_h( z_mo / ol(j,i) )          &
893                                            + psi_h( z0h(j,i) / ol(j,i) ) )
894                   ENDIF
895
896                ENDDO
897             ENDDO
898          ENDIF
899       ENDIF
900
901
902!
903!--    If required compute qr* and nr*
904       IF ( cloud_physics .AND. icloud_scheme == 0 .AND. precipitation )  THEN
905
906          !$OMP PARALLEL DO PRIVATE( k, z_mo )
907          ! ### WARNING: leads to severe compile error
908          !!!!!$acc kernels loop independent
909          DO  i = nxlg, nxrg
910             !!!!$acc loop independent
911             DO  j = nysg, nyng
912
913                k   = nzb_s_inner(j,i)
914                z_mo = zu(k+1) - zw(k)
915
916                   qrs(j,i) = kappa * ( qr(k+1,j,i) - qr(k,j,i) )              &
917                                     / ( LOG( z_mo / z0h(j,i) )                &
918                                         - psi_h( z_mo / ol(j,i) )             &
919                                         + psi_h( z0h(j,i) / ol(j,i) ) )
920
921                   nrs(j,i) = kappa * ( nr(k+1,j,i) - nr(k,j,i) )              &
922                                     / ( LOG( z_mo / z0h(j,i) )                &
923                                         - psi_h( z_mo / ol(j,i) )             &
924                                         + psi_h( z0h(j,i) / ol(j,i) ) )
925             ENDDO
926          ENDDO
927
928       ENDIF
929       !$acc end data
930
931    END SUBROUTINE calc_scaling_parameters
932
933
934
935!
936!-- Calculate surface fluxes usws, vsws, shf, qsws, (qrsws, nrsws)
937    SUBROUTINE calc_surface_fluxes
938
939       IMPLICIT NONE
940
941       REAL(wp) :: ol_mid !< Grid-interpolated L
942
943!
944!--    Compute u'w' for the total model domain.
945!--    First compute the corresponding component of u* and square it.
946       !$OMP PARALLEL DO PRIVATE( k, ol_mid, z_mo )
947       ! ### WARNING: leads to severe compile error
948       !!!!!$acc kernels loop
949       !$acc update host( ol, u, us, v )
950       DO  i = nxl, nxr
951          DO  j = nys, nyn
952
953             k   = nzb_u_inner(j,i)
954             z_mo = zu(k+1) - zw(k)
955!
956!--          Compute bulk Obukhov length for this point
957             ol_mid = 0.5_wp * ( ol(j,i-1) + ol(j,i) )
958
959             IF ( ol_mid == 0.0_wp )  THEN
960                ol_mid = MIN(ol(j,i-1), ol(j,i))
961             ENDIF
962
963             usws(j,i) = kappa * ( u(k+1,j,i) - u(k,j,i) )                     &
964                                 / ( LOG( z_mo / z0(j,i) )                     &
965                                     - psi_m( z_mo / ol_mid )                  &
966                                     + psi_m( z0(j,i) / ol_mid ) )
967
968             usws(j,i) = -usws(j,i) * 0.5_wp * ( us(j,i-1) + us(j,i) )
969          ENDDO
970       ENDDO
971
972!
973!--    Compute v'w' for the total model domain.
974!--    First compute the corresponding component of u* and square it.
975       !$OMP PARALLEL DO PRIVATE( k, ol_mid, z_mo )
976       ! ### WARNING: leads to severe compile error
977       !!!!!!$acc kernels loop
978       DO  i = nxl, nxr
979          DO  j = nys, nyn
980
981             k   = nzb_v_inner(j,i)
982             z_mo = zu(k+1) - zw(k)
983!
984!--          Compute bulk Obukhov length for this point
985             ol_mid = 0.5_wp * ( ol(j-1,i) + ol(j,i) )
986
987             IF ( ol_mid == 0.0_wp )  THEN
988                ol_mid = MIN(ol(j-1,i), ol(j-1,i))
989             ENDIF
990
991             vsws(j,i) = kappa * ( v(k+1,j,i) - v(k,j,i) )                     &
992                                 / ( LOG( z_mo / z0(j,i) )                     &
993                                     - psi_m( z_mo / ol_mid )                  &
994                                     + psi_m( z0(j,i) / ol_mid ) )
995
996             vsws(j,i) = -vsws(j,i) * 0.5_wp * ( us(j,i-1) + us(j,i) )
997
998          ENDDO
999       ENDDO
1000
1001!
1002!--    Exchange the boundaries for the momentum fluxes (is this still required?)
1003       !!!!!!!!$acc update host( usws, vsws )
1004       CALL exchange_horiz_2d( usws )
1005       CALL exchange_horiz_2d( vsws )
1006       !$acc update device( usws, vsws )
1007
1008!
1009!--    Compute the vertical kinematic heat flux
1010       IF ( .NOT. constant_heatflux .AND. ( simulated_time <=                &
1011            skip_time_do_lsm .OR. .NOT. land_surface ) )  THEN
1012          !$OMP PARALLEL DO
1013          !$acc kernels loop independent present( shf, ts, us )
1014          DO  i = nxlg, nxrg
1015             !$acc loop independent
1016             DO  j = nysg, nyng
1017                shf(j,i) = -ts(j,i) * us(j,i)
1018             ENDDO
1019          ENDDO
1020
1021       ENDIF
1022
1023!
1024!-- Compute the vertical water/scalar flux
1025       IF ( .NOT. constant_waterflux .AND. ( humidity .OR. passive_scalar )    &
1026            .AND. ( simulated_time <= skip_time_do_lsm .OR. .NOT.            &
1027            land_surface ) )  THEN
1028          !$OMP PARALLEL DO
1029          !$acc kernels loop independent present( qs, qsws, us )
1030          DO  i = nxlg, nxrg
1031             !$acc loop independent
1032             DO  j = nysg, nyng
1033                qsws(j,i) = -qs(j,i) * us(j,i)
1034             ENDDO
1035          ENDDO
1036
1037       ENDIF
1038
1039!
1040!--    Compute (turbulent) fluxes of rain water content and rain drop conc.
1041       IF ( cloud_physics  .AND.  icloud_scheme == 0  .AND.                    &
1042            precipitation )  THEN
1043          !$OMP PARALLEL DO
1044          !$acc kernels loop independent present( nrs, nrsws, qrs, qrsws, us )
1045          DO  i = nxlg, nxrg
1046             !$acc loop independent
1047             DO  j = nysg, nyng
1048                qrsws(j,i) = -qrs(j,i) * us(j,i)
1049                nrsws(j,i) = -nrs(j,i) * us(j,i)
1050             ENDDO
1051          ENDDO
1052       ENDIF
1053
1054!
1055!--    Bottom boundary condition for the TKE
1056       IF ( ibc_e_b == 2 )  THEN
1057          !$OMP PARALLEL DO
1058          !$acc kernels loop independent present( e, nzb_s_inner, us )
1059          DO  i = nxlg, nxrg
1060             !$acc loop independent
1061             DO  j = nysg, nyng
1062                e(nzb_s_inner(j,i)+1,j,i) = ( us(j,i) / 0.1_wp )**2
1063!
1064!--             As a test: cm = 0.4
1065!               e(nzb_s_inner(j,i)+1,j,i) = ( us(j,i) / 0.4_wp )**2
1066                e(nzb_s_inner(j,i),j,i)   = e(nzb_s_inner(j,i)+1,j,i)
1067             ENDDO
1068          ENDDO
1069       ENDIF
1070
1071    END SUBROUTINE calc_surface_fluxes
1072
1073
1074!
1075!-- Integrated stability function for momentum
1076    FUNCTION psi_m( zeta ) 
1077       
1078       USE kinds
1079
1080       IMPLICIT NONE
1081
1082       REAL(wp)            :: psi_m !< Integrated similarity function result
1083       REAL(wp)            :: zeta  !< Stability parameter z/L
1084       REAL(wp)            :: x     !< dummy variable
1085
1086       REAL(wp), PARAMETER :: a = 1.0_wp            !< constant
1087       REAL(wp), PARAMETER :: b = 0.66666666666_wp  !< constant
1088       REAL(wp), PARAMETER :: c = 5.0_wp            !< constant
1089       REAL(wp), PARAMETER :: d = 0.35_wp           !< constant
1090       REAL(wp), PARAMETER :: c_d_d = c / d         !< constant
1091       REAL(wp), PARAMETER :: bc_d_d = b * c / d    !< constant
1092
1093
1094       IF ( zeta < 0.0_wp )  THEN
1095          x = SQRT( SQRT(1.0_wp  - 16.0_wp * zeta ) )
1096          psi_m = pi * 0.5_wp - 2.0_wp * ATAN( x ) + LOG( ( 1.0_wp + x )**2    &
1097                  * ( 1.0_wp + x**2 ) * 0.125_wp )
1098       ELSE
1099
1100          psi_m = - b * ( zeta - c_d_d ) * EXP( -d * zeta ) - a * zeta         &
1101                   - bc_d_d
1102!
1103!--       Old version for stable conditions (only valid for z/L < 0.5)
1104!--       psi_m = - 5.0_wp * zeta
1105
1106       ENDIF
1107
1108    END FUNCTION psi_m
1109
1110
1111!
1112!-- Integrated stability function for heat and moisture
1113    FUNCTION psi_h( zeta ) 
1114       
1115       USE kinds
1116
1117       IMPLICIT NONE
1118
1119       REAL(wp)            :: psi_h !< Integrated similarity function result
1120       REAL(wp)            :: zeta  !< Stability parameter z/L
1121       REAL(wp)            :: x     !< dummy variable
1122
1123       REAL(wp), PARAMETER :: a = 1.0_wp            !< constant
1124       REAL(wp), PARAMETER :: b = 0.66666666666_wp  !< constant
1125       REAL(wp), PARAMETER :: c = 5.0_wp            !< constant
1126       REAL(wp), PARAMETER :: d = 0.35_wp           !< constant
1127       REAL(wp), PARAMETER :: c_d_d = c / d         !< constant
1128       REAL(wp), PARAMETER :: bc_d_d = b * c / d    !< constant
1129
1130
1131       IF ( zeta < 0.0_wp )  THEN
1132          x = SQRT(1.0_wp  - 16.0_wp * zeta )
1133          psi_h = 2.0_wp * LOG( (1.0_wp + x ) / 2.0_wp )
1134       ELSE
1135          psi_h = - b * ( zeta - c_d_d ) * EXP( -d * zeta ) - (1.0_wp          &
1136                  + 0.66666666666_wp * a * zeta )**1.5_wp - bc_d_d             &
1137                  + 1.0_wp
1138!
1139!--       Old version for stable conditions (only valid for z/L < 0.5)
1140!--       psi_h = - 5.0_wp * zeta
1141       ENDIF
1142
1143    END FUNCTION psi_h
1144
1145 END MODULE surface_layer_fluxes_mod
Note: See TracBrowser for help on using the repository browser.