/[lmdze]/trunk/Sources/dyn3d/Guide/Read_reanalyse/reanalyse2nat.f
ViewVC logotype

Annotation of /trunk/Sources/dyn3d/Guide/Read_reanalyse/reanalyse2nat.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 88 - (hide annotations)
Tue Mar 11 15:09:02 2014 UTC (10 years, 2 months ago) by guez
Original Path: trunk/dyn3d/Read_reanalyse/reanalyse2nat.f
File size: 3327 byte(s)
Removed useless argument mode of subroutine read_reanalyse.

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

  ViewVC Help
Powered by ViewVC 1.1.21