/[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 102 - (show annotations)
Tue Jul 15 13:43:24 2014 UTC (9 years, 10 months ago) by guez
Original Path: trunk/dyn3d/Read_reanalyse/nat2gcm.f
File size: 2457 byte(s)
Removed unused file "condsurf.f" (only useful for ocean slab).

day_step must be a multiple of 4 * iperiod if ok_guide.

Changed type of variable online of module conf_guide_m from integer to
logical. Value -1 was not useful, equivalent to not ok_guide.

Removed argument masse of procedure guide. masse is kept consistent
with ps throughout the run. masse need only be computed again just
after ps has been modified. In prodecure guide, replaced use of
remanent variable first by test on itau. Replaced test on variable
"test" by test on integer values.

In leapfrog, for the call to guide, replaced test on real values by
test on integer values.

Bug fix in tau2alpha: computation of dxdyv (following LMDZ revision 1040).

In procedure wrgrads, replaced badly chosen argument name "if" by i_f.

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

  ViewVC Help
Powered by ViewVC 1.1.21