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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 265 - (hide annotations)
Tue Mar 20 09:35:59 2018 UTC (6 years, 2 months ago) by guez
File size: 3068 byte(s)
Rename module dimens_m to dimensions.
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 173 subroutine reanalyse2nat(invert_y, psi, unc, vnc, tnc, qnc, pl, u, v, t, q, &
8     pk)
9 guez 37
10 guez 172 ! Inversion nord-sud de la grille et interpolation verticale sur
11     ! les niveaux du modèle.
12 guez 37
13 guez 173 USE comconst, ONLY: cpp, kappa
14     USE comgeom, ONLY: aireu_2d, airev_2d, aire_2d
15 guez 265 USE dimensions, ONLY: jjm, llm
16 guez 88 USE disvert_m, ONLY: ap, bp, preff
17     USE exner_hyb_m, ONLY: exner_hyb
18 guez 91 use massbar_m, only: massbar
19 guez 173 USE paramet_m, ONLY: iip1, jjp1, llmp1
20 guez 171 use pres2lev_m, only: pres2lev
21 guez 37
22 guez 173 logical, intent(in):: invert_y
23     real, intent(in):: psi(:, :) ! (iip1, jjp1)
24 guez 37
25 guez 173 real, intent(in):: unc(:, :, :) ! (iip1, jjp1, :)
26     real, intent(in):: vnc(:, :, :) ! (iip1, jjm, :)
27     real, intent(in):: tnc(:, :, :) ! (iip1, jjp1, :)
28     real, intent(in):: qnc(:, :, :) ! (iip1, jjp1, :)
29     real, intent(in):: pl(:)
30    
31     real, intent(out):: u(:, :, :) ! (iip1, jjp1, llm)
32     real, intent(out):: v(:, :, :) ! (iip1, jjm, llm)
33     real, intent(out):: t(:, :, :), q(:, :, :) ! (iip1, jjp1, llm)
34     real, intent(out):: pk(:, :, :) ! (iip1, jjp1, llm)
35    
36 guez 88 ! Local:
37 guez 37
38 guez 88 real zu(iip1, jjp1, llm), zv(iip1, jjm, llm)
39     real zt(iip1, jjp1, llm), zq(iip1, jjp1, llm)
40 guez 37
41 guez 88 real pext(iip1, jjp1, llm)
42     real pbarx(iip1, jjp1, llm), pbary(iip1, jjm, llm)
43     real plunc(iip1, jjp1, llm), plvnc(iip1, jjm, llm)
44     real plsnc(iip1, jjp1, llm)
45 guez 37
46 guez 88 real p(iip1, jjp1, llmp1)
47     real pks(iip1, jjp1)
48     real pls(iip1, jjp1, llm)
49 guez 178 real unskap
50 guez 37
51 guez 88 integer i, j, l
52 guez 37
53 guez 88 ! -----------------------------------------------------------------
54 guez 37
55 guez 88 ! calcul de la pression au milieu des couches
56     forall (l = 1: llm + 1) p(:, :, l) = ap(l) + bp(l) * psi
57 guez 137 CALL exner_hyb(psi, p, pks, pk)
58 guez 37
59 guez 88 ! Calcul de pls, pression au milieu des couches, en Pascals
60     unskap=1./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 173 call pres2lev(unc, zu, pl, plunc)
96     call pres2lev(vnc, zv, pl, plvnc )
97     call pres2lev(tnc, zt, pl, plsnc)
98     call pres2lev(qnc, zq, pl, plsnc)
99 guez 37
100 guez 173 if (invert_y) then
101     ! Inversion Nord/Sud
102     u=zu(:, jjp1:1:-1, :)
103     v=zv(:, jjm:1:-1, :)
104     t=zt(:, jjp1:1:-1, :)
105     q=zq(:, jjp1:1:-1, :)
106     else
107     u = zu
108     v = zv
109     t = zt
110     q = zq
111     end if
112 guez 37
113 guez 88 end subroutine reanalyse2nat
114 guez 37
115 guez 88 end module reanalyse2nat_m

  ViewVC Help
Powered by ViewVC 1.1.21