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

Contents of /trunk/dyn3d/Guide/Read_reanalyse/reanalyse2nat.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 172 - (show annotations)
Wed Sep 30 15:59:14 2015 UTC (8 years, 8 months ago) by guez
Original Path: trunk/Sources/dyn3d/Guide/Read_reanalyse/reanalyse2nat.f
File size: 3342 byte(s)
Just indented correctbid and nat2gcm.

The procedure read_reanalyse just reads the next time slab every time
it is called. No use keeping track of the time index in the calling
procedure, guide. It is simpler to do it in read_reanalyse. Also
simpler to read the number of vertical levels in read_reanalyse than
in guide, since we have already in read_reanalyse the input of
pressure levels. We then have to make the arrays containing reanalyses
static allocatable instead of automatic. Also only read pressure
levels at the first call of read_reanalyse instead of at every call.

masserea2 not used in guide. Remove it down the chain in
read_reanalyse and reanalyse2nat.

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

  ViewVC Help
Powered by ViewVC 1.1.21