SURFEX v8.1
General documentation of Surfex
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, DUO, U, CHS, DSO, D, DC, &
7  OHANDLE_SIC,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 !! S. Senesi 08/15 Add 2nd dimension name for SW bands to write_surf calls
38 !-------------------------------------------------------------------------------
39 !
40 !* 0. DECLARATIONS
41 ! ------------
42 !
44 USE modd_surf_atm_n, ONLY : surf_atm_t
47 !
48 #ifdef SFX_ARO
49 USE modd_io_surf_aro, ONLY : nblock
50 #endif
51 !
52 #ifdef SFX_OL
53 USE modd_io_surf_ol, ONLY : ldef
54 #endif
55 !
57 !
58 USE modd_surf_par, ONLY : xundef
59 !
60 USE modi_init_io_surf_n
62 USE modi_end_io_surf_n
63 !
64 USE yomhook ,ONLY : lhook, dr_hook
65 USE parkind1 ,ONLY : jprb
66 !
67 IMPLICIT NONE
68 !
69 !* 0.1 Declarations of arguments
70 ! -------------------------
71 !
72 !
73 TYPE(data_cover_t), INTENT(INOUT) :: DTCO
74 TYPE(diag_options_t), INTENT(INOUT) :: DUO
75 TYPE(surf_atm_t), INTENT(INOUT) :: U
76 TYPE(ch_seaflux_t), INTENT(INOUT) :: CHS
77 TYPE(diag_options_t), INTENT(INOUT) :: DSO
78 TYPE(diag_t), INTENT(INOUT) :: D
79 TYPE(diag_t), INTENT(INOUT) :: DC
80 LOGICAL, INTENT(IN) :: OHANDLE_SIC
81 !
82  CHARACTER(LEN=6), INTENT(IN) :: HPROGRAM ! program calling
83 !
84 !* 0.2 Declarations of local variables
85 ! -------------------------------
86 !
87 INTEGER :: IRESP ! IRESP : return-code if a problem appears
88  CHARACTER(LEN=12) :: YRECFM ! Name of the article to be read
89  CHARACTER(LEN=100):: YCOMMENT ! Comment string
90  CHARACTER(LEN=2) :: YNUM
91 !
92 LOGICAL :: GRESET
93 INTEGER :: JSV, JSW
94 LOGICAL :: GMISC
95 !
96 REAL(KIND=JPRB) :: ZHOOK_HANDLE
97 !
98 !-------------------------------------------------------------------------------
99 !
100 IF (lhook) CALL dr_hook('WRITE_DIAG_SEB_SEAFLUX_N',0,zhook_handle)
101 !
102 ! Initialisation for IO
103 !
104 greset=.true.
105 #ifdef SFX_ARO
106 greset=(nblock>0)
107 #endif
108 #ifdef SFX_OL
109 IF (ldef) greset = .false.
110 #endif
111 !
112  CALL init_io_surf_n(dtco, u,hprogram,'SEA ','SEAFLX','WRITE','SEAFLUX_DIAGNOSTICS.OUT.nc')
113 !
114 !
115 !* 1. Surface temperature :
116 ! ---------------------
117 !
118 gmisc=(dso%N2M>=1.OR.dso%LSURF_BUDGET.OR.dso%LSURF_BUDGETC)
119 !
120 IF (gmisc.AND.ohandle_sic) THEN
121  !
122  yrecfm='TS_SEA'
123  ycomment='X_Y_'//yrecfm//' (K)'
124  !
125  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XTS(:),iresp,hcomment=ycomment)
126  !
127  yrecfm='TSRAD_SEA'
128  ycomment='X_Y_'//yrecfm//' (K)'
129  !
130  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XTSRAD(:),iresp,hcomment=ycomment)
131  !
132 ENDIF
133 !
134 !* 2. Richardson number :
135 ! -----------------
136 !
137 IF (dso%N2M>=1) THEN
138  !
139  yrecfm='RI_SEA'
140  ycomment='X_Y_'//yrecfm
141  !
142  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XRI(:),iresp,hcomment=ycomment)
143  !
144 ENDIF
145  !
146  !* 3. Energy fluxes :
147  ! -------------
148  !
149 IF (dso%LSURF_BUDGET) THEN
150  !
151  yrecfm='RN_SEA'
152  ycomment='X_Y_'//yrecfm//' (W/m2)'
153  !
154  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XRN(:),iresp,hcomment=ycomment)
155  !
156  yrecfm='H_SEA'
157  ycomment='X_Y_'//yrecfm//' (W/m2)'
158  !
159  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XH(:),iresp,hcomment=ycomment)
160  !
161  yrecfm='LE_SEA'
162  ycomment='X_Y_'//yrecfm//' (W/m2)'
163  !
164  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XLE(:),iresp,hcomment=ycomment)
165  !
166  yrecfm='LEI_SEA'
167  ycomment='X_Y_'//yrecfm//' (W/m2)'
168  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XLEI(:),iresp,hcomment=ycomment)
169  !
170  yrecfm='GFLUX_SEA'
171  ycomment='X_Y_'//yrecfm//' (W/m2)'
172  !
173  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XGFLUX(:),iresp,hcomment=ycomment)
174  !
175  yrecfm='EVAP_SEA'
176  ycomment='X_Y_'//yrecfm//' (kg/m2/s)'
177  !
178  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XEVAP(:),iresp,hcomment=ycomment)
179  !
180  yrecfm='SUBL_SEA'
181  ycomment='X_Y_'//yrecfm//' (kg/m2/s)'
182  !
183  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XSUBL(:),iresp,hcomment=ycomment)
184  !
185  IF (dso%LRAD_BUDGET) THEN
186  !
187  yrecfm='SWD_SEA'
188  ycomment='X_Y_'//yrecfm//' (W/m2)'
189  !
190  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XSWD(:),iresp,hcomment=ycomment)
191  !
192  yrecfm='SWU_SEA'
193  ycomment='X_Y_'//yrecfm//' (W/m2)'
194  !
195  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XSWU(:),iresp,hcomment=ycomment)
196  !
197  yrecfm='LWD_SEA'
198  ycomment='X_Y_'//yrecfm//' (W/m2)'
199  !
200  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XLWD(:),iresp,hcomment=ycomment)
201  !
202  yrecfm='LWU_SEA'
203  ycomment='X_Y_'//yrecfm//' (W/m2)'
204  !
205  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XLWU(:),iresp,hcomment=ycomment)
206  !
207  IF (lallow_add_dim) THEN
208  !
209  yrecfm='SWD_SEA_'
210  ycomment='X_Y_'//yrecfm//' (W/m2)'
211  CALL write_surf(duo%CSELECT,&
212  hprogram,yrecfm,d%XSWBD(:,:),iresp,hcomment=ycomment, hnam_dim=yswband_dim_name)
213  !
214  yrecfm='SWU_SEA_'
215  ycomment='X_Y_'//yrecfm//' (W/m2)'
216  CALL write_surf(duo%CSELECT,&
217  hprogram,yrecfm,d%XSWBU(:,:),iresp,hcomment=ycomment, hnam_dim=yswband_dim_name)
218  !
219  ELSE
220  !
221  DO jsw=1, SIZE(d%XSWBD,2)
222  ynum=achar(48+jsw)
223  !
224  yrecfm='SWD_SEA_'//ynum
225  ycomment='X_Y_'//yrecfm//' (W/m2)'
226  !
227  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XSWBD(:,jsw),iresp,hcomment=ycomment)
228  !
229  yrecfm='SWU_SEA_'//ynum
230  ycomment='X_Y_'//yrecfm//' (W/m2)'
231  !
232  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XSWBU(:,jsw),iresp,hcomment=ycomment)
233  !
234  ENDDO
235  !
236  ENDIF
237  !
238  ENDIF
239  !
240  yrecfm='FMU_SEA'
241  ycomment='X_Y_'//yrecfm//' (kg/ms2)'
242  !
243  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XFMU(:),iresp,hcomment=ycomment)
244  !
245  yrecfm='FMV_SEA'
246  ycomment='X_Y_'//yrecfm//' (kg/ms2)'
247  !
248  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XFMV(:),iresp,hcomment=ycomment)
249  !
250 END IF
251 !
252 IF (dso%LSURF_BUDGET.OR.dso%LSURF_BUDGETC) THEN
253 !
254  yrecfm='TALB_SEA'
255  ycomment='total albedo over tile sea (-)'
256  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XALBT(:),iresp,hcomment=ycomment)
257 !
258 ENDIF
259 !
260 !* 4. transfer coefficients
261 ! ---------------------
262 !
263 IF (dso%LCOEF) THEN
264  !
265  yrecfm='CD_SEA'
266  ycomment='X_Y_'//yrecfm//' (W/s2)'
267  !
268  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XCD(:),iresp,hcomment=ycomment)
269  !
270  yrecfm='CH_SEA'
271  ycomment='X_Y_'//yrecfm//' (W/s)'
272  !
273  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XCH(:),iresp,hcomment=ycomment)
274  !
275  yrecfm='CE_SEA'
276  ycomment='X_Y_'//yrecfm//' (W/s/K)'
277  !
278  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XCE(:),iresp,hcomment=ycomment)
279  !
280  yrecfm='Z0_SEA'
281  ycomment='X_Y_'//yrecfm//' (M)'
282  !
283  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XZ0(:),iresp,hcomment=ycomment)
284  !
285  yrecfm='Z0H_SEA'
286  ycomment='X_Y_'//yrecfm//' (M)'
287  !
288  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XZ0H(:),iresp,hcomment=ycomment)
289  !
290 END IF
291 !
292 !
293 !* 5. Surface humidity
294 ! ----------------
295 !
296 IF (dso%LSURF_VARS) THEN
297  !
298  yrecfm='QS_SEA'
299  ycomment='X_Y_'//yrecfm//' (KG/KG)'
300  !
301  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XQS(:),iresp,hcomment=ycomment)
302  !
303 ENDIF
304 !
305 !
306 !* 6. parameters at 2 and 10 meters :
307 ! -----------------------------
308 !
309 IF (dso%N2M>=1) THEN
310  !
311  yrecfm='T2M_SEA'
312  ycomment='X_Y_'//yrecfm//' (K)'
313  !
314  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XT2M(:),iresp,hcomment=ycomment)
315  !
316  yrecfm='T2MMIN_SEA'
317  ycomment='X_Y_'//yrecfm//' (K)'
318  !
319  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XT2M_MIN(:),iresp,hcomment=ycomment)
320  IF(greset)d%XT2M_MIN(:)=xundef
321  !
322  yrecfm='T2MMAX_SEA'
323  ycomment='X_Y_'//yrecfm//' (K)'
324  !
325  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XT2M_MAX(:),iresp,hcomment=ycomment)
326  IF(greset)d%XT2M_MAX(:)=0.0
327  !
328  yrecfm='Q2M_SEA'
329  ycomment='X_Y_'//yrecfm//' (KG/KG)'
330  !
331  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XQ2M(:),iresp,hcomment=ycomment)
332  !
333  yrecfm='HU2M_SEA'
334  ycomment='X_Y_'//yrecfm//' (-)'
335  !
336  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XHU2M(:),iresp,hcomment=ycomment)
337  !
338  yrecfm='HU2MMIN_SEA'
339  ycomment='X_Y_'//yrecfm//' (-)'
340  !
341  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XHU2M_MIN(:),iresp,hcomment=ycomment)
342  IF(greset)d%XHU2M_MIN(:)=xundef
343  !
344  yrecfm='HU2MMAX_SEA'
345  ycomment='X_Y_'//yrecfm//' (-)'
346  !
347  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XHU2M_MAX(:),iresp,hcomment=ycomment)
348  IF(greset)d%XHU2M_MAX(:)=-xundef
349  !
350  yrecfm='ZON10M_SEA'
351  ycomment='X_Y_'//yrecfm//' (M/S)'
352  !
353  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XZON10M(:),iresp,hcomment=ycomment)
354  !
355  yrecfm='MER10M_SEA'
356  ycomment='X_Y_'//yrecfm//' (M/S)'
357  !
358  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XMER10M(:),iresp,hcomment=ycomment)
359  !
360  yrecfm='W10M_SEA'
361  ycomment='X_Y_'//yrecfm//' (M/S)'
362  !
363  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XWIND10M(:),iresp,hcomment=ycomment)
364  !
365  yrecfm='W10MMAX_SEA'
366  ycomment='X_Y_'//yrecfm//' (M/S)'
367  !
368  CALL write_surf(duo%CSELECT,hprogram,yrecfm,d%XWIND10M_MAX(:),iresp,hcomment=ycomment)
369  IF(greset)d%XWIND10M_MAX(:)=0.0
370  !
371 END IF
372 !
373 !
374 !* 7. chemical diagnostics:
375 ! --------------------
376 !
377 IF (chs%SVS%NBEQ>0 .AND. chs%CCH_DRY_DEP=="WES89 ") THEN
378  DO jsv = 1,SIZE(chs%CCH_NAMES,1)
379  yrecfm='DVSE'//trim(chs%CCH_NAMES(jsv))
380  WRITE(ycomment,'(A13,I3.3)')'(m/s) DV_SEA_',jsv
381  CALL write_surf(duo%CSELECT,hprogram,yrecfm,chs%XDEP(:,jsv),iresp,hcomment=ycomment)
382  END DO
383 ENDIF
384 !
385 IF (dso%LSURF_BUDGETC) THEN
386  !
387  CALL end_io_surf_n(hprogram)
388  CALL init_io_surf_n(dtco, u,hprogram,'SEA ','SEAFLX','WRITE','SEAFLUX_DIAG_CUMUL.OUT.nc')
389  !
390  yrecfm='RNC_SEA'
391  ycomment='X_Y_'//yrecfm//' (J/m2)'
392  !
393  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XRN(:),iresp,hcomment=ycomment)
394  !
395  yrecfm='HC_SEA'
396  ycomment='X_Y_'//yrecfm//' (J/m2)'
397  !
398  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XH(:),iresp,hcomment=ycomment)
399  !
400  yrecfm='LEC_SEA'
401  ycomment='X_Y_'//yrecfm//' (J/m2)'
402  !
403  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XLE(:),iresp,hcomment=ycomment)
404  !
405  yrecfm='LEIC_SEA'
406  ycomment='X_Y_'//yrecfm//' (W/m2)'
407  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XLEI(:),iresp,hcomment=ycomment)
408  !
409  yrecfm='GFLUXC_SEA'
410  ycomment='X_Y_'//yrecfm//' (J/m2)'
411  !
412  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XGFLUX(:),iresp,hcomment=ycomment)
413  !
414  yrecfm='EVAPC_SEA'
415  ycomment='X_Y_'//yrecfm//' (kg/m2)'
416  !
417  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XEVAP(:),iresp,hcomment=ycomment)
418  !
419  yrecfm='SUBLC_SEA'
420  ycomment='X_Y_'//yrecfm//' (kg/m2)'
421  !
422  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XSUBL(:),iresp,hcomment=ycomment)
423  !
424  IF (dso%LRAD_BUDGET .OR. (dso%LSURF_BUDGETC .AND. .NOT.duo%LRESET_BUDGETC)) THEN
425  !
426  yrecfm='SWDC_SEA'
427  ycomment='X_Y_'//yrecfm//' (J/m2)'
428  !
429  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XSWD(:),iresp,hcomment=ycomment)
430  !
431  yrecfm='SWUC_SEA'
432  ycomment='X_Y_'//yrecfm//' (J/m2)'
433  !
434  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XSWU(:),iresp,hcomment=ycomment)
435  !
436  yrecfm='LWDC_SEA'
437  ycomment='X_Y_'//yrecfm//' (J/m2)'
438  !
439  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XLWD(:),iresp,hcomment=ycomment)
440  !
441  yrecfm='LWUC_SEA'
442  ycomment='X_Y_'//yrecfm//' (J/m2)'
443  !
444  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XLWU(:),iresp,hcomment=ycomment)
445  !
446  ENDIF
447  !
448  yrecfm='FMUC_SEA'
449  ycomment='X_Y_'//yrecfm//' (kg/ms)'
450  !
451  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XFMU(:),iresp,hcomment=ycomment)
452  !
453  yrecfm='FMVC_SEA'
454  ycomment='X_Y_'//yrecfm//' (kg/ms)'
455  !
456  CALL write_surf(duo%CSELECT,hprogram,yrecfm,dc%XFMV(:),iresp,hcomment=ycomment)
457  !
458 END IF
459 !
460 !------------------------------------------------------------------------------
461 !
462 ! End of IO
463 !
464  CALL end_io_surf_n(hprogram)
465 !
466 IF (lhook) CALL dr_hook('WRITE_DIAG_SEB_SEAFLUX_N',1,zhook_handle)
467 !
468 !
469 END SUBROUTINE write_diag_seb_seaflux_n
static const char * trim(const char *name, int *n)
Definition: drhook.c:2383
real, parameter xundef
integer, parameter jprb
Definition: parkind1.F90:32
logical lallow_add_dim
Definition: modd_xios.F90:49
subroutine write_diag_seb_seaflux_n(DTCO, DUO, U, CHS, DSO, D, DC
character(len=30) yswband_dim_name
Definition: modd_xios.F90:69
subroutine end_io_surf_n(HPROGRAM)
Definition: end_io_surfn.F90:7
logical lhook
Definition: yomhook.F90:15
subroutine init_io_surf_n(DTCO, U, HPROGRAM, HMASK, HSCHEME, HACTION