--- trunk/phylmd/Interface_surf/interfoce_lim.f 2014/09/04 10:05:52 104 +++ trunk/Sources/phylmd/Interface_surf/interfoce_lim.f 2015/04/29 15:47:56 134 @@ -4,8 +4,8 @@ contains - SUBROUTINE interfoce_lim(itime, dtime, jour, klon, knon, knindex, & - debut, lmt_sst, pctsrf_new) + SUBROUTINE interfoce_lim(itime, dtime, jour, knindex, debut, lmt_sst, & + pctsrf_new) ! Cette routine sert d'interface entre le modèle atmosphérique et ! un fichier de conditions aux limites. @@ -13,51 +13,52 @@ ! Laurent FAIRHEAD, February 2000 use abort_gcm_m, only: abort_gcm + USE dimphy, ONLY: klon USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf USE netcdf, ONLY: nf90_noerr, nf90_nowrite use netcdf95, only: NF95_CLOSE, nf95_get_var, NF95_INQ_VARID, nf95_open + use nr_util, only: assert_eq integer, intent(IN):: itime ! numero du pas de temps courant real, intent(IN):: dtime ! pas de temps de la physique (en s) integer, intent(IN):: jour ! jour a lire dans l'annee - integer, intent(IN):: klon ! taille de la grille - integer, intent(IN):: knon ! nombre de points dans le domaine a traiter - integer, intent(in):: knindex(klon) + integer, intent(in):: knindex(:) ! (knon) ! index des points de la surface a traiter logical, intent(IN):: debut ! 1er appel a la physique (initialisation) - real, intent(out), dimension(knon):: lmt_sst - ! SST lues dans le fichier de CL + real, intent(out):: lmt_sst(:) ! (knon) + ! SST lues dans le fichier de conditions aux limites - real, intent(out), dimension(klon, nbsrf):: pctsrf_new + real, intent(out):: pctsrf_new(:, :) ! (klon, nbsrf) ! sous-maille fractionnelle ! Local: - integer ii + integer knon ! nombre de points dans le domaine a traiter INTEGER, save:: lmt_pas ! frequence de lecture des conditions limites ! (en pas de physique) - logical, save:: deja_lu ! pour indiquer que le jour a lire a deja - ! lu pour une surface precedente + logical, save:: deja_lu + ! pour indiquer que le jour à lire a déjà été lu pour une surface + ! précédente integer, save:: jour_lu ! Champs lus dans le fichier de conditions aux limites : - real, allocatable, save, dimension(:):: sst_lu - real, allocatable, save, dimension(:, :):: pct_tmp + real, allocatable, save:: sst_lu(:) + real, allocatable, save:: pct_tmp(:, :) - ! Pour Netcdf : - integer nid, nvarid - integer, dimension(2):: start, epais + integer ncid, varid ! pour NetCDF ! -------------------------------------------------- + knon = assert_eq(size(knindex), size(lmt_sst), "interfoce_lim knon") + if (debut .and. .not. allocated(sst_lu)) then - lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour + lmt_pas = nint(86400. / dtime) ! pour une lecture une fois par jour jour_lu = jour - 1 allocate(sst_lu(klon)) allocate(pct_tmp(klon, nbsrf)) @@ -67,44 +68,34 @@ ! Tester d'abord si c'est le moment de lire le fichier if (mod(itime - 1, lmt_pas) == 0 .and. .not. deja_lu) then - call NF95_OPEN ('limit.nc', NF90_NOWRITE, nid) - - ! La tranche de donnees a lire: - start(1) = 1 - start(2) = jour - epais(1) = klon - epais(2) = 1 + call NF95_OPEN ('limit.nc', NF90_NOWRITE, ncid) ! Fraction "ocean" - call NF95_INQ_VARID(nid, 'FOCE', nvarid) - call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_oce), start, epais) + call NF95_INQ_VARID(ncid, 'FOCE', varid) + call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_oce), start = (/1, jour/)) ! Fraction "glace de mer" - call NF95_INQ_VARID(nid, 'FSIC', nvarid) - call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_sic), start, epais) + call NF95_INQ_VARID(ncid, 'FSIC', varid) + call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_sic), start = (/1, jour/)) ! Fraction "terre" - call NF95_INQ_VARID(nid, 'FTER', nvarid) - call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_ter), start, epais) + call NF95_INQ_VARID(ncid, 'FTER', varid) + call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_ter), start = (/1, jour/)) ! Fraction "glacier terre" - call NF95_INQ_VARID(nid, 'FLIC', nvarid) - call NF95_GET_VAR(nid, nvarid, pct_tmp(:, is_lic), start, epais) + call NF95_INQ_VARID(ncid, 'FLIC', varid) + call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_lic), start = (/1, jour/)) - call NF95_INQ_VARID(nid, 'SST', nvarid) - call NF95_GET_VAR(nid, nvarid, sst_lu, start, epais) + call NF95_INQ_VARID(ncid, 'SST', varid) + call NF95_GET_VAR(ncid, varid, sst_lu, start = (/1, jour/)) - call NF95_CLOSE(nid) + call NF95_CLOSE(ncid) deja_lu = .true. jour_lu = jour endif ! Recopie des variables dans les champs de sortie - - do ii = 1, knon - lmt_sst(ii) = sst_lu(knindex(ii)) - enddo - + lmt_sst = sst_lu(knindex) pctsrf_new(:, is_oce) = pct_tmp(:, is_oce) pctsrf_new(:, is_sic) = pct_tmp(:, is_sic)