source: trunk/LMDZ.MARS/libf/phymars/write_output_mod.F90 @ 3183

Last change on this file since 3183 was 3169, checked in by emillour, 18 months ago

Mars PCM:
Minor fix for XIOS output of a scalar in write_output_mod.F90.
Also added more helpful error messages in xios_output_mod.F90.
EM

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