source: LMDZ6/trunk/libf/dyn3d/integrd.f90 @ 5435

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

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