/[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 114 - (hide annotations)
Fri Sep 19 11:41:35 2014 UTC (9 years, 8 months ago) by guez
Original Path: trunk/dyn3d/Guide/Read_reanalyse/reanalyse2nat.f
File size: 3360 byte(s)


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 guez 91 use massbar_m, only: massbar
20 guez 88 use massdair_m, only: massdair
21 guez 37
22 guez 88 integer nlevnc
23     real, intent(in):: psi(iip1, jjp1)
24     real unc(iip1, jjp1, nlevnc), vnc(iip1, jjm, nlevnc)
25     real tnc(iip1, jjp1, nlevnc)
26     real qnc(iip1, jjp1, nlevnc)
27     real pl(nlevnc)
28     real u(iip1, jjp1, llm), v(iip1, jjm, llm)
29     real t(iip1, jjp1, llm), q(iip1, jjp1, llm)
30     real masse(iip1, jjp1, llm)
31     real pk(iip1, jjp1, llm)
32 guez 37
33 guez 88 ! Local:
34 guez 37
35 guez 88 real zu(iip1, jjp1, llm), zv(iip1, jjm, llm)
36     real zt(iip1, jjp1, llm), zq(iip1, jjp1, llm)
37 guez 37
38 guez 88 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 guez 37
43 guez 88 real p(iip1, jjp1, llmp1)
44     real pks(iip1, jjp1)
45     real pkf(iip1, jjp1, llm)
46     real pls(iip1, jjp1, llm)
47     real prefkap, unskap
48 guez 37
49 guez 88 integer i, j, l
50 guez 37
51 guez 88 ! -----------------------------------------------------------------
52 guez 37
53 guez 88 ! calcul de la pression au milieu des couches
54     forall (l = 1: llm + 1) p(:, :, l) = ap(l) + bp(l) * psi
55     call massdair(p, masse)
56     CALL exner_hyb(psi, p, pks, pk, pkf)
57 guez 37
58 guez 88 ! Calcul de pls, pression au milieu des couches, en Pascals
59     unskap=1./kappa
60     prefkap = preff ** kappa
61     DO l = 1, llm
62     DO j=1, jjp1
63     DO i =1, iip1
64     pls(i, j, l) = preff * ( pk(i, j, l)/cpp) ** unskap
65     ENDDO
66     ENDDO
67     ENDDO
68 guez 37
69 guez 88 ! calcul des pressions pour les grilles u et v
70 guez 37
71 guez 88 do l=1, llm
72     do j=1, jjp1
73     do i=1, iip1
74     pext(i, j, l)=pls(i, j, l)*aire_2d(i, j)
75     enddo
76     enddo
77     enddo
78     call massbar(pext, pbarx, pbary )
79     do l=1, llm
80     do j=1, jjp1
81     do i=1, iip1
82     plunc(i, jjp1+1-j, l)=pbarx(i, j, l)/aireu_2d(i, j)
83     plsnc(i, jjp1+1-j, l)=pls(i, j, l)
84     enddo
85     enddo
86     enddo
87     do l=1, llm
88     do j=1, jjm
89     do i=1, iip1
90     plvnc(i, jjm+1-j, l)=pbary(i, j, l)/airev_2d(i, j)
91     enddo
92     enddo
93     enddo
94 guez 37
95 guez 88 call pres2lev(unc, zu, nlevnc, llm, pl, plunc, iip1, jjp1)
96     call pres2lev(vnc, zv, nlevnc, llm, pl, plvnc, iip1, jjm )
97     call pres2lev(tnc, zt, nlevnc, llm, pl, plsnc, iip1, jjp1)
98     call pres2lev(qnc, zq, nlevnc, llm, pl, plsnc, iip1, jjp1)
99 guez 37
100 guez 88 ! Inversion Nord/Sud
101     do l=1, llm
102     do j=1, jjp1
103     do i=1, iim
104     u(i, j, l)=zu(i, jjp1+1-j, l)
105     t(i, j, l)=zt(i, jjp1+1-j, l)
106     q(i, j, l)=zq(i, jjp1+1-j, l)
107     enddo
108     u(iip1, j, l)=u(1, j, l)
109     t(iip1, j, l)=t(1, j, l)
110     q(iip1, j, l)=q(1, j, l)
111     enddo
112     enddo
113 guez 37
114 guez 88 do l=1, llm
115     do j=1, jjm
116     do i=1, iim
117     v(i, j, l)=zv(i, jjm+1-j, l)
118     enddo
119     v(iip1, j, l)=v(1, j, l)
120     enddo
121     enddo
122 guez 37
123 guez 88 end subroutine reanalyse2nat
124 guez 37
125 guez 88 end module reanalyse2nat_m

  ViewVC Help
Powered by ViewVC 1.1.21