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

Annotation of /trunk/Sources/dyn3d/Guide/Read_reanalyse/read_reanalyse.f

Parent Directory Parent Directory | Revision Log Revision Log


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

1 guez 88 module read_reanalyse_m
2 guez 37
3 guez 88 IMPLICIT NONE
4 guez 3
5 guez 88 contains
6 guez 20
7 guez 88 subroutine read_reanalyse(timestep, psi, u, v, t, q, masse, nlevnc)
8 guez 3
9 guez 88 ! From LMDZ4/libf/dyn3d/read_reanalyse.F, version 1.3, 2005/04/15 12:31:21
10 guez 3
11 guez 88 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 guez 3
17 guez 88 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 guez 3
24 guez 88 ! Local:
25 guez 3
26 guez 88 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 guez 3
39 guez 88 ! -----------------------------------------------------------------
40 guez 3
41 guez 88 ! Initialisation de la lecture des fichiers
42 guez 3
43 guez 88 if (first) then
44     ncidpl=-99
45     print *, 'Intitialisation de read reanalsye'
46 guez 3
47 guez 88 ! 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 guez 3
55 guez 88 ! 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 guez 3
63 guez 88 ! 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 guez 3
71 guez 88 ! 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 guez 3
79 guez 88 ! 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 guez 3
91 guez 88 ! 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 guez 3
103 guez 88 ! lecture des champs u, v, T
104 guez 3
105 guez 88 ! dimensions pour les champs scalaires et le vent zonal
106 guez 3
107 guez 88 start(1)=1
108     start(2)=1
109     start(3)=1
110     start(4)=timestep
111 guez 3
112 guez 88 count(1)=iip1
113     count(2)=jjp1
114     count(3)=nlevnc
115     count(4)=1
116 guez 3
117 guez 88 ! mise a zero des tableaux
118 guez 3
119 guez 88 unc(:, :, :)=0.
120     vnc(:, :, :)=0.
121     tnc(:, :, :)=0.
122     Qnc(:, :, :)=0.
123 guez 3
124 guez 88 ! Vent zonal
125 guez 3
126 guez 88 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 guez 3
135 guez 88 ! Temperature
136 guez 3
137 guez 88 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 guez 3
146 guez 88 ! Humidite
147 guez 3
148 guez 88 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 guez 3
154 guez 88 count(2)=jjm
155     ! Vent meridien
156 guez 3
157 guez 88 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 guez 3
163 guez 88 start(3)=timestep
164     start(4)=0
165     count(2)=jjp1
166     count(3)=1
167     count(4)=0
168 guez 3
169 guez 88 ! Interpolation verticale sur les niveaux modele
170 guez 3
171 guez 88 call reanalyse2nat(nlevnc, psi, unc, vnc, tnc, Qnc, pl, u, v, t, Q, &
172     masse, pk)
173 guez 3
174 guez 88 call dump2d(iip1, jjm, v, 'V COUCHE APRES ')
175 guez 3
176 guez 88 ! Passage aux variables du modele (vents covariants, temperature
177     ! potentielle et humidite specifique)
178 guez 3
179 guez 88 call nat2gcm(u, v, t, Q, pk, u, v, t, Q)
180     print *, 'TIMESTEP ', timestep
181     first=.false.
182 guez 3
183 guez 88 end subroutine read_reanalyse
184 guez 3
185 guez 88 end module read_reanalyse_m

  ViewVC Help
Powered by ViewVC 1.1.21