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

Annotation of /trunk/Sources/phylmd/phyredem.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 12 - (hide annotations)
Mon Jul 21 16:05:07 2008 UTC (15 years, 10 months ago) by guez
Original Path: trunk/libf/phylmd/phyredem.f90
File size: 13041 byte(s)
-- Minor modification of input/output:

Created procedure "read_logic". Variables of module "logic" are read
by "read_logic" instead of "conf_gcm". Variable "offline" of module
"conf_gcm" is read from namelist instead of "*.def".

Deleted arguments "dtime", "co2_ppm_etat0", "solaire_etat0",
"tabcntr0" and local variables "radpas", "tab_cntrl" of
"phyetat0". "phyetat0" does not read "controle" in "startphy.nc" any
longer. "phyetat0" now reads global attribute "itau_phy" from
"startphy.nc". "phyredem" does not create variable "controle" in
"startphy.nc" any longer. "phyredem" now writes global attribute
"itau_phy" of "startphy.nc". Deleted argument "tabcntr0" of
"printflag". Removed diagnostic messages written by "printflag" for
comparison of the variable "controle" of "startphy.nc" and the
variables read from "*.def" or namelist input.

-- Removing unwanted functionality:

Removed variable "lunout" from module "iniprint", replaced everywhere
by standard output.

Removed case "ocean == 'couple'" in "clmain", "interfsurf_hq" and
"physiq". Removed procedure "interfoce_cpl".

-- Should not change anything at run time:

Automated creation of graphs in documentation. More documentation on
input files.

Converted Fortran files to free format: "phyredem.f90", "printflag.f90".

Split module "clesphy" into "clesphys" and "clesphys2".

Removed variables "conser", "leapf", "forward", "apphys", "apdiss" and
"statcl" from module "logic". Added arguments "conser" to "advect",
"leapf" to "integrd". Added local variables "forward", "leapf",
"apphys", "conser", "apdiss" in "leapfrog".

Added intent attributes.

Deleted arguments "dtime" of "phyredem", "pdtime" of "flxdtdq", "sh"
of "phytrac", "dt" of "yamada".

Deleted local variables "dtime", "co2_ppm_etat0", "solaire_etat0",
"length", "tabcntr0" in "physiq". Replaced all references to "dtime"
by references to "pdtphys".

