source: LMDZ6/branches/Amaury_dev/libf/dyn3d/integrd.F90 @ 5134

Last change on this file since 5134 was 5134, checked in by abarral, 4 months ago

Replace academic.h, alpale.h, comdissip.h, comdissipn.h, comdissnew.h by modules
Remove unused clesph0.h

  • 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.7 KB
Line 
1! $Id: integrd.F90 5134 2024-07-26 15:56:37Z abarral $
2
3SUBROUTINE integrd &
4        (nq, vcovm1, ucovm1, tetam1, psm1, massem1, &
5        dv, du, dteta, dq, dp, vcov, ucov, teta, q, ps, masse, phis & !,finvmaold
6        )
7
8  USE control_mod, ONLY: planet_type
9  USE comconst_mod, ONLY: pi
10  USE logic_mod, ONLY: leapf
11  USE comvert_mod, ONLY: ap, bp
12  USE temps_mod, ONLY: dt
13  USE lmdz_iniprint, ONLY: lunout, prt_level
14  USE lmdz_ssum_scopy, ONLY: scopy, ssum
15
16
17  IMPLICIT NONE
18
19
20  !=======================================================================
21  !
22  !   Auteur:  P. Le Van
23  !   -------
24  !
25  !   objet:
26  !   ------
27  !
28  !   Incrementation des tendances dynamiques
29  !
30  !=======================================================================
31  !-----------------------------------------------------------------------
32  !   Declarations:
33  !   -------------
34
35  INCLUDE "dimensions.h"
36  INCLUDE "paramet.h"
37  INCLUDE "comgeom.h"
38
39  !   Arguments:
40  !   ----------
41
42  INTEGER, INTENT(IN) :: nq ! number of tracers to handle in this routine
43  REAL, INTENT(INOUT) :: vcov(ip1jm, llm) ! covariant meridional wind
44  REAL, INTENT(INOUT) :: ucov(ip1jmp1, llm) ! covariant zonal wind
45  REAL, INTENT(INOUT) :: teta(ip1jmp1, llm) ! potential temperature
46  REAL, INTENT(INOUT) :: q(ip1jmp1, llm, nq) ! advected tracers
47  REAL, INTENT(INOUT) :: ps(ip1jmp1) ! surface pressure
48  REAL, INTENT(INOUT) :: masse(ip1jmp1, llm) ! atmospheric mass
49  REAL, INTENT(IN) :: phis(ip1jmp1) ! ground geopotential !!! unused
50  ! values at previous time step
51  REAL, INTENT(INOUT) :: vcovm1(ip1jm, llm)
52  REAL, INTENT(INOUT) :: ucovm1(ip1jmp1, llm)
53  REAL, INTENT(INOUT) :: tetam1(ip1jmp1, llm)
54  REAL, INTENT(INOUT) :: psm1(ip1jmp1)
55  REAL, INTENT(INOUT) :: massem1(ip1jmp1, llm)
56  ! the tendencies to add
57  REAL, INTENT(IN) :: dv(ip1jm, llm)
58  REAL, INTENT(IN) :: du(ip1jmp1, llm)
59  REAL, INTENT(IN) :: dteta(ip1jmp1, llm)
60  REAL, INTENT(IN) :: dp(ip1jmp1)
61  REAL, INTENT(IN) :: dq(ip1jmp1, llm, nq) !!! unused
62  ! REAL,INTENT(OUT) :: finvmaold(ip1jmp1,llm) !!! unused
63
64  !   Local:
65  !   ------
66
67  REAL :: vscr(ip1jm), uscr(ip1jmp1), hscr(ip1jmp1), pscr(ip1jmp1)
68  REAL :: massescr(ip1jmp1, llm)
69  ! REAL finvmasse(ip1jmp1,llm)
70  REAL :: p(ip1jmp1, llmp1)
71  REAL :: tpn, tps, tppn(iim), tpps(iim)
72  REAL :: qpn, qps, qppn(iim), qpps(iim)
73  REAL :: deltap(ip1jmp1, llm)
74
75  INTEGER :: l, ij, iq, i, j
76
77  !-----------------------------------------------------------------------
78
79  DO  l = 1, llm
80    DO  ij = 1, iip1
81      ucov(ij, l) = 0.
82      ucov(ij + ip1jm, l) = 0.
83      uscr(ij) = 0.
84      uscr(ij + ip1jm) = 0.
85    ENDDO
86  ENDDO
87
88
89  !    ............    integration  de       ps         ..............
90
91  CALL SCOPY(ip1jmp1 * llm, masse, 1, massescr, 1)
92
93  DO ij = 1, ip1jmp1
94    pscr (ij) = ps(ij)
95    ps (ij) = psm1(ij) + dt * dp(ij)
96  ENDDO
97  !
98  DO ij = 1, ip1jmp1
99    IF(ps(ij)<0.) THEN
100      WRITE(lunout, *) "integrd: negative surface pressure ", ps(ij)
101      WRITE(lunout, *) " at node ij =", ij
102      ! since ij=j+(i-1)*jjp1 , we have
103      j = modulo(ij, jjp1)
104      i = 1 + (ij - j) / jjp1
105      WRITE(lunout, *) " lon = ", rlonv(i) * 180. / pi, " deg", &
106              " lat = ", rlatu(j) * 180. / pi, " deg"
107      CALL abort_gcm("integrd", "", 1)
108    ENDIF
109  ENDDO
110  !
111  DO  ij = 1, iim
112    tppn(ij) = aire(ij) * ps(ij)
113    tpps(ij) = aire(ij + ip1jm) * ps(ij + ip1jm)
114  ENDDO
115  tpn = SSUM(iim, tppn, 1) / apoln
116  tps = SSUM(iim, tpps, 1) / apols
117  DO ij = 1, iip1
118    ps(ij) = tpn
119    ps(ij + ip1jm) = tps
120  ENDDO
121  !
122  !  ... Calcul  de la nouvelle masse d'air au dernier temps integre t+1 ...
123  !
124  CALL pression (ip1jmp1, ap, bp, ps, p)
125  CALL massdair (p, masse)
126
127  ! Ehouarn : we don't use/need finvmaold and finvmasse,
128  ! so might as well not compute them
129  ! CALL   SCOPY( ijp1llm  , masse, 1, finvmasse,  1      )
130  ! CALL filtreg( finvmasse, jjp1, llm, -2, 2, .TRUE., 1  )
131  !
132
133  !    ............   integration  de  ucov, vcov,  h     ..............
134
135  DO l = 1, llm
136
137    DO ij = iip2, ip1jm
138      uscr(ij) = ucov(ij, l)
139      ucov(ij, l) = ucovm1(ij, l) + dt * du(ij, l)
140    ENDDO
141
142    DO ij = 1, ip1jm
143      vscr(ij) = vcov(ij, l)
144      vcov(ij, l) = vcovm1(ij, l) + dt * dv(ij, l)
145    ENDDO
146
147    DO ij = 1, ip1jmp1
148      hscr(ij) = teta(ij, l)
149      teta (ij, l) = tetam1(ij, l) * massem1(ij, l) / masse(ij, l) &
150              + dt * dteta(ij, l) / masse(ij, l)
151    ENDDO
152
153    !   ....  Calcul de la valeur moyenne, unique  aux poles pour  teta    ......
154    !
155    !
156    DO  ij = 1, iim
157      tppn(ij) = aire(ij) * teta(ij, l)
158      tpps(ij) = aire(ij + ip1jm) * teta(ij + ip1jm, l)
159    ENDDO
160    tpn = SSUM(iim, tppn, 1) / apoln
161    tps = SSUM(iim, tpps, 1) / apols
162
163    DO ij = 1, iip1
164      teta(ij, l) = tpn
165      teta(ij + ip1jm, l) = tps
166    ENDDO
167    !
168
169    IF(leapf)  THEN
170      CALL SCOPY (ip1jmp1, uscr(1), 1, ucovm1(1, l), 1)
171      CALL SCOPY (ip1jm, vscr(1), 1, vcovm1(1, l), 1)
172      CALL SCOPY (ip1jmp1, hscr(1), 1, tetam1(1, l), 1)
173    END IF
174
175  ENDDO ! of DO l = 1,llm
176
177
178  !
179  !   .......  integration de   q   ......
180  !
181  !$$$      IF( iadv(1).NE.3.AND.iadv(2).NE.3 )    THEN
182  !$$$c
183  !$$$       IF( forward .OR.  leapf )  THEN
184  !$$$        DO iq = 1,2
185  !$$$        DO  l = 1,llm
186  !$$$        DO ij = 1,ip1jmp1
187  !$$$        q(ij,l,iq) = ( q(ij,l,iq)*finvmaold(ij,l) + dtvr *dq(ij,l,iq) )/
188  !$$$     $                            finvmasse(ij,l)
189  !$$$        ENDDO
190  !$$$        ENDDO
191  !$$$        ENDDO
192  !$$$       ELSE
193  !$$$         DO iq = 1,2
194  !$$$         DO  l = 1,llm
195  !$$$         DO ij = 1,ip1jmp1
196  !$$$         q( ij,l,iq ) = q( ij,l,iq ) * finvmaold(ij,l) / finvmasse(ij,l)
197  !$$$         ENDDO
198  !$$$         ENDDO
199  !$$$         ENDDO
200  !$$$
201  !$$$       END IF
202  !$$$c
203  !$$$      ENDIF
204
205  IF (planet_type=="earth") THEN
206    ! Earth-specific treatment of first 2 tracers (water)
207    DO l = 1, llm
208      DO ij = 1, ip1jmp1
209        deltap(ij, l) = p(ij, l) - p(ij, l + 1)
210      ENDDO
211    ENDDO
212
213    CALL qminimum(q, nq, deltap)
214
215    !
216    !    .....  Calcul de la valeur moyenne, unique  aux poles pour  q .....
217    !
218
219    DO iq = 1, nq
220      DO l = 1, llm
221
222        DO ij = 1, iim
223          qppn(ij) = aire(ij) * q(ij, l, iq)
224          qpps(ij) = aire(ij + ip1jm) * q(ij + ip1jm, l, iq)
225        ENDDO
226        qpn = SSUM(iim, qppn, 1) / apoln
227        qps = SSUM(iim, qpps, 1) / apols
228
229        DO ij = 1, iip1
230          q(ij, l, iq) = qpn
231          q(ij + ip1jm, l, iq) = qps
232        ENDDO
233
234      ENDDO
235    ENDDO
236
237    ! Ehouarn: forget about finvmaold
238    ! CALL  SCOPY( ijp1llm , finvmasse, 1, finvmaold, 1 )
239
240  ENDIF ! of if (planet_type.EQ."earth")
241  !
242  !
243  ! .....   FIN  de l'integration  de   q    .......
244
245  !    .................................................................
246
247  IF(leapf)  THEN
248    CALL SCOPY (ip1jmp1, pscr, 1, psm1, 1)
249    CALL SCOPY (ip1jmp1 * llm, massescr, 1, massem1, 1)
250  END IF
251
252END SUBROUTINE integrd
Note: See TracBrowser for help on using the repository browser.