/[lmdze]/trunk/dyn3d/Read_reanalyse/read_reanalyse.f
ViewVC logotype

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

  ViewVC Help
Powered by ViewVC 1.1.21