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

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

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

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

Legend:
Removed from v.3  
changed lines
  Added in v.156

  ViewVC Help
Powered by ViewVC 1.1.21