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