SURFEX  V8_0
Surfex V8_0 release
 All Classes Files Functions Variables
write_diag_seb_seafluxn.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 write_diag_seb_seaflux_n (DTCO, DGU, U, CHS, DGS, S, &
7  hprogram)
8 ! #################################
9 !
10 !!**** *WRITE_DIAG_SEB_SEAFLUX_n* - write the SEAFLUX diagnostic fields
11 !!
12 !! PURPOSE
13 !! -------
14 !!
15 !!
16 !!** METHOD
17 !! ------
18 !!
19 !! REFERENCE
20 !! ---------
21 !!
22 !!
23 !! AUTHOR
24 !! ------
25 !! V. Masson *Meteo France*
26 !!
27 !! MODIFICATIONS
28 !! -------------
29 !! Original 01/2004
30 !! Modified 01/2006 : sea flux parameterization.
31 !! Modified 08/2009 : cumulated diag
32 !! B. Decharme 06/2013 : Add evap and sublimation diag
33 !! Delete LPROVAR_TO_DIAG here
34 !! S.Senesi 01/2014 : add diags on seaice
35 !! S. Belamari 06/2014 : Introduce GRESET to avoid errors due to NBLOCK=0
36 !! when coupled with ARPEGE/ALADIN/AROME
37 !! B. Decharme 02/2016 : NBLOCK instead of LCOUNTW for compilation in AAA
38 !-------------------------------------------------------------------------------
39 !
40 !* 0. DECLARATIONS
41 ! ------------
42 !
43 !
46 USE modd_surf_atm_n, ONLY : surf_atm_t
49 USE modd_seaflux_n, ONLY : seaflux_t
50 !
51 !
52 #ifdef SFX_ARO
53 USE modd_io_surf_aro, ONLY : nblock
54 #endif
55 !
56 !
57 USE modd_surf_par, ONLY : xundef
58 !
59 !
60 !
61 !
62 !
63 USE modi_init_io_surf_n
65 USE modi_end_io_surf_n
66 !
67 USE yomhook ,ONLY : lhook, dr_hook
68 USE parkind1 ,ONLY : jprb
69 !
70 IMPLICIT NONE
71 !
72 !* 0.1 Declarations of arguments
73 ! -------------------------
74 !
75 !
76 TYPE(data_cover_t), INTENT(INOUT) :: dtco
77 TYPE(diag_surf_atm_t), INTENT(INOUT) :: dgu
78 TYPE(surf_atm_t), INTENT(INOUT) :: u
79 TYPE(ch_seaflux_t), INTENT(INOUT) :: chs
80 TYPE(diag_seaflux_t), INTENT(INOUT) :: dgs
81 TYPE(seaflux_t), INTENT(INOUT) :: s
82 !
83  CHARACTER(LEN=6), INTENT(IN) :: hprogram ! program calling
84 !
85 !* 0.2 Declarations of local variables
86 ! -------------------------------
87 !
88 INTEGER :: iresp ! IRESP : return-code if a problem appears
89  CHARACTER(LEN=12) :: yrecfm ! Name of the article to be read
90  CHARACTER(LEN=100):: ycomment ! Comment string
91  CHARACTER(LEN=2) :: ynum
92 !
93 LOGICAL :: greset
94 INTEGER :: jsv, jsw
95 LOGICAL :: gmisc
96 !
97 REAL(KIND=JPRB) :: zhook_handle
98 !
99 !-------------------------------------------------------------------------------
100 !
101 IF (lhook) CALL dr_hook('WRITE_DIAG_SEB_SEAFLUX_N',0,zhook_handle)
102 !
103 ! Initialisation for IO
104 !
105 greset=.true.
106 #ifdef SFX_ARO
107 greset=(nblock>0)
108 #endif
109 !
110  CALL init_io_surf_n(dtco, dgu, u, &
111  hprogram,'SEA ','SEAFLX','WRITE')
112 !
113 !
114 !* 1. Surface temperature :
115 ! ---------------------
116 !
117 gmisc=(dgs%N2M>=1.OR.dgs%LSURF_BUDGET.OR.dgs%LSURF_BUDGETC)
118 !
119 IF (gmisc.AND.s%LHANDLE_SIC) THEN
120  !
121  yrecfm='TS_SEA'
122  ycomment='X_Y_'//yrecfm//' (K)'
123  !
124  CALL write_surf(dgu, u, &
125  hprogram,yrecfm,dgs%XTS(:),iresp,hcomment=ycomment)
126  !
127  yrecfm='TSRAD_SEA'
128  ycomment='X_Y_'//yrecfm//' (K)'
129  !
130  CALL write_surf(dgu, u, &
131  hprogram,yrecfm,dgs%XTSRAD(:),iresp,hcomment=ycomment)
132  !
133 ENDIF
134 !
135 !* 2. Richardson number :
136 ! -----------------
137 !
138 IF (dgs%N2M>=1) THEN
139  !
140  yrecfm='RI_SEA'
141  ycomment='X_Y_'//yrecfm
142  !
143  CALL write_surf(dgu, u, &
144  hprogram,yrecfm,dgs%XRI(:),iresp,hcomment=ycomment)
145  !
146 ENDIF
147  !
148  !* 3. Energy fluxes :
149  ! -------------
150  !
151 IF (dgs%LSURF_BUDGET) THEN
152  !
153  yrecfm='RN_SEA'
154  ycomment='X_Y_'//yrecfm//' (W/m2)'
155  !
156  CALL write_surf(dgu, u, &
157  hprogram,yrecfm,dgs%XRN(:),iresp,hcomment=ycomment)
158  !
159  yrecfm='H_SEA'
160  ycomment='X_Y_'//yrecfm//' (W/m2)'
161  !
162  CALL write_surf(dgu, u, &
163  hprogram,yrecfm,dgs%XH(:),iresp,hcomment=ycomment)
164  !
165  yrecfm='LE_SEA'
166  ycomment='X_Y_'//yrecfm//' (W/m2)'
167  !
168  CALL write_surf(dgu, u, &
169  hprogram,yrecfm,dgs%XLE(:),iresp,hcomment=ycomment)
170  !
171  yrecfm='LEI_SEA'
172  ycomment='X_Y_'//yrecfm//' (W/m2)'
173  CALL write_surf(dgu, u, &
174  hprogram,yrecfm,dgs%XLE_ICE(:),iresp,hcomment=ycomment)
175  !
176  yrecfm='GFLUX_SEA'
177  ycomment='X_Y_'//yrecfm//' (W/m2)'
178  !
179  CALL write_surf(dgu, u, &
180  hprogram,yrecfm,dgs%XGFLUX(:),iresp,hcomment=ycomment)
181  !
182  yrecfm='EVAP_SEA'
183  ycomment='X_Y_'//yrecfm//' (kg/m2/s)'
184  !
185  CALL write_surf(dgu, u, &
186  hprogram,yrecfm,dgs%XEVAP(:),iresp,hcomment=ycomment)
187  !
188  yrecfm='SUBL_SEA'
189  ycomment='X_Y_'//yrecfm//' (kg/m2/s)'
190  !
191  CALL write_surf(dgu, u, &
192  hprogram,yrecfm,dgs%XSUBL(:),iresp,hcomment=ycomment)
193  !
194  IF (dgs%LRAD_BUDGET) THEN
195  !
196  yrecfm='SWD_SEA'
197  ycomment='X_Y_'//yrecfm//' (W/m2)'
198  !
199  CALL write_surf(dgu, u, &
200  hprogram,yrecfm,dgs%XSWD(:),iresp,hcomment=ycomment)
201  !
202  yrecfm='SWU_SEA'
203  ycomment='X_Y_'//yrecfm//' (W/m2)'
204  !
205  CALL write_surf(dgu, u, &
206  hprogram,yrecfm,dgs%XSWU(:),iresp,hcomment=ycomment)
207  !
208  yrecfm='LWD_SEA'
209  ycomment='X_Y_'//yrecfm//' (W/m2)'
210  !
211  CALL write_surf(dgu, u, &
212  hprogram,yrecfm,dgs%XLWD(:),iresp,hcomment=ycomment)
213  !
214  yrecfm='LWU_SEA'
215  ycomment='X_Y_'//yrecfm//' (W/m2)'
216  !
217  CALL write_surf(dgu, u, &
218  hprogram,yrecfm,dgs%XLWU(:),iresp,hcomment=ycomment)
219  !
220  DO jsw=1, SIZE(dgs%XSWBD,2)
221  ynum=achar(48+jsw)
222  !
223  yrecfm='SWD_SEA_'//ynum
224  ycomment='X_Y_'//yrecfm//' (W/m2)'
225  !
226  CALL write_surf(dgu, u, &
227  hprogram,yrecfm,dgs%XSWBD(:,jsw),iresp,hcomment=ycomment)
228  !
229  yrecfm='SWU_SEA_'//ynum
230  ycomment='X_Y_'//yrecfm//' (W/m2)'
231  !
232  CALL write_surf(dgu, u, &
233  hprogram,yrecfm,dgs%XSWBU(:,jsw),iresp,hcomment=ycomment)
234  !
235  ENDDO
236  !
237  ENDIF
238  !
239  yrecfm='FMU_SEA'
240  ycomment='X_Y_'//yrecfm//' (kg/ms2)'
241  !
242  CALL write_surf(dgu, u, &
243  hprogram,yrecfm,dgs%XFMU(:),iresp,hcomment=ycomment)
244  !
245  yrecfm='FMV_SEA'
246  ycomment='X_Y_'//yrecfm//' (kg/ms2)'
247  !
248  CALL write_surf(dgu, u, &
249  hprogram,yrecfm,dgs%XFMV(:),iresp,hcomment=ycomment)
250  !
251 END IF
252 !
253 IF (dgs%LSURF_BUDGETC) THEN
254  !
255  yrecfm='RNC_SEA'
256  ycomment='X_Y_'//yrecfm//' (J/m2)'
257  !
258  CALL write_surf(dgu, u, &
259  hprogram,yrecfm,dgs%XRNC(:),iresp,hcomment=ycomment)
260  !
261  yrecfm='HC_SEA'
262  ycomment='X_Y_'//yrecfm//' (J/m2)'
263  !
264  CALL write_surf(dgu, u, &
265  hprogram,yrecfm,dgs%XHC(:),iresp,hcomment=ycomment)
266  !
267  yrecfm='LEC_SEA'
268  ycomment='X_Y_'//yrecfm//' (J/m2)'
269  !
270  CALL write_surf(dgu, u, &
271  hprogram,yrecfm,dgs%XLEC(:),iresp,hcomment=ycomment)
272  !
273  yrecfm='LEIC_SEA'
274  ycomment='X_Y_'//yrecfm//' (W/m2)'
275  CALL write_surf(dgu, u, &
276  hprogram,yrecfm,dgs%XLEC_ICE(:),iresp,hcomment=ycomment)
277  !
278  yrecfm='GFLUXC_SEA'
279  ycomment='X_Y_'//yrecfm//' (J/m2)'
280  !
281  CALL write_surf(dgu, u, &
282  hprogram,yrecfm,dgs%XGFLUXC(:),iresp,hcomment=ycomment)
283  !
284  yrecfm='EVAPC_SEA'
285  ycomment='X_Y_'//yrecfm//' (kg/m2)'
286  !
287  CALL write_surf(dgu, u, &
288  hprogram,yrecfm,dgs%XEVAPC(:),iresp,hcomment=ycomment)
289  !
290  yrecfm='SUBLC_SEA'
291  ycomment='X_Y_'//yrecfm//' (kg/m2)'
292  !
293  CALL write_surf(dgu, u, &
294  hprogram,yrecfm,dgs%XSUBLC(:),iresp,hcomment=ycomment)
295  !
296  IF (dgs%LRAD_BUDGET .OR. (dgs%LSURF_BUDGETC .AND. .NOT.dgu%LRESET_BUDGETC)) THEN
297  !
298  yrecfm='SWDC_SEA'
299  ycomment='X_Y_'//yrecfm//' (J/m2)'
300  !
301  CALL write_surf(dgu, u, &
302  hprogram,yrecfm,dgs%XSWDC(:),iresp,hcomment=ycomment)
303  !
304  yrecfm='SWUC_SEA'
305  ycomment='X_Y_'//yrecfm//' (J/m2)'
306  !
307  CALL write_surf(dgu, u, &
308  hprogram,yrecfm,dgs%XSWUC(:),iresp,hcomment=ycomment)
309  !
310  yrecfm='LWDC_SEA'
311  ycomment='X_Y_'//yrecfm//' (J/m2)'
312  !
313  CALL write_surf(dgu, u, &
314  hprogram,yrecfm,dgs%XLWDC(:),iresp,hcomment=ycomment)
315  !
316  yrecfm='LWUC_SEA'
317  ycomment='X_Y_'//yrecfm//' (J/m2)'
318  !
319  CALL write_surf(dgu, u, &
320  hprogram,yrecfm,dgs%XLWUC(:),iresp,hcomment=ycomment)
321  !
322  ENDIF
323  !
324  yrecfm='FMUC_SEA'
325  ycomment='X_Y_'//yrecfm//' (kg/ms)'
326  !
327  CALL write_surf(dgu, u, &
328  hprogram,yrecfm,dgs%XFMUC(:),iresp,hcomment=ycomment)
329  !
330  yrecfm='FMVC_SEA'
331  ycomment='X_Y_'//yrecfm//' (kg/ms)'
332  !
333  CALL write_surf(dgu, u, &
334  hprogram,yrecfm,dgs%XFMVC(:),iresp,hcomment=ycomment)
335  !
336 END IF
337 !
338 IF (dgs%LSURF_BUDGET.OR.dgs%LSURF_BUDGETC) THEN
339 !
340  yrecfm='TALB_SEA'
341  ycomment='total albedo over tile sea (-)'
342  CALL write_surf(dgu, u, &
343  hprogram,yrecfm,dgs%XALBT(:),iresp,hcomment=ycomment)
344 !
345 ENDIF
346 !
347 !* 4. transfer coefficients
348 ! ---------------------
349 !
350 IF (dgs%LCOEF) THEN
351  !
352  yrecfm='CD_SEA'
353  ycomment='X_Y_'//yrecfm//' (W/s2)'
354  !
355  CALL write_surf(dgu, u, &
356  hprogram,yrecfm,dgs%XCD(:),iresp,hcomment=ycomment)
357  !
358  yrecfm='CH_SEA'
359  ycomment='X_Y_'//yrecfm//' (W/s)'
360  !
361  CALL write_surf(dgu, u, &
362  hprogram,yrecfm,dgs%XCH(:),iresp,hcomment=ycomment)
363  !
364  yrecfm='CE_SEA'
365  ycomment='X_Y_'//yrecfm//' (W/s/K)'
366  !
367  CALL write_surf(dgu, u, &
368  hprogram,yrecfm,dgs%XCE(:),iresp,hcomment=ycomment)
369  !
370  yrecfm='Z0_SEA'
371  ycomment='X_Y_'//yrecfm//' (M)'
372  !
373  CALL write_surf(dgu, u, &
374  hprogram,yrecfm,dgs%XZ0(:),iresp,hcomment=ycomment)
375  !
376  yrecfm='Z0H_SEA'
377  ycomment='X_Y_'//yrecfm//' (M)'
378  !
379  CALL write_surf(dgu, u, &
380  hprogram,yrecfm,dgs%XZ0H(:),iresp,hcomment=ycomment)
381  !
382 END IF
383 !
384 !
385 !* 5. Surface humidity
386 ! ----------------
387 !
388 IF (dgs%LSURF_VARS) THEN
389  !
390  yrecfm='QS_SEA'
391  ycomment='X_Y_'//yrecfm//' (KG/KG)'
392  !
393  CALL write_surf(dgu, u, &
394  hprogram,yrecfm,dgs%XQS(:),iresp,hcomment=ycomment)
395  !
396 ENDIF
397 !
398 !
399 !* 6. parameters at 2 and 10 meters :
400 ! -----------------------------
401 !
402 IF (dgs%N2M>=1) THEN
403  !
404  yrecfm='T2M_SEA'
405  ycomment='X_Y_'//yrecfm//' (K)'
406  !
407  CALL write_surf(dgu, u, &
408  hprogram,yrecfm,dgs%XT2M(:),iresp,hcomment=ycomment)
409  !
410  yrecfm='T2MMIN_SEA'
411  ycomment='X_Y_'//yrecfm//' (K)'
412  !
413  CALL write_surf(dgu, u, &
414  hprogram,yrecfm,dgs%XT2M_MIN(:),iresp,hcomment=ycomment)
415  IF(greset)dgs%XT2M_MIN(:)=xundef
416  !
417  yrecfm='T2MMAX_SEA'
418  ycomment='X_Y_'//yrecfm//' (K)'
419  !
420  CALL write_surf(dgu, u, &
421  hprogram,yrecfm,dgs%XT2M_MAX(:),iresp,hcomment=ycomment)
422  IF(greset)dgs%XT2M_MAX(:)=0.0
423  !
424  yrecfm='Q2M_SEA'
425  ycomment='X_Y_'//yrecfm//' (KG/KG)'
426  !
427  CALL write_surf(dgu, u, &
428  hprogram,yrecfm,dgs%XQ2M(:),iresp,hcomment=ycomment)
429  !
430  yrecfm='HU2M_SEA'
431  ycomment='X_Y_'//yrecfm//' (-)'
432  !
433  CALL write_surf(dgu, u, &
434  hprogram,yrecfm,dgs%XHU2M(:),iresp,hcomment=ycomment)
435  !
436  yrecfm='HU2MMIN_SEA'
437  ycomment='X_Y_'//yrecfm//' (-)'
438  !
439  CALL write_surf(dgu, u, &
440  hprogram,yrecfm,dgs%XHU2M_MIN(:),iresp,hcomment=ycomment)
441  IF(greset)dgs%XHU2M_MIN(:)=xundef
442  !
443  yrecfm='HU2MMAX_SEA'
444  ycomment='X_Y_'//yrecfm//' (-)'
445  !
446  CALL write_surf(dgu, u, &
447  hprogram,yrecfm,dgs%XHU2M_MAX(:),iresp,hcomment=ycomment)
448  IF(greset)dgs%XHU2M_MAX(:)=-xundef
449  !
450  yrecfm='ZON10M_SEA'
451  ycomment='X_Y_'//yrecfm//' (M/S)'
452  !
453  CALL write_surf(dgu, u, &
454  hprogram,yrecfm,dgs%XZON10M(:),iresp,hcomment=ycomment)
455  !
456  yrecfm='MER10M_SEA'
457  ycomment='X_Y_'//yrecfm//' (M/S)'
458  !
459  CALL write_surf(dgu, u, &
460  hprogram,yrecfm,dgs%XMER10M(:),iresp,hcomment=ycomment)
461  !
462  yrecfm='W10M_SEA'
463  ycomment='X_Y_'//yrecfm//' (M/S)'
464  !
465  CALL write_surf(dgu, u, &
466  hprogram,yrecfm,dgs%XWIND10M(:),iresp,hcomment=ycomment)
467  !
468  yrecfm='W10MMAX_SEA'
469  ycomment='X_Y_'//yrecfm//' (M/S)'
470  !
471  CALL write_surf(dgu, u, &
472  hprogram,yrecfm,dgs%XWIND10M_MAX(:),iresp,hcomment=ycomment)
473  IF(greset)dgs%XWIND10M_MAX(:)=0.0
474  !
475 END IF
476 !
477 !
478 !* 7. chemical diagnostics:
479 ! --------------------
480 !
481 IF (chs%SVS%NBEQ>0 .AND. chs%CCH_DRY_DEP=="WES89 ") THEN
482  DO jsv = 1,SIZE(chs%CCH_NAMES,1)
483  yrecfm='DV_SEA_'//trim(chs%CCH_NAMES(jsv))
484  WRITE(ycomment,'(A13,I3.3)')'(m/s) DV_SEA_',jsv
485  CALL write_surf(dgu, u, &
486  hprogram,yrecfm,chs%XDEP(:,jsv),iresp,hcomment=ycomment)
487  END DO
488 ENDIF
489 !
490 !------------------------------------------------------------------------------
491 !
492 ! End of IO
493 !
494  CALL end_io_surf_n(hprogram)
495 !
496 IF (lhook) CALL dr_hook('WRITE_DIAG_SEB_SEAFLUX_N',1,zhook_handle)
497 !
498 !
499 END SUBROUTINE write_diag_seb_seaflux_n
subroutine init_io_surf_n(DTCO, DGU, U, HPROGRAM, HMASK, HSCHEME, HACTION)
subroutine end_io_surf_n(HPROGRAM)
Definition: end_io_surfn.F90:6
subroutine write_diag_seb_seaflux_n(DTCO, DGU, U, CHS, DGS, S, HPROGRAM)