source: LMDZ6/branches/Amaury_dev/libf/dyn3d/top_bound.F90 @ 5105

Last change on this file since 5105 was 5103, checked in by abarral, 5 months ago

Handle CPP_INLANDSIS in lmdz_cppkeys_wrapper.F90
Remove obsolete key wrgrads_thermcell, _ADV_HALO, _ADV_HALLO, isminmax
Remove redundant uses of CPPKEY_INCA (thanks acozic)
Remove obsolete misc/write_field.F90
Remove unused ioipsl_* wrappers
Remove calls to WriteField_u with wrong signature
Convert .F -> .[fF]90
(lint) uppercase fortran operators
[note: 1d and iso still broken - working on it]

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 6.0 KB
Line 
1! $Id: top_bound.F90 5103 2024-07-23 13:29:36Z abarral $
2
3SUBROUTINE top_bound(vcov, ucov, teta, masse, dt)
4
5  USE comconst_mod, ONLY: iflag_top_bound, mode_top_bound, &
6          tau_top_bound
7  USE comvert_mod, ONLY: presnivs, preff, scaleheight
8
9  IMPLICIT NONE
10  !
11  include "dimensions.h"
12  include "paramet.h"
13  include "comgeom2.h"
14
15
16  ! ..  DISSIPATION LINEAIRE A HAUT NIVEAU, RUN MESO,
17  ! F. LOTT DEC. 2006
18  !                             (  10/12/06  )
19
20  !=======================================================================
21  !
22  !   Auteur:  F. LOTT
23  !   -------
24  !
25  !   Objet:
26  !   ------
27  !
28  !   Dissipation linéaire (ex top_bound de la physique)
29  !
30  !=======================================================================
31
32  ! top_bound sponge layer model:
33  ! Quenching is modeled as: A(t)=Am+A0*exp(-lambda*t)
34  ! where Am is the zonal average of the field (or zero), and lambda the inverse
35  ! of the characteristic quenching/relaxation time scale
36  ! Thus, assuming Am to be time-independent, field at time t+dt is given by:
37  ! A(t+dt)=A(t)-(A(t)-Am)*(1-exp(-lambda*t))
38  ! Moreover lambda can be a function of model level (see below), and relaxation
39  ! can be toward the average zonal field or just zero (see below).
40
41  ! NB: top_bound sponge is only called from leapfrog if ok_strato=.TRUE.
42
43  ! sponge parameters: (loaded/set in conf_gcm.F ; stored in comconst_mod)
44  !    iflag_top_bound=0 for no sponge
45  !    iflag_top_bound=1 for sponge over 4 topmost layers
46  !    iflag_top_bound=2 for sponge from top to ~1% of top layer pressure
47  !    mode_top_bound=0: no relaxation
48  !    mode_top_bound=1: u and v relax towards 0
49  !    mode_top_bound=2: u and v relax towards their zonal mean
50  !    mode_top_bound=3: u,v and pot. temp. relax towards their zonal mean
51  !    tau_top_bound : inverse of charactericstic relaxation time scale at
52  !                   the topmost layer (Hz)
53
54  include "comdissipn.h"
55  include "iniprint.h"
56
57  !   Arguments:
58  !   ----------
59
60  real, intent(inout) :: ucov(iip1, jjp1, llm) ! covariant zonal wind
61  real, intent(inout) :: vcov(iip1, jjm, llm) ! covariant meridional wind
62  real, intent(inout) :: teta(iip1, jjp1, llm) ! potential temperature
63  real, intent(in) :: masse(iip1, jjp1, llm) ! mass of atmosphere
64  real, intent(in) :: dt ! time step (s) of sponge model
65
66  !   Local:
67  !   ------
68
69  REAL :: massebx(iip1, jjp1, llm), masseby(iip1, jjm, llm), zm
70  REAL :: uzon(jjp1, llm), vzon(jjm, llm), tzon(jjp1, llm)
71
72  integer :: i
73  REAL, SAVE :: rdamp(llm) ! quenching coefficient
74  real, save :: lambda(llm) ! inverse or quenching time scale (Hz)
75
76  LOGICAL, SAVE :: first = .TRUE.
77
78  INTEGER :: j, l
79
80  if (iflag_top_bound==0) return
81
82  if (first) then
83    if (iflag_top_bound==1) then
84      ! sponge quenching over the topmost 4 atmospheric layers
85      lambda(:) = 0.
86      lambda(llm) = tau_top_bound
87      lambda(llm - 1) = tau_top_bound / 2.
88      lambda(llm - 2) = tau_top_bound / 4.
89      lambda(llm - 3) = tau_top_bound / 8.
90    else if (iflag_top_bound==2) then
91      ! sponge quenching over topmost layers down to pressures which are
92      ! higher than 100 times the topmost layer pressure
93      lambda(:) = tau_top_bound &
94              * max(presnivs(llm) / presnivs(:) - 0.01, 0.)
95    endif
96
97    ! quenching coefficient rdamp(:)
98    ! rdamp(:)=dt*lambda(:) ! Explicit Euler approx.
99    rdamp(:) = 1. - exp(-lambda(:) * dt)
100
101    write(lunout, *)'TOP_BOUND mode', mode_top_bound
102    write(lunout, *)'Sponge layer coefficients'
103    write(lunout, *)'p (Pa)  z(km)  tau(s)   1./tau (Hz)'
104    do l = 1, llm
105      if (rdamp(l)/=0.) then
106        write(lunout, '(6(1pe12.4,1x))') &
107                presnivs(l), log(preff / presnivs(l)) * scaleheight, &
108                1. / lambda(l), lambda(l)
109      endif
110    enddo
111    first = .FALSE.
112  endif ! of if (first)
113
114  CALL massbar(masse, massebx, masseby)
115
116  ! ! compute zonal average of vcov and u
117  if (mode_top_bound>=2) then
118    do l = 1, llm
119      do j = 1, jjm
120        vzon(j, l) = 0.
121        zm = 0.
122        do i = 1, iim
123          ! NB: we can work using vcov zonal mean rather than v since the
124          ! cv coefficient (which relates the two) only varies with latitudes
125          vzon(j, l) = vzon(j, l) + vcov(i, j, l) * masseby(i, j, l)
126          zm = zm + masseby(i, j, l)
127        enddo
128        vzon(j, l) = vzon(j, l) / zm
129      enddo
130    enddo
131
132    do l = 1, llm
133      do j = 2, jjm ! excluding poles
134        uzon(j, l) = 0.
135        zm = 0.
136        do i = 1, iim
137          uzon(j, l) = uzon(j, l) + massebx(i, j, l) * ucov(i, j, l) / cu(i, j)
138          zm = zm + massebx(i, j, l)
139        enddo
140        uzon(j, l) = uzon(j, l) / zm
141      enddo
142    enddo
143  else ! ucov and vcov will relax towards 0
144    vzon(:, :) = 0.
145    uzon(:, :) = 0.
146  endif ! of if (mode_top_bound.ge.2)
147
148  ! ! compute zonal average of potential temperature, if necessary
149  if (mode_top_bound>=3) then
150    do l = 1, llm
151      do j = 2, jjm ! excluding poles
152        zm = 0.
153        tzon(j, l) = 0.
154        do i = 1, iim
155          tzon(j, l) = tzon(j, l) + teta(i, j, l) * masse(i, j, l)
156          zm = zm + masse(i, j, l)
157        enddo
158        tzon(j, l) = tzon(j, l) / zm
159      enddo
160    enddo
161  endif ! of if (mode_top_bound.ge.3)
162
163  if (mode_top_bound>=1) then
164    ! ! Apply sponge quenching on vcov:
165    do l = 1, llm
166      do i = 1, iip1
167        do j = 1, jjm
168          vcov(i, j, l) = vcov(i, j, l) &
169                  - rdamp(l) * (vcov(i, j, l) - vzon(j, l))
170        enddo
171      enddo
172    enddo
173
174    ! ! Apply sponge quenching on ucov:
175    do l = 1, llm
176      do i = 1, iip1
177        do j = 2, jjm ! excluding poles
178          ucov(i, j, l) = ucov(i, j, l) &
179                  - rdamp(l) * (ucov(i, j, l) - cu(i, j) * uzon(j, l))
180        enddo
181      enddo
182    enddo
183  endif ! of if (mode_top_bound.ge.1)
184
185  if (mode_top_bound>=3) then
186    ! ! Apply sponge quenching on teta:
187    do l = 1, llm
188      do i = 1, iip1
189        do j = 2, jjm ! excluding poles
190          teta(i, j, l) = teta(i, j, l) &
191                  - rdamp(l) * (teta(i, j, l) - tzon(j, l))
192        enddo
193      enddo
194    enddo
195  endif ! of if (mode_top_bound.ge.3)
196
197END SUBROUTINE top_bound
Note: See TracBrowser for help on using the repository browser.