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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 173 - (show annotations)
Tue Oct 6 15:57:02 2015 UTC (8 years, 7 months ago) by guez
File size: 3144 byte(s)
correctbid did nothing. (Not used either in LMDZ since revision 1170.)

Avoid aliasing in arguments of nat2gcm: use a single set of arguments
with intent inout. Argument q of nat2gcm was not used.

pres2lev now accepts po in any monotonic order. So the input files for
nudging can now have the pressure coordinate in any order. Also, we
read the latitude coordinate from the input files for nudging and we
invert order if necessary so the input files for nudging can now have
the latitude coordinate in any order.

In pre2lev, no need for lmomx: use automatic arrays.

Removed variable ncep of module conf_guide_m. Instead, we find out
what the pressure coordinate is with find_coord.


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

  ViewVC Help
Powered by ViewVC 1.1.21