SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
writesurf_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 writesurf_ocean_n (DGU, U, &
7  o, or, &
8  hprogram)
9 ! ########################################
10 !
11 !!**** *WRITE_OCEAN_n* - writes OCEAN fields
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 03/2007
37 !! Modified 07/2012, P. Le Moigne : CMO1D phasing
38 !! 11/2014 (David BARBARY) : Write oceanic level
39 !-------------------------------------------------------------------------------
40 !
41 !* 0. DECLARATIONS
42 ! ------------
43 !
44 !
45 !
46 !
48 USE modd_surf_atm_n, ONLY : surf_atm_t
49 !
50 USE modd_ocean_n, ONLY : ocean_t
51 USE modd_ocean_rel_n, ONLY : ocean_rel_t
52 !
54 !
56 !
57 !
58 USE yomhook ,ONLY : lhook, dr_hook
59 USE parkind1 ,ONLY : jprb
60 !
61 IMPLICIT NONE
62 !
63 !* 0.1 Declarations of arguments
64 ! -------------------------
65 !
66 !
67 !
68 TYPE(diag_surf_atm_t), INTENT(INOUT) :: dgu
69 TYPE(surf_atm_t), INTENT(INOUT) :: u
70 !
71 TYPE(ocean_t), INTENT(INOUT) :: o
72 TYPE(ocean_rel_t), INTENT(INOUT) :: or
73 !
74  CHARACTER(LEN=6), INTENT(IN) :: hprogram ! program calling
75 
76 !
77 !* 0.2 Declarations of local variables
78 ! -------------------------------
79 !
80 INTEGER :: iresp ! IRESP : return-code if a problem appears
81  CHARACTER(LEN=12) :: yrecfm ! Name of the article to be read
82  CHARACTER(LEN=4 ) :: ylvl
83  CHARACTER(LEN=100):: ycomment ! Comment string
84  CHARACTER(LEN=14) :: yform ! Writing format
85 !
86 INTEGER :: jlevel ! loop counter on oceanic levels
87 !
88 REAL, DIMENSION(SIZE(O%XSEAT,1)) :: zwork ! 1D array to write data in file
89 REAL(KIND=JPRB) :: zhook_handle
90 !
91 !-------------------------------------------------------------------------------
92 IF (lhook) CALL dr_hook('WRITESURF_OCEAN_N',0,zhook_handle)
93 !
94 !* flag to define if OCEAN is used
95 !
96 yrecfm='SEA_OCEAN'
97 ycomment='flag to use OCEAN model'
98  CALL write_surf(dgu, u, &
99  hprogram,yrecfm,o%LMERCATOR,iresp,hcomment=ycomment)
100 !
101 IF (.NOT. o%LMERCATOR .AND. lhook) CALL dr_hook('WRITESURF_OCEAN_N',1,zhook_handle)
102 IF (.NOT. o%LMERCATOR) RETURN
103 !
104 ! Write Oceanic Level
105 yrecfm='SEA_NBLEVEL'
106 ycomment='Number of OCEAN levels'
107  CALL write_surf(dgu, u, &
108  hprogram,yrecfm,nockmax,iresp,hcomment=ycomment)
109 
110 DO jlevel=nockmin+1,nockmax
111  WRITE(ylvl,'(I4)') jlevel
112  yrecfm='LEVL_OC'//adjustl(ylvl(:len_trim(ylvl)))
113  yform='(A21,I1.1,A4)'
114  IF (jlevel >= 10) yform='(A21,I2.2,A4)'
115  WRITE(ycomment,fmt=yform) 'Depth of OCEAN level ',jlevel,' (m)'
116  CALL write_surf(dgu, u, &
117  hprogram,yrecfm,xzhoc(jlevel),iresp,hcomment=ycomment)
118 END DO
119 !
120 ! Relaxation time
121 yrecfm='TAU_REL_OC'
122 ycomment='Relaxation time of ocean model (s)'
123  CALL write_surf(dgu, u, &
124  hprogram,yrecfm,or%XTAU_REL,iresp,hcomment=ycomment)
125 !
126 yrecfm='LREL_CUR_OC'
127 ycomment='FLAG FOR RELAXATION ON CURRENT'
128  CALL write_surf(dgu, u, &
129  hprogram,yrecfm,or%LREL_CUR,iresp,hcomment=ycomment)
130 !
131 yrecfm='LREL_TS_OC'
132 ycomment='FLAG FOR RELAXATION ON T,S'
133  CALL write_surf(dgu, u, &
134  hprogram,yrecfm,or%LREL_TS,iresp,hcomment=ycomment)
135 !
136 yrecfm='LFLX_NULL_OC'
137 ycomment='FLAG FOR ZERO FLUX '
138  CALL write_surf(dgu, u, &
139  hprogram,yrecfm,or%LFLUX_NULL,iresp,hcomment=ycomment)
140 
141 yrecfm='LFLX_CORR_OC'
142 ycomment='FLAG FOR FLUX CORRECTION '
143  CALL write_surf(dgu, u, &
144  hprogram,yrecfm,or%LFLX_CORR,iresp,hcomment=ycomment)
145 
146 yrecfm='CORR_FLX_OC'
147 ycomment='FLUX CORRECTION COEFF (W/M2/K)'
148  CALL write_surf(dgu, u, &
149  hprogram,yrecfm,or%XQCORR,iresp,hcomment=ycomment)
150 !
151 yrecfm='LDIAPYC_OC'
152 ycomment='FLAG FOR DIAPYCNAL MIXING '
153  CALL write_surf(dgu, u, &
154  hprogram,yrecfm,or%LDIAPYCNAL,iresp,hcomment=ycomment)
155 !
156 !-------------------------------------------------------------------------------
157 !
158 !* 3. Prognostic fields:
159 ! -----------------
160 !
161 ! tendance surface liee au flux non solaire
162 !
163  jlevel=1
164  WRITE(ylvl,'(I4)') jlevel
165  yrecfm='DTFNSOL'//adjustl(ylvl(:len_trim(ylvl)))
166  yform='(A11,I1.1,A5)'
167  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
168 !RJ: extended ascii should be avoided in I/O
169  WRITE(ycomment,fmt=yform) 'X_Y_DTFNSOL',jlevel,' (°C)'
170  zwork=o%XDTFNSOL(:)
171  CALL write_surf(dgu, u, &
172  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
173 !
174 !* oceanic temperature
175 !
176 DO jlevel=nockmin+1,nockmax
177  WRITE(ylvl,'(I4)') jlevel
178  yrecfm='TEMP_OC'//adjustl(ylvl(:len_trim(ylvl)))
179  yform='(A11,I1.1,A5)'
180  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
181 !RJ: extended ascii should be avoided in I/O
182  WRITE(ycomment,fmt=yform) 'X_Y_TEMP_OC',jlevel,' (°C)'
183  zwork=o%XSEAT(:,jlevel)
184  CALL write_surf(dgu, u, &
185  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
186 END DO
187 !
188 ! Budget terms lie au flux solaire
189 DO jlevel=nockmin+1,nockmax
190  WRITE(ylvl,'(I4)') jlevel
191  yrecfm='DTFSOL'//adjustl(ylvl(:len_trim(ylvl)))
192  yform='(A11,I1.1,A5)'
193  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
194 !RJ: extended ascii should be avoided in I/O
195  WRITE(ycomment,fmt=yform) 'X_Y_DTFSOL',jlevel,' (°C/s)'
196  zwork=o%XDTFSOL(:,jlevel)
197  CALL write_surf(dgu, u, &
198  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
199 END DO
200 
201 !
202 !* relaxation profile for oceanic temperature
203 DO jlevel=nockmin+1,nockmax
204  WRITE(ylvl,'(I4)') jlevel
205  yrecfm='T_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
206  yform='(A11,I1.1,A5)'
207  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
208 !RJ: extended ascii should be avoided in I/O
209  WRITE(ycomment,fmt=yform) 'X_Y_T_OC_REL',jlevel,' (°C)'
210  zwork=or%XSEAT_REL(:,jlevel)
211  CALL write_surf(dgu, u, &
212  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
213 END DO
214 !
215 !* oceanic salinity
216 !
217 DO jlevel=nockmin+1,nockmax
218  WRITE(ylvl,'(I4)') jlevel
219  yrecfm='SALT_OC'//adjustl(ylvl(:len_trim(ylvl)))
220  yform='(A11,I1.1,A5)'
221  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
222  WRITE(ycomment,fmt=yform) 'X_Y_SALT_OC',jlevel,'(psu)'
223  zwork=o%XSEAS(:,jlevel)
224  CALL write_surf(dgu, u, &
225  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
226 END DO
227 !
228 !* oceanic salinity profile of relxation
229 !
230 DO jlevel=nockmin+1,nockmax
231  WRITE(ylvl,'(I4)') jlevel
232  yrecfm='S_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
233  yform='(A11,I1.1,A5)'
234  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
235  WRITE(ycomment,fmt=yform) 'X_Y_S_OC_REL',jlevel,'(psu)'
236  zwork=or%XSEAS_REL(:,jlevel)
237  CALL write_surf(dgu, u, &
238  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
239 END DO
240 !
241 !* oceanic zonal current
242 !
243 DO jlevel=nockmin+1,nockmax
244  WRITE(ylvl,'(I4)') jlevel
245  yrecfm='U_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
246  yform='(A11,I1.1,A5)'
247  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
248  WRITE(ycomment,fmt=yform) 'X_Y_U_OC_REL',jlevel,' M/S'
249  zwork=or%XSEAU_REL(:,jlevel)
250  CALL write_surf(dgu, u, &
251  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
252 END DO
253 !
254 DO jlevel=nockmin+1,nockmax
255  WRITE(ylvl,'(I4)') jlevel
256  yrecfm='UCUR_OC'//adjustl(ylvl(:len_trim(ylvl)))
257  yform='(A11,I1.1,A5)'
258  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
259  WRITE(ycomment,fmt=yform) 'X_Y_UCUR_OC',jlevel,' (m/s)'
260  zwork=o%XSEAU(:,jlevel)
261  CALL write_surf(dgu, u, &
262  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
263 END DO
264 !
265 !* oceanic meridian current
266 !
267 DO jlevel=nockmin+1,nockmax
268  WRITE(ylvl,'(I4)') jlevel
269  yrecfm='V_OC_REL'//adjustl(ylvl(:len_trim(ylvl)))
270  yform='(A11,I1.1,A5)'
271  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
272  WRITE(ycomment,fmt=yform) 'X_Y_V_OC_REL',jlevel,' M/S'
273  zwork=or%XSEAV_REL(:,jlevel)
274  CALL write_surf(dgu, u, &
275  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
276 END DO
277 !
278 DO jlevel=nockmin+1,nockmax
279  WRITE(ylvl,'(I4)') jlevel
280  yrecfm='VCUR_OC'//adjustl(ylvl(:len_trim(ylvl)))
281  yform='(A11,I1.1,A5)'
282  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
283  WRITE(ycomment,fmt=yform) 'X_Y_VCUR_OC',jlevel,'(m/s)'
284  zwork=o%XSEAV(:,jlevel)
285  CALL write_surf(dgu, u, &
286  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
287 END DO
288 !
289 !* oceanic turbulent kinetic energy
290 !
291 DO jlevel=nockmin+1,nockmax
292  WRITE(ylvl,'(I4)') jlevel
293  yrecfm='TKE_OC'//adjustl(ylvl(:len_trim(ylvl)))
294  yform='(A11,I1.1,A5)'
295  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
296  WRITE(ycomment,fmt=yform) 'X_Y_TKE_OC',jlevel,' (J)'
297  zwork=o%XSEAE(:,jlevel)
298  CALL write_surf(dgu, u, &
299  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
300 END DO
301 !
302 !* oceanic mixing coefficient
303 !
304 DO jlevel=nockmin+1,nockmax
305  WRITE(ylvl,'(I4)') jlevel
306  yrecfm='KMEL_OC'//adjustl(ylvl(:len_trim(ylvl)))
307  yform='(A11,I1.1,A5)'
308  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
309  WRITE(ycomment,fmt=yform) 'X_Y_KMEL_OC',jlevel,' (m2/s2)'
310  zwork=o%XKMEL(:,jlevel)
311  CALL write_surf(dgu, u, &
312  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
313 END DO
314 
315 !
316 !-------------------------------------------------------------------------------
317 !
318 !* 4. Semi-prognostic fields:
319 ! ----------------------
320 !* bathymetry indice
321 !
322 DO jlevel=nockmin+1,nockmax
323  WRITE(ylvl,'(I4)') jlevel
324  yrecfm='SEAINDBATH'//adjustl(ylvl(:len_trim(ylvl)))
325  yform='(A11,I1.1,A5)'
326  IF (jlevel >= 10) yform='(A11,I2.2,A5)'
327  WRITE(ycomment,fmt=yform) 'X_Y_SEAINDBATH',jlevel,' (J)'
328  zwork=o%XSEABATH(:,jlevel)
329  CALL write_surf(dgu, u, &
330  hprogram,yrecfm,zwork,iresp,hcomment=ycomment)
331 END DO
332 !
333 !-------------------------------------------------------------------------------
334 !Sea Surface Salinity
335 !
336 yrecfm='SSS_OC'
337 ycomment='SSS_OC'
338  CALL write_surf(dgu, u, &
339  hprogram,yrecfm,o%XSEAS(:,nockmin),iresp,hcomment=ycomment)
340 !-------------------------------------------------------------------------------
341 !Sea Surface Salinity
342 !
343 yrecfm='SEA_HMO'
344 ycomment='X_Y_'
345  CALL write_surf(dgu, u, &
346  hprogram,yrecfm,o%XSEAHMO(:),iresp,hcomment=ycomment)
347 !-------------------------------------------------------------------------------
348 IF (lhook) CALL dr_hook('WRITESURF_OCEAN_N',1,zhook_handle)
349 !
350 END SUBROUTINE writesurf_ocean_n
subroutine writesurf_ocean_n(DGU, U, O, OR, HPROGRAM)