/[lmdze]/trunk/phylmd/Interface_surf/limit_read_sst.f90
ViewVC logotype

Contents of /trunk/phylmd/Interface_surf/limit_read_sst.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 328 - (show annotations)
Thu Jun 13 14:40:06 2019 UTC (5 years ago) by guez
File size: 1315 byte(s)
Change all `.f` suffixes to `.f90`. (The opposite was done in revision
82.)  Because of change of philosopy in GNUmakefile: we already had a
rewritten rule for `.f`, so it does not make the makefile longer to
replace it by a rule for `.f90`. And it spares us options of
makedepf90 and of the compiler. Also we prepare the way for a simpler
`CMakeLists.txt`.

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

  ViewVC Help
Powered by ViewVC 1.1.21