1 guez 12 SUBROUTINE phyredem(fichnom,radpas,rlat,rlon,pctsrf,tsol,tsoil, &
2     tslab,seaice,qsurf,qsol,snow,albedo,alblw,evap,rain_fall,snow_fall, &
3     solsw,sollw,fder,radsol,frugs,agesno,zmea,zstd,zsig,zgam,zthe,zpic, &
4     zval,rugsrel,t_ancien,q_ancien,rnebcon,ratqs,clwcon,run_off_lic_0)
5 guez 3
6 guez 12 ! From phylmd/phyredem.F,v 1.3 2005/05/25 13:10:09
7     ! Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818
8     ! Objet: Ecriture de l'etat de démarrage ou redémarrage pour la physique
9 guez 3
10 guez 12 USE indicesol, ONLY : is_lic, is_oce, is_sic, is_ter, nbsrf
11     USE dimphy, ONLY : klev, klon, zmasq
12     USE dimsoil, ONLY : nsoilmx
13     USE temps, ONLY : itau_phy
14     USE netcdf95, ONLY : nf95_create, nf95_put_att
15     USE netcdf, ONLY : nf90_clobber, nf90_global
16 guez 3
17 guez 12 IMPLICIT NONE
18 guez 3
19 guez 12 INCLUDE 'netcdf.inc'
20    
21     CHARACTER(len=*) fichnom
22     INTEGER :: radpas
23     REAL, INTENT (IN) :: rlat(klon), rlon(klon)
24     REAL :: tsol(klon,nbsrf)
25     REAL :: tsoil(klon,nsoilmx,nbsrf)
26     !IM "slab" ocean
27     REAL :: tslab(klon), seaice(klon) !IM "slab" ocean
28     REAL :: qsurf(klon,nbsrf)
29     REAL :: qsol(klon)
30     REAL :: snow(klon,nbsrf)
31     REAL :: albedo(klon,nbsrf)
32     !IM BEG
33     REAL :: alblw(klon,nbsrf)
34     !IM END
35     REAL :: evap(klon,nbsrf)
36     REAL :: rain_fall(klon)
37     REAL :: snow_fall(klon)
38     REAL :: solsw(klon)
39     REAL :: sollw(klon)
40     REAL :: fder(klon)
41     REAL :: radsol(klon)
42     REAL :: frugs(klon,nbsrf)
43     REAL :: agesno(klon,nbsrf)
44     REAL :: zmea(klon)
45     REAL :: zstd(klon)
46     REAL :: zsig(klon)
47     REAL :: zgam(klon)
48     REAL :: zthe(klon)
49     REAL :: zpic(klon)
50     REAL :: zval(klon)
51     REAL :: rugsrel(klon)
52     REAL :: pctsrf(klon,nbsrf)
53     REAL :: t_ancien(klon,klev), q_ancien(klon,klev)
54     REAL :: clwcon(klon,klev), rnebcon(klon,klev), ratqs(klon,klev)
55     REAL :: run_off_lic_0(klon)
56    
57     INTEGER :: nid, nvarid, idim2, idim3
58     INTEGER :: ierr
59    
60     INTEGER :: isoil, nsrf
61     CHARACTER (7) :: str7
62     CHARACTER (2) :: str2
63    
64     !------------------------------------------------------------
65    
66     PRINT *, 'Call sequence information: phyredem'
67     CALL nf95_create(fichnom,nf90_clobber,nid)
68    
69     call nf95_put_att(nid, nf90_global, 'title', 'Fichier redémarrage physique')
70     call nf95_put_att(nid, nf90_global, "itau_phy", itau_phy)
71    
72     ierr = nf_def_dim(nid,'points_physiques',klon,idim2)
73     ierr = nf_def_dim(nid,'horizon_vertical',klon*klev,idim3)
74    
75     ierr = nf_def_var(nid,'longitude',nf_float,1,idim2,nvarid)
76     ierr = nf_put_att_text(nid,nvarid,'title',32, &
77     'Longitudes de la grille physique')
78     ierr = nf_enddef(nid)
79     ierr = nf_put_var_real(nid,nvarid,rlon)
80    
81     ierr = nf_redef(nid)
82     ierr = nf_def_var(nid,'latitude',nf_float,1,idim2,nvarid)
83     ierr = nf_put_att_text(nid,nvarid,'title',31, &
84     'Latitudes de la grille physique')
85     ierr = nf_enddef(nid)
86     ierr = nf_put_var_real(nid,nvarid,rlat)
87    
88     ! PB ajout du masque terre/mer
89    
90     ierr = nf_redef(nid)
91     ierr = nf_def_var(nid,'masque',nf_float,1,idim2,nvarid)
92     ierr = nf_put_att_text(nid,nvarid,'title',16,'masque terre mer')
93     ierr = nf_enddef(nid)
94     ierr = nf_put_var_real(nid,nvarid,zmasq)
95     ! BP ajout des fraction de chaque sous-surface
96    
97     ! 1. fraction de terre
98    
99     ierr = nf_redef(nid)
100     ierr = nf_def_var(nid,'FTER',nf_float,1,idim2,nvarid)
101     ierr = nf_put_att_text(nid,nvarid,'title',21,'fraction de continent')
102     ierr = nf_enddef(nid)
103     ierr = nf_put_var_real(nid,nvarid,pctsrf(1:klon,is_ter))
104    
105     ! 2. Fraction de glace de terre
106    
107     ierr = nf_redef(nid)
108     ierr = nf_def_var(nid,'FLIC',nf_float,1,idim2,nvarid)
109     ierr = nf_put_att_text(nid,nvarid,'title',24,'fraction glace de terre')
110     ierr = nf_enddef(nid)
111     ierr = nf_put_var_real(nid,nvarid,pctsrf(1:klon,is_lic))
112    
113     ! 3. fraction ocean
114    
115     ierr = nf_redef(nid)
116     ierr = nf_def_var(nid,'FOCE',nf_float,1,idim2,nvarid)
117     ierr = nf_put_att_text(nid,nvarid,'title',14,'fraction ocean')
118     ierr = nf_enddef(nid)
119     ierr = nf_put_var_real(nid,nvarid,pctsrf(1:klon,is_oce))
120    
121     ! 4. Fraction glace de mer
122    
123     ierr = nf_redef(nid)
124     ierr = nf_def_var(nid,'FSIC',nf_float,1,idim2,nvarid)
125     ierr = nf_put_att_text(nid,nvarid,'title',18,'fraction glace mer')
126     ierr = nf_enddef(nid)
127     ierr = nf_put_var_real(nid,nvarid,pctsrf(1:klon,is_sic))
128    
129    
130    
131     DO nsrf = 1, nbsrf
132     IF (nsrf<=99) THEN
133     WRITE (str2,'(i2.2)') nsrf
134     ierr = nf_redef(nid)
135     ierr = nf_def_var(nid,'TS'//str2,nf_float,1,idim2,nvarid)
136     ierr = nf_put_att_text(nid,nvarid,'title',28, &
137     'Temperature de surface No.'//str2)
138     ierr = nf_enddef(nid)
139     ELSE
140     PRINT *, 'Trop de sous-mailles'
141     STOP 1
142     END IF
143     ierr = nf_put_var_real(nid,nvarid,tsol(1,nsrf))
144     END DO
145    
146     DO nsrf = 1, nbsrf
147     DO isoil = 1, nsoilmx
148     IF (isoil<=99 .AND. nsrf<=99) THEN
149     WRITE (str7,'(i2.2,"srf",i2.2)') isoil, nsrf
150     ierr = nf_redef(nid)
151     ierr = nf_def_var(nid,'Tsoil'//str7,nf_float,1,idim2,nvarid)
152     ierr = nf_put_att_text(nid,nvarid,'title',29, &
153     'Temperature du sol No.'//str7)
154     ierr = nf_enddef(nid)
155 guez 3 ELSE
156 guez 12 PRINT *, 'Trop de couches'
157     STOP 1
158     END IF
159     ierr = nf_put_var_real(nid,nvarid,tsoil(1,isoil,nsrf))
160     END DO
161     END DO
162 guez 3
163 guez 12 !IM "slab" ocean
164     ierr = nf_redef(nid)
165     ierr = nf_def_var(nid,'TSLAB',nf_float,1,idim2,nvarid)
166     ierr = nf_put_att_text(nid,nvarid,'title',33, &
167     'Ecart de la SST (pour slab-ocean)')
168     ierr = nf_enddef(nid)
169     ierr = nf_put_var_real(nid,nvarid,tslab)
170    
171     ierr = nf_redef(nid)
172     ierr = nf_def_var(nid,'SEAICE',nf_float,1,idim2,nvarid)
173     ierr = nf_put_att_text(nid,nvarid,'title',33, &
174     'Glace de mer kg/m2 (pour slab-ocean)')
175     ierr = nf_enddef(nid)
176     ierr = nf_put_var_real(nid,nvarid,seaice)
177    
178     DO nsrf = 1, nbsrf
179     IF (nsrf<=99) THEN
180     WRITE (str2,'(i2.2)') nsrf
181     ierr = nf_redef(nid)
182     ierr = nf_def_var(nid,'QS'//str2,nf_float,1,idim2,nvarid)
183     ierr = nf_put_att_text(nid,nvarid,'title',25, &
184     'Humidite de surface No.'//str2)
185     ierr = nf_enddef(nid)
186     ELSE
187     PRINT *, 'Trop de sous-mailles'
188     STOP 1
189     END IF
190     ierr = nf_put_var_real(nid,nvarid,qsurf(1,nsrf))
191     END DO
192    
193     ierr = nf_redef(nid)
194     ierr = nf_def_var(nid,'QSOL',nf_float,1,idim2,nvarid)
195     ierr = nf_put_att_text(nid,nvarid,'title',20,'Eau dans le sol (mm)')
196     ierr = nf_enddef(nid)
197     ierr = nf_put_var_real(nid,nvarid,qsol)
198    
199     DO nsrf = 1, nbsrf
200     IF (nsrf<=99) THEN
201     WRITE (str2,'(i2.2)') nsrf
202     ierr = nf_redef(nid)
203     ierr = nf_def_var(nid,'ALBE'//str2,nf_float,1,idim2,nvarid)
204     ierr = nf_put_att_text(nid,nvarid,'title',23, &
205     'albedo de surface No.'//str2)
206     ierr = nf_enddef(nid)
207     ELSE
208     PRINT *, 'Trop de sous-mailles'
209     STOP 1
210     END IF
211     ierr = nf_put_var_real(nid,nvarid,albedo(1,nsrf))
212     END DO
213    
214     !IM BEG albedo LW
215     DO nsrf = 1, nbsrf
216     IF (nsrf<=99) THEN
217     WRITE (str2,'(i2.2)') nsrf
218     ierr = nf_redef(nid)
219     ierr = nf_def_var(nid,'ALBLW'//str2,nf_float,1,idim2,nvarid)
220     ierr = nf_put_att_text(nid,nvarid,'title',23, &
221     'albedo LW de surface No.'//str2)
222     ierr = nf_enddef(nid)
223     ELSE
224     PRINT *, 'Trop de sous-mailles'
225     STOP 1
226     END IF
227     ierr = nf_put_var_real(nid,nvarid,alblw(1,nsrf))
228     END DO
229     !IM END albedo LW
230    
231     DO nsrf = 1, nbsrf
232     IF (nsrf<=99) THEN
233     WRITE (str2,'(i2.2)') nsrf
234     ierr = nf_redef(nid)
235     ierr = nf_def_var(nid,'EVAP'//str2,nf_float,1,idim2,nvarid)
236     ierr = nf_put_att_text(nid,nvarid,'title',28, &
237     'Evaporation de surface No.'//str2)
238     ierr = nf_enddef(nid)
239     ELSE
240     PRINT *, 'Trop de sous-mailles'
241     STOP 1
242     END IF
243     ierr = nf_put_var_real(nid,nvarid,evap(1,nsrf))
244     END DO
245    
246    
247     DO nsrf = 1, nbsrf
248     IF (nsrf<=99) THEN
249     WRITE (str2,'(i2.2)') nsrf
250     ierr = nf_redef(nid)
251     ierr = nf_def_var(nid,'SNOW'//str2,nf_float,1,idim2,nvarid)
252     ierr = nf_put_att_text(nid,nvarid,'title',22, &
253     'Neige de surface No.'//str2)
254     ierr = nf_enddef(nid)
255     ELSE
256     PRINT *, 'Trop de sous-mailles'
257     STOP 1
258     END IF
259     ierr = nf_put_var_real(nid,nvarid,snow(1,nsrf))
260     END DO
261    
262    
263     ierr = nf_redef(nid)
264     ierr = nf_def_var(nid,'RADS',nf_float,1,idim2,nvarid)
265     ierr = nf_put_att_text(nid,nvarid,'title',28, &
266     'Rayonnement net a la surface')
267     ierr = nf_enddef(nid)
268     ierr = nf_put_var_real(nid,nvarid,radsol)
269    
270     ierr = nf_redef(nid)
271     ierr = nf_def_var(nid,'solsw',nf_float,1,idim2,nvarid)
272     ierr = nf_put_att_text(nid,nvarid,'title',32, &
273     'Rayonnement solaire a la surface')
274     ierr = nf_enddef(nid)
275     ierr = nf_put_var_real(nid,nvarid,solsw)
276    
277     ierr = nf_redef(nid)
278     ierr = nf_def_var(nid,'sollw',nf_float,1,idim2,nvarid)
279     ierr = nf_put_att_text(nid,nvarid,'title',27, &
280     'Rayonnement IF a la surface')
281     ierr = nf_enddef(nid)
282     ierr = nf_put_var_real(nid,nvarid,sollw)
283    
284     ierr = nf_redef(nid)
285     ierr = nf_def_var(nid,'fder',nf_float,1,idim2,nvarid)
286     ierr = nf_put_att_text(nid,nvarid,'title',14,'Derive de flux')
287     ierr = nf_enddef(nid)
288     ierr = nf_put_var_real(nid,nvarid,fder)
289    
290     ierr = nf_redef(nid)
291     ierr = nf_def_var(nid,'rain_f',nf_float,1,idim2,nvarid)
292     ierr = nf_put_att_text(nid,nvarid,'title',21,'precipitation liquide')
293     ierr = nf_enddef(nid)
294     ierr = nf_put_var_real(nid,nvarid,rain_fall)
295    
296     ierr = nf_redef(nid)
297     ierr = nf_def_var(nid,'snow_f',nf_float,1,idim2,nvarid)
298     ierr = nf_put_att_text(nid,nvarid,'title',20,'precipitation solide')
299     ierr = nf_enddef(nid)
300     ierr = nf_put_var_real(nid,nvarid,snow_fall)
301    
302     DO nsrf = 1, nbsrf
303     IF (nsrf<=99) THEN
304     WRITE (str2,'(i2.2)') nsrf
305     ierr = nf_redef(nid)
306     ierr = nf_def_var(nid,'RUG'//str2,nf_float,1,idim2,nvarid)
307     ierr = nf_put_att_text(nid,nvarid,'title',23, &
308     'rugosite de surface No.'//str2)
309     ierr = nf_enddef(nid)
310     ELSE
311     PRINT *, 'Trop de sous-mailles'
312     STOP 1
313     END IF
314     ierr = nf_put_var_real(nid,nvarid,frugs(1,nsrf))
315     END DO
316    
317     DO nsrf = 1, nbsrf
318     IF (nsrf<=99) THEN
319     WRITE (str2,'(i2.2)') nsrf
320     ierr = nf_redef(nid)
321     ierr = nf_def_var(nid,'AGESNO'//str2,nf_float,1,idim2,nvarid)
322     ierr = nf_put_att_text(nid,nvarid,'title',15, &
323     'Age de la neige surface No.'//str2)
324     ierr = nf_enddef(nid)
325     ELSE
326     PRINT *, 'Trop de sous-mailles'
327     STOP 1
328     END IF
329     ierr = nf_put_var_real(nid,nvarid,agesno(1,nsrf))
330     END DO
331    
332     ierr = nf_redef(nid)
333     ierr = nf_def_var(nid,'ZMEA',nf_float,1,idim2,nvarid)
334     ierr = nf_enddef(nid)
335     ierr = nf_put_var_real(nid,nvarid,zmea)
336    
337     ierr = nf_redef(nid)
338     ierr = nf_def_var(nid,'ZSTD',nf_float,1,idim2,nvarid)
339     ierr = nf_enddef(nid)
340     ierr = nf_put_var_real(nid,nvarid,zstd)
341     ierr = nf_redef(nid)
342     ierr = nf_def_var(nid,'ZSIG',nf_float,1,idim2,nvarid)
343     ierr = nf_enddef(nid)
344     ierr = nf_put_var_real(nid,nvarid,zsig)
345     ierr = nf_redef(nid)
346     ierr = nf_def_var(nid,'ZGAM',nf_float,1,idim2,nvarid)
347     ierr = nf_enddef(nid)
348     ierr = nf_put_var_real(nid,nvarid,zgam)
349     ierr = nf_redef(nid)
350     ierr = nf_def_var(nid,'ZTHE',nf_float,1,idim2,nvarid)
351     ierr = nf_enddef(nid)
352     ierr = nf_put_var_real(nid,nvarid,zthe)
353     ierr = nf_redef(nid)
354     ierr = nf_def_var(nid,'ZPIC',nf_float,1,idim2,nvarid)
355     ierr = nf_enddef(nid)
356     ierr = nf_put_var_real(nid,nvarid,zpic)
357     ierr = nf_redef(nid)
358     ierr = nf_def_var(nid,'ZVAL',nf_float,1,idim2,nvarid)
359     ierr = nf_enddef(nid)
360     ierr = nf_put_var_real(nid,nvarid,zval)
361     ierr = nf_redef(nid)
362     ierr = nf_def_var(nid,'RUGSREL',nf_float,1,idim2,nvarid)
363     ierr = nf_enddef(nid)
364     ierr = nf_put_var_real(nid,nvarid,rugsrel)
365    
366     ierr = nf_redef(nid)
367     ierr = nf_def_var(nid,'TANCIEN',nf_float,1,idim3,nvarid)
368     ierr = nf_enddef(nid)
369     ierr = nf_put_var_real(nid,nvarid,t_ancien)
370    
371     ierr = nf_redef(nid)
372     ierr = nf_def_var(nid,'QANCIEN',nf_float,1,idim3,nvarid)
373     ierr = nf_enddef(nid)
374     ierr = nf_put_var_real(nid,nvarid,q_ancien)
375    
376     ierr = nf_redef(nid)
377     ierr = nf_def_var(nid,'RUGMER',nf_float,1,idim2,nvarid)
378     ierr = nf_put_att_text(nid,nvarid,'title',28, &
379     'Longueur de rugosite sur mer')
380     ierr = nf_enddef(nid)
381     ierr = nf_put_var_real(nid,nvarid,frugs(1,is_oce))
382    
383     ierr = nf_redef(nid)
384     ierr = nf_def_var(nid,'CLWCON',nf_float,1,idim2,nvarid)
385     ierr = nf_put_att_text(nid,nvarid,'title',28,'Eau liquide convective')
386     ierr = nf_enddef(nid)
387     ierr = nf_put_var_real(nid,nvarid,clwcon)
388    
389     ierr = nf_redef(nid)
390     ierr = nf_def_var(nid,'RNEBCON',nf_float,1,idim2,nvarid)
391     ierr = nf_put_att_text(nid,nvarid,'title',28,'Nebulosite convective')
392     ierr = nf_enddef(nid)
393     ierr = nf_put_var_real(nid,nvarid,rnebcon)
394    
395     ierr = nf_redef(nid)
396     ierr = nf_def_var(nid,'RATQS',nf_float,1,idim2,nvarid)
397     ierr = nf_put_att_text(nid,nvarid,'title',5,'Ratqs')
398     ierr = nf_enddef(nid)
399     ierr = nf_put_var_real(nid,nvarid,ratqs)
400    
401     ! run_off_lic_0
402    
403     ierr = nf_redef(nid)
404     ierr = nf_def_var(nid,'RUNOFFLIC0',nf_float,1,idim2,nvarid)
405     ierr = nf_put_att_text(nid,nvarid,'title',10,'Runofflic0')
406     ierr = nf_enddef(nid)
407     ierr = nf_put_var_real(nid,nvarid,run_off_lic_0)
408    
409    
410     ierr = nf_close(nid)
411    
412     END SUBROUTINE phyredem

  ViewVC Help
Powered by ViewVC 1.1.21