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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 10 - (show annotations)
Fri Apr 18 14:45:53 2008 UTC (16 years, 1 month ago) by guez
File size: 9137 byte(s)
Added NetCDF directory "/home/guez/include" in "g95.mk" and
"nag_tools.mk".

Added some "intent" attributes in "PVtheta", "advtrac", "caladvtrac",
"calfis", "diagedyn", "dissip", "vlspltqs", "aeropt", "ajsec",
"calltherm", "clmain", "cltrac", "cltracrn", "concvl", "conema3",
"conflx", "fisrtilp", "newmicro", "nuage", "diagcld1", "diagcld2",
"drag_noro", "lift_noro", "SUGWD", "physiq", "phytrac", "radlwsw", "thermcell".

Removed the case "ierr == 0" in "abort_gcm"; moved call to "histclo"
and messages for end of run from "abort_gcm" to "gcm"; replaced call
to "abort_gcm" in "leapfrog" by exit from outer loop.

In "calfis": removed argument "pp" and variable "unskap"; changed
"pksurcp" from scalar to rank 2; use "pressure_var"; rewrote
computation of "zplev", "zplay", "ztfi", "pcvgt" using "dyn_phy";
added computation of "pls".

Removed unused variable in "dynredem0".

In "exner_hyb": changed "dellta" from scalar to rank 1; replaced call
to "ssum" by call to "sum"; removed variables "xpn" and "xps";
replaced some loops by array expressions.

In "leapfrog": use "pressure_var"; deleted variables "p", "longcles".

Converted common blocks "YOECUMF", "YOEGWD" to modules.

Removed argument "pplay" in "cvltr", "diagetpq", "nflxtr".

Created module "raddimlw" from include file "raddimlw.h".

Corrected call to "new_unit" in "test_disvert".

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

  ViewVC Help
Powered by ViewVC 1.1.21