/[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 178 - (hide annotations)
Fri Mar 11 18:47:26 2016 UTC (8 years, 2 months ago) by guez
File size: 3228 byte(s)
Moved variables date0, deltat, datasz_max, ncvar_ids, point, buff_pos,
buffer, regular from module histcom_var to modules where they are
defined.

Removed procedure ioipslmpp, useless for a sequential program.

Added argument datasz_max to histwrite_real (to avoid circular
dependency with histwrite).

Removed useless variables and computations everywhere.

Changed real litteral constants from default kind to double precision
in lwb, lwu, lwvn, sw1s, swtt, swtt1, swu.

Removed unused arguments: paer of sw, sw1s, sw2s, swclr; pcldsw of
sw1s, sw2s; pdsig, prayl of swr; co2_ppm of clmain, clqh; tsol of
transp_lay; nsrf of screenp; kcrit and kknu of gwstress; pstd of
orosetup.

Added output of relative humidity.

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 175 ! lecture conditions limites
11 guez 104 ! Cette routine sert d'interface entre le modèle atmosphérique et
12     ! un fichier de conditions aux limites.
13 guez 54
14 guez 104 ! Laurent FAIRHEAD, February 2000
15 guez 54
16 guez 106 USE dimphy, ONLY: klon
17 guez 104 USE indicesol, ONLY: is_lic, is_oce, is_sic, is_ter, nbsrf
18 guez 178 USE netcdf, ONLY: nf90_nowrite
19 guez 104 use netcdf95, only: NF95_CLOSE, nf95_get_var, NF95_INQ_VARID, nf95_open
20 guez 178 use nr_util, only: assert
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     INTEGER, save:: lmt_pas ! frequence de lecture des conditions limites
40 guez 54 ! (en pas de physique)
41 guez 104
42 guez 106 logical, save:: deja_lu
43     ! pour indiquer que le jour à lire a déjà été lu pour une surface
44     ! précédente
45 guez 54
46 guez 104 integer, save:: jour_lu
47 guez 54
48 guez 104 ! Champs lus dans le fichier de conditions aux limites :
49 guez 106 real, allocatable, save:: sst_lu(:)
50     real, allocatable, save:: pct_tmp(:, :)
51 guez 54
52 guez 106 integer ncid, varid ! pour NetCDF
53 guez 104
54 guez 54 ! --------------------------------------------------
55    
56 guez 178 call assert(size(knindex) == size(lmt_sst), "interfoce_lim knon")
57 guez 106
58 guez 54 if (debut .and. .not. allocated(sst_lu)) then
59 guez 106 lmt_pas = nint(86400. / dtime) ! pour une lecture une fois par jour
60 guez 54 jour_lu = jour - 1
61     allocate(sst_lu(klon))
62     allocate(pct_tmp(klon, nbsrf))
63     endif
64    
65     if ((jour - jour_lu) /= 0) deja_lu = .false.
66    
67     ! Tester d'abord si c'est le moment de lire le fichier
68 guez 104 if (mod(itime - 1, lmt_pas) == 0 .and. .not. deja_lu) then
69 guez 106 call NF95_OPEN ('limit.nc', NF90_NOWRITE, ncid)
70 guez 54
71 guez 104 ! Fraction "ocean"
72 guez 106 call NF95_INQ_VARID(ncid, 'FOCE', varid)
73     call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_oce), start = (/1, jour/))
74 guez 54
75 guez 104 ! Fraction "glace de mer"
76 guez 106 call NF95_INQ_VARID(ncid, 'FSIC', varid)
77     call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_sic), start = (/1, jour/))
78 guez 54
79 guez 104 ! Fraction "terre"
80 guez 106 call NF95_INQ_VARID(ncid, 'FTER', varid)
81     call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_ter), start = (/1, jour/))
82 guez 54
83 guez 104 ! Fraction "glacier terre"
84 guez 106 call NF95_INQ_VARID(ncid, 'FLIC', varid)
85     call NF95_GET_VAR(ncid, varid, pct_tmp(:, is_lic), start = (/1, jour/))
86 guez 54
87 guez 106 call NF95_INQ_VARID(ncid, 'SST', varid)
88     call NF95_GET_VAR(ncid, varid, sst_lu, start = (/1, jour/))
89 guez 54
90 guez 106 call NF95_CLOSE(ncid)
91 guez 54 deja_lu = .true.
92     jour_lu = jour
93     endif
94    
95     ! Recopie des variables dans les champs de sortie
96 guez 106 lmt_sst = sst_lu(knindex)
97 guez 54 pctsrf_new(:, is_oce) = pct_tmp(:, is_oce)
98     pctsrf_new(:, is_sic) = pct_tmp(:, is_sic)
99    
100     END SUBROUTINE interfoce_lim
101    
102     end module interfoce_lim_m

  ViewVC Help
Powered by ViewVC 1.1.21