/[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

revision 156 by guez, Thu Jul 16 17:39:10 2015 UTC revision 207 by guez, Thu Sep 1 10:30:53 2016 UTC
# Line 4  module phyredem_m Line 4  module phyredem_m
4    
5  contains  contains
6    
7    SUBROUTINE phyredem(fichnom, pctsrf, tsol, tsoil, tslab, seaice, qsurf, &    SUBROUTINE phyredem(pctsrf, ftsol, ftsoil, qsurf, qsol, snow, albedo, evap, &
8         qsol, snow, albedo, evap, rain_fall, snow_fall, solsw, sollw, &         rain_fall, snow_fall, solsw, sollw, fder, radsol, frugs, agesno, zmea, &
9         fder, radsol, frugs, agesno, zmea, zstd, zsig, zgam, zthe, zpic, zval, &         zstd, zsig, zgam, zthe, zpic, zval, t_ancien, q_ancien, rnebcon, &
10         t_ancien, q_ancien, rnebcon, ratqs, clwcon, run_off_lic_0, sig1, w01)         ratqs, clwcon, run_off_lic_0, sig1, w01)
11    
12      ! From phylmd/phyredem.F, version 1.3, 2005/05/25 13:10:09      ! From phylmd/phyredem.F, version 1.3, 2005/05/25 13:10:09
13      ! Author: Z. X. Li (LMD/CNRS)      ! Author: Z. X. Li (LMD/CNRS)
# Line 17  contains Line 17  contains
17      ! pour la physique      ! pour la physique
18    
19      USE dimphy, ONLY: klev, klon, zmasq      USE dimphy, ONLY: klev, klon, zmasq
     USE dimsoil, ONLY: nsoilmx  
20      USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf      USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf
21      USE netcdf, ONLY: nf90_clobber, nf90_global, nf90_float      USE netcdf95, ONLY: nf95_inq_varid, nf95_put_var, nf95_close
22      USE netcdf95, ONLY: nf95_create, nf95_put_att, nf95_def_dim, &      use phyredem0_m, only: ncid_restartphy
          nf95_def_var, nf95_enddef, nf95_put_var, nf95_close  
     use phyetat0_m, only: rlat, rlon  
     USE temps, ONLY: itau_phy  
23    
     CHARACTER(len=*), INTENT(IN):: fichnom  
24      REAL, INTENT(IN):: pctsrf(:, :) ! (klon, nbsrf)      REAL, INTENT(IN):: pctsrf(:, :) ! (klon, nbsrf)
25      REAL, INTENT(IN):: tsol(:, :) ! (klon, nbsrf)      REAL, INTENT(IN):: ftsol(:, :) ! (klon, nbsrf)
26      REAL, INTENT(IN):: tsoil(:, :, :) ! (klon, nsoilmx, nbsrf)      REAL, INTENT(IN):: ftsoil(:, :, :) ! (klon, nsoilmx, nbsrf)
     REAL, INTENT(IN):: tslab(:), seaice(:) ! (klon) slab ocean  
27      REAL, INTENT(IN):: qsurf(:, :) ! (klon, nbsrf)      REAL, INTENT(IN):: qsurf(:, :) ! (klon, nbsrf)
28    
29      REAL, intent(in):: qsol(:) ! (klon)      REAL, intent(in):: qsol(:) ! (klon)
# Line 63  contains Line 57  contains
57      ! vertical velocity within adiabatic updraft      ! vertical velocity within adiabatic updraft
58    
59      ! Local:      ! Local:
60        integer varid
     INTEGER ncid, idim2, idim3, dimid_nbsrf, dimid_nsoilmx  
     integer varid, varid_run_off_lic_0, varid_sig1, varid_w01, varid_rlon  
     integer varid_rlat, varid_zmasq, varid_fter, varid_flic, varid_foce  
     integer varid_fsic, varid_ts, varid_tsoil, varid_tslab, varid_SEAICE  
     integer varid_qs, varid_QSOL, varid_albe, varid_evap, varid_snow  
     integer varid_rads, varid_solsw, varid_sollw, varid_fder, varid_rain_f  
     integer varid_snow_f, varid_rug, varid_agesno, varid_zmea, varid_zstd  
     integer varid_zsig, varid_zgam, varid_zthe, varid_zpic, varid_zval  
     integer varid_tancien, varid_qancien, varid_rugmer, varid_clwcon  
     integer varid_rnebcon, varid_ratqs  
61    
62      !------------------------------------------------------------      !------------------------------------------------------------
63    
64      PRINT *, 'Call sequence information: phyredem'      PRINT *, 'Call sequence information: phyredem'
     CALL nf95_create(fichnom, nf90_clobber, ncid)  
65    
66      call nf95_put_att(ncid, nf90_global, 'title', &      call nf95_inq_varid(ncid_restartphy, "masque", varid)
67           'start file for the physics code')      call nf95_put_var(ncid_restartphy, varid, zmasq)
68      call nf95_put_att(ncid, nf90_global, "itau_phy", itau_phy)  
69        call nf95_inq_varid(ncid_restartphy, "FTER", varid)
70      call nf95_def_dim(ncid, 'points_physiques', klon, idim2)      call nf95_put_var(ncid_restartphy, varid, pctsrf(:, is_ter))
71      call nf95_def_dim(ncid, 'klev', klev, idim3)  
72      call nf95_def_dim(ncid, 'nbsrf', nbsrf, dimid_nbsrf)      call nf95_inq_varid(ncid_restartphy, "FLIC", varid)
73      call nf95_def_dim(ncid, 'nsoilmx', nsoilmx, dimid_nsoilmx)      call nf95_put_var(ncid_restartphy, varid, pctsrf(:, is_lic))
74    
75      call nf95_def_var(ncid, 'longitude', nf90_float, idim2, varid_rlon)      call nf95_inq_varid(ncid_restartphy, "FOCE", varid)
76      call nf95_def_var(ncid, 'latitude', nf90_float, idim2, varid_rlat)      call nf95_put_var(ncid_restartphy, varid, pctsrf(:, is_oce))
77    
78      call nf95_def_var(ncid, 'masque', nf90_float, idim2, varid_zmasq)      call nf95_inq_varid(ncid_restartphy, "FSIC", varid)
79      call nf95_put_att(ncid, varid_zmasq, 'title', 'masque terre mer')      call nf95_put_var(ncid_restartphy, varid, pctsrf(:, is_sic))
80    
81      ! Fractions de chaque sous-surface      call nf95_inq_varid(ncid_restartphy, "TS", varid)
82        call nf95_put_var(ncid_restartphy, varid, ftsol)
83      call nf95_def_var(ncid, 'FTER', nf90_float, idim2, varid_fter)  
84      call nf95_put_att(ncid, varid_fter, 'title', 'fraction de continent')      call nf95_inq_varid(ncid_restartphy, "Tsoil", varid)
85        call nf95_put_var(ncid_restartphy, varid, ftsoil)
86      call nf95_def_var(ncid, 'FLIC', nf90_float, idim2, varid_flic)  
87      call nf95_put_att(ncid, varid_flic, 'title', 'fraction glace de terre')      call nf95_inq_varid(ncid_restartphy, "QS", varid)
88        call nf95_put_var(ncid_restartphy, varid, qsurf)
89      call nf95_def_var(ncid, 'FOCE', nf90_float, idim2, varid_foce)  
90      call nf95_put_att(ncid, varid_foce, 'title', 'fraction ocean')      call nf95_inq_varid(ncid_restartphy, "QSOL", varid)
91        call nf95_put_var(ncid_restartphy, varid, qsol)
92      call nf95_def_var(ncid, 'FSIC', nf90_float, idim2, varid_fsic)  
93      call nf95_put_att(ncid, varid_fsic, 'title', 'fraction glace mer')      call nf95_inq_varid(ncid_restartphy, "ALBE", varid)
94        call nf95_put_var(ncid_restartphy, varid, albedo)
95      call nf95_def_var(ncid, 'TS', nf90_float, (/idim2, dimid_nbsrf/), varid_ts)  
96      call nf95_put_att(ncid, varid_ts, 'title', 'surface temperature')      call nf95_inq_varid(ncid_restartphy, "EVAP", varid)
97        call nf95_put_var(ncid_restartphy, varid, evap)
98      call nf95_def_var(ncid, 'Tsoil', nf90_float, (/idim2, dimid_nsoilmx, &  
99           dimid_nbsrf/), varid_tsoil)      call nf95_inq_varid(ncid_restartphy, "SNOW", varid)
100      call nf95_put_att(ncid, varid_tsoil, 'title', 'soil temperature')      call nf95_put_var(ncid_restartphy, varid, snow)
101    
102      call nf95_def_var(ncid, 'TSLAB', nf90_float, idim2, varid_tslab)      call nf95_inq_varid(ncid_restartphy, "RADS", varid)
103      call nf95_put_att(ncid, varid_tslab, 'title', &      call nf95_put_var(ncid_restartphy, varid, radsol)
104           'Ecart de la SST (pour slab-ocean)')  
105        call nf95_inq_varid(ncid_restartphy, "solsw", varid)
106      call nf95_def_var(ncid, 'SEAICE', nf90_float, idim2, varid_SEAICE)      call nf95_put_var(ncid_restartphy, varid, solsw)
107      call nf95_put_att(ncid, varid_SEAICE, 'title', &  
108           'Glace de mer kg/m2 (pour slab-ocean)')      call nf95_inq_varid(ncid_restartphy, "sollw", varid)
109        call nf95_put_var(ncid_restartphy, varid, sollw)
110      call nf95_def_var(ncid, 'QS', nf90_float, (/idim2, dimid_nbsrf/), varid_qs)  
111      call nf95_put_att(ncid, varid, 'title', 'Humidite de surface')      call nf95_inq_varid(ncid_restartphy, "fder", varid)
112        call nf95_put_var(ncid_restartphy, varid, fder)
113      call nf95_def_var(ncid, 'QSOL', nf90_float, idim2, varid_QSOL)  
114      call nf95_put_att(ncid, varid_QSOL, 'title', 'Eau dans le sol (mm)')      call nf95_inq_varid(ncid_restartphy, "rain_f", varid)
115        call nf95_put_var(ncid_restartphy, varid, rain_fall)
116      call nf95_def_var(ncid, 'ALBE', nf90_float, (/idim2, dimid_nbsrf/), &  
117           varid_albe)      call nf95_inq_varid(ncid_restartphy, "snow_f", varid)
118      call nf95_put_att(ncid, varid_albe, 'title', 'albedo de surface')      call nf95_put_var(ncid_restartphy, varid, snow_fall)
119    
120      call nf95_def_var(ncid, 'EVAP', nf90_float, (/idim2, dimid_nbsrf/), &      call nf95_inq_varid(ncid_restartphy, "RUG", varid)
121           varid_evap)      call nf95_put_var(ncid_restartphy, varid, frugs)
122      call nf95_put_att(ncid, varid_evap, 'title', 'Evaporation de surface')  
123        call nf95_inq_varid(ncid_restartphy, "AGESNO", varid)
124      call nf95_def_var(ncid, 'SNOW', nf90_float, (/idim2, dimid_nbsrf/), &      call nf95_put_var(ncid_restartphy, varid, agesno)
125           varid_snow)  
126      call nf95_put_att(ncid, varid_snow, 'title', 'Neige de surface')      call nf95_inq_varid(ncid_restartphy, "ZMEA", varid)
127        call nf95_put_var(ncid_restartphy, varid, zmea)
128      call nf95_def_var(ncid, 'RADS', nf90_float, idim2, varid_rads)  
129      call nf95_put_att(ncid, varid_rads, 'title', &      call nf95_inq_varid(ncid_restartphy, "ZSTD", varid)
130           'Rayonnement net a la surface')      call nf95_put_var(ncid_restartphy, varid, zstd)
131    
132      call nf95_def_var(ncid, 'solsw', nf90_float, idim2, varid_solsw)      call nf95_inq_varid(ncid_restartphy, "ZSIG", varid)
133      call nf95_put_att(ncid, varid_solsw, 'title', &      call nf95_put_var(ncid_restartphy, varid, zsig)
134           'Rayonnement solaire a la surface')  
135        call nf95_inq_varid(ncid_restartphy, "ZGAM", varid)
136      call nf95_def_var(ncid, 'sollw', nf90_float, idim2, varid_sollw)      call nf95_put_var(ncid_restartphy, varid, zgam)
137      call nf95_put_att(ncid, varid_sollw, 'title', &  
138           'Rayonnement IF a la surface')      call nf95_inq_varid(ncid_restartphy, "ZTHE", varid)
139        call nf95_put_var(ncid_restartphy, varid, zthe)
140      call nf95_def_var(ncid, 'fder', nf90_float, idim2, varid_fder)  
141      call nf95_put_att(ncid, varid_fder, 'title', 'Derive de flux')      call nf95_inq_varid(ncid_restartphy, "ZPIC", varid)
142        call nf95_put_var(ncid_restartphy, varid, zpic)
143      call nf95_def_var(ncid, 'rain_f', nf90_float, idim2, varid_rain_f)  
144      call nf95_put_att(ncid, varid_rain_f, 'title', 'precipitation liquide')      call nf95_inq_varid(ncid_restartphy, "ZVAL", varid)
145        call nf95_put_var(ncid_restartphy, varid, zval)
146      call nf95_def_var(ncid, 'snow_f', nf90_float, idim2, varid_snow_f)  
147      call nf95_put_att(ncid, varid_snow_f, 'title', 'precipitation solide')      call nf95_inq_varid(ncid_restartphy, "TANCIEN", varid)
148        call nf95_put_var(ncid_restartphy, varid, t_ancien)
149      call nf95_def_var(ncid, 'RUG', nf90_float, (/idim2, dimid_nbsrf/), &  
150           varid_rug)      call nf95_inq_varid(ncid_restartphy, "QANCIEN", varid)
151      call nf95_put_att(ncid, varid_rug, 'title', 'rugosite de surface')      call nf95_put_var(ncid_restartphy, varid, q_ancien)
152    
153      call nf95_def_var(ncid, 'AGESNO', nf90_float, (/idim2, dimid_nbsrf/), &      call nf95_inq_varid(ncid_restartphy, "RUGMER", varid)
154           varid_agesno)      call nf95_put_var(ncid_restartphy, varid, frugs(:, is_oce))
155      call nf95_put_att(ncid, varid_agesno, 'title', 'Age de la neige surface')  
156        call nf95_inq_varid(ncid_restartphy, "CLWCON", varid)
157      call nf95_def_var(ncid, 'ZMEA', nf90_float, idim2, varid_zmea)      call nf95_put_var(ncid_restartphy, varid, clwcon(:, 1))
158      call nf95_def_var(ncid, 'ZSTD', nf90_float, idim2, varid_zstd)  
159      call nf95_def_var(ncid, 'ZSIG', nf90_float, idim2, varid_zsig)      call nf95_inq_varid(ncid_restartphy, "RNEBCON", varid)
160      call nf95_def_var(ncid, 'ZGAM', nf90_float, idim2, varid_zgam)      call nf95_put_var(ncid_restartphy, varid, rnebcon(:, 1))
161      call nf95_def_var(ncid, 'ZTHE', nf90_float, idim2, varid_zthe)  
162      call nf95_def_var(ncid, 'ZPIC', nf90_float, idim2, varid_zpic)      call nf95_inq_varid(ncid_restartphy, "RATQS", varid)
163      call nf95_def_var(ncid, 'ZVAL', nf90_float, idim2, varid_zval)      call nf95_put_var(ncid_restartphy, varid, ratqs(:, 1))
164      call nf95_def_var(ncid, 'TANCIEN', nf90_float, (/idim2, idim3/), &  
165           varid_tancien)      call nf95_inq_varid(ncid_restartphy, "RUNOFFLIC0", varid)
166      call nf95_def_var(ncid, 'QANCIEN', nf90_float, (/idim2, idim3/), &      call nf95_put_var(ncid_restartphy, varid, run_off_lic_0)
167           varid_qancien)  
168        call nf95_inq_varid(ncid_restartphy, "sig1", varid)
169      call nf95_def_var(ncid, 'RUGMER', nf90_float, idim2, varid_rugmer)      call nf95_put_var(ncid_restartphy, varid, sig1)
170      call nf95_put_att(ncid, varid_rugmer, 'title', &  
171           'Longueur de rugosite sur mer')      call nf95_inq_varid(ncid_restartphy, "w01", varid)
172        call nf95_put_var(ncid_restartphy, varid, w01)
     call nf95_def_var(ncid, 'CLWCON', nf90_float, idim2, varid_clwcon)  
     call nf95_put_att(ncid, varid_clwcon, 'title', 'Eau liquide convective')  
   
     call nf95_def_var(ncid, 'RNEBCON', nf90_float, idim2, varid_rnebcon)  
     call nf95_put_att(ncid, varid_rnebcon, 'title', 'Nebulosite convective')  
   
     call nf95_def_var(ncid, 'RATQS', nf90_float, idim2, varid_ratqs)  
     call nf95_put_att(ncid, varid_ratqs, 'title', 'Ratqs')  
   
     call nf95_def_var(ncid, 'RUNOFFLIC0', nf90_float, idim2, &  
          varid_run_off_lic_0)  
     call nf95_put_att(ncid, varid_run_off_lic_0, 'title', 'Runofflic0')  
   
     call nf95_def_var(ncid, 'sig1', nf90_float, (/idim2, idim3/), varid_sig1)  
     call nf95_put_att(ncid, varid_sig1, 'long_name', &  
          'section adiabatic updraft')  
   
     call nf95_def_var(ncid, 'w01', nf90_float, (/idim2, idim3/), varid_w01)  
     call nf95_put_att(ncid, varid_w01, 'long_name', &  
          'vertical velocity within adiabatic updraft')  
   
     call nf95_enddef(ncid)  
   
     call nf95_put_var(ncid, varid_rlon, rlon)  
     call nf95_put_var(ncid, varid_rlat, rlat)  
     call nf95_put_var(ncid, varid_zmasq, zmasq)  
     call nf95_put_var(ncid, varid_fter, pctsrf(:, is_ter))  
     call nf95_put_var(ncid, varid_flic, pctsrf(:, is_lic))  
     call nf95_put_var(ncid, varid_foce, pctsrf(:, is_oce))  
     call nf95_put_var(ncid, varid_fsic, pctsrf(:, is_sic))  
     call nf95_put_var(ncid, varid_ts, tsol)  
     call nf95_put_var(ncid, varid_tsoil, tsoil)  
     call nf95_put_var(ncid, varid_tslab, tslab)  
     call nf95_put_var(ncid, varid_SEAICE, seaice)  
     call nf95_put_var(ncid, varid_qs, qsurf)  
     call nf95_put_var(ncid, varid_QSOL, qsol)  
     call nf95_put_var(ncid, varid_albe, albedo)  
     call nf95_put_var(ncid, varid_evap, evap)  
     call nf95_put_var(ncid, varid_snow, snow)  
     call nf95_put_var(ncid, varid_rads, radsol)  
     call nf95_put_var(ncid, varid_solsw, solsw)  
     call nf95_put_var(ncid, varid_sollw, sollw)  
     call nf95_put_var(ncid, varid_fder, fder)  
     call nf95_put_var(ncid, varid_rain_f, rain_fall)  
     call nf95_put_var(ncid, varid_snow_f, snow_fall)  
     call nf95_put_var(ncid, varid_rug, frugs)  
     call nf95_put_var(ncid, varid_agesno, agesno)  
     call nf95_put_var(ncid, varid_zmea, zmea)  
     call nf95_put_var(ncid, varid_zstd, zstd)  
     call nf95_put_var(ncid, varid_zsig, zsig)  
     call nf95_put_var(ncid, varid_zgam, zgam)  
     call nf95_put_var(ncid, varid_zthe, zthe)  
     call nf95_put_var(ncid, varid_zpic, zpic)  
     call nf95_put_var(ncid, varid_zval, zval)  
     call nf95_put_var(ncid, varid_tancien, t_ancien)  
     call nf95_put_var(ncid, varid_qancien, q_ancien)  
     call nf95_put_var(ncid, varid_rugmer, frugs(:, is_oce))  
     call nf95_put_var(ncid, varid_clwcon, clwcon(:, 1))  
     call nf95_put_var(ncid, varid_rnebcon, rnebcon(:, 1))  
     call nf95_put_var(ncid, varid_ratqs, ratqs(:, 1))  
     call nf95_put_var(ncid, varid_run_off_lic_0, run_off_lic_0)  
     call nf95_put_var(ncid, varid_sig1, sig1)  
     call nf95_put_var(ncid, varid_w01, w01)  
173    
174      call nf95_close(ncid)      call nf95_close(ncid_restartphy)
175    
176    END SUBROUTINE phyredem    END SUBROUTINE phyredem
177    

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

  ViewVC Help
Powered by ViewVC 1.1.21