/[lmdze]/trunk/libf/dyn3d/Read_reanalyse/read_reanalyse.f90
ViewVC logotype

Contents of /trunk/libf/dyn3d/Read_reanalyse/read_reanalyse.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 44 - (show annotations)
Wed Apr 13 12:29:18 2011 UTC (13 years, 1 month ago) by guez
File size: 6955 byte(s)
Removed argument "pdteta" of "calfis", because it was not used.

Created module "conf_guide_m", containing procedure
"conf_guide". Moved module variables from "guide_m" to "conf_guide_m".

In module "getparam", removed "ini_getparam" and "fin_getparam" from
generic interface "getpar".

Created module variables in "tau2alpha_m" to replace common "comdxdy".

1 subroutine read_reanalyse(timestep,psi &
2 ,u,v,t,q,masse,ps,mode,nlevnc)
3
4 !
5 ! $Header: /home/cvsroot/LMDZ4/libf/dyn3d/read_reanalyse.F,v 1.3 2005/04/15 12:31:21 lmdzadmin Exp $
6 !
7 !
8 !
9 ! mode=0 variables naturelles
10 ! mode=1 variabels GCM
11
12 ! -----------------------------------------------------------------
13 ! Declarations
14 ! -----------------------------------------------------------------
15 use dimens_m
16 use paramet_m
17 use comvert
18 use comgeom
19 use conf_guide_m
20 use netcdf
21
22 IMPLICIT NONE
23
24 ! common
25 ! ------
26
27 include "netcdf.inc"
28
29
30 ! arguments
31 ! ---------
32 integer nlevnc
33 integer timestep,mode,l
34
35 real psi(iip1,jjp1)
36 real u(iip1,jjp1,llm),v(iip1,jjm,llm)
37 real t(iip1,jjp1,llm),ps(iip1,jjp1),q(iip1,jjp1,llm)
38 real masse(iip1,jjp1,llm),pk(iip1,jjp1,llm)
39
40
41 ! local
42 ! -----
43 integer ncidu,varidu,ncidv,varidv,ncidt,varidt,ncidps,varidps
44 integer ncidpl
45 integer varidpl,ncidQ,varidQ
46 save ncidu,varidu,ncidv,varidv,ncidt,varidt,ncidps,varidps
47 save ncidpl
48 save varidpl,ncidQ,varidQ
49
50 real unc(iip1,jjp1,nlevnc),vnc(iip1,jjm,nlevnc)
51 real tnc(iip1,jjp1,nlevnc),psnc(iip1,jjp1)
52 real Qnc(iip1,jjp1,nlevnc)
53 real pl(nlevnc)
54
55 integer start(4),count(4),status
56
57 real rcode
58 logical first
59 save first
60
61 data first/.true./
62
63
64
65 ! -----------------------------------------------------------------
66 ! Initialisation de la lecture des fichiers
67 ! -----------------------------------------------------------------
68 if (first) then
69 ncidpl=-99
70 print*,'Intitialisation de read reanalsye'
71
72 ! Vent zonal
73 if (guide_u) then
74 rcode=nf90_open('u.nc',nf90_nowrite,ncidu)
75 rcode = nf90_inq_varid(ncidu, 'UWND', varidu)
76 print*,'ncidu,varidu',ncidu,varidu
77 if (ncidpl.eq.-99) ncidpl=ncidu
78 endif
79
80 ! Vent meridien
81 if (guide_v) then
82 rcode=nf90_open('v.nc',nf90_nowrite,ncidv)
83 rcode = nf90_inq_varid(ncidv, 'VWND', varidv)
84 print*,'ncidv,varidv',ncidv,varidv
85 if (ncidpl.eq.-99) ncidpl=ncidv
86 endif
87
88 ! Temperature
89 if (guide_T) then
90 rcode=nf90_open('T.nc',nf90_nowrite,ncidt)
91 rcode = nf90_inq_varid(ncidt, 'AIR', varidt)
92 print*,'ncidt,varidt',ncidt,varidt
93 if (ncidpl.eq.-99) ncidpl=ncidt
94 endif
95
96 ! Humidite
97 if (guide_Q) then
98 rcode=nf90_open('hur.nc',nf90_nowrite,ncidQ)
99 rcode = nf90_inq_varid(ncidQ, 'RH', varidQ)
100 print*,'ncidQ,varidQ',ncidQ,varidQ
101 if (ncidpl.eq.-99) ncidpl=ncidQ
102 endif
103
104 ! Pression de surface
105 if (guide_P) then
106 rcode=nf90_open('ps.nc',nf90_nowrite,ncidps)
107 rcode = nf90_inq_varid(ncidps, 'SP', varidps)
108 print*,'ncidps,varidps',ncidps,varidps
109 endif
110
111 ! Coordonnee verticale
112 if (ncep) then
113 print*,'Vous etes entrain de lire des donnees NCEP'
114 rcode = nf90_inq_varid(ncidpl, 'LEVEL', varidpl)
115 else
116 print*,'Vous etes entrain de lire des donnees ECMWF'
117 rcode = nf90_inq_varid(ncidpl, 'PRESSURE', varidpl)
118 endif
119 print*,'ncidu,varidpl',ncidu,varidpl
120 endif
121 print*,'ok1'
122
123 ! Niveaux de pression
124 print*,'WARNING!!! Il n y a pas de test de coherence'
125 print*,'sur le nombre de niveaux verticaux dans le fichier nc'
126 status=NF_GET_VARA_REAL(ncidpl,varidpl,1,nlevnc,pl)
127 ! passage en pascal
128 pl(:)=100.*pl(:)
129 if (first) then
130 do l=1,nlevnc
131 print*,'PL(',l,')=',pl(l)
132 enddo
133 endif
134
135 ! -----------------------------------------------------------------
136 ! lecture des champs u, v, T, ps
137 ! -----------------------------------------------------------------
138
139 ! dimensions pour les champs scalaires et le vent zonal
140 ! -----------------------------------------------------
141
142 start(1)=1
143 start(2)=1
144 start(3)=1
145 start(4)=timestep
146
147 count(1)=iip1
148 count(2)=jjp1
149 count(3)=nlevnc
150 count(4)=1
151
152 ! mise a zero des tableaux
153 ! ------------------------
154 unc(:,:,:)=0.
155 vnc(:,:,:)=0.
156 tnc(:,:,:)=0.
157 Qnc(:,:,:)=0.
158
159 ! Vent zonal
160 ! ----------
161
162 if (guide_u) then
163 print*,'avant la lecture de UNCEP nd de niv:',nlevnc
164 status=NF_GET_VARA_REAL(ncidu,varidu,start,count,unc)
165 ! call dump2d(iip1,jjp1,unc,'VENT NCEP ')
166 ! call dump2d(iip1,40,unc(1,1,nlevnc),'VENT NCEP ')
167 print*,'WARNING!!! Correction bidon pour palier a un '
168 print*,'probleme dans la creation des fichiers nc'
169 call correctbid(iim,jjp1*nlevnc,unc)
170 call dump2d(iip1,jjp1,unc,'UNC COUCHE 1 ')
171 endif
172
173 ! Temperature
174 ! -----------
175
176 print*,'ncidt=',ncidt,'varidt=',varidt,'start=',start
177 print*,'count=',count
178 if (guide_T) then
179 status=NF_GET_VARA_REAL(ncidt,varidt,start,count,tnc)
180 call dump2d(iip1,jjp1,tnc,'TNC COUCHE 1 AAA ')
181 call correctbid(iim,jjp1*nlevnc,tnc)
182 call dump2d(iip1,jjp1,tnc,'TNC COUCHE 1 BBB ')
183 endif
184
185 ! Humidite
186 ! --------
187
188 if (guide_Q) then
189 status=NF_GET_VARA_REAL(ncidQ,varidQ,start,count,Qnc)
190 call correctbid(iim,jjp1*nlevnc,Qnc)
191 call dump2d(iip1,jjp1,Qnc,'QNC COUCHE 1 ')
192 endif
193
194 count(2)=jjm
195 ! Vent meridien
196 ! -------------
197
198 if (guide_v) then
199 status=NF_GET_VARA_REAL(ncidv,varidv,start,count,vnc)
200 call correctbid(iim,jjm*nlevnc,vnc)
201 call dump2d(iip1,jjm,vnc,'VNC COUCHE 1 ')
202 endif
203
204 start(3)=timestep
205 start(4)=0
206 count(2)=jjp1
207 count(3)=1
208 count(4)=0
209
210 ! Pression de surface
211 ! -------------------
212
213 if (guide_P) then
214 status=NF_GET_VARA_REAL(ncidps,varidps,start,count,psnc)
215 call dump2d(iip1,jjp1,psnc,'PSNC COUCHE 1 ')
216 call correctbid(iim,jjp1,psnc)
217 endif
218
219
220
221 ! -----------------------------------------------------------------
222 ! Interpollation verticale sur les niveaux modele
223 ! -----------------------------------------------------------------
224 call reanalyse2nat(nlevnc,psi,unc,vnc,tnc,Qnc,psnc,pl,u,v,t,Q &
225 ,ps,masse,pk)
226
227 call dump2d(iip1,jjm,v,'V COUCHE APRES ')
228
229
230 ! -----------------------------------------------------------------
231 ! Passage aux variables du modele (vents covariants, temperature
232 ! potentielle et humidite specifique)
233 ! -----------------------------------------------------------------
234 call nat2gcm(u,v,t,Q,pk,u,v,t,Q)
235 print*,'TIMESTEP ',timestep
236 if(mode.ne.1) stop'mode pas egal 0'
237 ! call dump2d(iip1,jjm,v,'VCOV COUCHE 1 ')
238
239 ! Lignes introduites a une epoque pour un probleme oublie...
240 ! do l=1,llm
241 ! do i=1,iip1
242 ! v(i,1,l)=0.
243 ! v(i,jjm,l)=0.
244 ! enddo
245 ! enddo
246 first=.false.
247
248 return
249 end

  ViewVC Help
Powered by ViewVC 1.1.21