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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 15 - (show annotations)
Fri Aug 1 15:24:12 2008 UTC (15 years, 9 months ago) by guez
File size: 9136 byte(s)
-- Minor modification of input/output:

Added variable "Sigma_O3_Royer" to "histday.nc". "ecrit_day" is not
modified in "physiq". Removed variables "pyu1", "pyv1", "ftsol1",
"ftsol2", "ftsol3", "ftsol4", "psrf1", "psrf2", "psrf3", "psrf4"
"mfu", "mfd", "en_u", "en_d", "de_d", "de_u", "coefh" from
"histrac.nc".

Variable "raz_date" of module "conf_gcm_m" has logical type instead of
integer type.

-- Should not change any result at run time:

Modified calls to "IOIPSL_Lionel" procedures because the interfaces of
these procedures have been simplified.

Changed name of variable in module "start_init_orog_m": "masque" to
"mask".

Created a module containing procedure "phyredem".

Removed arguments "punjours", "pdayref" and "ptimestep" of procedure
"iniphysiq".

Renamed procedure "gr_phy_write" to "gr_phy_write_2d". Created
procedure "gr_phy_write_3d".

Removed procedures "ini_undefstd", "moy_undefSTD", "calcul_STDlev",
"calcul_divers".

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

  ViewVC Help
Powered by ViewVC 1.1.21