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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 66 - (hide annotations)
Thu Sep 20 13:00:41 2012 UTC (11 years, 8 months ago) by guez
Original Path: trunk/libf/dyn3d/Read_reanalyse/nat2gcm.f90
File size: 2450 byte(s)
Changed name of module "comvert" to "disvert_m". Changed constant
1. to 0.3 in vertical sampling "strato".

1 guez 37
2     !===========================================================================
3     subroutine nat2gcm(u,v,t,rh,pk,ucov,vcov,teta,q)
4     !===========================================================================
5    
6     use dimens_m
7     use paramet_m
8     use comconst
9 guez 66 use disvert_m
10 guez 37 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 ps(iip1,jjp1)
19    
20     real ucov(iip1,jjp1,llm),vcov(iip1,jjm,llm)
21     real teta(iip1,jjp1,llm),q(iip1,jjp1,llm)
22    
23     real pres(iip1,jjp1,llm),qsat(iip1,jjp1,llm)
24    
25     real unskap
26    
27     integer i,j,l
28    
29    
30     print*,'Entree dans nat2gcm'
31     ! ucov(:,:,:)=0.
32     ! do l=1,llm
33     ! ucov(:,2:jjm,l)=u(:,2:jjm,l)*cu_2d(:,2:jjm)
34     ! enddo
35     ! ucov(iip1,:,:)=ucov(1,:,:)
36    
37     ! teta(:,:,:)=t(:,:,:)*cpp/pk(:,:,:)
38     ! teta(iip1,:,:)=teta(1,:,:)
39    
40     ! calcul de ucov et de la temperature potentielle
41     do l=1,llm
42     do j=1,jjp1
43     do i=1,iim
44     ucov(i,j,l)=u(i,j,l)*cu_2d(i,j)
45     teta(i,j,l)=t(i,j,l)*cpp/pk(i,j,l)
46     enddo
47     ucov(iip1,j,l)=ucov(1,j,l)
48     teta(iip1,j,l)=teta(1,j,l)
49     enddo
50     do i=1,iip1
51     ucov(i,1,l)=0.
52     ucov(i,jjp1,l)=0.
53     teta(i,1,l)=teta(1,1,l)
54     teta(i,jjp1,l)=teta(1,jjp1,l)
55     enddo
56     enddo
57    
58     ! calcul de ucov
59     do l=1,llm
60     do j=1,jjm
61     do i=1,iim
62     vcov(i,j,l)=v(i,j,l)*cv_2d(i,j)
63     enddo
64     vcov(iip1,j,l)=vcov(1,j,l)
65     enddo
66     enddo
67    
68     ! call dump2d(iip1,jjp1,teta,'TETA EN BAS ')
69     ! call dump2d(iip1,jjp1,teta(1,1,llm),'TETA EN HAUT ')
70    
71     ! Humidite relative -> specifique
72     ! -------------------------------
73     if (1.eq.0) then
74     ! FINALEMENT ON GUIDE EN HUMIDITE RELATIVE
75     print*,'calcul de unskap'
76     unskap = 1./ kappa
77     print*,'calcul de pres'
78     pres(:,:,:)=preff*(pk(:,:,:)/cpp)**unskap
79     print*,'calcul de qsat'
80     qsat = q_sat(t, pres)
81     print*,'calcul de q'
82     ! ATTENTION : humidites relatives en %
83     rh(:,:,:)=max(rh(:,:,:)*0.01,1.e-6)
84     q(:,:,:)=qsat(:,:,:)*rh(:,:,:)
85     print*,'calcul de q OK'
86    
87     call dump2d(iip1,jjp1,pres,'PRESSION PREMIERE COUCHE ')
88     call dump2d(iip1,jjp1,q,'HUMIDITE SPECIFIQUE COUCHE 1 ')
89     endif
90    
91    
92     return
93     end

  ViewVC Help
Powered by ViewVC 1.1.21