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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 154 - (show annotations)
Tue Jul 7 17:49:23 2015 UTC (8 years, 11 months ago) by guez
Original Path: trunk/Sources/phylmd/Interface_surf/interfsur_lim.f
File size: 2533 byte(s)
Removed argument dtphys of physiq. Use it directly from comconst in
physiq instead.

Donwgraded variables eignfnu, eignfnv of module inifgn_m to dummy
arguments of SUBROUTINE inifgn. They were not used elsewhere than in
the calling procedure inifilr. Renamed argument dv of inifgn to eignval_v.

Made alboc and alboc_cd independent of the size of arguments. Now we
can call them only at indices knindex in interfsurf_hq, where we need
them. Fixed a bug in alboc_cd: rmu0 was modified, and the
corresponding actual argument in interfsurf_hq is an intent(in)
argument of interfsurf_hq.

Variables of size knon instead of klon in interfsur_lim and interfsurf_hq.

Removed argument alb_new of interfsurf_hq because it was the same than
alblw. Simplified test on cycle_diurne, following LMDZ.

Moved tests on nbapp_rad from physiq to read_clesphys2. No need for
separate counter itaprad, we can use itap. Define lmt_pas and radpas
from integer input parameters instead of real-type computed values.

1 module interfsur_lim_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE interfsur_lim(itime, dtime, jour, knindex, debut, alb_new, z0_new)
8
9 ! Cette routine sert d'interface entre le mod\`ele atmosph\'erique et
10 ! un fichier de conditions aux limites.
11
12 ! Laurent FAIRHEAD, February 2000
13
14 use abort_gcm_m, only: abort_gcm
15 USE dimphy, ONLY: klon
16 use netcdf, only: NF90_NOWRITE
17 use netcdf95, only: NF95_close, NF95_GET_VAR, NF95_INQ_VARID, NF95_OPEN
18
19 integer, intent(IN):: itime ! numero du pas de temps courant
20 real, intent(IN):: dtime ! pas de temps de la physique (en s)
21 integer, intent(IN):: jour ! jour a lire dans l'annee
22
23 integer, intent(in):: knindex(:) ! (knon)
24 ! index des points de la surface \`a traiter
25
26 logical, intent(IN):: debut ! premier appel \`a la physique (initialisation)
27 real, intent(out):: alb_new(:) ! (knon) albedo lu
28 real, intent(out):: z0_new(:) ! (klon) longueur de rugosit\'e lue
29
30 ! Local:
31
32 integer knon ! nombre de points dans le domaine a traiter
33
34 integer, save:: lmt_pas ! frequence de lecture des conditions limites
35 ! (en pas de physique)
36
37 logical, save:: deja_lu_sur
38 ! jour \`a lire d\'ej\`a lu pour une surface pr\'ec\'edente
39
40 integer, save:: jour_lu_sur
41
42 ! Champs lus dans le fichier de conditions aux limites :
43 real, allocatable, save:: alb_lu(:), rug_lu(:)
44
45 integer ncid, varid
46
47 !------------------------------------------------------------
48
49 knon = size(knindex)
50
51 if (debut) then
52 lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
53 jour_lu_sur = jour - 1
54 allocate(alb_lu(klon))
55 allocate(rug_lu(klon))
56 endif
57
58 if (jour - jour_lu_sur /= 0) deja_lu_sur = .false.
59
60 ! Tester d'abord si c'est le moment de lire le fichier
61 if (mod(itime - 1, lmt_pas) == 0 .and. .not. deja_lu_sur) then
62 call NF95_OPEN('limit.nc', NF90_NOWRITE, ncid)
63
64 ! Lecture Albedo
65 call NF95_INQ_VARID(ncid, 'ALB', varid)
66 call NF95_GET_VAR(ncid, varid, alb_lu, start=(/1, jour/))
67
68 ! Lecture rugosit\'e
69 call NF95_INQ_VARID(ncid, 'RUG', varid)
70 call NF95_GET_VAR(ncid, varid, rug_lu, start=(/1, jour/))
71
72 call NF95_CLOSE(ncid)
73 deja_lu_sur = .true.
74 jour_lu_sur = jour
75 endif
76
77 ! Recopie des variables dans les champs de sortie
78 alb_new = alb_lu(knindex(:knon))
79 z0_new(:knon) = rug_lu(knindex(:knon))
80 z0_new(knon + 1:) = 999999.
81
82 END SUBROUTINE interfsur_lim
83
84 end module interfsur_lim_m

  ViewVC Help
Powered by ViewVC 1.1.21