/[lmdze]/trunk/libf/dyn3d/dynredem0.f90
ViewVC logotype

Contents of /trunk/libf/dyn3d/dynredem0.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 57 - (show annotations)
Mon Jan 30 12:54:02 2012 UTC (12 years, 3 months ago) by guez
File size: 9215 byte(s)
Write used namelists to file "" instead of standard output.

Avoid aliasing in "inidissip" in calls to "divgrad2", "divgrad",
"gradiv2", "gradiv", "nxgraro2" and "nxgrarot". Add a degenerate
dimension to arrays so they have rank 3, like the dummy arguments in
"divgrad2", "divgrad", "gradiv2", "gradiv", "nxgraro2" and "nxgrarot".

Extract the initialization part from "bilan_dyn" and make a separate
procedure, "init_dynzon", from it.

Move variables from modules "iniprint" and "logic" to module
"conf_gcm_m".

Promote internal procedures of "fxy" to private procedures of module
"fxy_m".

Extracted documentation from "inigeom". Removed useless "save"
attributes. Removed useless intermediate variables. Extracted
processing of poles from loop on latitudes. Write coordinates to file
"longitude_latitude.txt" instead of standard output.

Do not use ozone tracer for radiative transfer.

1 MODULE dynredem0_m
2
3 IMPLICIT NONE
4
5 CONTAINS
6
7 SUBROUTINE dynredem0(fichnom, iday_end, phis)
8
9 ! From dyn3d/dynredem.F, version 1.2 2004/06/22 11:45:30
10 ! Ecriture du fichier de redémarrage au format NetCDF (initialisation)
11
12 USE comconst, ONLY : cpp, daysec, dtvr, g, kappa, omeg, rad
13 USE comvert, ONLY : ap, bp, nivsig, nivsigs, pa, preff, presnivs
14 USE comgeom, ONLY : aire_2d, cu_2d, cv_2d, rlatu, rlatv, rlonu, rlonv
15 USE dimens_m, ONLY : iim, jjm, llm, nqmx
16 USE ener, ONLY : ang0, etot0, ptot0, stot0, ztot0
17 USE calendar, ONLY : ju2ymds, ymds2ju
18 USE iniadvtrac_m, ONLY : tname, ttext
19 use conf_gcm_m, ONLY : fxyhypb, ysinus
20 USE netcdf95, ONLY : nf95_close, nf95_create, nf95_def_dim, &
21 nf95_def_var, nf95_enddef, nf95_inq_varid, nf95_put_att, &
22 nf95_put_var
23 USE netcdf, ONLY : nf90_clobber, nf90_float, nf90_global, &
24 nf90_unlimited
25 USE paramet_m, ONLY : iip1, jjp1, llmp1
26 USE serre, ONLY : clat, clon, dzoomx, dzoomy, grossismx, grossismy, &
27 taux, tauy
28 USE temps, ONLY : annee_ref, day_ref
29
30 CHARACTER (len=*), INTENT (IN) :: fichnom
31 INTEGER, INTENT (IN) :: iday_end
32 REAL, INTENT (IN) :: phis(:, :)
33
34 ! Local:
35
36 INTEGER :: iq, l
37 INTEGER, PARAMETER:: length = 100
38 REAL :: tab_cntrl(length) ! tableau des parametres du run
39
40 ! Variables locales pour NetCDF:
41
42 INTEGER :: dims2(2), dims3(3), dims4(4)
43 INTEGER :: idim_index
44 INTEGER :: idim_rlonu, idim_rlonv, idim_rlatu, idim_rlatv
45 INTEGER :: idim_s, idim_sig
46 INTEGER :: idim_tim
47 INTEGER :: ncid, varid
48
49 REAL :: zjulian, hours
50 INTEGER :: yyears0, jjour0, mmois0
51 CHARACTER (len=30) :: unites
52
53 !-----------------------------------------------------------------------
54
55 PRINT *, 'Call sequence information: dynredem0'
56
57 CALL ymds2ju(annee_ref, 1, iday_end, 0., zjulian)
58 CALL ju2ymds(zjulian, yyears0, mmois0, jjour0, hours)
59
60 DO l = 1, length
61 tab_cntrl(l) = 0.
62 END DO
63 tab_cntrl(1) = real(iim)
64 tab_cntrl(2) = real(jjm)
65 tab_cntrl(3) = real(llm)
66 tab_cntrl(4) = real(day_ref)
67 tab_cntrl(5) = real(annee_ref)
68 tab_cntrl(6) = rad
69 tab_cntrl(7) = omeg
70 tab_cntrl(8) = g
71 tab_cntrl(9) = cpp
72 tab_cntrl(10) = kappa
73 tab_cntrl(11) = daysec
74 tab_cntrl(12) = dtvr
75 tab_cntrl(13) = etot0
76 tab_cntrl(14) = ptot0
77 tab_cntrl(15) = ztot0
78 tab_cntrl(16) = stot0
79 tab_cntrl(17) = ang0
80 tab_cntrl(18) = pa
81 tab_cntrl(19) = preff
82
83 ! Paramètres pour le zoom :
84
85 tab_cntrl(20) = clon
86 tab_cntrl(21) = clat
87 tab_cntrl(22) = grossismx
88 tab_cntrl(23) = grossismy
89
90 IF (fxyhypb) THEN
91 tab_cntrl(24) = 1.
92 tab_cntrl(25) = dzoomx
93 tab_cntrl(26) = dzoomy
94 tab_cntrl(27) = 0.
95 tab_cntrl(28) = taux
96 tab_cntrl(29) = tauy
97 ELSE
98 tab_cntrl(24) = 0.
99 tab_cntrl(25) = dzoomx
100 tab_cntrl(26) = dzoomy
101 tab_cntrl(27) = 0.
102 tab_cntrl(28) = 0.
103 tab_cntrl(29) = 0.
104 IF (ysinus) tab_cntrl(27) = 1.
105 END IF
106
107 tab_cntrl(30) = real(iday_end)
108
109 CALL nf95_create(fichnom, nf90_clobber, ncid)
110 CALL nf95_put_att(ncid, nf90_global, 'title', &
111 'Fichier de démarrage dynamique')
112
113 ! Definir les dimensions du fichiers:
114
115 CALL nf95_def_dim(ncid, 'index', length, idim_index)
116 CALL nf95_def_dim(ncid, 'rlonu', iip1, idim_rlonu)
117 CALL nf95_def_dim(ncid, 'rlatu', jjp1, idim_rlatu)
118 CALL nf95_def_dim(ncid, 'rlonv', iip1, idim_rlonv)
119 CALL nf95_def_dim(ncid, 'rlatv', jjm, idim_rlatv)
120 CALL nf95_def_dim(ncid, 'sigs', llm, idim_s)
121 CALL nf95_def_dim(ncid, 'sig', llmp1, idim_sig)
122 CALL nf95_def_dim(ncid, 'temps', nf90_unlimited, idim_tim)
123
124 ! Definir et enregistrer certains champs invariants:
125
126 CALL nf95_def_var(ncid, 'controle', nf90_float, idim_index, varid)
127 CALL nf95_put_att(ncid, varid, 'title', 'Parametres de controle')
128
129 CALL nf95_def_var(ncid, 'rlonu', nf90_float, idim_rlonu, varid)
130 CALL nf95_put_att(ncid, varid, 'title', 'Longitudes des points U')
131
132 CALL nf95_def_var(ncid, 'rlatu', nf90_float, idim_rlatu, varid)
133 CALL nf95_put_att(ncid, varid, 'title', 'Latitudes des points U')
134
135 CALL nf95_def_var(ncid, 'rlonv', nf90_float, idim_rlonv, varid)
136 CALL nf95_put_att(ncid, varid, 'title', 'Longitudes des points V')
137
138 CALL nf95_def_var(ncid, 'rlatv', nf90_float, idim_rlatv, varid)
139 CALL nf95_put_att(ncid, varid, 'title', 'Latitudes des points V')
140
141 CALL nf95_def_var(ncid, 'nivsigs', nf90_float, idim_s, varid)
142 CALL nf95_put_att(ncid, varid, 'title', 'Numero naturel des couches s')
143
144 CALL nf95_def_var(ncid, 'nivsig', nf90_float, idim_sig, varid)
145 CALL nf95_put_att(ncid, varid, 'title', &
146 'Numero naturel des couches sigma')
147
148 CALL nf95_def_var(ncid, 'ap', nf90_float, idim_sig, varid)
149 CALL nf95_put_att(ncid, varid, 'title', 'Coefficient A pour hybride')
150
151 CALL nf95_def_var(ncid, 'bp', nf90_float, idim_sig, varid)
152 CALL nf95_put_att(ncid, varid, 'title', 'Coefficient B pour hybride')
153
154 CALL nf95_def_var(ncid, 'presnivs', nf90_float, idim_s, varid)
155
156 ! Coefficients de passage cov. <-> contra. <--> naturel
157
158 dims2(1) = idim_rlonu
159 dims2(2) = idim_rlatu
160 CALL nf95_def_var(ncid, 'cu', nf90_float, dims2, varid)
161 CALL nf95_put_att(ncid, varid, 'title', 'Coefficient de passage pour U')
162
163 dims2(1) = idim_rlonv
164 dims2(2) = idim_rlatv
165 CALL nf95_def_var(ncid, 'cv', nf90_float, dims2, varid)
166 CALL nf95_put_att(ncid, varid, 'title', 'Coefficient de passage pour V')
167
168 ! Aire de chaque maille:
169
170 dims2(1) = idim_rlonv
171 dims2(2) = idim_rlatu
172 CALL nf95_def_var(ncid, 'aire', nf90_float, dims2, varid)
173 CALL nf95_put_att(ncid, varid, 'title', 'Aires de chaque maille')
174
175 ! Geopentiel au sol:
176
177 dims2(1) = idim_rlonv
178 dims2(2) = idim_rlatu
179 CALL nf95_def_var(ncid, 'phisinit', nf90_float, dims2, varid)
180 CALL nf95_put_att(ncid, varid, 'title', 'Geopotentiel au sol')
181
182 ! Definir les variables pour pouvoir les enregistrer plus tard:
183
184 CALL nf95_def_var(ncid, 'temps', nf90_float, idim_tim, varid)
185 CALL nf95_put_att(ncid, varid, 'title', 'Temps de simulation')
186 WRITE (unites, 200) yyears0, mmois0, jjour0
187 200 FORMAT ('days since ', I4, '-', I2.2, '-', I2.2, ' 00:00:00')
188 CALL nf95_put_att(ncid, varid, 'units', unites)
189
190
191 dims4(1) = idim_rlonu
192 dims4(2) = idim_rlatu
193 dims4(3) = idim_s
194 dims4(4) = idim_tim
195 CALL nf95_def_var(ncid, 'ucov', nf90_float, dims4, varid)
196 CALL nf95_put_att(ncid, varid, 'title', 'Vitesse U')
197
198 dims4(1) = idim_rlonv
199 dims4(2) = idim_rlatv
200 dims4(3) = idim_s
201 dims4(4) = idim_tim
202 CALL nf95_def_var(ncid, 'vcov', nf90_float, dims4, varid)
203 CALL nf95_put_att(ncid, varid, 'title', 'Vitesse V')
204
205 dims4(1) = idim_rlonv
206 dims4(2) = idim_rlatu
207 dims4(3) = idim_s
208 dims4(4) = idim_tim
209 CALL nf95_def_var(ncid, 'teta', nf90_float, dims4, varid)
210 CALL nf95_put_att(ncid, varid, 'title', 'Temperature')
211
212 dims4(1) = idim_rlonv
213 dims4(2) = idim_rlatu
214 dims4(3) = idim_s
215 dims4(4) = idim_tim
216 DO iq = 1, nqmx
217 CALL nf95_def_var(ncid, tname(iq), nf90_float, dims4, varid)
218 CALL nf95_put_att(ncid, varid, 'title', ttext(iq))
219 END DO
220
221 dims4(1) = idim_rlonv
222 dims4(2) = idim_rlatu
223 dims4(3) = idim_s
224 dims4(4) = idim_tim
225 CALL nf95_def_var(ncid, 'masse', nf90_float, dims4, varid)
226 CALL nf95_put_att(ncid, varid, 'title', 'C est quoi ?')
227
228 dims3(1) = idim_rlonv
229 dims3(2) = idim_rlatu
230 dims3(3) = idim_tim
231 CALL nf95_def_var(ncid, 'ps', nf90_float, dims3, varid)
232 CALL nf95_put_att(ncid, varid, 'title', 'Pression au sol')
233
234 CALL nf95_enddef(ncid)
235
236 CALL nf95_inq_varid(ncid, 'controle', varid)
237 CALL nf95_put_var(ncid, varid, tab_cntrl)
238
239 CALL nf95_inq_varid(ncid, 'rlonu', varid)
240 CALL nf95_put_var(ncid, varid, rlonu)
241
242 CALL nf95_inq_varid(ncid, 'rlatu', varid)
243 CALL nf95_put_var(ncid, varid, rlatu)
244
245 CALL nf95_inq_varid(ncid, 'rlonv', varid)
246 CALL nf95_put_var(ncid, varid, rlonv)
247
248 CALL nf95_inq_varid(ncid, 'rlatv', varid)
249 CALL nf95_put_var(ncid, varid, rlatv)
250
251 CALL nf95_inq_varid(ncid, 'nivsigs', varid)
252 CALL nf95_put_var(ncid, varid, nivsigs)
253
254 CALL nf95_inq_varid(ncid, 'nivsig', varid)
255 CALL nf95_put_var(ncid, varid, nivsig)
256
257 CALL nf95_inq_varid(ncid, 'ap', varid)
258 CALL nf95_put_var(ncid, varid, ap)
259
260 CALL nf95_inq_varid(ncid, 'bp', varid)
261 CALL nf95_put_var(ncid, varid, bp)
262
263 CALL nf95_inq_varid(ncid, 'presnivs', varid)
264 CALL nf95_put_var(ncid, varid, presnivs)
265
266 CALL nf95_inq_varid(ncid, 'cu', varid)
267 CALL nf95_put_var(ncid, varid, cu_2d)
268
269 CALL nf95_inq_varid(ncid, 'cv', varid)
270 CALL nf95_put_var(ncid, varid, cv_2d)
271
272 CALL nf95_inq_varid(ncid, 'aire', varid)
273 CALL nf95_put_var(ncid, varid, aire_2d)
274
275 CALL nf95_inq_varid(ncid, 'phisinit', varid)
276 CALL nf95_put_var(ncid, varid, phis)
277
278 CALL nf95_close(ncid) ! fermer le fichier
279
280 PRINT *, 'iim, jjm, llm, iday_end', iim, jjm, llm, iday_end
281 PRINT *, 'rad, omeg, g, cpp, kappa', rad, omeg, g, cpp, kappa
282
283 END SUBROUTINE dynredem0
284
285 END MODULE dynredem0_m

  ViewVC Help
Powered by ViewVC 1.1.21