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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (hide annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years, 1 month ago) by guez
File size: 6893 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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

  ViewVC Help
Powered by ViewVC 1.1.21