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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 301 - (hide annotations)
Thu Aug 2 17:23:07 2018 UTC (5 years, 9 months ago) by guez
File size: 1875 byte(s)
Move the call to conf_interface up to physiq, so there is no need to
test first call inside pbl_surface for this.

run_off_lic in fonte_neige was computed but not used. Pass it up to
pbl_surface so we can output it (following LMDZ).

Simplify the logic in interfsur_lim so we do not need debut.

Remove the tests on the order of surface types in interfsurf_hq. Just
add comments in indicesol.

1 guez 54 module interfsur_lim_m
2    
3     implicit none
4    
5     contains
6    
7 guez 301 SUBROUTINE interfsur_lim(jour, knindex, albedo, z0_new)
8 guez 54
9 guez 150 ! Cette routine sert d'interface entre le mod\`ele atmosph\'erique et
10 guez 54 ! un fichier de conditions aux limites.
11    
12 guez 104 ! Laurent FAIRHEAD, February 2000
13 guez 54
14 guez 301 use conf_gcm_m, only: lmt_pas
15 guez 98 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 guez 191 use time_phylmdz, only: itap
19 guez 54
20 guez 98 integer, intent(IN):: jour ! jour a lire dans l'annee
21 guez 54
22 guez 106 integer, intent(in):: knindex(:) ! (knon)
23 guez 150 ! index des points de la surface \`a traiter
24 guez 98
25 guez 156 real, intent(out):: albedo(:) ! (knon) albedo lu
26 guez 282 real, intent(out):: z0_new(:) ! (knon) longueur de rugosit\'e lue
27 guez 98
28     ! Local:
29    
30     logical, save:: deja_lu_sur
31 guez 150 ! jour \`a lire d\'ej\`a lu pour une surface pr\'ec\'edente
32 guez 54
33 guez 301 integer:: jour_lu_sur = - 1
34 guez 54
35 guez 98 ! Champs lus dans le fichier de conditions aux limites :
36 guez 301 real, save:: alb_lu(klon), rug_lu(klon)
37 guez 98
38     integer ncid, varid
39    
40 guez 54 !------------------------------------------------------------
41    
42 guez 98 if (jour - jour_lu_sur /= 0) deja_lu_sur = .false.
43 guez 54
44     ! Tester d'abord si c'est le moment de lire le fichier
45 guez 191 if (mod(itap - 1, lmt_pas) == 0 .and. .not. deja_lu_sur) then
46 guez 98 call NF95_OPEN('limit.nc', NF90_NOWRITE, ncid)
47 guez 54
48     ! Lecture Albedo
49 guez 98 call NF95_INQ_VARID(ncid, 'ALB', varid)
50     call NF95_GET_VAR(ncid, varid, alb_lu, start=(/1, jour/))
51 guez 54
52 guez 150 ! Lecture rugosit\'e
53 guez 98 call NF95_INQ_VARID(ncid, 'RUG', varid)
54     call NF95_GET_VAR(ncid, varid, rug_lu, start=(/1, jour/))
55 guez 54
56 guez 98 call NF95_CLOSE(ncid)
57 guez 54 deja_lu_sur = .true.
58     jour_lu_sur = jour
59     endif
60    
61     ! Recopie des variables dans les champs de sortie
62 guez 156 albedo = alb_lu(knindex)
63 guez 282 z0_new = rug_lu(knindex)
64 guez 54
65     END SUBROUTINE interfsur_lim
66    
67     end module interfsur_lim_m

  ViewVC Help
Powered by ViewVC 1.1.21