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

Last change on this file since 2282 was 2277, checked in by kanani, 7 years ago

code documentation and cleanup

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