/[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 85 - (hide annotations)
Thu Mar 6 17:35:22 2014 UTC (10 years, 2 months ago) by guez
Original Path: trunk/dyn3d/Read_reanalyse/reanalyse2nat.f
File size: 3210 byte(s)
Removed option to guide surface pressure because it was not
functional: psrea1 was not defined in procedure guide. Removed local
variables psrea1 and psrea2 of procedure guide. ps becomes an
"intent(in)" argument in guide. Removed case guide_p in guide. Removed
variable guide_p of module conf_guide_m. Removed case guide_p and
argument ps in read_reanalyse. Removed case guide_p and argument ps in
reanalyse2nat.

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

  ViewVC Help
Powered by ViewVC 1.1.21