/[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 175 - (show annotations)
Fri Feb 5 16:02:34 2016 UTC (8 years, 3 months ago) by guez
File size: 3346 byte(s)
Added argument itau_phy to ini_histins, phyetat0, phytrac and
phyredem0. Removed variable itau_phy of module temps. Avoiding side
effect in etat0 and phyetat0. The procedures ini_histins, phyetat0,
phytrac and phyredem0 are all called by physiq so there is no
cascading variable penalty.

In procedure inifilr, made the condition on colat0 weaker to allow for
rounding error.

Removed arguments flux_o, flux_g and t_slab of clmain, flux_o and
flux_g of clqh and interfsurf_hq, tslab and seaice of phyetat0 and
phyredem. NetCDF variables TSLAB and SEAICE no longer in
restartphy.nc. All these variables were related to the not-implemented
slab ocean. seaice and tslab were just set to 0 in phyetat0 and never
used nor changed. flux_o and flux_g were computed in clmain but never
used in physiq.

Removed argument swnet of clqh. Was used only to compute a local
variable, swdown, which was not used.

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

  ViewVC Help
Powered by ViewVC 1.1.21