source: trunk/LMDZ.GENERIC/libf/phystd/write_output_mod.F90 @ 3523

Last change on this file since 3523 was 3522, checked in by afalco, 9 days ago

Generic PCM: imported write_output() subroutine from Mars PCM to write both diagfi and XIOS.

Still some bug on a few variables from generic condensation, to be checked.

AF

File size: 9.5 KB
Line 
1MODULE write_output_mod
2
3implicit none
4
5private
6
7INTERFACE write_output
8    MODULE PROCEDURE write_output_d0, write_output_d1, write_output_d2, &
9                     write_output_i0, write_output_i1, write_output_i2, &
10                     write_output_l0, write_output_l1, write_output_l2
11END INTERFACE write_output
12
13public write_output
14
15!-----------------------------------------------------------------------
16contains
17!-----------------------------------------------------------------------
18
19SUBROUTINE write_output_d0(field_name,title,units,field)
20! For a scalar
21
22#ifdef CPP_XIOS
23    ! use xios_output_mod, only: xios_is_active_field
24    use xios_output_mod, only: send_xios_field
25#endif
26
27implicit none
28
29include "dimensions.h"
30
31integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
32character(*), intent(in) :: field_name
33character(*), intent(in) :: title
34character(*), intent(in) :: units
35real,         intent(in) :: field
36
37logical :: is_active ! For XIOS, should this field be sent or not
38
39call writediagfi(ngrid,field_name,title,units,0,(/field/))
40#ifdef CPP_XIOS
41!is_active=xios_is_active_field(field_name)
42is_active=.true.
43    ! only send the field to xios if the user asked for it
44    if (is_active) call send_xios_field(field_name,field)
45#endif
46
47END SUBROUTINE write_output_d0
48
49!-----------------------------------------------------------------------
50
51SUBROUTINE write_output_d1(field_name,title,units,field)
52! For a surface field
53
54#ifdef CPP_XIOS
55    ! use xios_output_mod, only: xios_is_active_field
56    use xios_output_mod, only: send_xios_field
57#endif
58
59implicit none
60
61include "dimensions.h"
62
63integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
64character(*),       intent(in) :: field_name
65character(*),       intent(in) :: title
66character(*),       intent(in) :: units
67real, dimension(:), intent(in) :: field
68
69logical :: is_active ! For XIOS, should this field be sent or not
70
71call writediagfi(ngrid,field_name,title,units,2,field)
72#ifdef CPP_XIOS
73!is_active=xios_is_active_field(field_name)
74is_active=.true.
75    ! only send the field to xios if the user asked for it
76    if (is_active) call send_xios_field(field_name,field)
77#endif
78
79END SUBROUTINE write_output_d1
80
81!-----------------------------------------------------------------------
82
83SUBROUTINE write_output_d2(field_name,title,units,field)
84! For a "3D" horizontal-vertical field
85
86#ifdef CPP_XIOS
87    ! use xios_output_mod, only: xios_is_active_field
88    use xios_output_mod, only: send_xios_field
89#endif
90
91use comsoil_h,         only: nsoilmx
92use writediagsoil_mod, only: writediagsoil
93
94implicit none
95
96include "dimensions.h"
97
98integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
99character(*),         intent(in) :: field_name
100character(*),         intent(in) :: title
101character(*),         intent(in) :: units
102real, dimension(:,:), intent(in) :: field
103
104logical :: is_active ! For XIOS, should this field be sent or not
105
106if (size(field,2) == nsoilmx) then
107    call writediagsoil(ngrid,field_name,title,units,3,field)
108else
109    call writediagfi(ngrid,field_name,title,units,3,field)
110endif
111
112#ifdef CPP_XIOS
113!is_active=xios_is_active_field(field_name)
114is_active=.true.
115    ! only send the field to xios if the user asked for it
116    if (is_active) call send_xios_field(field_name,field)
117#endif
118
119END SUBROUTINE write_output_d2
120
121!-----------------------------------------------------------------------
122
123SUBROUTINE write_output_i0(field_name,title,units,field)
124! For a scalar
125
126#ifdef CPP_XIOS
127    ! use xios_output_mod, only: xios_is_active_field
128    use xios_output_mod, only: send_xios_field
129#endif
130
131implicit none
132
133include "dimensions.h"
134
135integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
136character(*), intent(in) :: field_name
137character(*), intent(in) :: title
138character(*), intent(in) :: units
139integer,      intent(in) :: field
140
141logical :: is_active ! For XIOS, should this field be sent or not
142
143call writediagfi(ngrid,field_name,title,units,0,(/real(field)/))
144#ifdef CPP_XIOS
145!is_active=xios_is_active_field(field_name)
146is_active=.true.
147    ! only send the field to xios if the user asked for it
148    if (is_active) call send_xios_field(field_name,(/real(field)/))
149#endif
150
151END SUBROUTINE write_output_i0
152
153!-----------------------------------------------------------------------
154
155SUBROUTINE write_output_i1(field_name,title,units,field)
156! For a surface field
157
158#ifdef CPP_XIOS
159    ! use xios_output_mod, only: xios_is_active_field
160    use xios_output_mod, only: send_xios_field
161#endif
162
163implicit none
164
165include "dimensions.h"
166
167integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
168character(*),          intent(in) :: field_name
169character(*),          intent(in) :: title
170character(*),          intent(in) :: units
171integer, dimension(:), intent(in) :: field
172
173logical :: is_active ! For XIOS, should this field be sent or not
174
175call writediagfi(ngrid,field_name,title,units,2,real(field))
176#ifdef CPP_XIOS
177!is_active=xios_is_active_field(field_name)
178is_active=.true.
179    ! only send the field to xios if the user asked for it
180    if (is_active) call send_xios_field(field_name,real(field))
181#endif
182
183END SUBROUTINE write_output_i1
184
185!-----------------------------------------------------------------------
186
187SUBROUTINE write_output_i2(field_name,title,units,field)
188! For a "3D" horizontal-vertical field
189
190#ifdef CPP_XIOS
191    ! use xios_output_mod, only: xios_is_active_field
192    use xios_output_mod, only: send_xios_field
193#endif
194
195use comsoil_h,         only: nsoilmx
196use writediagsoil_mod, only: writediagsoil
197
198implicit none
199
200include "dimensions.h"
201
202integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
203character(*),            intent(in) :: field_name
204character(*),            intent(in) :: title
205character(*),            intent(in) :: units
206integer, dimension(:,:), intent(in) :: field
207
208logical :: is_active ! For XIOS, should this field be sent or not
209
210if (size(field,2) == nsoilmx) then
211    call writediagsoil(ngrid,field_name,title,units,3,real(field))
212else
213    call writediagfi(ngrid,field_name,title,units,3,real(field))
214endif
215#ifdef CPP_XIOS
216!is_active=xios_is_active_field(field_name)
217is_active=.true.
218    ! only send the field to xios if the user asked for it
219    if (is_active) call send_xios_field(field_name,real(field))
220#endif
221
222END SUBROUTINE write_output_i2
223
224!-----------------------------------------------------------------------
225
226SUBROUTINE write_output_l0(field_name,title,units,field)
227! For a scalar
228
229#ifdef CPP_XIOS
230    ! use xios_output_mod, only: xios_is_active_field
231    use xios_output_mod, only: send_xios_field
232#endif
233
234implicit none
235
236include "dimensions.h"
237
238integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
239character(*), intent(in) :: field_name
240character(*), intent(in) :: title
241character(*), intent(in) :: units
242logical,      intent(in) :: field
243! Local argument used to convert logical to real array
244real, dimension(1) :: field_real
245logical :: is_active ! For XIOS, should this field be sent or not
246
247field_real = 0.
248if (field) field_real = 1.
249
250call writediagfi(ngrid,field_name,title,units,0,field_real)
251#ifdef CPP_XIOS
252!is_active=xios_is_active_field(field_name)
253is_active=.true.
254    ! only send the field to xios if the user asked for it
255    if (is_active) call send_xios_field(field_name,field_real)
256#endif
257
258END SUBROUTINE write_output_l0
259
260!-----------------------------------------------------------------------
261
262SUBROUTINE write_output_l1(field_name,title,units,field)
263! For a surface field
264
265#ifdef CPP_XIOS
266    ! use xios_output_mod, only: xios_is_active_field
267    use xios_output_mod, only: send_xios_field
268#endif
269
270implicit none
271
272include "dimensions.h"
273
274integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
275character(*),          intent(in) :: field_name
276character(*),          intent(in) :: title
277character(*),          intent(in) :: units
278logical, dimension(:), intent(in) :: field
279! Local argument used to convert logical to real
280real, dimension(ngrid) :: field_real
281logical :: is_active ! For XIOS, should this field be sent or not
282
283field_real = 0.
284where (field) field_real = 1.
285
286call writediagfi(ngrid,field_name,title,units,2,field_real)
287#ifdef CPP_XIOS
288!is_active=xios_is_active_field(field_name)
289is_active=.true.
290    ! only send the field to xios if the user asked for it
291    if (is_active) call send_xios_field(field_name,field_real)
292#endif
293
294END SUBROUTINE write_output_l1
295
296!-----------------------------------------------------------------------
297
298SUBROUTINE write_output_l2(field_name,title,units,field)
299! For a "3D" horizontal-vertical field
300
301#ifdef CPP_XIOS
302    ! use xios_output_mod, only: xios_is_active_field
303    use xios_output_mod, only: send_xios_field
304#endif
305
306use comsoil_h,         only: nsoilmx
307use writediagsoil_mod, only: writediagsoil
308
309implicit none
310
311include "dimensions.h"
312
313integer, parameter :: ngrid = 2 + (jjm - 1)*iim - 1/jjm
314character(*),            intent(in) :: field_name
315character(*),            intent(in) :: title
316character(*),            intent(in) :: units
317logical, dimension(:,:), intent(in) :: field
318! Local argument used to convert logical to real
319real, allocatable, dimension(:,:) :: field_real
320logical :: is_active ! For XIOS, should this field be sent or not
321
322allocate(field_real(size(field,1),size(field,2)))
323field_real = 0.
324where (field) field_real = 1.
325
326if (size(field,2) == nsoilmx) then
327    call writediagsoil(ngrid,field_name,title,units,3,field_real)
328else
329    call writediagfi(ngrid,field_name,title,units,3,field_real)
330endif
331
332#ifdef CPP_XIOS
333! is_active=xios_is_active_field(field_name)
334is_active=.true.
335    ! only send the field to xios if the user asked for it
336    if (is_active) call send_xios_field(field_name,field_real)
337#endif
338
339deallocate(field_real)
340
341END SUBROUTINE write_output_l2
342
343END MODULE write_output_mod
Note: See TracBrowser for help on using the repository browser.