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