/[lmdze]/trunk/phylmd/phyredem.f
ViewVC logotype

Contents of /trunk/phylmd/phyredem.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 156 - (show annotations)
Thu Jul 16 17:39:10 2015 UTC (8 years, 10 months ago) by guez
Original Path: trunk/Sources/phylmd/phyredem.f
File size: 11153 byte(s)
In procedure cltracrn, no need for local variable zx_trs, use directly
local_trs.

In (re)startphy.nc, agglomerate variables for different surface types
into a single variable with an added dimension.

In phyredem, bring together all definitions, do not use redef.

1 module phyredem_m
2
3 IMPLICIT NONE
4
5 contains
6
7 SUBROUTINE phyredem(fichnom, pctsrf, tsol, tsoil, tslab, seaice, qsurf, &
8 qsol, snow, albedo, evap, rain_fall, snow_fall, solsw, sollw, &
9 fder, radsol, frugs, agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, &
10 t_ancien, q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, w01)
11
12 ! From phylmd/phyredem.F, version 1.3, 2005/05/25 13:10:09
13 ! Author: Z. X. Li (LMD/CNRS)
14 ! Date: 1993/08/18
15
16 ! Objet : \'ecriture de l'\'etat de d\'emarrage ou red\'emarrage
17 ! pour la physique
18
19 USE dimphy, ONLY: klev, klon, zmasq
20 USE dimsoil, ONLY: nsoilmx
21 USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf
22 USE netcdf, ONLY: nf90_clobber, nf90_global, nf90_float
23 USE netcdf95, ONLY: nf95_create, nf95_put_att, nf95_def_dim, &
24 nf95_def_var, nf95_enddef, nf95_put_var, nf95_close
25 use phyetat0_m, only: rlat, rlon
26 USE temps, ONLY: itau_phy
27
28 CHARACTER(len=*), INTENT(IN):: fichnom
29 REAL, INTENT(IN):: pctsrf(:, :) ! (klon, nbsrf)
30 REAL, INTENT(IN):: tsol(:, :) ! (klon, nbsrf)
31 REAL, INTENT(IN):: tsoil(:, :, :) ! (klon, nsoilmx, nbsrf)
32 REAL, INTENT(IN):: tslab(:), seaice(:) ! (klon) slab ocean
33 REAL, INTENT(IN):: qsurf(:, :) ! (klon, nbsrf)
34
35 REAL, intent(in):: qsol(:) ! (klon)
36 ! column-density of water in soil, in kg m-2
37
38 REAL, INTENT(IN):: snow(klon, nbsrf)
39 REAL, INTENT(IN):: albedo(klon, nbsrf)
40 REAL, INTENT(IN):: evap(klon, nbsrf)
41 REAL, INTENT(IN):: rain_fall(klon)
42 REAL, INTENT(IN):: snow_fall(klon)
43 REAL, INTENT(IN):: solsw(klon)
44 REAL, INTENT(IN):: sollw(klon)
45 REAL, INTENT(IN):: fder(klon)
46 REAL, INTENT(IN):: radsol(klon)
47 REAL, INTENT(IN):: frugs(klon, nbsrf)
48 REAL, INTENT(IN):: agesno(klon, nbsrf)
49 REAL, INTENT(IN):: zmea(klon)
50 REAL, intent(in):: zstd(klon)
51 REAL, intent(in):: zsig(klon)
52 REAL, intent(in):: zgam(klon)
53 REAL, intent(in):: zthe(klon)
54 REAL, intent(in):: zpic(klon)
55 REAL, intent(in):: zval(klon)
56 REAL, intent(in):: t_ancien(klon, klev), q_ancien(klon, klev)
57 REAL, intent(in):: rnebcon(klon, klev), ratqs(klon, klev)
58 REAL, intent(in):: clwcon(klon, klev)
59 REAL, intent(in):: run_off_lic_0(klon)
60 real, intent(in):: sig1(klon, klev) ! section adiabatic updraft
61
62 real, intent(in):: w01(klon, klev)
63 ! vertical velocity within adiabatic updraft
64
65 ! Local:
66
67 INTEGER ncid, idim2, idim3, dimid_nbsrf, dimid_nsoilmx
68 integer varid, varid_run_off_lic_0, varid_sig1, varid_w01, varid_rlon
69 integer varid_rlat, varid_zmasq, varid_fter, varid_flic, varid_foce
70 integer varid_fsic, varid_ts, varid_tsoil, varid_tslab, varid_SEAICE
71 integer varid_qs, varid_QSOL, varid_albe, varid_evap, varid_snow
72 integer varid_rads, varid_solsw, varid_sollw, varid_fder, varid_rain_f
73 integer varid_snow_f, varid_rug, varid_agesno, varid_zmea, varid_zstd
74 integer varid_zsig, varid_zgam, varid_zthe, varid_zpic, varid_zval
75 integer varid_tancien, varid_qancien, varid_rugmer, varid_clwcon
76 integer varid_rnebcon, varid_ratqs
77
78 !------------------------------------------------------------
79
80 PRINT *, 'Call sequence information: phyredem'
81 CALL nf95_create(fichnom, nf90_clobber, ncid)
82
83 call nf95_put_att(ncid, nf90_global, 'title', &
84 'start file for the physics code')
85 call nf95_put_att(ncid, nf90_global, "itau_phy", itau_phy)
86
87 call nf95_def_dim(ncid, 'points_physiques', klon, idim2)
88 call nf95_def_dim(ncid, 'klev', klev, idim3)
89 call nf95_def_dim(ncid, 'nbsrf', nbsrf, dimid_nbsrf)
90 call nf95_def_dim(ncid, 'nsoilmx', nsoilmx, dimid_nsoilmx)
91
92 call nf95_def_var(ncid, 'longitude', nf90_float, idim2, varid_rlon)
93 call nf95_def_var(ncid, 'latitude', nf90_float, idim2, varid_rlat)
94
95 call nf95_def_var(ncid, 'masque', nf90_float, idim2, varid_zmasq)
96 call nf95_put_att(ncid, varid_zmasq, 'title', 'masque terre mer')
97
98 ! Fractions de chaque sous-surface
99
100 call nf95_def_var(ncid, 'FTER', nf90_float, idim2, varid_fter)
101 call nf95_put_att(ncid, varid_fter, 'title', 'fraction de continent')
102
103 call nf95_def_var(ncid, 'FLIC', nf90_float, idim2, varid_flic)
104 call nf95_put_att(ncid, varid_flic, 'title', 'fraction glace de terre')
105
106 call nf95_def_var(ncid, 'FOCE', nf90_float, idim2, varid_foce)
107 call nf95_put_att(ncid, varid_foce, 'title', 'fraction ocean')
108
109 call nf95_def_var(ncid, 'FSIC', nf90_float, idim2, varid_fsic)
110 call nf95_put_att(ncid, varid_fsic, 'title', 'fraction glace mer')
111
112 call nf95_def_var(ncid, 'TS', nf90_float, (/idim2, dimid_nbsrf/), varid_ts)
113 call nf95_put_att(ncid, varid_ts, 'title', 'surface temperature')
114
115 call nf95_def_var(ncid, 'Tsoil', nf90_float, (/idim2, dimid_nsoilmx, &
116 dimid_nbsrf/), varid_tsoil)
117 call nf95_put_att(ncid, varid_tsoil, 'title', 'soil temperature')
118
119 call nf95_def_var(ncid, 'TSLAB', nf90_float, idim2, varid_tslab)
120 call nf95_put_att(ncid, varid_tslab, 'title', &
121 'Ecart de la SST (pour slab-ocean)')
122
123 call nf95_def_var(ncid, 'SEAICE', nf90_float, idim2, varid_SEAICE)
124 call nf95_put_att(ncid, varid_SEAICE, 'title', &
125 'Glace de mer kg/m2 (pour slab-ocean)')
126
127 call nf95_def_var(ncid, 'QS', nf90_float, (/idim2, dimid_nbsrf/), varid_qs)
128 call nf95_put_att(ncid, varid, 'title', 'Humidite de surface')
129
130 call nf95_def_var(ncid, 'QSOL', nf90_float, idim2, varid_QSOL)
131 call nf95_put_att(ncid, varid_QSOL, 'title', 'Eau dans le sol (mm)')
132
133 call nf95_def_var(ncid, 'ALBE', nf90_float, (/idim2, dimid_nbsrf/), &
134 varid_albe)
135 call nf95_put_att(ncid, varid_albe, 'title', 'albedo de surface')
136
137 call nf95_def_var(ncid, 'EVAP', nf90_float, (/idim2, dimid_nbsrf/), &
138 varid_evap)
139 call nf95_put_att(ncid, varid_evap, 'title', 'Evaporation de surface')
140
141 call nf95_def_var(ncid, 'SNOW', nf90_float, (/idim2, dimid_nbsrf/), &
142 varid_snow)
143 call nf95_put_att(ncid, varid_snow, 'title', 'Neige de surface')
144
145 call nf95_def_var(ncid, 'RADS', nf90_float, idim2, varid_rads)
146 call nf95_put_att(ncid, varid_rads, 'title', &
147 'Rayonnement net a la surface')
148
149 call nf95_def_var(ncid, 'solsw', nf90_float, idim2, varid_solsw)
150 call nf95_put_att(ncid, varid_solsw, 'title', &
151 'Rayonnement solaire a la surface')
152
153 call nf95_def_var(ncid, 'sollw', nf90_float, idim2, varid_sollw)
154 call nf95_put_att(ncid, varid_sollw, 'title', &
155 'Rayonnement IF a la surface')
156
157 call nf95_def_var(ncid, 'fder', nf90_float, idim2, varid_fder)
158 call nf95_put_att(ncid, varid_fder, 'title', 'Derive de flux')
159
160 call nf95_def_var(ncid, 'rain_f', nf90_float, idim2, varid_rain_f)
161 call nf95_put_att(ncid, varid_rain_f, 'title', 'precipitation liquide')
162
163 call nf95_def_var(ncid, 'snow_f', nf90_float, idim2, varid_snow_f)
164 call nf95_put_att(ncid, varid_snow_f, 'title', 'precipitation solide')
165
166 call nf95_def_var(ncid, 'RUG', nf90_float, (/idim2, dimid_nbsrf/), &
167 varid_rug)
168 call nf95_put_att(ncid, varid_rug, 'title', 'rugosite de surface')
169
170 call nf95_def_var(ncid, 'AGESNO', nf90_float, (/idim2, dimid_nbsrf/), &
171 varid_agesno)
172 call nf95_put_att(ncid, varid_agesno, 'title', 'Age de la neige surface')
173
174 call nf95_def_var(ncid, 'ZMEA', nf90_float, idim2, varid_zmea)
175 call nf95_def_var(ncid, 'ZSTD', nf90_float, idim2, varid_zstd)
176 call nf95_def_var(ncid, 'ZSIG', nf90_float, idim2, varid_zsig)
177 call nf95_def_var(ncid, 'ZGAM', nf90_float, idim2, varid_zgam)
178 call nf95_def_var(ncid, 'ZTHE', nf90_float, idim2, varid_zthe)
179 call nf95_def_var(ncid, 'ZPIC', nf90_float, idim2, varid_zpic)
180 call nf95_def_var(ncid, 'ZVAL', nf90_float, idim2, varid_zval)
181 call nf95_def_var(ncid, 'TANCIEN', nf90_float, (/idim2, idim3/), &
182 varid_tancien)
183 call nf95_def_var(ncid, 'QANCIEN', nf90_float, (/idim2, idim3/), &
184 varid_qancien)
185
186 call nf95_def_var(ncid, 'RUGMER', nf90_float, idim2, varid_rugmer)
187 call nf95_put_att(ncid, varid_rugmer, 'title', &
188 'Longueur de rugosite sur mer')
189
190 call nf95_def_var(ncid, 'CLWCON', nf90_float, idim2, varid_clwcon)
191 call nf95_put_att(ncid, varid_clwcon, 'title', 'Eau liquide convective')
192
193 call nf95_def_var(ncid, 'RNEBCON', nf90_float, idim2, varid_rnebcon)
194 call nf95_put_att(ncid, varid_rnebcon, 'title', 'Nebulosite convective')
195
196 call nf95_def_var(ncid, 'RATQS', nf90_float, idim2, varid_ratqs)
197 call nf95_put_att(ncid, varid_ratqs, 'title', 'Ratqs')
198
199 call nf95_def_var(ncid, 'RUNOFFLIC0', nf90_float, idim2, &
200 varid_run_off_lic_0)
201 call nf95_put_att(ncid, varid_run_off_lic_0, 'title', 'Runofflic0')
202
203 call nf95_def_var(ncid, 'sig1', nf90_float, (/idim2, idim3/), varid_sig1)
204 call nf95_put_att(ncid, varid_sig1, 'long_name', &
205 'section adiabatic updraft')
206
207 call nf95_def_var(ncid, 'w01', nf90_float, (/idim2, idim3/), varid_w01)
208 call nf95_put_att(ncid, varid_w01, 'long_name', &
209 'vertical velocity within adiabatic updraft')
210
211 call nf95_enddef(ncid)
212
213 call nf95_put_var(ncid, varid_rlon, rlon)
214 call nf95_put_var(ncid, varid_rlat, rlat)
215 call nf95_put_var(ncid, varid_zmasq, zmasq)
216 call nf95_put_var(ncid, varid_fter, pctsrf(:, is_ter))
217 call nf95_put_var(ncid, varid_flic, pctsrf(:, is_lic))
218 call nf95_put_var(ncid, varid_foce, pctsrf(:, is_oce))
219 call nf95_put_var(ncid, varid_fsic, pctsrf(:, is_sic))
220 call nf95_put_var(ncid, varid_ts, tsol)
221 call nf95_put_var(ncid, varid_tsoil, tsoil)
222 call nf95_put_var(ncid, varid_tslab, tslab)
223 call nf95_put_var(ncid, varid_SEAICE, seaice)
224 call nf95_put_var(ncid, varid_qs, qsurf)
225 call nf95_put_var(ncid, varid_QSOL, qsol)
226 call nf95_put_var(ncid, varid_albe, albedo)
227 call nf95_put_var(ncid, varid_evap, evap)
228 call nf95_put_var(ncid, varid_snow, snow)
229 call nf95_put_var(ncid, varid_rads, radsol)
230 call nf95_put_var(ncid, varid_solsw, solsw)
231 call nf95_put_var(ncid, varid_sollw, sollw)
232 call nf95_put_var(ncid, varid_fder, fder)
233 call nf95_put_var(ncid, varid_rain_f, rain_fall)
234 call nf95_put_var(ncid, varid_snow_f, snow_fall)
235 call nf95_put_var(ncid, varid_rug, frugs)
236 call nf95_put_var(ncid, varid_agesno, agesno)
237 call nf95_put_var(ncid, varid_zmea, zmea)
238 call nf95_put_var(ncid, varid_zstd, zstd)
239 call nf95_put_var(ncid, varid_zsig, zsig)
240 call nf95_put_var(ncid, varid_zgam, zgam)
241 call nf95_put_var(ncid, varid_zthe, zthe)
242 call nf95_put_var(ncid, varid_zpic, zpic)
243 call nf95_put_var(ncid, varid_zval, zval)
244 call nf95_put_var(ncid, varid_tancien, t_ancien)
245 call nf95_put_var(ncid, varid_qancien, q_ancien)
246 call nf95_put_var(ncid, varid_rugmer, frugs(:, is_oce))
247 call nf95_put_var(ncid, varid_clwcon, clwcon(:, 1))
248 call nf95_put_var(ncid, varid_rnebcon, rnebcon(:, 1))
249 call nf95_put_var(ncid, varid_ratqs, ratqs(:, 1))
250 call nf95_put_var(ncid, varid_run_off_lic_0, run_off_lic_0)
251 call nf95_put_var(ncid, varid_sig1, sig1)
252 call nf95_put_var(ncid, varid_w01, w01)
253
254 call nf95_close(ncid)
255
256 END SUBROUTINE phyredem
257
258 end module phyredem_m

  ViewVC Help
Powered by ViewVC 1.1.21