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

Last change on this file since 5136 was 5136, checked in by abarral, 8 weeks ago

Put comgeom.h, comgeom2.h into modules

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