/[lmdze]/trunk/Sources/phylmd/Interface_surf/read_sst.f
ViewVC logotype

Annotation of /trunk/Sources/phylmd/Interface_surf/read_sst.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 221 - (hide annotations)
Thu Apr 20 14:44:47 2017 UTC (7 years, 2 months ago) by guez
File size: 1275 byte(s)
clcdrag is no longer used in LMDZ. Replaced by cdrag in LMDZ. In cdrag
in LMDZ, zxli is a symbolic constant, false. So removed case zxli true
in LMDZE.

read_sst is called zero (if no ocean point on the whole planet) time or
once per call of physiq. If mod(itap - 1, lmt_pas) == 0 then we have
advanced in time of lmt_pas and deja_lu is necessarily false.

qsat[sl] and dqsat[sl] were never called.

Added output of qsurf in histins, following LMDZ.

Last dummy argument dtime of phystokenc is always the same as first
dummy argument pdtphys, removed dtime.

Removed make rules for nag_xref95, since it does not exist any longer.

1 guez 202 module read_sst_m
2    
3     implicit none
4    
5     contains
6    
7 guez 221 SUBROUTINE read_sst(julien, knindex, lmt_sst)
8 guez 202
9     ! From interfoce_lim
10    
11 guez 221 use conf_gcm_m, only: lmt_pas
12 guez 202 USE dimphy, ONLY: klon
13     USE netcdf, ONLY: nf90_nowrite
14     use netcdf95, only: NF95_CLOSE, nf95_get_var, NF95_INQ_VARID, nf95_open
15     use nr_util, only: assert
16     use time_phylmdz, only: itap
17    
18 guez 221 integer, intent(IN):: julien ! jour a lire dans l'annee
19 guez 202
20     integer, intent(in):: knindex(:) ! (knon)
21     ! index des points de la surface a traiter
22    
23     real, intent(out):: lmt_sst(:) ! (knon)
24     ! SST lues dans le fichier de conditions aux limites
25    
26     ! Local:
27    
28     ! Champ lu dans le fichier de conditions aux limites :
29 guez 221 real, save:: sst_lu(klon)
30 guez 202
31     integer ncid, varid ! pour NetCDF
32    
33     ! --------------------------------------------------
34    
35     call assert(size(knindex) == size(lmt_sst), "read_sst knon")
36    
37     ! Tester d'abord si c'est le moment de lire le fichier
38 guez 221 if (mod(itap - 1, lmt_pas) == 0) then
39 guez 202 call NF95_OPEN ('limit.nc', NF90_NOWRITE, ncid)
40    
41     call NF95_INQ_VARID(ncid, 'SST', varid)
42 guez 221 call NF95_GET_VAR(ncid, varid, sst_lu, start = (/1, julien/))
43 guez 202
44     call NF95_CLOSE(ncid)
45     endif
46    
47     lmt_sst = sst_lu(knindex)
48    
49     END SUBROUTINE read_sst
50    
51     end module read_sst_m

  ViewVC Help
Powered by ViewVC 1.1.21