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