source: LMDZ6/trunk/libf/dyn3dmem/caldyn_loc.F90 @ 5250

Last change on this file since 5250 was 5246, checked in by abarral, 29 hours ago

Convert fixed-form to free-form sources .F -> .{f,F}90
(WIP: some .F remain, will be handled in subsequent commits)

  • 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
File size: 6.0 KB
Line 
1!
2! $Id: $
3!
4#undef DEBUG_IO
5!#define DEBUG_IO
6
7SUBROUTINE caldyn_loc &
8        (itau,ucov,vcov,teta,ps,masse,pk,pkf,phis , &
9        phi,conser,du,dv,dteta,dp,w,pbaru,pbarv,time )
10  USE parallel_lmdz
11  USE Write_Field_loc
12  USE caldyn_mod, ONLY: vcont, ucont, ang, p, massebx, masseby, &
13        vorpot, ecin, bern, massebxy, convm
14  USE comvert_mod, ONLY: ap, bp
15
16  IMPLICIT NONE
17
18  !=======================================================================
19  !
20  !  Auteur :  P. Le Van
21  !
22  !   Objet:
23  !   ------
24  !
25  !   Calcul des tendances dynamiques.
26  !
27  ! Modif 04/93 F.Forget
28  !=======================================================================
29
30  !-----------------------------------------------------------------------
31  !   0. Declarations:
32  !   ----------------
33
34  include "dimensions.h"
35  include "paramet.h"
36  include "comgeom.h"
37
38  !   Arguments:
39  !   ----------
40
41  LOGICAL,INTENT(IN) :: conser ! triggers printing some diagnostics ! not used
42  INTEGER,INTENT(IN) :: itau ! time step index ! not used
43  REAL,INTENT(IN) :: vcov(ijb_v:ije_v,llm) ! covariant meridional wind
44  REAL,INTENT(IN) :: ucov(ijb_u:ije_u,llm) ! covariant zonal wind
45  REAL,INTENT(IN) :: teta(ijb_u:ije_u,llm) ! potential temperature
46  REAL,INTENT(IN) :: ps(ijb_u:ije_u) ! surface pressure
47  REAL,INTENT(IN) :: phis(ijb_u:ije_u) ! geopotential at the surface
48  REAL,INTENT(IN) :: pk(iip1,jjb_u:jje_u,llm) ! Exner at mid-layer
49  REAL,INTENT(IN) :: pkf(ijb_u:ije_u,llm) ! filtered Exner
50  REAL,INTENT(IN) :: phi(ijb_u:ije_u,llm) ! geopotential
51  REAL,INTENT(OUT) :: masse(ijb_u:ije_u,llm) ! air mass
52  REAL,INTENT(OUT) :: dv(ijb_v:ije_v,llm) ! tendency on vcov
53  REAL,INTENT(OUT) :: du(ijb_u:ije_u,llm) ! tendency on ucov
54  REAL,INTENT(OUT) :: dteta(ijb_u:ije_u,llm) ! tenddency on teta
55  REAL,INTENT(OUT) :: dp(ijb_u:ije_u) ! tendency on ps
56  REAL,INTENT(OUT) :: w(ijb_u:ije_u,llm) ! vertical velocity
57  REAL,INTENT(OUT) :: pbaru(ijb_u:ije_u,llm) ! mass flux in the zonal direction
58  REAL,INTENT(OUT) :: pbarv(ijb_v:ije_v,llm) ! mass flux in the meridional direction
59  REAL,INTENT(IN) :: time ! current time
60
61  !   Local:
62  !   ------
63
64  INTEGER :: ij,l,ijb,ije,ierr
65
66
67  !-----------------------------------------------------------------------
68  !   Compute dynamical tendencies:
69  !--------------------------------
70
71  ! ! compute contravariant winds ucont() and vcont
72  CALL covcont_loc  ( llm    , ucov    , vcov , ucont, vcont     )
73  ! ! compute pressure p()
74  CALL pression_loc ( ip1jmp1, ap      , bp   ,  ps  , p         )
75  !ym      CALL psextbar (   ps   , psexbarxy                          )
76!$OMP BARRIER
77  ! ! compute mass in each atmospheric mesh: masse()
78  CALL massdair_loc (    p   , masse                             )
79  ! ! compute X and Y-averages of mass, massebx() and masseby()
80  CALL massbar_loc  (   masse, massebx , masseby                 )
81  ! ! compute XY-average of mass, massebxy()
82  call massbarxy_loc(   masse, massebxy                          )
83  ! ! compute mass fluxes pbaru() and pbarv()
84  CALL flumass_loc  ( massebx, masseby,vcont,ucont,pbaru,pbarv   )
85  ! ! compute dteta() , horizontal converging flux of theta
86  CALL dteta1_loc   (   teta , pbaru   , pbarv, dteta            )
87  ! ! compute convm(), horizontal converging flux of mass
88  CALL convmas1_loc  (   pbaru, pbarv   , convm                  )
89!$OMP BARRIER
90  CALL convmas2_loc  (   convm                      )
91!$OMP BARRIER
92#ifdef DEBUG_IO
93  call WriteField_u('ucont',ucont)
94  call WriteField_v('vcont',vcont)
95  call WriteField_u('p',p)
96  call WriteField_u('masse',masse)
97  call WriteField_u('massebx',massebx)
98  call WriteField_v('masseby',masseby)
99  call WriteField_v('massebxy',massebxy)
100  call WriteField_u('pbaru',pbaru)
101  call WriteField_v('pbarv',pbarv)
102  call WriteField_u('dteta',dteta)
103  call WriteField_u('convm',convm)
104#endif
105
106!$OMP BARRIER
107!$OMP MASTER
108  ijb=ij_begin
109  ije=ij_end
110  ! ! compute pressure variation due to mass convergence
111  DO ij =ijb, ije
112     dp( ij ) = convm( ij,1 ) / airesurg( ij )
113  ENDDO
114!$OMP END MASTER
115!$OMP BARRIER
116
117  ! ! compute vertical velocity w()
118  CALL vitvert_loc ( convm  , w                                )
119  ! ! compute potential vorticity vorpot()
120  CALL tourpot_loc ( vcov   , ucov  , massebxy  , vorpot       )
121  ! ! compute rotation induced du() and dv()
122  CALL dudv1_loc   ( vorpot , pbaru , pbarv     , du     , dv  )
123
124#ifdef DEBUG_IO
125  call WriteField_u('w',w)
126  call WriteField_v('vorpot',vorpot)
127  call WriteField_u('du',du)
128  call WriteField_v('dv',dv)
129#endif
130
131  ! ! compute kinetic energy ecin()
132  CALL enercin_loc ( vcov   , ucov  , vcont   , ucont  , ecin  )
133  ! ! compute Bernouilli function bern()
134  CALL bernoui_loc ( ip1jmp1, llm   , phi       , ecin   , bern)
135  ! ! compute and add du() and dv() contributions from Bernouilli and pressure
136  CALL dudv2_loc   ( teta   , pkf   , bern      , du     , dv  )
137
138#ifdef DEBUG_IO
139  call WriteField_u('ecin',ecin)
140  call WriteField_u('bern',bern)
141  call WriteField_u('du',du)
142  call WriteField_v('dv',dv)
143  call WriteField_u('pkf',pkf)
144#endif
145
146  ijb=ij_begin-iip1
147  ije=ij_end+iip1
148
149  if (pole_nord) ijb=ij_begin
150  if (pole_sud) ije=ij_end
151
152!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
153  DO l=1,llm
154     DO ij=ijb,ije
155        ang(ij,l) = ucov(ij,l) + constang(ij)
156    ENDDO
157  ENDDO
158!$OMP END DO
159
160  ! ! compute vertical advection contributions to du(), dv() and dteta()
161  CALL advect_new_loc(ang,vcov,teta,w,massebx,masseby,du,dv,dteta)
162
163  !  WARNING probleme de peridocite de dv sur les PC/linux. Pb d'arrondi
164       ! probablement. Observe sur le code compile avec pgf90 3.0-1
165  ijb=ij_begin
166  ije=ij_end
167  if (pole_sud) ije=ij_end-iip1
168
169!$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
170  DO l = 1, llm
171     DO ij = ijb, ije, iip1
172       IF( dv(ij,l).NE.dv(ij+iim,l) )  THEN
173      ! PRINT *,'!!!ATTENTION!!! probleme de periodicite sur vcov',
174  !    ,   ' dans caldyn'
175      ! PRINT *,' l,  ij = ', l, ij, ij+iim,dv(ij+iim,l),dv(ij,l)
176      dv(ij+iim,l) = dv(ij,l)
177      endif
178     enddo
179  enddo
180!$OMP END DO NOWAIT
181
182  ! Ehouarn: NB: output of control variables not implemented...
183
184  RETURN
185END SUBROUTINE caldyn_loc
Note: See TracBrowser for help on using the repository browser.