/[lmdze]/trunk/libf/phylmd/Interface_surf/interfoce_lim.f90
ViewVC logotype

Contents of /trunk/libf/phylmd/Interface_surf/interfoce_lim.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 54 - (show annotations)
Tue Dec 6 15:07:04 2011 UTC (12 years, 5 months ago) by guez
File size: 6953 byte(s)
Removed Numerical Recipes procedure "ran1". Replaced calls to "ran1"
in "inidissip" by calls to intrinsic procedures.

Split file "interface_surf.f90" into a file with a module containing
only variables, "interface_surf", and single-procedure files. Gathered
files into directory "Interface_surf".

Added argument "cdivu" to "gradiv" and "gradiv2", "cdivh" to
"divgrad2" and "divgrad", and "crot" to "nxgraro2" and
"nxgrarot". "dissip" now uses variables "cdivu", "cdivh" and "crot"
from module "inidissip_m", so it can pass them to "gradiv2",
etc. Thanks to this modification, we avoid a circular dependency
betwwen "inidissip.f90" and "gradiv2.f90", etc. The value -1. used by
"gradiv2", for instance, during computation of eigenvalues is not the
value "cdivu" computed by "inidissip".

Extracted procedure "start_inter_3d" from module "startdyn", to its
own module.

In "inidissip", unrolled loop on "ii". I find it clearer now.

Moved variables "matriceun", "matriceus", "matricevn", "matricevs",
"matrinvn" and "matrinvs" from module "parafilt" to module
"inifilr_m". Moved variables "jfiltnu", "jfiltnv", "jfiltsu",
"jfiltsv" from module "coefils" to module "inifilr_m".

