source: palm/trunk/SOURCE/nudging_mod.f90 @ 2256

Last change on this file since 2256 was 2233, checked in by suehring, 8 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 21.5 KB
RevLine 
[1850]1!> @file nudging_mod.f90
[2000]2!------------------------------------------------------------------------------!
[1239]3! This file is part of PALM.
4!
[2000]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.
[1239]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!
[2101]17! Copyright 1997-2017 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1239]19!
20! Current revisions:
21! ------------------
[1383]22!
[2233]23!
[1383]24! Former revisions:
25! -----------------
26! $Id: nudging_mod.f90 2233 2017-05-30 18:08:54Z suehring $
27!
[2233]28! 2232 2017-05-30 17:47:52Z suehring
29! Adjustments to new topography concept
30!
[2001]31! 2000 2016-08-20 18:09:15Z knoop
32! Forced header and separation lines into 80 columns
33!
[1851]34! 1850 2016-04-08 13:29:27Z maronga
35! Module renamed
36!
37!
[1758]38! 1757 2016-02-22 15:49:32Z maronga
39! Bugfix: allow for using higher vertical resolution in nudging file than grid
40! spacing in the LES model
41!
[1683]42! 1682 2015-10-07 23:56:08Z knoop
43! Code annotations made doxygen readable
44!
[1399]45! 1398 2014-05-07 11:15:00Z heinze
46! Subroutine nudge_ref is extended to set u_init and v_init to the current
47! nudging profiles
48!
[1383]49! 1382 2014-04-30 12:15:41Z boeske
[1382]50! Changed the weighting factor that is used in the summation of nudging
51! tendencies for profile data output from weight_pres to weight_substep,
52! added Neumann boundary conditions for profile data output of nudging terms at
53! nzt+1
[1366]54!
[1381]55! 1380 2014-04-28 12:40:45Z heinze
56! Subroutine nudge_ref added to account for proper upper scalar boundary
57! conditions in case of nudging
58!
[1366]59! 1365 2014-04-22 15:03:56Z boeske
[1365]60! Variable t renamed nt, variable currtnudge renamed tmp_tnudge,
61! summation of nudging tendencies for data output added
62! +sums_ls_l, tmp_tend
63! Added new subroutine calc_tnudge, which calculates the current nudging time
64! scale at each time step
[1354]65!
[1356]66! 1355 2014-04-10 10:21:29Z heinze
67! Error message specified.
68!
[1354]69! 1353 2014-04-08 15:21:23Z heinze
70! REAL constants provided with KIND-attribute
71!
[1321]72! 1320 2014-03-20 08:40:49Z raasch
[1320]73! ONLY-attribute added to USE-statements,
74! kind-parameters added to all INTEGER and REAL declaration statements,
75! kinds are defined in new module kinds,
76! old module precision_kind is removed,
77! revision history before 2012 removed,
78! comment fields (!:) to be used for variable explanations added to
79! all variable declaration statements
[1239]80!
[1319]81! 1318 2014-03-17 13:35:16Z raasch
82! module interfaces removed
83!
[1269]84! 1268 2013-12-12 09:47:53Z heinze
85! bugfix: argument of calc_mean_profile corrected
86!
[1252]87! 1251 2013-11-07 08:14:30Z heinze
88! bugfix: calculate dtm and dtp also in vector version
89!
[1250]90! 1249 2013-11-06 10:45:47Z heinze
91! remove call of user module
92! reformatting
93!
[1242]94! 1241 2013-10-30 11:36:58Z heinze
95! Initial revision
[1239]96!
97! Description:
98! ------------
[1682]99!> Nudges u, v, pt and q to given profiles on a relaxation timescale tnudge.
100!> Profiles are read in from NUDGIN_DATA. Code is based on Neggers et al. (2012)
101!> and also part of DALES and UCLA-LES.
[1239]102!--------------------------------------------------------------------------------!
[1682]103 MODULE nudge_mod
104 
[1239]105
106    PRIVATE
[1380]107    PUBLIC init_nudge, calc_tnudge, nudge, nudge_ref
[1239]108    SAVE
109
110    INTERFACE nudge
111       MODULE PROCEDURE nudge
112       MODULE PROCEDURE nudge_ij
113    END INTERFACE nudge
114
115 CONTAINS
116
[1682]117!------------------------------------------------------------------------------!
118! Description:
119! ------------
120!> @todo Missing subroutine description.
121!------------------------------------------------------------------------------!
[1239]122    SUBROUTINE init_nudge
123
[1320]124       USE arrays_3d,                                                          &
[1365]125           ONLY:  ptnudge, qnudge, timenudge, tmp_tnudge, tnudge, unudge,      &
126                  vnudge, wnudge, zu
[1239]127
[1320]128       USE control_parameters,                                                 &
129           ONLY:  dt_3d, lptnudge, lqnudge, lunudge, lvnudge, lwnudge,         &
130                   message_string, ntnudge
131
132       USE indices,                                                            &
133           ONLY:  nzb, nzt
134
135       USE kinds
136
[1239]137       IMPLICIT NONE
138
139
[1682]140       INTEGER(iwp) ::  finput = 90  !<
141       INTEGER(iwp) ::  ierrn        !<
142       INTEGER(iwp) ::  k            !<
143       INTEGER(iwp) ::  nt            !<
[1239]144
[1682]145       CHARACTER(1) ::  hash     !<
[1320]146
[1682]147       REAL(wp) ::  highheight   !<
148       REAL(wp) ::  highqnudge   !<
149       REAL(wp) ::  highptnudge  !<
150       REAL(wp) ::  highunudge   !<
151       REAL(wp) ::  highvnudge   !<
152       REAL(wp) ::  highwnudge   !<
153       REAL(wp) ::  hightnudge   !<
[1320]154
[1682]155       REAL(wp) ::  lowheight    !<
156       REAL(wp) ::  lowqnudge    !<
157       REAL(wp) ::  lowptnudge   !<
158       REAL(wp) ::  lowunudge    !<
159       REAL(wp) ::  lowvnudge    !<
160       REAL(wp) ::  lowwnudge    !<
161       REAL(wp) ::  lowtnudge    !<
[1320]162
[1682]163       REAL(wp) ::  fac          !<
[1320]164
[1239]165       ALLOCATE( ptnudge(nzb:nzt+1,1:ntnudge), qnudge(nzb:nzt+1,1:ntnudge), &
166                 tnudge(nzb:nzt+1,1:ntnudge), unudge(nzb:nzt+1,1:ntnudge),  &
167                 vnudge(nzb:nzt+1,1:ntnudge), wnudge(nzb:nzt+1,1:ntnudge)  )
168
[1365]169       ALLOCATE( tmp_tnudge(nzb:nzt) )
170
[1239]171       ALLOCATE( timenudge(0:ntnudge) )
172
[1353]173       ptnudge = 0.0_wp; qnudge = 0.0_wp; tnudge = 0.0_wp; unudge = 0.0_wp
174       vnudge = 0.0_wp; wnudge = 0.0_wp; timenudge = 0.0_wp
[1365]175!
176!--    Initialize array tmp_nudge with a current nudging time scale of 6 hours
177       tmp_tnudge = 21600.0_wp
[1239]178
[1365]179       nt = 0
[1249]180       OPEN ( finput, FILE='NUDGING_DATA', STATUS='OLD', &
181              FORM='FORMATTED', IOSTAT=ierrn )
[1239]182
[1249]183       IF ( ierrn /= 0 )  THEN
[1239]184          message_string = 'file NUDGING_DATA does not exist'
185          CALL message( 'nudging', 'PA0365', 1, 2, 0, 6, 0 )
186       ENDIF
187
188       ierrn = 0
189
190 rloop:DO
[1365]191          nt = nt + 1
[1239]192          hash = "#"
[1320]193          ierrn = 1 ! not zero
[1239]194!
195!--       Search for the next line consisting of "# time",
196!--       from there onwards the profiles will be read
197          DO WHILE ( .NOT. ( hash == "#" .AND. ierrn == 0 ) ) 
198         
[1365]199            READ ( finput, *, IOSTAT=ierrn ) hash, timenudge(nt)
[1249]200            IF ( ierrn < 0 )  EXIT rloop
[1239]201
202          ENDDO
203
204          ierrn = 0
[1249]205          READ ( finput, *, IOSTAT=ierrn ) lowheight, lowtnudge, lowunudge,   &
206                                           lowvnudge, lowwnudge , lowptnudge, &
207                                           lowqnudge
[1239]208
[1249]209          IF ( ierrn /= 0 )  THEN
[1239]210             message_string = 'errors in file NUDGING_DATA'
211             CALL message( 'nudging', 'PA0366', 1, 2, 0, 6, 0 )
212          ENDIF
213
214          ierrn = 0
[1249]215          READ ( finput, *, IOSTAT=ierrn ) highheight, hightnudge, highunudge,   &
216                                           highvnudge, highwnudge , highptnudge, &
217                                           highqnudge
[1239]218
[1249]219          IF ( ierrn /= 0 )  THEN
[1239]220             message_string = 'errors in file NUDGING_DATA'
221             CALL message( 'nudging', 'PA0366', 1, 2, 0, 6, 0 )
222          ENDIF
223
224          DO  k = nzb, nzt+1
[1757]225             DO WHILE ( highheight < zu(k) )
[1239]226                lowheight  = highheight
227                lowtnudge  = hightnudge
228                lowunudge  = highunudge
229                lowvnudge  = highvnudge
230                lowwnudge  = highwnudge
231                lowptnudge = highptnudge
232                lowqnudge  = highqnudge
233 
234                ierrn = 0
[1249]235                READ ( finput, *, IOSTAT=ierrn )  highheight , hightnudge , &
236                                                  highunudge , highvnudge , &
237                                                  highwnudge , highptnudge, &
238                                                  highqnudge
239                IF (ierrn /= 0 )  THEN
[1355]240                   WRITE( message_string, * ) 'zu(nzt+1) = ', zu(nzt+1), 'm is ',&
241                        'higher than the maximum height in NUDING_DATA which ',  &
242                        'is ', lowheight, 'm. Interpolation on PALM ',           &
243                        'grid is not possible.'
244                   CALL message( 'nudging', 'PA0364', 1, 2, 0, 6, 0 )
[1239]245                ENDIF
[1757]246             ENDDO
[1239]247
248!
249!--          Interpolation of prescribed profiles in space
250
[1249]251             fac = ( highheight - zu(k) ) / ( highheight - lowheight )
[1239]252
[1365]253             tnudge(k,nt)  = fac * lowtnudge  + ( 1.0_wp - fac ) * hightnudge
254             unudge(k,nt)  = fac * lowunudge  + ( 1.0_wp - fac ) * highunudge
255             vnudge(k,nt)  = fac * lowvnudge  + ( 1.0_wp - fac ) * highvnudge
256             wnudge(k,nt)  = fac * lowwnudge  + ( 1.0_wp - fac ) * highwnudge
257             ptnudge(k,nt) = fac * lowptnudge + ( 1.0_wp - fac ) * highptnudge
258             qnudge(k,nt)  = fac * lowqnudge  + ( 1.0_wp - fac ) * highqnudge
[1239]259          ENDDO
260
261       ENDDO rloop
262
[1249]263       CLOSE ( finput )
[1239]264
265!
266!--    Prevent nudging if nudging profiles exhibt too small values
[1241]267!--    not used so far
[1353]268       lptnudge  = ANY( ABS( ptnudge ) > 1.0e-8_wp )
269       lqnudge   = ANY( ABS( qnudge )  > 1.0e-8_wp )
270       lunudge   = ANY( ABS( unudge )  > 1.0e-8_wp )
271       lvnudge   = ANY( ABS( vnudge )  > 1.0e-8_wp )
272       lwnudge   = ANY( ABS( wnudge )  > 1.0e-8_wp )
[1239]273
274    END SUBROUTINE init_nudge
275
[1365]276
[1682]277!------------------------------------------------------------------------------!
278! Description:
279! ------------
280!> @todo Missing subroutine description.
281!------------------------------------------------------------------------------!
[1365]282    SUBROUTINE calc_tnudge ( time )
283
284       USE arrays_3d,                                                          &
285           ONLY:  timenudge, tmp_tnudge, tnudge
286
287       USE control_parameters,                                                 &
288           ONLY:  dt_3d 
289
290       USE indices,                                                            &
291           ONLY:  nzb, nzt
292
293       USE kinds
294
295       IMPLICIT NONE
296
297
[1682]298       REAL(wp) ::  dtm         !<
299       REAL(wp) ::  dtp         !<
300       REAL(wp) ::  time        !<
[1365]301
[1682]302       INTEGER(iwp) ::  k   !<
303       INTEGER(iwp) ::  nt  !<
[1365]304
305       nt = 1
306       DO WHILE ( time > timenudge(nt) )
307         nt = nt+1
308       ENDDO
309       IF ( time /= timenudge(1) ) THEN
310         nt = nt-1
311       ENDIF
312
313       dtm = ( time - timenudge(nt) ) / ( timenudge(nt+1) - timenudge(nt) )
314       dtp = ( timenudge(nt+1) - time ) / ( timenudge(nt+1) - timenudge(nt) )
315
316       DO  k = nzb, nzt
317          tmp_tnudge(k) = MAX( dt_3d, tnudge(k,nt) * dtp + tnudge(k,nt+1) * dtm )
318       ENDDO
319
320    END SUBROUTINE calc_tnudge
321
[1239]322!------------------------------------------------------------------------------!
[1682]323! Description:
324! ------------
325!> Call for all grid points
[1239]326!------------------------------------------------------------------------------!
327    SUBROUTINE nudge ( time, prog_var )
328
[1320]329       USE arrays_3d,                                                          &
[1365]330           ONLY:  pt, ptnudge, q, qnudge, tend, timenudge, tmp_tnudge, tnudge, &
331                  u, unudge, v, vnudge
[1239]332
[1320]333       USE control_parameters,                                                 &
[1365]334           ONLY:  dt_3d, intermediate_timestep_count, message_string
[1320]335
336       USE indices,                                                            &
[2232]337           ONLY:  nxl, nxr, nys, nyn, nzb, nzt, wall_flags_0
[1320]338
[1365]339       USE kinds
[1320]340
341       USE statistics,                                                         &
[1382]342           ONLY:  hom, sums_ls_l, weight_substep
[1320]343
[1239]344       IMPLICIT NONE
345
[1682]346       CHARACTER (LEN=*) ::  prog_var  !<
[1239]347
[1682]348       REAL(wp) ::  tmp_tend    !<
349       REAL(wp) ::  dtm         !<
350       REAL(wp) ::  dtp         !<
351       REAL(wp) ::  time        !<
[1239]352
[1682]353       INTEGER(iwp) ::  i  !<
354       INTEGER(iwp) ::  j  !<
355       INTEGER(iwp) ::  k  !<
356       INTEGER(iwp) ::  nt  !<
[1239]357
358
[1365]359       nt = 1
360       DO WHILE ( time > timenudge(nt) )
361         nt = nt+1
[1251]362       ENDDO
363       IF ( time /= timenudge(1) ) THEN
[1365]364         nt = nt-1
[1251]365       ENDIF
366
[1365]367       dtm = ( time - timenudge(nt) ) / ( timenudge(nt+1) - timenudge(nt) )
368       dtp = ( timenudge(nt+1) - time ) / ( timenudge(nt+1) - timenudge(nt) )
[1251]369
[1239]370       SELECT CASE ( prog_var )
371
372          CASE ( 'u' )
373
374             DO  i = nxl, nxr
375                DO  j = nys, nyn
[1382]376
[2232]377                   DO  k = nzb+1, nzt
[1239]378
[1365]379                      tmp_tend = - ( hom(k,1,1,0) - ( unudge(k,nt) * dtp +     &
380                                     unudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]381
[2232]382                      tend(k,j,i) = tend(k,j,i) + tmp_tend *                   &
383                                        MERGE( 1.0_wp, 0.0_wp,                 &
384                                               BTEST( wall_flags_0(k,j,i), 1 ) )
[1239]385
[1365]386                      sums_ls_l(k,6) = sums_ls_l(k,6) + tmp_tend *             &
[1382]387                                     weight_substep(intermediate_timestep_count)
[1239]388                   ENDDO
[1382]389                 
390                   sums_ls_l(nzt+1,6) = sums_ls_l(nzt,6)
391 
[1239]392                ENDDO
393            ENDDO
394
395          CASE ( 'v' )
396
397             DO  i = nxl, nxr
398                DO  j = nys, nyn
[1382]399
[2232]400                   DO  k = nzb+1, nzt
[1239]401
[1365]402                      tmp_tend = - ( hom(k,1,2,0) - ( vnudge(k,nt) * dtp +     &
403                                     vnudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]404
[2232]405                      tend(k,j,i) = tend(k,j,i) + tmp_tend *                   &
406                                        MERGE( 1.0_wp, 0.0_wp,                 &
407                                               BTEST( wall_flags_0(k,j,i), 2 ) )
[1239]408
[1365]409                      sums_ls_l(k,7) = sums_ls_l(k,7) + tmp_tend *             &
[1382]410                                     weight_substep(intermediate_timestep_count)
[1239]411                   ENDDO
[1382]412                 
413                   sums_ls_l(nzt+1,7) = sums_ls_l(nzt,7)
414
[1239]415                ENDDO
416            ENDDO
417
418          CASE ( 'pt' )
419
420             DO  i = nxl, nxr
421                DO  j = nys, nyn
[1382]422
[2232]423                   DO  k = nzb+1, nzt
[1239]424
[1365]425                      tmp_tend = - ( hom(k,1,4,0) - ( ptnudge(k,nt) * dtp +    &
426                                     ptnudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]427
[2232]428                      tend(k,j,i) = tend(k,j,i) + tmp_tend *                   &
429                                        MERGE( 1.0_wp, 0.0_wp,                 &
430                                               BTEST( wall_flags_0(k,j,i), 0 ) )
[1239]431
[1365]432                      sums_ls_l(k,4) = sums_ls_l(k,4) + tmp_tend *             &
[1382]433                                     weight_substep(intermediate_timestep_count)
[1239]434                   ENDDO
[1382]435
436                   sums_ls_l(nzt+1,4) = sums_ls_l(nzt,4)
437
[1239]438                ENDDO
439            ENDDO
440
441          CASE ( 'q' )
442
443             DO  i = nxl, nxr
444                DO  j = nys, nyn
[1382]445
[2232]446                   DO  k = nzb+1, nzt
[1239]447
[1365]448                      tmp_tend = - ( hom(k,1,41,0) - ( qnudge(k,nt) * dtp +    &
449                                     qnudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]450
[2232]451                      tend(k,j,i) = tend(k,j,i) + tmp_tend *                   &
452                                        MERGE( 1.0_wp, 0.0_wp,                 &
453                                               BTEST( wall_flags_0(k,j,i), 0 ) )
[1239]454
[1365]455                      sums_ls_l(k,5) = sums_ls_l(k,5) + tmp_tend *             &
[1382]456                                     weight_substep(intermediate_timestep_count)
[1239]457                   ENDDO
[1382]458                 
459                   sums_ls_l(nzt+1,5) = sums_ls_l(nzt,5)
460
[1239]461                ENDDO
462            ENDDO
463
464          CASE DEFAULT
465             message_string = 'unknown prognostic variable "' // prog_var // '"'
466             CALL message( 'nudge', 'PA0367', 1, 2, 0, 6, 0 )
467
468       END SELECT
469
470    END SUBROUTINE nudge
471
472
473!------------------------------------------------------------------------------!
[1682]474! Description:
475! ------------
476!> Call for grid point i,j
[1239]477!------------------------------------------------------------------------------!
478
479    SUBROUTINE nudge_ij( i, j, time, prog_var )
480
[1320]481       USE arrays_3d,                                                          &
[1365]482           ONLY:  pt, ptnudge, q, qnudge, tend, timenudge, tmp_tnudge, tnudge, &
483                  u, unudge, v, vnudge
[1239]484
[1320]485       USE control_parameters,                                                 &
[1365]486           ONLY:  dt_3d, intermediate_timestep_count, message_string
[1320]487
488       USE indices,                                                            &
[2232]489           ONLY:  nxl, nxr, nys, nyn, nzb, nzt, wall_flags_0
[1320]490
[1365]491       USE kinds
[1320]492
493       USE statistics,                                                         &
[1382]494           ONLY:  hom, sums_ls_l, weight_substep
[1320]495
[1239]496       IMPLICIT NONE
497
498
[1682]499       CHARACTER (LEN=*) ::  prog_var  !<
[1239]500
[1682]501       REAL(wp) ::  tmp_tend    !<
502       REAL(wp) ::  dtm         !<
503       REAL(wp) ::  dtp         !<
504       REAL(wp) ::  time        !<
[1239]505
[1682]506       INTEGER(iwp) ::  i  !<
507       INTEGER(iwp) ::  j  !<
508       INTEGER(iwp) ::  k  !<
509       INTEGER(iwp) ::  nt  !<
[1239]510
[1320]511
[1365]512       nt = 1
513       DO WHILE ( time > timenudge(nt) )
514         nt = nt+1
[1239]515       ENDDO
[1249]516       IF ( time /= timenudge(1) )  THEN
[1365]517         nt = nt-1
[1239]518       ENDIF
519
[1365]520       dtm = ( time - timenudge(nt) ) / ( timenudge(nt+1) - timenudge(nt) )
521       dtp = ( timenudge(nt+1) - time ) / ( timenudge(nt+1) - timenudge(nt) )
[1239]522
523       SELECT CASE ( prog_var )
524
525          CASE ( 'u' )
526
[2232]527             DO  k = nzb+1, nzt
[1239]528
[1365]529                tmp_tend = - ( hom(k,1,1,0) - ( unudge(k,nt) * dtp +           &
530                               unudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]531
[2232]532                tend(k,j,i) = tend(k,j,i) + tmp_tend *                         &
533                                        MERGE( 1.0_wp, 0.0_wp,                 &
534                                               BTEST( wall_flags_0(k,j,i), 1 ) )
[1365]535
536                sums_ls_l(k,6) = sums_ls_l(k,6) + tmp_tend                     &
[1382]537                                 * weight_substep(intermediate_timestep_count)
[1239]538             ENDDO
539
[1382]540             sums_ls_l(nzt+1,6) = sums_ls_l(nzt,6)
541
[1239]542          CASE ( 'v' )
543
[2232]544             DO  k = nzb+1, nzt
[1239]545
[1365]546                tmp_tend = - ( hom(k,1,2,0) - ( vnudge(k,nt) * dtp +           &
547                               vnudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]548
[2232]549                tend(k,j,i) = tend(k,j,i) + tmp_tend *                         &
550                                        MERGE( 1.0_wp, 0.0_wp,                 &
551                                               BTEST( wall_flags_0(k,j,i), 2 ) )
[1365]552
553                sums_ls_l(k,7) = sums_ls_l(k,7) + tmp_tend                     &
[1382]554                                 * weight_substep(intermediate_timestep_count)
[1239]555             ENDDO
556
[1382]557             sums_ls_l(nzt+1,7) = sums_ls_l(nzt,7)
[1239]558
559          CASE ( 'pt' )
560
[2232]561             DO  k = nzb+1, nzt
[1239]562
[1365]563                tmp_tend = - ( hom(k,1,4,0) - ( ptnudge(k,nt) * dtp +          &
564                               ptnudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]565
[2232]566                tend(k,j,i) = tend(k,j,i) + tmp_tend *                         &
567                                        MERGE( 1.0_wp, 0.0_wp,                 &
568                                               BTEST( wall_flags_0(k,j,i), 0 ) )
[1365]569
570                sums_ls_l(k,4) = sums_ls_l(k,4) + tmp_tend                     &
[1382]571                                 * weight_substep(intermediate_timestep_count)
[1239]572             ENDDO
573
[1382]574             sums_ls_l(nzt+1,4) = sums_ls_l(nzt,4)
[1239]575
[1382]576
[1239]577          CASE ( 'q' )
578
[2232]579             DO  k = nzb+1, nzt
[1239]580
[1365]581                tmp_tend = - ( hom(k,1,41,0) - ( qnudge(k,nt) * dtp +          &
582                               qnudge(k,nt+1) * dtm ) ) / tmp_tnudge(k)
[1239]583
[2232]584                tend(k,j,i) = tend(k,j,i) + tmp_tend *                         &
585                                        MERGE( 1.0_wp, 0.0_wp,                 &
586                                               BTEST( wall_flags_0(k,j,i), 0 ) )
[1365]587
588                sums_ls_l(k,5) = sums_ls_l(k,5) + tmp_tend                     &
[1382]589                                 * weight_substep(intermediate_timestep_count)
[1239]590             ENDDO
591
[1382]592             sums_ls_l(nzt+1,5) = sums_ls_l(nzt,5)
593
[1239]594          CASE DEFAULT
595             message_string = 'unknown prognostic variable "' // prog_var // '"'
596             CALL message( 'nudge', 'PA0367', 1, 2, 0, 6, 0 )
597
598       END SELECT
599
600
601    END SUBROUTINE nudge_ij
602
[1380]603
[1682]604!------------------------------------------------------------------------------!
605! Description:
606! ------------
607!> @todo Missing subroutine description.
608!------------------------------------------------------------------------------!
[1380]609    SUBROUTINE nudge_ref ( time )
610
611       USE arrays_3d,                                                          &
[1398]612           ONLY:  time_vert, ptnudge, pt_init, qnudge, q_init, unudge, u_init, &
613                  vnudge, v_init
[1380]614
615       USE kinds
616
617
618       IMPLICIT NONE
619
[1682]620       INTEGER(iwp) ::  nt                    !<
[1380]621
[1682]622       REAL(wp)             ::  fac           !<
623       REAL(wp), INTENT(in) ::  time          !<
[1380]624
625!
626!--    Interpolation in time of NUDGING_DATA for pt_init and q_init. This is
627!--    needed for correct upper boundary conditions for pt and q and in case that
628!      large scale subsidence as well as scalar Rayleigh-damping are used
629       nt = 1
630       DO WHILE ( time > time_vert(nt) )
631          nt = nt + 1
632       ENDDO
633       IF ( time /= time_vert(nt) )  THEN
634        nt = nt - 1
635       ENDIF
636
637       fac = ( time-time_vert(nt) ) / ( time_vert(nt+1)-time_vert(nt) )
638
639       pt_init = ptnudge(:,nt) + fac * ( ptnudge(:,nt+1) - ptnudge(:,nt) )
640       q_init  = qnudge(:,nt) + fac * ( qnudge(:,nt+1) - qnudge(:,nt) )
[1398]641       u_init  = unudge(:,nt) + fac * ( unudge(:,nt+1) - unudge(:,nt) )
642       v_init  = vnudge(:,nt) + fac * ( vnudge(:,nt+1) - vnudge(:,nt) )
[1380]643
644    END SUBROUTINE nudge_ref
645
[1239]646 END MODULE nudge_mod
Note: See TracBrowser for help on using the repository browser.