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

Contents of /trunk/Sources/phylmd/Interface_surf/interfoce_lim.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 191 - (show annotations)
Mon May 9 19:56:28 2016 UTC (8 years ago) by guez
File size: 3179 byte(s)
Extracted the call to read_comdissnew out of conf_gcm.

Made ok_instan a variable of module clesphys, itau_phy a variable of
module phyetat0_m, nid_ins a variable of module ini_histins_m, itap a
variable of new module time_phylmdz, so that histwrite_phy can be
called from any procedure without the need to cascade those variables
into that procedure. Made itau_w a variable of module time_phylmdz so
that it is computed only once per time step of physics.

Extracted variables of module clesphys which were in namelist
conf_phys_nml into their own namelist, clesphys_nml, and created
procedure read_clesphys reading clesphys_nml, to avoid side effect.

No need for double precision in procedure getso4fromfile. Assume there
is a single variable for the whole year in the NetCDF file instead of
one variable per month.

Created generic procedure histwrite_phy and removed procedure
write_histins, following LMDZ. histwrite_phy has only two arguments,
can be called from anywhere, and should manage the logic of writing or
not writing into various history files with various operations. So the
test on ok_instan goes inside histwrite_phy.

Test for raz_date in phyetat0 instead of physiq to avoid side effect.

Created procedure increment_itap to avoid side effect.

Removed unnecessary differences between procedures readsulfate and
readsulfate_pi.

1 module interfoce_lim_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE interfoce_lim(dtime, jour, knindex, debut, lmt_sst, pctsrf_new)
8
9 ! lecture conditions limites
10 ! Cette routine sert d'interface entre le modèle atmosphérique et
11 ! un fichier de conditions aux limites.
12
13 ! Laurent FAIRHEAD, February 2000
14
15 USE dimphy, ONLY: klon
16 USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf
17 USE netcdf, ONLY: nf90_nowrite
18 use netcdf95, only: NF95_CLOSE, nf95_get_var, NF95_INQ_VARID, nf95_open
19 use nr_util, only: assert
20 use time_phylmdz, only: itap
21
22 real, intent(IN):: dtime ! pas de temps de la physique (en s)
23 integer, intent(IN):: jour ! jour a lire dans l'annee
24
25 integer, intent(in):: knindex(:) ! (knon)
26 ! index des points de la surface a traiter
27
28 logical, intent(IN):: debut ! 1er appel a la physique (initialisation)
29
30 real, intent(out):: lmt_sst(:) ! (knon)
31 ! SST lues dans le fichier de conditions aux limites
32
33 real, intent(out):: pctsrf_new(:, :) ! (klon, nbsrf)
34 ! sous-maille fractionnelle
35
36 ! Local:
37
38 INTEGER, save:: lmt_pas ! frequence de lecture des conditions limites
39 ! (en pas de physique)
40
41 logical, save:: deja_lu
42 ! pour indiquer que le jour à lire a déjà été lu pour une surface
43 ! précédente
44
45 integer, save:: jour_lu
46
47 ! Champs lus dans le fichier de conditions aux limites :
48 real, allocatable, save:: sst_lu(:)
49 real, allocatable, save:: pct_tmp(:, :)
50
51 integer ncid, varid ! pour NetCDF
52
53 ! --------------------------------------------------
54
55 call assert(size(knindex) == size(lmt_sst), "interfoce_lim knon")
56
57 if (debut .and. .not. allocated(sst_lu)) then
58 lmt_pas = nint(86400. / dtime) ! pour une lecture une fois par jour
59 jour_lu = jour - 1
60 allocate(sst_lu(klon))
61 allocate(pct_tmp(klon, nbsrf))
62 endif
63
64 if ((jour - jour_lu) /= 0) deja_lu = .false.
65
66 ! Tester d'abord si c'est le moment de lire le fichier
67 if (mod(itap - 1, lmt_pas) == 0 .and. .not. deja_lu) then
68 call NF95_OPEN ('limit.nc', NF90_NOWRITE, ncid)
69
70 ! Fraction "ocean"
71 call NF95_INQ_VARID(ncid, 'FOCE', varid)
72 call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_oce), start = (/1, jour/))
73
74 ! Fraction "glace de mer"
75 call NF95_INQ_VARID(ncid, 'FSIC', varid)
76 call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_sic), start = (/1, jour/))
77
78 ! Fraction "terre"
79 call NF95_INQ_VARID(ncid, 'FTER', varid)
80 call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_ter), start = (/1, jour/))
81
82 ! Fraction "glacier terre"
83 call NF95_INQ_VARID(ncid, 'FLIC', varid)
84 call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_lic), start = (/1, jour/))
85
86 call NF95_INQ_VARID(ncid, 'SST', varid)
87 call NF95_GET_VAR(ncid, varid, sst_lu, start = (/1, jour/))
88
89 call NF95_CLOSE(ncid)
90 deja_lu = .true.
91 jour_lu = jour
92 endif
93
94 ! Recopie des variables dans les champs de sortie
95 lmt_sst = sst_lu(knindex)
96 pctsrf_new(:, is_oce) = pct_tmp(:, is_oce)
97 pctsrf_new(:, is_sic) = pct_tmp(:, is_sic)
98
99 END SUBROUTINE interfoce_lim
100
101 end module interfoce_lim_m

  ViewVC Help
Powered by ViewVC 1.1.21