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 |
! 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 dimphy, ONLY: klon |
17 |
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 |
use nr_util, only: assert_eq |
21 |
|
22 |
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 |
|
26 |
integer, intent(in):: knindex(:) ! (knon) |
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):: lmt_sst(:) ! (knon) |
32 |
! SST lues dans le fichier de conditions aux limites |
33 |
|
34 |
real, intent(out):: pctsrf_new(:, :) ! (klon, nbsrf) |
35 |
! sous-maille fractionnelle |
36 |
|
37 |
! Local: |
38 |
|
39 |
integer knon ! nombre de points dans le domaine a traiter |
40 |
|
41 |
INTEGER, save:: lmt_pas ! frequence de lecture des conditions limites |
42 |
! (en pas de physique) |
43 |
|
44 |
logical, save:: deja_lu |
45 |
! pour indiquer que le jour à lire a déjà été lu pour une surface |
46 |
! précédente |
47 |
|
48 |
integer, save:: jour_lu |
49 |
|
50 |
! Champs lus dans le fichier de conditions aux limites : |
51 |
real, allocatable, save:: sst_lu(:) |
52 |
real, allocatable, save:: pct_tmp(:, :) |
53 |
|
54 |
integer ncid, varid ! pour NetCDF |
55 |
|
56 |
! -------------------------------------------------- |
57 |
|
58 |
knon = assert_eq(size(knindex), size(lmt_sst), "interfoce_lim knon") |
59 |
|
60 |
if (debut .and. .not. allocated(sst_lu)) then |
61 |
lmt_pas = nint(86400. / dtime) ! pour une lecture une fois par jour |
62 |
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 |
if (mod(itime - 1, lmt_pas) == 0 .and. .not. deja_lu) then |
71 |
call NF95_OPEN ('limit.nc', NF90_NOWRITE, ncid) |
72 |
|
73 |
! Fraction "ocean" |
74 |
call NF95_INQ_VARID(ncid, 'FOCE', varid) |
75 |
call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_oce), start = (/1, jour/)) |
76 |
|
77 |
! Fraction "glace de mer" |
78 |
call NF95_INQ_VARID(ncid, 'FSIC', varid) |
79 |
call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_sic), start = (/1, jour/)) |
80 |
|
81 |
! Fraction "terre" |
82 |
call NF95_INQ_VARID(ncid, 'FTER', varid) |
83 |
call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_ter), start = (/1, jour/)) |
84 |
|
85 |
! Fraction "glacier terre" |
86 |
call NF95_INQ_VARID(ncid, 'FLIC', varid) |
87 |
call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_lic), start = (/1, jour/)) |
88 |
|
89 |
call NF95_INQ_VARID(ncid, 'SST', varid) |
90 |
call NF95_GET_VAR(ncid, varid, sst_lu, start = (/1, jour/)) |
91 |
|
92 |
call NF95_CLOSE(ncid) |
93 |
deja_lu = .true. |
94 |
jour_lu = jour |
95 |
endif |
96 |
|
97 |
! Recopie des variables dans les champs de sortie |
98 |
lmt_sst = sst_lu(knindex) |
99 |
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 |