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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 104 - (show annotations)
Thu Sep 4 10:05:52 2014 UTC (9 years, 9 months ago) by guez
File size: 3444 byte(s)
Removed procedure sortvarc0. Called sortvarc with an additional
argument resetvarc instead. (Following LMDZ.) Moved current time
computations and some printing statements from sortvarc to
caldyn. Could then remove arguments itau and time_0 of sortvarc, and
could remove "use dynetat0". Better to keep "dynetat0.f" as a gcm-only
file.

Moved some variables from module ener to module sortvarc.

Split file "mathelp.f" into single-procedure files.

Removed unused argument nadv of adaptdt. Removed dimension arguments
of bernoui.

Removed unused argument nisurf of interfoce_lim. Changed the size of
argument lmt_sst of interfoce_lim from klon to knon. Removed case when
newlmt is false.

dynredem1 is called only once in each run, either ce0l or gcm. So
variable nb in call to nf95_put_var was always 1. Removed variable nb.

Removed dimension arguments of calcul_fluxs. Removed unused arguments
precip_rain, precip_snow, snow of calcul_fluxs. Changed the size of
all the arrays in calcul_fluxs from klon to knon.

Removed dimension arguments of fonte_neige. Changed the size of all
the arrays in fonte_neige from klon to knon.

Changed the size of arguments tsurf and tsurf_new of interfsurf_hq
from klon to knon. Changed the size of argument ptsrf of soil from
klon to knon.

1 module interfoce_lim_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE interfoce_lim(itime, dtime, jour, klon, knon, knindex, &
8 debut, lmt_sst, pctsrf_new)
9
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 abort_gcm_m, only: abort_gcm
16 USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf
17 USE netcdf, ONLY: nf90_noerr, nf90_nowrite
18 use netcdf95, only: NF95_CLOSE, nf95_get_var, NF95_INQ_VARID, nf95_open
19
20 integer, intent(IN):: itime ! numero du pas de temps courant
21 real, intent(IN):: dtime ! pas de temps de la physique (en s)
22 integer, intent(IN):: jour ! jour a lire dans l'annee
23 integer, intent(IN):: klon ! taille de la grille
24 integer, intent(IN):: knon ! nombre de points dans le domaine a traiter
25
26 integer, intent(in):: knindex(klon)
27 ! index des points de la surface a traiter
28
29 logical, intent(IN):: debut ! 1er appel a la physique (initialisation)
30
31 real, intent(out), dimension(knon):: lmt_sst
32 ! SST lues dans le fichier de CL
33
34 real, intent(out), dimension(klon, nbsrf):: pctsrf_new
35 ! sous-maille fractionnelle
36
37 ! Local:
38
39 integer ii
40
41 INTEGER, save:: lmt_pas ! frequence de lecture des conditions limites
42 ! (en pas de physique)
43
44 logical, save:: deja_lu ! pour indiquer que le jour a lire a deja
45 ! lu pour une surface precedente
46
47 integer, save:: jour_lu
48
49 ! Champs lus dans le fichier de conditions aux limites :
50 real, allocatable, save, dimension(:):: sst_lu
51 real, allocatable, save, dimension(:, :):: pct_tmp
52
53 ! Pour Netcdf :
54 integer nid, nvarid
55 integer, dimension(2):: start, epais
56
57 ! --------------------------------------------------
58
59 if (debut .and. .not. allocated(sst_lu)) then
60 lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
61 jour_lu = jour - 1
62 allocate(sst_lu(klon))
63 allocate(pct_tmp(klon, nbsrf))
64 endif
65
66 if ((jour - jour_lu) /= 0) deja_lu = .false.
67
68 ! Tester d'abord si c'est le moment de lire le fichier
69 if (mod(itime - 1, lmt_pas) == 0 .and. .not. deja_lu) then
70 call NF95_OPEN ('limit.nc', NF90_NOWRITE, nid)
71
72 ! La tranche de donnees a lire:
73 start(1) = 1
74 start(2) = jour
75 epais(1) = klon
76 epais(2) = 1
77
78 ! Fraction "ocean"
79 call NF95_INQ_VARID(nid, 'FOCE', nvarid)
80 call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_oce), start, epais)
81
82 ! Fraction "glace de mer"
83 call NF95_INQ_VARID(nid, 'FSIC', nvarid)
84 call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_sic), start, epais)
85
86 ! Fraction "terre"
87 call NF95_INQ_VARID(nid, 'FTER', nvarid)
88 call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_ter), start, epais)
89
90 ! Fraction "glacier terre"
91 call NF95_INQ_VARID(nid, 'FLIC', nvarid)
92 call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_lic), start, epais)
93
94 call NF95_INQ_VARID(nid, 'SST', nvarid)
95 call NF95_GET_VAR(nid, nvarid, sst_lu, start, epais)
96
97 call NF95_CLOSE(nid)
98 deja_lu = .true.
99 jour_lu = jour
100 endif
101
102 ! Recopie des variables dans les champs de sortie
103
104 do ii = 1, knon
105 lmt_sst(ii) = sst_lu(knindex(ii))
106 enddo
107
108 pctsrf_new(:, is_oce) = pct_tmp(:, is_oce)
109 pctsrf_new(:, is_sic) = pct_tmp(:, is_sic)
110
111 END SUBROUTINE interfoce_lim
112
113 end module interfoce_lim_m

  ViewVC Help
Powered by ViewVC 1.1.21