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