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