1 module interfoce_lim_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE interfoce_lim(itime, dtime, jour, &
8 klon, nisurf, knon, knindex, &
9 debut, &
10 lmt_sst, pctsrf_new)
11
12 ! Cette routine sert d'interface entre le modele atmospherique et
13 ! un fichier de conditions aux limites
14
15 ! L. Fairhead 02/2000
16
17 use abort_gcm_m, only: abort_gcm
18 use indicesol
19
20 integer, intent(IN) :: itime ! numero du pas de temps courant
21 real , intent(IN) :: dtime ! pas de temps de la physique (en s)
22 integer, intent(IN) :: jour ! jour a lire dans l'annee
23 integer, intent(IN) :: nisurf ! index de la surface a traiter (1 = sol continental)
24 integer, intent(IN) :: knon ! nombre de points dans le domaine a traiter
25 integer, intent(IN) :: klon ! taille de la grille
26 integer, dimension(klon), intent(in) :: knindex ! index des points de la surface a traiter
27 logical, intent(IN) :: debut ! logical: 1er appel a la physique (initialisation)
28
29 ! Parametres de sortie
30 ! output:
31 ! lmt_sst SST lues dans le fichier de CL
32 ! pctsrf_new sous-maille fractionnelle
33 real, intent(out), dimension(klon) :: lmt_sst
34 real, intent(out), dimension(klon, nbsrf) :: pctsrf_new
35
36 ! Variables locales
37 integer :: ii
38 INTEGER, save :: lmt_pas ! frequence de lecture des conditions limites
39 ! (en pas de physique)
40 logical, save :: deja_lu ! pour indiquer que le jour a lire a deja
41 ! lu pour une surface precedente
42 integer, save :: jour_lu
43 integer :: ierr
44 character (len = 20) :: modname = 'interfoce_lim'
45 character (len = 80) :: abort_message
46 logical, save :: newlmt = .TRUE.
47 logical, save :: check = .FALSE.
48 ! Champs lus dans le fichier de CL
49 real, allocatable , save, dimension(:) :: sst_lu, rug_lu, nat_lu
50 real, allocatable , save, dimension(:, :) :: pct_tmp
51
52 ! quelques variables pour netcdf
53
54 include "netcdf.inc"
55 integer :: nid, nvarid
56 integer, dimension(2) :: start, epais
57
58 ! --------------------------------------------------
59
60 if (debut .and. .not. allocated(sst_lu)) then
61 lmt_pas = nint(86400./dtime * 1.0) ! pour une lecture une fois par jour
62 jour_lu = jour - 1
63 allocate(sst_lu(klon))
64 allocate(nat_lu(klon))
65 allocate(pct_tmp(klon, nbsrf))
66 endif
67
68 if ((jour - jour_lu) /= 0) deja_lu = .false.
69
70 if (check) write(*, *)modname, ' :: jour, jour_lu, deja_lu', jour, jour_lu, &
71 deja_lu
72 if (check) write(*, *)modname, ' :: itime, lmt_pas ', itime, lmt_pas, dtime
73
74 ! Tester d'abord si c'est le moment de lire le fichier
75 if (mod(itime-1, lmt_pas) == 0 .and. .not. deja_lu) then
76
77 ! Ouverture du fichier
78
79 ierr = NF_OPEN ('limit.nc', NF_NOWRITE, nid)
80 if (ierr.NE.NF_NOERR) then
81 abort_message &
82 = 'Pb d''ouverture du fichier de conditions aux limites'
83 call abort_gcm(modname, abort_message, 1)
84 endif
85
86 ! La tranche de donnees a lire:
87
88 start(1) = 1
89 start(2) = jour
90 epais(1) = klon
91 epais(2) = 1
92
93 if (newlmt) then
94
95 ! Fraction "ocean"
96
97 ierr = NF_INQ_VARID(nid, 'FOCE', nvarid)
98 if (ierr /= NF_NOERR) then
99 abort_message = 'Le champ <FOCE> est absent'
100 call abort_gcm(modname, abort_message, 1)
101 endif
102 ierr = NF_GET_VARA_REAL(nid, nvarid, start, epais, pct_tmp(1, is_oce))
103 if (ierr /= NF_NOERR) then
104 abort_message = 'Lecture echouee pour <FOCE>'
105 call abort_gcm(modname, abort_message, 1)
106 endif
107
108 ! Fraction "glace de mer"
109
110 ierr = NF_INQ_VARID(nid, 'FSIC', nvarid)
111 if (ierr /= NF_NOERR) then
112 abort_message = 'Le champ <FSIC> est absent'
113 call abort_gcm(modname, abort_message, 1)
114 endif
115 ierr = NF_GET_VARA_REAL(nid, nvarid, start, epais, pct_tmp(1, is_sic))
116 if (ierr /= NF_NOERR) then
117 abort_message = 'Lecture echouee pour <FSIC>'
118 call abort_gcm(modname, abort_message, 1)
119 endif
120
121 ! Fraction "terre"
122
123 ierr = NF_INQ_VARID(nid, 'FTER', nvarid)
124 if (ierr /= NF_NOERR) then
125 abort_message = 'Le champ <FTER> est absent'
126 call abort_gcm(modname, abort_message, 1)
127 endif
128 ierr = NF_GET_VARA_REAL(nid, nvarid, start, epais, pct_tmp(1, is_ter))
129 if (ierr /= NF_NOERR) then
130 abort_message = 'Lecture echouee pour <FTER>'
131 call abort_gcm(modname, abort_message, 1)
132 endif
133
134 ! Fraction "glacier terre"
135
136 ierr = NF_INQ_VARID(nid, 'FLIC', nvarid)
137 if (ierr /= NF_NOERR) then
138 abort_message = 'Le champ <FLIC> est absent'
139 call abort_gcm(modname, abort_message, 1)
140 endif
141 ierr = NF_GET_VARA_REAL(nid, nvarid, start, epais, pct_tmp(1, is_lic))
142 if (ierr /= NF_NOERR) then
143 abort_message = 'Lecture echouee pour <FLIC>'
144 call abort_gcm(modname, abort_message, 1)
145 endif
146
147 else ! on en est toujours a rnatur
148
149 ierr = NF_INQ_VARID(nid, 'NAT', nvarid)
150 if (ierr /= NF_NOERR) then
151 abort_message = 'Le champ <NAT> est absent'
152 call abort_gcm(modname, abort_message, 1)
153 endif
154 ierr = NF_GET_VARA_REAL(nid, nvarid, start, epais, nat_lu)
155 if (ierr /= NF_NOERR) then
156 abort_message = 'Lecture echouee pour <NAT>'
157 call abort_gcm(modname, abort_message, 1)
158 endif
159
160 ! Remplissage des fractions de surface
161 ! nat = 0, 1, 2, 3 pour ocean, terre, glacier, seaice
162
163 pct_tmp = 0.0
164 do ii = 1, klon
165 pct_tmp(ii, nint(nat_lu(ii)) + 1) = 1.
166 enddo
167
168
169 ! On se retrouve avec ocean en 1 et terre en 2 alors qu'on veut le contraire
170
171 pctsrf_new = pct_tmp
172 pctsrf_new (:, 2)= pct_tmp (:, 1)
173 pctsrf_new (:, 1)= pct_tmp (:, 2)
174 pct_tmp = pctsrf_new
175 endif ! fin test sur newlmt
176
177 ! Lecture SST
178
179 ierr = NF_INQ_VARID(nid, 'SST', nvarid)
180 if (ierr /= NF_NOERR) then
181 abort_message = 'Le champ <SST> est absent'
182 call abort_gcm(modname, abort_message, 1)
183 endif
184 ierr = NF_GET_VARA_REAL(nid, nvarid, start, epais, sst_lu)
185 if (ierr /= NF_NOERR) then
186 abort_message = 'Lecture echouee pour <SST>'
187 call abort_gcm(modname, abort_message, 1)
188 endif
189
190
191 ! Fin de lecture
192
193 ierr = NF_CLOSE(nid)
194 deja_lu = .true.
195 jour_lu = jour
196 endif
197
198 ! Recopie des variables dans les champs de sortie
199
200 lmt_sst = 999999999.
201 do ii = 1, knon
202 lmt_sst(ii) = sst_lu(knindex(ii))
203 enddo
204
205 pctsrf_new(:, is_oce) = pct_tmp(:, is_oce)
206 pctsrf_new(:, is_sic) = pct_tmp(:, is_sic)
207
208 END SUBROUTINE interfoce_lim
209
210 end module interfoce_lim_m

  ViewVC Help
Powered by ViewVC 1.1.21