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

Contents of /trunk/dyn3d/Read_reanalyse/read_reanalyse.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 88 - (show annotations)
Tue Mar 11 15:09:02 2014 UTC (10 years, 2 months ago) by guez
File size: 5234 byte(s)
Removed useless argument mode of subroutine read_reanalyse.

1 module read_reanalyse_m
2
3 IMPLICIT NONE
4
5 contains
6
7 subroutine read_reanalyse(timestep, psi, u, v, t, q, masse, nlevnc)
8
9 ! From LMDZ4/libf/dyn3d/read_reanalyse.F, version 1.3, 2005/04/15 12:31:21
10
11 USE conf_guide_m, ONLY: guide_q, guide_t, guide_u, guide_v, ncep
12 USE dimens_m, ONLY: iim, jjm, llm
13 USE netcdf, ONLY: nf90_get_var, nf90_inq_varid, nf90_nowrite, nf90_open
14 USE paramet_m, ONLY: iip1, jjp1
15 use reanalyse2nat_m, only: reanalyse2nat
16
17 integer timestep
18 real, intent(in):: psi(iip1, jjp1)
19 real u(iip1, jjp1, llm), v(iip1, jjm, llm)
20 real t(iip1, jjp1, llm), q(iip1, jjp1, llm)
21 real masse(iip1, jjp1, llm)
22 integer nlevnc
23
24 ! Local:
25
26 integer l
27 real pk(iip1, jjp1, llm)
28 integer, save:: ncidu, varidu, ncidv, varidv, ncidt, varidt
29 integer, save:: ncidpl
30 integer, save:: varidpl, ncidQ, varidQ
31 real unc(iip1, jjp1, nlevnc), vnc(iip1, jjm, nlevnc)
32 real tnc(iip1, jjp1, nlevnc)
33 real Qnc(iip1, jjp1, nlevnc)
34 real pl(nlevnc)
35 integer start(4), count(4), status
36 real rcode
37 logical:: first = .true.
38
39 ! -----------------------------------------------------------------
40
41 ! Initialisation de la lecture des fichiers
42
43 if (first) then
44 ncidpl=-99
45 print *, 'Intitialisation de read reanalsye'
46
47 ! Vent zonal
48 if (guide_u) then
49 rcode=nf90_open('u.nc', nf90_nowrite, ncidu)
50 rcode = nf90_inq_varid(ncidu, 'UWND', varidu)
51 print *, 'ncidu, varidu', ncidu, varidu
52 if (ncidpl.eq.-99) ncidpl=ncidu
53 endif
54
55 ! Vent meridien
56 if (guide_v) then
57 rcode=nf90_open('v.nc', nf90_nowrite, ncidv)
58 rcode = nf90_inq_varid(ncidv, 'VWND', varidv)
59 print *, 'ncidv, varidv', ncidv, varidv
60 if (ncidpl.eq.-99) ncidpl=ncidv
61 endif
62
63 ! Temperature
64 if (guide_T) then
65 rcode=nf90_open('T.nc', nf90_nowrite, ncidt)
66 rcode = nf90_inq_varid(ncidt, 'AIR', varidt)
67 print *, 'ncidt, varidt', ncidt, varidt
68 if (ncidpl.eq.-99) ncidpl=ncidt
69 endif
70
71 ! Humidite
72 if (guide_Q) then
73 rcode=nf90_open('hur.nc', nf90_nowrite, ncidQ)
74 rcode = nf90_inq_varid(ncidQ, 'RH', varidQ)
75 print *, 'ncidQ, varidQ', ncidQ, varidQ
76 if (ncidpl.eq.-99) ncidpl=ncidQ
77 endif
78
79 ! Coordonnee verticale
80 if (ncep) then
81 print *, 'Vous etes entrain de lire des donnees NCEP'
82 rcode = nf90_inq_varid(ncidpl, 'LEVEL', varidpl)
83 else
84 print *, 'Vous etes entrain de lire des donnees ECMWF'
85 rcode = nf90_inq_varid(ncidpl, 'PRESSURE', varidpl)
86 endif
87 print *, 'ncidu, varidpl', ncidu, varidpl
88 endif
89 print *, 'ok1'
90
91 ! Niveaux de pression
92 print *, 'WARNING!!! Il n y a pas de test de coherence'
93 print *, 'sur le nombre de niveaux verticaux dans le fichier nc'
94 status=NF90_GET_VAR(ncidpl, varidpl, pl)
95 ! passage en pascal
96 pl(:)=100.*pl(:)
97 if (first) then
98 do l=1, nlevnc
99 print *, 'PL(', l, ')=', pl(l)
100 enddo
101 endif
102
103 ! lecture des champs u, v, T
104
105 ! dimensions pour les champs scalaires et le vent zonal
106
107 start(1)=1
108 start(2)=1
109 start(3)=1
110 start(4)=timestep
111
112 count(1)=iip1
113 count(2)=jjp1
114 count(3)=nlevnc
115 count(4)=1
116
117 ! mise a zero des tableaux
118
119 unc(:, :, :)=0.
120 vnc(:, :, :)=0.
121 tnc(:, :, :)=0.
122 Qnc(:, :, :)=0.
123
124 ! Vent zonal
125
126 if (guide_u) then
127 print *, 'avant la lecture de UNCEP nd de niv:', nlevnc
128 status=NF90_GET_VAR(ncidu, varidu, unc, start, count)
129 print *, 'WARNING!!! Correction bidon pour palier a un '
130 print *, 'probleme dans la creation des fichiers nc'
131 call correctbid(iim, jjp1*nlevnc, unc)
132 call dump2d(iip1, jjp1, unc, 'UNC COUCHE 1 ')
133 endif
134
135 ! Temperature
136
137 print *, 'ncidt=', ncidt, 'varidt=', varidt, 'start=', start
138 print *, 'count=', count
139 if (guide_T) then
140 status=NF90_GET_VAR(ncidt, varidt, tnc, start, count)
141 call dump2d(iip1, jjp1, tnc, 'TNC COUCHE 1 AAA ')
142 call correctbid(iim, jjp1*nlevnc, tnc)
143 call dump2d(iip1, jjp1, tnc, 'TNC COUCHE 1 BBB ')
144 endif
145
146 ! Humidite
147
148 if (guide_Q) then
149 status=NF90_GET_VAR(ncidQ, varidQ, Qnc, start, count)
150 call correctbid(iim, jjp1*nlevnc, Qnc)
151 call dump2d(iip1, jjp1, Qnc, 'QNC COUCHE 1 ')
152 endif
153
154 count(2)=jjm
155 ! Vent meridien
156
157 if (guide_v) then
158 status=NF90_GET_VAR(ncidv, varidv, vnc, start, count)
159 call correctbid(iim, jjm*nlevnc, vnc)
160 call dump2d(iip1, jjm, vnc, 'VNC COUCHE 1 ')
161 endif
162
163 start(3)=timestep
164 start(4)=0
165 count(2)=jjp1
166 count(3)=1
167 count(4)=0
168
169 ! Interpolation verticale sur les niveaux modele
170
171 call reanalyse2nat(nlevnc, psi, unc, vnc, tnc, Qnc, pl, u, v, t, Q, &
172 masse, pk)
173
174 call dump2d(iip1, jjm, v, 'V COUCHE APRES ')
175
176 ! Passage aux variables du modele (vents covariants, temperature
177 ! potentielle et humidite specifique)
178
179 call nat2gcm(u, v, t, Q, pk, u, v, t, Q)
180 print *, 'TIMESTEP ', timestep
181 first=.false.
182
183 end subroutine read_reanalyse
184
185 end module read_reanalyse_m

  ViewVC Help
Powered by ViewVC 1.1.21