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

Contents of /trunk/dyn3d/Guide/Read_reanalyse/nat2gcm.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/nat2gcm.f
File size: 1982 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 subroutine nat2gcm(u,v,t,rh,pk,ucov,vcov,teta,q)
2
3 ! Passage aux variables du modele (vents covariants, temperature
4 ! potentielle et humidite specifique)
5
6 use dimens_m
7 use paramet_m
8 use comconst
9 use disvert_m
10 use comgeom
11 use q_sat_m, only: q_sat
12 use guide_m
13 implicit none
14
15
16 real u(iip1,jjp1,llm),v(iip1,jjm,llm)
17 real t(iip1,jjp1,llm),pk(iip1,jjp1,llm),rh(iip1,jjp1,llm)
18 real ucov(iip1,jjp1,llm),vcov(iip1,jjm,llm)
19 real teta(iip1,jjp1,llm),q(iip1,jjp1,llm)
20
21 real pres(iip1,jjp1,llm),qsat(iip1,jjp1,llm)
22
23 real unskap
24
25 integer i,j,l
26
27
28 print*,'Entree dans nat2gcm'
29 ! ucov(:,:,:)=0.
30 ! do l=1,llm
31 ! ucov(:,2:jjm,l)=u(:,2:jjm,l)*cu_2d(:,2:jjm)
32 ! enddo
33 ! ucov(iip1,:,:)=ucov(1,:,:)
34
35 ! teta(:,:,:)=t(:,:,:)*cpp/pk(:,:,:)
36 ! teta(iip1,:,:)=teta(1,:,:)
37
38 ! calcul de ucov et de la temperature potentielle
39 do l=1,llm
40 do j=1,jjp1
41 do i=1,iim
42 ucov(i,j,l)=u(i,j,l)*cu_2d(i,j)
43 teta(i,j,l)=t(i,j,l)*cpp/pk(i,j,l)
44 enddo
45 ucov(iip1,j,l)=ucov(1,j,l)
46 teta(iip1,j,l)=teta(1,j,l)
47 enddo
48 do i=1,iip1
49 ucov(i,1,l)=0.
50 ucov(i,jjp1,l)=0.
51 teta(i,1,l)=teta(1,1,l)
52 teta(i,jjp1,l)=teta(1,jjp1,l)
53 enddo
54 enddo
55
56 ! calcul de ucov
57 do l=1,llm
58 do j=1,jjm
59 do i=1,iim
60 vcov(i,j,l)=v(i,j,l)*cv_2d(i,j)
61 enddo
62 vcov(iip1,j,l)=vcov(1,j,l)
63 enddo
64 enddo
65
66 ! Humidite relative -> specifique
67 ! -------------------------------
68 if (1.eq.0) then
69 ! FINALEMENT ON GUIDE EN HUMIDITE RELATIVE
70 print*,'calcul de unskap'
71 unskap = 1./ kappa
72 print*,'calcul de pres'
73 pres(:,:,:)=preff*(pk(:,:,:)/cpp)**unskap
74 print*,'calcul de qsat'
75 qsat = q_sat(t, pres)
76 print*,'calcul de q'
77 ! ATTENTION : humidites relatives en %
78 rh(:,:,:)=max(rh(:,:,:)*0.01,1.e-6)
79 q(:,:,:)=qsat(:,:,:)*rh(:,:,:)
80 print*,'calcul de q OK'
81
82 endif
83
84 end subroutine nat2gcm

  ViewVC Help
Powered by ViewVC 1.1.21