/[lmdze]/trunk/libf/dyn3d/Read_reanalyse/reanalyse2nat.f90
ViewVC logotype

Contents of /trunk/libf/dyn3d/Read_reanalyse/reanalyse2nat.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 44 - (show annotations)
Wed Apr 13 12:29:18 2011 UTC (13 years ago) by guez
File size: 4036 byte(s)
Removed argument "pdteta" of "calfis", because it was not used.

Created module "conf_guide_m", containing procedure
"conf_guide". Moved module variables from "guide_m" to "conf_guide_m".

In module "getparam", removed "ini_getparam" and "fin_getparam" from
generic interface "getpar".

Created module variables in "tau2alpha_m" to replace common "comdxdy".

1
2
3 !===========================================================================
4 subroutine reanalyse2nat(nlevnc,psi &
5 ,unc,vnc,tnc,qnc,psnc,pl,u,v,t,q &
6 ,ps,masse,pk)
7 !===========================================================================
8
9 ! -----------------------------------------------------------------
10 ! Inversion Nord/sud de la grille + interpollation sur les niveaux
11 ! verticaux du modele.
12 ! -----------------------------------------------------------------
13
14 use dimens_m
15 use paramet_m
16 use comconst
17 use comvert
18 use comgeom
19 use exner_hyb_m, only: exner_hyb
20 use conf_guide_m
21
22 implicit none
23
24
25 integer nlevnc
26 real psi(iip1,jjp1)
27 real u(iip1,jjp1,llm),v(iip1,jjm,llm)
28 real t(iip1,jjp1,llm),ps(iip1,jjp1),q(iip1,jjp1,llm)
29
30 real pl(nlevnc)
31 real unc(iip1,jjp1,nlevnc),vnc(iip1,jjm,nlevnc)
32 real tnc(iip1,jjp1,nlevnc),psnc(iip1,jjp1)
33 real qnc(iip1,jjp1,nlevnc)
34
35 real zu(iip1,jjp1,llm),zv(iip1,jjm,llm)
36 real zt(iip1,jjp1,llm),zq(iip1,jjp1,llm)
37
38 real pext(iip1,jjp1,llm)
39 real pbarx(iip1,jjp1,llm),pbary(iip1,jjm,llm)
40 real plunc(iip1,jjp1,llm),plvnc(iip1,jjm,llm)
41 real plsnc(iip1,jjp1,llm)
42
43 real p(iip1,jjp1,llmp1),pk(iip1,jjp1,llm),pks(iip1,jjp1)
44 real pkf(iip1,jjp1,llm)
45 real masse(iip1,jjp1,llm),pls(iip1,jjp1,llm)
46 real prefkap,unskap
47
48
49 integer i,j,l
50
51
52 ! -----------------------------------------------------------------
53 ! calcul de la pression au milieu des couches.
54 ! -----------------------------------------------------------------
55
56 forall (l = 1: llm + 1) p(:, :, l) = ap(l) + bp(l) * psi
57 call massdair(p,masse)
58 CALL exner_hyb(psi,p,pks,pk,pkf)
59
60 ! .... Calcul de pls , pression au milieu des couches ,en Pascals
61 unskap=1./kappa
62 prefkap = preff ** kappa
63 ! PRINT *,' Pref kappa unskap ',preff,kappa,unskap
64 DO l = 1, llm
65 DO j=1,jjp1
66 DO i =1, iip1
67 pls(i,j,l) = preff * ( pk(i,j,l)/cpp) ** unskap
68 ENDDO
69 ENDDO
70 ENDDO
71
72
73 ! -----------------------------------------------------------------
74 ! calcul des pressions pour les grilles u et v
75 ! -----------------------------------------------------------------
76
77 do l=1,llm
78 do j=1,jjp1
79 do i=1,iip1
80 pext(i,j,l)=pls(i,j,l)*aire_2d(i,j)
81 enddo
82 enddo
83 enddo
84 call massbar(pext, pbarx, pbary )
85 do l=1,llm
86 do j=1,jjp1
87 do i=1,iip1
88 plunc(i,jjp1+1-j,l)=pbarx(i,j,l)/aireu_2d(i,j)
89 plsnc(i,jjp1+1-j,l)=pls(i,j,l)
90 enddo
91 enddo
92 enddo
93 do l=1,llm
94 do j=1,jjm
95 do i=1,iip1
96 plvnc(i,jjm+1-j,l)=pbary(i,j,l)/airev_2d(i,j)
97 enddo
98 enddo
99 enddo
100
101 ! -----------------------------------------------------------------
102
103 if (guide_P) then
104 do j=1,jjp1
105 do i=1,iim
106 ps(i,j)=psnc(i,jjp1+1-j)
107 enddo
108 ps(iip1,j)=ps(1,j)
109 enddo
110 endif
111
112
113 ! -----------------------------------------------------------------
114 call pres2lev(unc,zu,nlevnc,llm,pl,plunc,iip1,jjp1)
115 call pres2lev(vnc,zv,nlevnc,llm,pl,plvnc,iip1,jjm )
116 call pres2lev(tnc,zt,nlevnc,llm,pl,plsnc,iip1,jjp1)
117 call pres2lev(qnc,zq,nlevnc,llm,pl,plsnc,iip1,jjp1)
118
119 ! call dump2d(iip1,jjp1,ps,'PS ')
120 ! call dump2d(iip1,jjp1,psu,'PS ')
121 ! call dump2d(iip1,jjm,psv,'PS ')
122 ! Inversion Nord/Sud
123 do l=1,llm
124 do j=1,jjp1
125 do i=1,iim
126 u(i,j,l)=zu(i,jjp1+1-j,l)
127 t(i,j,l)=zt(i,jjp1+1-j,l)
128 q(i,j,l)=zq(i,jjp1+1-j,l)
129 enddo
130 u(iip1,j,l)=u(1,j,l)
131 t(iip1,j,l)=t(1,j,l)
132 q(iip1,j,l)=q(1,j,l)
133 enddo
134 enddo
135
136 do l=1,llm
137 do j=1,jjm
138 do i=1,iim
139 v(i,j,l)=zv(i,jjm+1-j,l)
140 enddo
141 v(iip1,j,l)=v(1,j,l)
142 enddo
143 enddo
144
145 return
146 end

  ViewVC Help
Powered by ViewVC 1.1.21