LMDZ
flumass_p.F
Go to the documentation of this file.
1  SUBROUTINE flumass_p(massebx,masseby, vcont, ucont, pbaru, pbarv)
3  IMPLICIT NONE
4 
5 c=======================================================================
6 c
7 c Auteurs: P. Le Van, F. Hourdin .
8 c -------
9 c
10 c Objet:
11 c ------
12 c
13 c *********************************************************************
14 c .... calcul du flux de masse aux niveaux s ......
15 c *********************************************************************
16 c massebx,masseby,vcont et ucont sont des argum. d'entree pour le s-pg .
17 c pbaru et pbarv sont des argum.de sortie pour le s-pg .
18 c
19 c=======================================================================
20 
21 
22 #include "dimensions.h"
23 #include "paramet.h"
24 #include "comgeom.h"
25 
26  REAL massebx( ip1jmp1,llm ),masseby( ip1jm,llm ) ,
27  * vcont( ip1jm,llm ),ucont( ip1jmp1,llm ),pbaru( ip1jmp1,llm ),
28  * pbarv( ip1jm,llm )
29 
30  REAL apbarun( iip1 ),apbarus( iip1 )
31 
32  REAL sairen,saireun,saires,saireus,ctn,cts,ctn0,cts0
33  INTEGER l,ij,i
34  INTEGER ijb,ije
35 
36  EXTERNAL ssum
37  REAL SSUM
38 
39 c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
40  DO 5 l = 1,llm
41 
42  ijb=ij_begin
43  ije=ij_end+iip1
44 
45  if (pole_nord) ijb=ij_begin+iip1
46  if (pole_sud) ije=ij_end-iip1
47 
48  DO 1 ij = ijb,ije
49  pbaru( ij,l ) = massebx( ij,l ) * ucont( ij,l )
50  1 CONTINUE
51 
52  ijb=ij_begin-iip1
53  ije=ij_end+iip1
54 
55  if (pole_nord) ijb=ij_begin
56  if (pole_sud) ije=ij_end-iip1
57 
58  DO 3 ij = ijb,ije
59  pbarv( ij,l ) = masseby( ij,l ) * vcont( ij,l )
60  3 CONTINUE
61 
62  5 CONTINUE
63 c$OMP END DO NOWAIT
64 c ................................................................
65 c calcul de la composante du flux de masse en x aux poles .......
66 c ................................................................
67 c par la resolution d'1 systeme de 2 equations .
68 
69 c la premiere equat.decrivant le calcul de la divergence en 1 point i
70 c du pole,ce calcul etant itere de i=1 a i=im .
71 c c.a.d ,
72 c ( ( 0.5*pbaru(i)-0.5*pbaru(i-1) - pbarv(i))/aire(i) =
73 c - somme de ( pbarv(n) )/aire pole
74 
75 c l'autre equat.specifiant que la moyenne du flux de masse au pole est =0.
76 c c.a.d somme de pbaru(n)*aire locale(n) = 0.
77 
78 c on en revient ainsi a determiner la constante additive commune aux pbaru
79 c qui representait pbaru(0,j,l) dans l'equat.du calcul de la diverg.au pt
80 c i=1 .
81 c i variant de 1 a im
82 c n variant de 1 a im
83 
84  IF (pole_nord) THEN
85 
86  sairen = ssum( iim, aire( 1 ), 1 )
87  saireun= ssum( iim, aireu( 1 ), 1 )
88 
89 c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
90  DO l = 1,llm
91 
92  ctn = ssum( iim, pbarv( 1 ,l), 1 )/ sairen
93 
94  pbaru(1,l)=pbarv(1,l) - ctn * aire(1)
95 
96  DO i = 2,iim
97  pbaru(i,l) = pbaru(i- 1,l ) +
98  * pbarv(i,l) - ctn * aire(i )
99  ENDDO
100 
101  DO i = 1,iim
102  apbarun(i) = aireu( i ) * pbaru( i , l)
103  ENDDO
104 
105  ctn0 = -ssum( iim,apbarun,1 )/saireun
106 
107  DO i = 1,iim
108  pbaru( i , l) = 2. * ( pbaru( i , l) + ctn0 )
109  ENDDO
110 
111  pbaru( iip1 ,l ) = pbaru( 1 ,l )
112 
113  ENDDO
114 c$OMP END DO NOWAIT
115 
116  ENDIF
117 
118 
119  IF (pole_sud) THEN
120 
121  saires = ssum( iim, aire( ip1jm+1 ), 1 )
122  saireus= ssum( iim, aireu( ip1jm+1 ), 1 )
123 
124 c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
125  DO l = 1,llm
126 
127  cts = ssum( iim, pbarv(ip1jmi1+ 1,l), 1 )/ saires
128  pbaru(ip1jm+1,l)= - pbarv(ip1jmi1+1,l) + cts * aire(ip1jm+1)
129 
130  DO i = 2,iim
131  pbaru(i+ ip1jm,l) = pbaru(i+ip1jm-1,l) -
132  * pbarv(i+ip1jmi1,l)+cts*aire(i+ip1jm)
133  ENDDO
134 
135  DO i = 1,iim
136  apbarus(i) = aireu(i +ip1jm) * pbaru(i +ip1jm, l)
137  ENDDO
138 
139  cts0 = -ssum( iim,apbarus,1 )/saireus
140 
141  DO i = 1,iim
142  pbaru(i+ ip1jm, l) = 2. * ( pbaru(i +ip1jm, l) + cts0 )
143  ENDDO
144 
145  pbaru( ip1jmp1,l ) = pbaru( ip1jm +1,l )
146 
147  ENDDO
148 c$OMP END DO NOWAIT
149  ENDIF
150 
151  RETURN
152  END
!$Header llmm1 INTEGER ip1jmi1
Definition: paramet.h:14
!$Header llmm1 INTEGER ip1jmp1
Definition: paramet.h:14
integer, save ij_end
logical, save pole_sud
!$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
!$Header!CDK comgeom COMMON comgeom aire
Definition: comgeom.h:25
!$Header llmm1 INTEGER ip1jm
Definition: paramet.h:14
!$Header!CDK comgeom COMMON comgeom aireu
Definition: comgeom.h:25
logical, save pole_nord
subroutine flumass_p(massebx, masseby, vcont, ucont, pbaru, pbarv)
Definition: flumass_p.F:2
integer, save ij_begin
c c zjulian c cym CALL iim cym klev iim
Definition: ini_bilKP_ave.h:24
real function ssum(n, sx, incx)
Definition: cray.F:27