SURFEX v8.1
General documentation of Surfex
read_oceann.F90
Go to the documentation of this file.
1 !SFX_LIC Copyright 1994-2014 CNRS, Meteo-France and Universite Paul Sabatier
2 !SFX_LIC This is part of the SURFEX software governed by the CeCILL-C licence
3 !SFX_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
4 !SFX_LIC for details. version 1.
5 ! #########
6  SUBROUTINE read_ocean_n (DTCO, O, OR, U, HPROGRAM)
7 ! #########################################
8 !
9 !!**** *READ_OCEAN_n* - read oceanic variables
10 !!
11 !!
12 !! PURPOSE
13 !! -------
14 !!
15 !!** METHOD
16 !! ------
17 !!
18 !! EXTERNAL
19 !! --------
20 !!
21 !!
22 !! IMPLICIT ARGUMENTS
23 !! ------------------
24 !!
25 !! REFERENCE
26 !! ---------
27 !!
28 !!
29 !! AUTHOR
30 !! ------
31 !! C. Lebeaupin Brossier *Meteo France*
32 !!
33 !! MODIFICATIONS
34 !! -------------
35 !! Original 04/2007
36 !! Modofied 07/2012, P. Le Moigne : CMO1D phasing
37 !-------------------------------------------------------------------------------
38 !
39 !* 0. DECLARATIONS
40 ! ------------
41 !
42 !
43 !
44 !
46 USE modd_ocean_n, ONLY : ocean_t
47 USE modd_ocean_rel_n, ONLY : ocean_rel_t
48 USE modd_surf_atm_n, ONLY : surf_atm_t
49 !
50 USE modd_surf_par, ONLY : xundef
52 !
53 !
55 USE modi_ocean_mercatorvergrid
56 !
57 !
58 USE yomhook ,ONLY : lhook, dr_hook
59 USE parkind1 ,ONLY : jprb
60 !
61 USE modi_get_type_dim_n
62 !
63 IMPLICIT NONE
64 !
65 !* 0.1 Declarations of arguments
66 ! -------------------------
67 !
68 !
69 TYPE(data_cover_t), INTENT(INOUT) :: DTCO
70 TYPE(ocean_t), INTENT(INOUT) :: O
71 TYPE(ocean_rel_t), INTENT(INOUT) :: OR
72 TYPE(surf_atm_t), INTENT(INOUT) :: U
73 !
74  CHARACTER(LEN=6), INTENT(IN) :: HPROGRAM ! calling program
75 !
76 !* 0.2 Declarations of local variables
77 ! -------------------------------
78 !
79 INTEGER :: ILU ! 1D physical dimension
80 !
81 INTEGER :: IRESP ! Error code after redding
82 !
83  CHARACTER(LEN=4) :: YLVL
84 !
85  CHARACTER(LEN=12) :: YRECFM ! Name of the article to be read
86  CHARACTER(LEN=14) :: YFORM ! Writing format
87 REAL, DIMENSION(:),ALLOCATABLE :: ZWORK ! 1D array to write data in file
88 !
89 INTEGER :: JLEVEL ! loop counter on oceanic levels
90 INTEGER :: J ! loop counter on sea grid points
91 !
92 INTEGER :: IVERSION ! surface version
93 REAL(KIND=JPRB) :: ZHOOK_HANDLE
94 !
95 !-------------------------------------------------------------------------------
96 IF (lhook) CALL dr_hook('READ_OCEAN_N',0,zhook_handle)
97 o%NOCTCOUNT=0
98 !
99 yrecfm='VERSION'
100  CALL read_surf(hprogram,yrecfm,iversion,iresp)
101 !
102 !* flag to use or not Ocean model
103 !
104 IF (iversion<=3) THEN
105  o%LMERCATOR=.false.
106 ELSE
107  yrecfm='SEA_OCEAN'
108  CALL read_surf(hprogram,yrecfm,o%LMERCATOR,iresp)
109 ENDIF
110 !
111 IF (.NOT. o%LMERCATOR) THEN
112  ALLOCATE(o%XSEAT(0,0))
113  ALLOCATE(o%XSEAS(0,0))
114  ALLOCATE(o%XSEAU(0,0))
115  ALLOCATE(o%XSEAV(0,0))
116  ALLOCATE(or%XSEAT_REL(0,0))
117  ALLOCATE(or%XSEAS_REL(0,0))
118  ALLOCATE(or%XSEAU_REL(0,0))
119  ALLOCATE(or%XSEAV_REL(0,0))
120  !
121  ALLOCATE(o%XSEAE(0,0))
122  ALLOCATE(o%XSEABATH(0,0))
123  ALLOCATE(o%XSEAHMO(0))
124  ALLOCATE(o%XLE (0,0))
125  ALLOCATE(o%XLK (0,0))
126  ALLOCATE(o%XKMEL (0,0))
127  ALLOCATE(o%XKMELM (0,0))
128  ALLOCATE(o%XSEATEND (0))
129  ALLOCATE(o%XDTFSOL(0,0))
130  ALLOCATE(o%XDTFNSOL(0))
131  IF (lhook) CALL dr_hook('READ_OCEAN_N',1,zhook_handle)
132  RETURN
133 ENDIF
134 !
135 !-------------------------------------------------------------------------------
136 !
137 nockmin = 0
138 yrecfm='SEA_NBLEVEL'
139  CALL read_surf(hprogram,yrecfm,nockmax,iresp)
140 !
141 ALLOCATE(xzhoc(nockmin:nockmax))
142 xzhoc(nockmin) = 0.
143 ! Read vertical grid
144 DO jlevel = nockmin+1,nockmax
145  WRITE(ylvl,'(I4)') jlevel
146  yrecfm='LEVL_OC'//adjustl(ylvl(:len_trim(ylvl)))
147  CALL read_surf(hprogram,yrecfm,xzhoc(jlevel),iresp)
148 END DO
149 !
151 !
152 ! Relaxation time and logical
153 yrecfm='TAU_REL_OC'
154  CALL read_surf(hprogram,yrecfm,or%XTAU_REL,iresp)
155 !
156 yrecfm='LREL_CUR_OC'
157  CALL read_surf(hprogram,yrecfm,or%LREL_CUR,iresp)
158 
159 yrecfm='LREL_TS_OC'
160  CALL read_surf(hprogram,yrecfm,or%LREL_TS,iresp)
161 yrecfm='LFLX_NULL_OC'
162  CALL read_surf(hprogram,yrecfm,or%LFLUX_NULL,iresp)
163 yrecfm='LFLX_CORR_OC'
164  CALL read_surf(hprogram,yrecfm,or%LFLX_CORR,iresp)
165 yrecfm='CORR_FLX_OC'
166  CALL read_surf(hprogram,yrecfm,or%XQCORR,iresp)
167 yrecfm='LDIAPYC_OC'
168  CALL read_surf(hprogram,yrecfm,or%LDIAPYCNAL,iresp)
169 !
170 !* 1D physical dimension
171 !
172 yrecfm='SIZE_SEA'
173  CALL get_type_dim_n(dtco, u, 'SEA ',ilu)
174 !
175 !* 2. Prognostic fields:
176 ! -----------------
177 !
178 ALLOCATE(zwork(ilu))
179 !* oceanic temperature
180 !
181 ALLOCATE(o%XSEAT(ilu,nockmin:nockmax))
182 !
183 DO jlevel=nockmin+1,nockmax
184  WRITE(ylvl,'(I4)') jlevel
185  yrecfm='TEMP_OC'//adjustl(ylvl(:len_trim(ylvl)))
186  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
187  o%XSEAT(:,jlevel)=zwork(:)
188 END DO
189 o%XSEAT(:,nockmin)=o%XSEAT(:,nockmin+1)
190 !
191 !* relaxation profile for oceanic temperature
192 !
193 ALLOCATE(or%XSEAT_REL(ilu,nockmin:nockmax))
194 !
195 DO jlevel=nockmin+1,nockmax
196  WRITE(ylvl,'(I4)') jlevel
197  yrecfm='T_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
198  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
199  or%XSEAT_REL(:,jlevel)=zwork(:)
200 END DO
201 or%XSEAT_REL(:,nockmin)=or%XSEAT_REL(:,nockmin+1)
202 !
203 !* oceanic salinity
204 !
205 ALLOCATE(o%XSEAS(ilu,nockmin:nockmax))
206 !
207 DO jlevel=nockmin+1,nockmax
208  WRITE(ylvl,'(I4)') jlevel
209  yrecfm='SALT_OC'//adjustl(ylvl(:len_trim(ylvl)))
210  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
211  o%XSEAS(:,jlevel)=zwork(:)
212 END DO
213 o%XSEAS(:,nockmin)=o%XSEAS(:,nockmin+1)
214 !
215 !* oceanic salinity profile of relaxation
216 !
217 ALLOCATE(or%XSEAS_REL(ilu,nockmin:nockmax))
218 !
219 DO jlevel=nockmin+1,nockmax
220  WRITE(ylvl,'(I4)') jlevel
221  yrecfm='S_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
222  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
223  or%XSEAS_REL(:,jlevel)=zwork(:)
224 END DO
225 or%XSEAS_REL(:,nockmin)=or%XSEAS_REL(:,nockmin+1)
226 !
227 !* oceanic current
228 !
229 ALLOCATE(or%XSEAU_REL(ilu,nockmin:nockmax))
230 ALLOCATE(or%XSEAV_REL(ilu,nockmin:nockmax))
231 !
232 DO jlevel=nockmin+1,nockmax
233  WRITE(ylvl,'(I4)') jlevel
234  yrecfm='U_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
235  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
236  or%XSEAU_REL(:,jlevel)=zwork(:)
237 END DO
238 or%XSEAU_REL(:,nockmin)=or%XSEAU_REL(:,nockmin+1)
239 !
240 DO jlevel=nockmin+1,nockmax
241  WRITE(ylvl,'(I4)') jlevel
242  yrecfm='V_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
243  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
244  or%XSEAV_REL(:,jlevel)=zwork(:)
245 END DO
246 or%XSEAV_REL(:,nockmin)=or%XSEAV_REL(:,nockmin+1)
247 !
248 ALLOCATE(o%XSEAU(ilu,nockmin:nockmax))
249 ALLOCATE(o%XSEAV(ilu,nockmin:nockmax))
250 !
251 DO jlevel=nockmin+1,nockmax
252  WRITE(ylvl,'(I4)') jlevel
253  yrecfm='UCUR_OC'//adjustl(ylvl(:len_trim(ylvl)))
254  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
255  o%XSEAU(:,jlevel)=zwork(:)
256 END DO
257 DO jlevel=nockmin+1,nockmax
258  WRITE(ylvl,'(I4)') jlevel
259  yrecfm='VCUR_OC'//adjustl(ylvl(:len_trim(ylvl)))
260  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
261  o%XSEAV(:,jlevel)=zwork(:)
262 END DO
263 o%XSEAU(:,nockmin)=o%XSEAU(:,nockmin+1)
264 o%XSEAV(:,nockmin)=o%XSEAV(:,nockmin+1)
265 !
266 !* oceanic turbulent kinetic energy
267 !
268 ALLOCATE(o%XSEAE(ilu,nockmin:nockmax))
269 !
270 DO jlevel=nockmin+1,nockmax
271  WRITE(ylvl,'(I4)') jlevel
272  yrecfm='TKE_OC'//adjustl(ylvl(:len_trim(ylvl)))
273  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
274  o%XSEAE(:,jlevel)=zwork(:)
275 END DO
276 o%XSEAE(:,nockmin)=o%XSEAE(:,nockmin+1)
277 !
278 !
279 !-------------------------------------------------------------------------------
280 !
281 !* 4. Semi-prognostic fields:
282 ! ----------------------
283 !
284 !* bathymetry indice
285 !
286 ALLOCATE(o%XSEABATH(ilu,nockmin:nockmax))
287 !
288 DO jlevel=nockmin+1,nockmax
289  WRITE(ylvl,'(I4)') jlevel
290  yrecfm='SEAINDBATH'//adjustl(ylvl(:len_trim(ylvl)))
291  CALL read_surf(hprogram,yrecfm,zwork(:),iresp)
292  o%XSEABATH(:,jlevel)=zwork(:)
293 END DO
294 o%XSEABATH(:,nockmin)=1.
295 !
296 !-------------------------------------------------------------------------------
297 !Complete undefined values for the oceanic 1D model convergence
298 DO j=1,ilu
299  DO jlevel=nockmin+2,nockmax
300  IF (o%XSEABATH(j,jlevel)==0.) THEN
301  o%XSEAT(j,jlevel)=o%XSEAT(j,jlevel-1)
302  o%XSEAS(j,jlevel)=o%XSEAS(j,jlevel-1)
303  o%XSEAU(j,jlevel)=o%XSEAU(j,jlevel-1)
304  o%XSEAV(j,jlevel)=o%XSEAV(j,jlevel-1)
305  o%XSEAE(j,jlevel)=o%XSEAE(j,jlevel-1)
306  !
307  or%XSEAT_REL(j,jlevel)=or%XSEAT_REL(j,jlevel-1)
308  or%XSEAS_REL(j,jlevel)=or%XSEAS_REL(j,jlevel-1)
309  or%XSEAU_REL(j,jlevel)=or%XSEAU_REL(j,jlevel-1)
310  or%XSEAV_REL(j,jlevel)=or%XSEAV_REL(j,jlevel-1)
311  ENDIF
312  ENDDO
313 ENDDO
314 !
315 DEALLOCATE(zwork)
316 !-------------------------------------------------------------------------------
317 ALLOCATE(o%XSEAHMO(ilu))
318 yrecfm='SEA_HMO'
319  CALL read_surf(hprogram,yrecfm,o%XSEAHMO(:),iresp)
320 !
321 !-------------------------------------------------------------------------------
322 ALLOCATE(o%XLE (SIZE(o%XSEAT,1),nockmin:nockmax))
323 ALLOCATE(o%XLK (SIZE(o%XSEAT,1),nockmin:nockmax))
324 ALLOCATE(o%XKMEL (SIZE(o%XSEAT,1),nockmin:nockmax))
325 ALLOCATE(o%XKMELM (SIZE(o%XSEAT,1),nockmin:nockmax))
326 o%XLE(:,:) =xundef
327 o%XLK(:,:) =xundef
328 o%XKMEL(:,:) =xundef
329 o%XKMELM(:,:) =xundef
330 !
331 ALLOCATE(o%XSEATEND (SIZE(o%XSEAT,1)))
332 o%XSEATEND(:) =xundef
333 !
334 ALLOCATE(o%XDTFSOL(ilu,nockmin:nockmax))
335 ALLOCATE(o%XDTFNSOL(ilu))
336 !
337 o%XDTFSOL(:,:) = xundef
338 o%XDTFNSOL(:) = xundef
339 !
340 IF (lhook) CALL dr_hook('READ_OCEAN_N',1,zhook_handle)
341 !
342 !------------------------------------------------------------------------------
343 END SUBROUTINE read_ocean_n
subroutine get_type_dim_n(DTCO, U, HTYPE, KDIM)
subroutine read_ocean_n(DTCO, O, OR, U, HPROGRAM)
Definition: read_oceann.F90:7
real, dimension(:), pointer xzhoc
real, parameter xundef
integer, parameter jprb
Definition: parkind1.F90:32
logical lhook
Definition: yomhook.F90:15
integer, save nockmax
subroutine ocean_mercatorvergrid
integer, save nockmin