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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21