/[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 108 - (show annotations)
Tue Sep 16 14:00:41 2014 UTC (9 years, 7 months ago) by guez
File size: 4412 byte(s)
Imported writefield from LMDZ. Close at the end of gcm the files which
were created by writefiled (not done in LMDZ).

Removed procedures for the output of Grads files. Removed calls to
dump2d. In guide, replaced calls to wrgrads by calls to writefield.

In vlspltqs, removed redundant programming of saturation
pressure. Call foeew from module FCTTRE instead.

Bug fix in interpre: size of w exceeding size of correponding actual
argument wg in advtrac.

In leapfrog, call guide until the end of the run, instead of six hours
before the end.

Bug fix in readsulfate_preind: type of arguments.

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 if (ncidpl.eq.-99) ncidpl=ncidu
52 endif
53
54 ! Vent meridien
55 if (guide_v) then
56 rcode=nf90_open('v.nc', nf90_nowrite, ncidv)
57 rcode = nf90_inq_varid(ncidv, 'VWND', varidv)
58 if (ncidpl.eq.-99) ncidpl=ncidv
59 endif
60
61 ! Temperature
62 if (guide_T) then
63 rcode=nf90_open('T.nc', nf90_nowrite, ncidt)
64 rcode = nf90_inq_varid(ncidt, 'AIR', varidt)
65 if (ncidpl.eq.-99) ncidpl=ncidt
66 endif
67
68 ! Humidite
69 if (guide_Q) then
70 rcode=nf90_open('hur.nc', nf90_nowrite, ncidQ)
71 rcode = nf90_inq_varid(ncidQ, 'RH', varidQ)
72 if (ncidpl.eq.-99) ncidpl=ncidQ
73 endif
74
75 ! Coordonnee verticale
76 if (ncep) then
77 print *, 'Vous etes entrain de lire des donnees NCEP'
78 rcode = nf90_inq_varid(ncidpl, 'LEVEL', varidpl)
79 else
80 print *, 'Vous etes entrain de lire des donnees ECMWF'
81 rcode = nf90_inq_varid(ncidpl, 'PRESSURE', varidpl)
82 endif
83 endif
84
85 ! Niveaux de pression
86
87 ! Warning: il n y a pas de test de coherence sur le nombre de
88 ! niveaux verticaux dans le fichier nc'
89 status=NF90_GET_VAR(ncidpl, varidpl, pl)
90 ! passage en pascal
91 pl(:)=100.*pl(:)
92 if (first) then
93 do l=1, nlevnc
94 print *, 'PL(', l, ')=', pl(l)
95 enddo
96 endif
97
98 ! lecture des champs u, v, T
99
100 ! dimensions pour les champs scalaires et le vent zonal
101
102 start(1)=1
103 start(2)=1
104 start(3)=1
105 start(4)=timestep
106
107 count(1)=iip1
108 count(2)=jjp1
109 count(3)=nlevnc
110 count(4)=1
111
112 ! mise a zero des tableaux
113
114 unc(:, :, :)=0.
115 vnc(:, :, :)=0.
116 tnc(:, :, :)=0.
117 Qnc(:, :, :)=0.
118
119 ! Vent zonal
120
121 if (guide_u) then
122 status=NF90_GET_VAR(ncidu, varidu, unc, start, count)
123 ! Warning Correction bidon pour palier a un probleme dans la
124 ! creation des fichiers nc
125 call correctbid(iim, jjp1*nlevnc, unc)
126 endif
127
128 ! Temperature
129
130 if (guide_T) then
131 status=NF90_GET_VAR(ncidt, varidt, tnc, start, count)
132 call correctbid(iim, jjp1*nlevnc, tnc)
133 endif
134
135 ! Humidite
136
137 if (guide_Q) then
138 status=NF90_GET_VAR(ncidQ, varidQ, Qnc, start, count)
139 call correctbid(iim, jjp1*nlevnc, Qnc)
140 endif
141
142 count(2)=jjm
143 ! Vent meridien
144
145 if (guide_v) then
146 status=NF90_GET_VAR(ncidv, varidv, vnc, start, count)
147 call correctbid(iim, jjm*nlevnc, vnc)
148 endif
149
150 start(3)=timestep
151 start(4)=0
152 count(2)=jjp1
153 count(3)=1
154 count(4)=0
155
156 ! Interpolation verticale sur les niveaux modele
157
158 call reanalyse2nat(nlevnc, psi, unc, vnc, tnc, Qnc, pl, u, v, t, Q, &
159 masse, pk)
160
161 ! Passage aux variables du modele (vents covariants, temperature
162 ! potentielle et humidite specifique)
163
164 call nat2gcm(u, v, t, Q, pk, u, v, t, Q)
165 first=.false.
166
167 end subroutine read_reanalyse
168
169 end module read_reanalyse_m

  ViewVC Help
Powered by ViewVC 1.1.21