1 |
|
|
! |
2 |
|
|
! $Header$ |
3 |
|
|
! |
4 |
|
|
SUBROUTINE interfoce_lim(itime, dtime, jour, & |
5 |
|
|
knon, knindex, & |
6 |
|
|
debut, & |
7 |
|
|
lmt_sst_p, pctsrf_new_p) |
8 |
|
|
|
9 |
|
|
USE mod_grid_phy_lmdz |
10 |
|
|
USE mod_phys_lmdz_para |
11 |
|
|
USE indice_sol_mod |
12 |
|
|
|
13 |
|
|
IMPLICIT NONE |
14 |
|
|
|
15 |
|
|
INCLUDE "netcdf.inc" |
16 |
|
|
|
17 |
|
|
! Cette routine sert d'interface entre le modele atmospherique et un fichier |
18 |
|
|
! de conditions aux limites |
19 |
|
|
! |
20 |
|
|
! L. Fairhead 02/2000 |
21 |
|
|
! |
22 |
|
|
! input: |
23 |
|
|
! itime numero du pas de temps courant |
24 |
|
|
! dtime pas de temps de la physique (en s) |
25 |
|
|
! jour jour a lire dans l'annee |
26 |
|
|
! nisurf index de la surface a traiter (1 = sol continental) |
27 |
|
|
! knon nombre de points dans le domaine a traiter |
28 |
|
|
! knindex index des points de la surface a traiter |
29 |
|
|
! klon taille de la grille |
30 |
|
|
! debut logical: 1er appel a la physique (initialisation) |
31 |
|
|
! |
32 |
|
|
! output: |
33 |
|
|
! lmt_sst_p SST lues dans le fichier de CL |
34 |
|
|
! pctsrf_new-p sous-maille fractionnelle |
35 |
|
|
! |
36 |
|
|
|
37 |
|
|
|
38 |
|
|
! Parametres d'entree |
39 |
|
|
!**************************************************************************************** |
40 |
|
|
INTEGER, INTENT(IN) :: itime |
41 |
|
|
INTEGER, INTENT(IN) :: jour |
42 |
|
|
INTEGER, INTENT(IN) :: knon |
43 |
|
|
INTEGER, DIMENSION(klon_loc), INTENT(IN) :: knindex |
44 |
|
|
REAL , INTENT(IN) :: dtime |
45 |
|
|
LOGICAL, INTENT(IN) :: debut |
46 |
|
|
|
47 |
|
|
! Parametres de sortie |
48 |
|
|
!**************************************************************************************** |
49 |
|
|
REAL, INTENT(OUT), DIMENSION(klon_loc) :: lmt_sst_p |
50 |
|
|
REAL, INTENT(OUT), DIMENSION(klon_loc,nbsrf) :: pctsrf_new_p |
51 |
|
|
|
52 |
|
|
|
53 |
|
|
! Variables locales avec l'attribut SAVE |
54 |
|
|
!**************************************************************************************** |
55 |
|
|
! frequence de lecture des conditions limites (en pas de physique) |
56 |
|
|
INTEGER,SAVE :: lmt_pas |
57 |
|
|
!$OMP THREADPRIVATE(lmt_pas) |
58 |
|
|
! pour indiquer que le jour a lire est deja lu pour une surface precedente |
59 |
|
|
LOGICAL,SAVE :: deja_lu |
60 |
|
|
!$OMP THREADPRIVATE(deja_lu) |
61 |
|
|
INTEGER,SAVE :: jour_lu |
62 |
|
|
!$OMP THREADPRIVATE(jour_lu) |
63 |
|
|
CHARACTER (len = 20),SAVE :: fich ='limit.nc' |
64 |
|
|
!$OMP THREADPRIVATE(fich) |
65 |
|
|
LOGICAL, SAVE :: newlmt = .TRUE. |
66 |
|
|
!$OMP THREADPRIVATE(newlmt) |
67 |
|
|
LOGICAL, SAVE :: check = .FALSE. |
68 |
|
|
!$OMP THREADPRIVATE(check) |
69 |
|
|
REAL, ALLOCATABLE , SAVE, DIMENSION(:) :: sst_lu_p |
70 |
|
|
!$OMP THREADPRIVATE(sst_lu_p) |
71 |
|
|
REAL, ALLOCATABLE , SAVE, DIMENSION(:,:) :: pct_tmp_p |
72 |
|
|
!$OMP THREADPRIVATE(pct_tmp_p) |
73 |
|
|
|
74 |
|
|
! Variables locales |
75 |
|
|
!**************************************************************************************** |
76 |
|
|
INTEGER :: nid, nvarid |
77 |
|
|
INTEGER :: ii |
78 |
|
|
INTEGER :: ierr |
79 |
|
|
INTEGER, DIMENSION(2) :: start, epais |
80 |
|
|
CHARACTER (len = 20) :: modname = 'interfoce_lim' |
81 |
|
|
CHARACTER (len = 80) :: abort_message |
82 |
|
|
REAL, DIMENSION(klon_glo,nbsrf) :: pctsrf_new |
83 |
|
|
REAL, DIMENSION(klon_glo,nbsrf) :: pct_tmp |
84 |
|
|
REAL, DIMENSION(klon_glo) :: sst_lu |
85 |
|
|
REAL, DIMENSION(klon_glo) :: nat_lu |
86 |
|
|
! |
87 |
|
|
! Fin declaration |
88 |
|
|
!**************************************************************************************** |
89 |
|
|
|
90 |
|
|
!**************************************************************************************** |
91 |
|
|
! Start calculation |
92 |
|
|
! |
93 |
|
|
!**************************************************************************************** |
94 |
|
|
IF (debut .AND. .NOT. ALLOCATED(sst_lu_p)) THEN |
95 |
|
|
lmt_pas = NINT(86400./dtime * 1.0) ! pour une lecture une fois par jour |
96 |
|
|
jour_lu = jour - 1 |
97 |
|
|
ALLOCATE(sst_lu_p(klon_loc)) |
98 |
|
|
ALLOCATE(pct_tmp_p(klon_loc,nbsrf)) |
99 |
|
|
ENDIF |
100 |
|
|
|
101 |
|
|
IF ((jour - jour_lu) /= 0) deja_lu = .FALSE. |
102 |
|
|
|
103 |
|
|
IF (check) WRITE(*,*) modname, ' :: jour, jour_lu, deja_lu', jour, jour_lu, deja_lu |
104 |
|
|
IF (check) WRITE(*,*) modname, ' :: itime, lmt_pas ', itime, lmt_pas,dtime |
105 |
|
|
|
106 |
|
|
!**************************************************************************************** |
107 |
|
|
! Ouverture et lecture du fichier pour le master process si c'est le bon moment |
108 |
|
|
! |
109 |
|
|
!**************************************************************************************** |
110 |
|
|
! Tester d'abord si c'est le moment de lire le fichier |
111 |
|
|
IF (MOD(itime-1, lmt_pas) == 0 .AND. .NOT. deja_lu) THEN |
112 |
|
|
|
113 |
|
|
!$OMP MASTER |
114 |
|
|
IF (is_mpi_root) THEN |
115 |
|
|
|
116 |
|
|
fich = TRIM(fich) |
117 |
|
|
ierr = NF_OPEN (fich, NF_NOWRITE,nid) |
118 |
|
|
IF (ierr.NE.NF_NOERR) THEN |
119 |
|
|
abort_message = 'Pb d''ouverture du fichier de conditions aux limites' |
120 |
|
|
CALL abort_physic(modname,abort_message,1) |
121 |
|
|
ENDIF |
122 |
|
|
|
123 |
|
|
! La tranche de donnees a lire: |
124 |
|
|
|
125 |
|
|
start(1) = 1 |
126 |
|
|
start(2) = jour |
127 |
|
|
epais(1) = klon_glo |
128 |
|
|
epais(2) = 1 |
129 |
|
|
|
130 |
|
|
IF (newlmt) THEN |
131 |
|
|
! |
132 |
|
|
! Fraction "ocean" |
133 |
|
|
! |
134 |
|
|
ierr = NF_INQ_VARID(nid, 'FOCE', nvarid) |
135 |
|
|
IF (ierr /= NF_NOERR) THEN |
136 |
|
|
abort_message = 'Le champ <FOCE> est absent' |
137 |
|
|
CALL abort_physic(modname,abort_message,1) |
138 |
|
|
ENDIF |
139 |
|
|
#ifdef NC_DOUBLE |
140 |
|
|
ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_oce)) |
141 |
|
|
#else |
142 |
|
|
ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_oce)) |
143 |
|
|
#endif |
144 |
|
|
IF (ierr /= NF_NOERR) THEN |
145 |
|
|
abort_message = 'Lecture echouee pour <FOCE>' |
146 |
|
|
CALL abort_physic(modname,abort_message,1) |
147 |
|
|
ENDIF |
148 |
|
|
! |
149 |
|
|
! Fraction "glace de mer" |
150 |
|
|
! |
151 |
|
|
ierr = NF_INQ_VARID(nid, 'FSIC', nvarid) |
152 |
|
|
IF (ierr /= NF_NOERR) THEN |
153 |
|
|
abort_message = 'Le champ <FSIC> est absent' |
154 |
|
|
CALL abort_physic(modname,abort_message,1) |
155 |
|
|
ENDIF |
156 |
|
|
#ifdef NC_DOUBLE |
157 |
|
|
ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_sic)) |
158 |
|
|
#else |
159 |
|
|
ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_sic)) |
160 |
|
|
#endif |
161 |
|
|
IF (ierr /= NF_NOERR) THEN |
162 |
|
|
abort_message = 'Lecture echouee pour <FSIC>' |
163 |
|
|
CALL abort_physic(modname,abort_message,1) |
164 |
|
|
ENDIF |
165 |
|
|
! |
166 |
|
|
! Fraction "terre" |
167 |
|
|
! |
168 |
|
|
ierr = NF_INQ_VARID(nid, 'FTER', nvarid) |
169 |
|
|
IF (ierr /= NF_NOERR) THEN |
170 |
|
|
abort_message = 'Le champ <FTER> est absent' |
171 |
|
|
CALL abort_physic(modname,abort_message,1) |
172 |
|
|
ENDIF |
173 |
|
|
#ifdef NC_DOUBLE |
174 |
|
|
ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_ter)) |
175 |
|
|
#else |
176 |
|
|
ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_ter)) |
177 |
|
|
#endif |
178 |
|
|
IF (ierr /= NF_NOERR) THEN |
179 |
|
|
abort_message = 'Lecture echouee pour <FTER>' |
180 |
|
|
CALL abort_physic(modname,abort_message,1) |
181 |
|
|
ENDIF |
182 |
|
|
! |
183 |
|
|
! Fraction "glacier terre" |
184 |
|
|
! |
185 |
|
|
ierr = NF_INQ_VARID(nid, 'FLIC', nvarid) |
186 |
|
|
IF (ierr /= NF_NOERR) THEN |
187 |
|
|
abort_message = 'Le champ <FLIC> est absent' |
188 |
|
|
CALL abort_physic(modname,abort_message,1) |
189 |
|
|
ENDIF |
190 |
|
|
#ifdef NC_DOUBLE |
191 |
|
|
ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais,pct_tmp(1,is_lic)) |
192 |
|
|
#else |
193 |
|
|
ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais,pct_tmp(1,is_lic)) |
194 |
|
|
#endif |
195 |
|
|
IF (ierr /= NF_NOERR) THEN |
196 |
|
|
abort_message = 'Lecture echouee pour <FLIC>' |
197 |
|
|
CALL abort_physic(modname,abort_message,1) |
198 |
|
|
ENDIF |
199 |
|
|
! |
200 |
|
|
ELSE ! on en est toujours a rnatur |
201 |
|
|
! |
202 |
|
|
ierr = NF_INQ_VARID(nid, 'NAT', nvarid) |
203 |
|
|
IF (ierr /= NF_NOERR) THEN |
204 |
|
|
abort_message = 'Le champ <NAT> est absent' |
205 |
|
|
CALL abort_physic(modname,abort_message,1) |
206 |
|
|
ENDIF |
207 |
|
|
#ifdef NC_DOUBLE |
208 |
|
|
ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, nat_lu) |
209 |
|
|
#else |
210 |
|
|
ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, nat_lu) |
211 |
|
|
#endif |
212 |
|
|
IF (ierr /= NF_NOERR) THEN |
213 |
|
|
abort_message = 'Lecture echouee pour <NAT>' |
214 |
|
|
CALL abort_physic(modname,abort_message,1) |
215 |
|
|
ENDIF |
216 |
|
|
! |
217 |
|
|
! Remplissage des fractions de surface |
218 |
|
|
! nat = 0, 1, 2, 3 pour ocean, terre, glacier, seaice |
219 |
|
|
! |
220 |
|
|
pct_tmp = 0.0 |
221 |
|
|
DO ii = 1, klon_glo |
222 |
|
|
pct_tmp(ii,NINT(nat_lu(ii)) + 1) = 1. |
223 |
|
|
ENDDO |
224 |
|
|
|
225 |
|
|
! |
226 |
|
|
! On se retrouve avec ocean en 1 et terre en 2 alors qu'on veut le contraire |
227 |
|
|
! |
228 |
|
|
pctsrf_new = pct_tmp |
229 |
|
|
pctsrf_new (:,2)= pct_tmp (:,1) |
230 |
|
|
pctsrf_new (:,1)= pct_tmp (:,2) |
231 |
|
|
pct_tmp = pctsrf_new |
232 |
|
|
ENDIF ! fin test sur newlmt |
233 |
|
|
! |
234 |
|
|
! Lecture SST |
235 |
|
|
! |
236 |
|
|
ierr = NF_INQ_VARID(nid, 'SST', nvarid) |
237 |
|
|
IF (ierr /= NF_NOERR) THEN |
238 |
|
|
abort_message = 'Le champ <SST> est absent' |
239 |
|
|
CALL abort_physic(modname,abort_message,1) |
240 |
|
|
ENDIF |
241 |
|
|
#ifdef NC_DOUBLE |
242 |
|
|
ierr = NF_GET_VARA_DOUBLE(nid,nvarid,start,epais, sst_lu) |
243 |
|
|
#else |
244 |
|
|
ierr = NF_GET_VARA_REAL(nid,nvarid,start,epais, sst_lu) |
245 |
|
|
#endif |
246 |
|
|
IF (ierr /= NF_NOERR) THEN |
247 |
|
|
abort_message = 'Lecture echouee pour <SST>' |
248 |
|
|
CALL abort_physic(modname,abort_message,1) |
249 |
|
|
ENDIF |
250 |
|
|
|
251 |
|
|
!**************************************************************************************** |
252 |
|
|
! Fin de lecture, fermeture de fichier |
253 |
|
|
! |
254 |
|
|
!**************************************************************************************** |
255 |
|
|
ierr = NF_CLOSE(nid) |
256 |
|
|
ENDIF ! is_mpi_root |
257 |
|
|
|
258 |
|
|
!$OMP END MASTER |
259 |
|
|
!$OMP BARRIER |
260 |
|
|
|
261 |
|
|
|
262 |
|
|
!**************************************************************************************** |
263 |
|
|
! Distribue les variables sur tous les processus |
264 |
|
|
! |
265 |
|
|
!**************************************************************************************** |
266 |
|
|
CALL Scatter(sst_lu,sst_lu_p) |
267 |
|
|
CALL Scatter(pct_tmp(:,is_oce),pct_tmp_p(:,is_oce)) |
268 |
|
|
CALL Scatter(pct_tmp(:,is_sic),pct_tmp_p(:,is_sic)) |
269 |
|
|
deja_lu = .TRUE. |
270 |
|
|
jour_lu = jour |
271 |
|
|
ENDIF |
272 |
|
|
|
273 |
|
|
!**************************************************************************************** |
274 |
|
|
! Recopie des variables dans les champs de sortie |
275 |
|
|
! |
276 |
|
|
!**************************************************************************************** |
277 |
|
|
lmt_sst_p = 999999999. |
278 |
|
|
|
279 |
|
|
DO ii = 1, knon |
280 |
|
|
lmt_sst_p(ii) = sst_lu_p(knindex(ii)) |
281 |
|
|
ENDDO |
282 |
|
|
|
283 |
|
|
DO ii=1,klon_loc |
284 |
|
|
pctsrf_new_p(ii,is_oce)=pct_tmp_p(ii,is_oce) |
285 |
|
|
pctsrf_new_p(ii,is_sic)=pct_tmp_p(ii,is_sic) |
286 |
|
|
ENDDO |
287 |
|
|
|
288 |
|
|
|
289 |
|
|
END SUBROUTINE interfoce_lim |