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

Last change on this file since 3188 was 3169, checked in by emillour, 2 years 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
RevLine 
[2932]1MODULE write_output_mod
[3092]2
[3160]3implicit none
[3092]4
[3160]5private
[2970]6
[3092]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
[2932]12
[3160]13public write_output
[3092]14
[3160]15!-----------------------------------------------------------------------
16contains
17!-----------------------------------------------------------------------
[2932]18
[3160]19SUBROUTINE write_output_d0(field_name,title,units,field)
[3161]20! For a scalar
[3160]21
[2934]22#ifdef CPP_XIOS
[3160]23    use xios_output_mod, only: xios_is_active_field
24    use xios_output_mod, only: send_xios_field
[2934]25#endif
[3092]26
[3160]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/))
[3092]38#ifdef CPP_XIOS
[3055]39    ! only send the field to xios if the user asked for it
[3169]40    if (xios_is_active_field(field_name)) call send_xios_field(field_name,field)
[2932]41#endif
[3092]42
[3160]43END SUBROUTINE write_output_d0
[2932]44
[3160]45!-----------------------------------------------------------------------
[3092]46
[3160]47SUBROUTINE write_output_d1(field_name,title,units,field)
48! For a surface field
49
[2934]50#ifdef CPP_XIOS
[3161]51    use xios_output_mod, only: xios_is_active_field
52    use xios_output_mod, only: send_xios_field
[2934]53#endif
[3092]54
[3160]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)
[3092]66#ifdef CPP_XIOS
[3055]67    ! only send the field to xios if the user asked for it
[3160]68    if (xios_is_active_field(field_name)) call send_xios_field(field_name,field)
[2932]69#endif
[3092]70
[3160]71END SUBROUTINE write_output_d1
[2932]72
[3160]73!-----------------------------------------------------------------------
[3092]74
[3160]75SUBROUTINE write_output_d2(field_name,title,units,field)
76! For a "3D" horizontal-vertical field
77
[2934]78#ifdef CPP_XIOS
[3160]79    use xios_output_mod, only: xios_is_active_field
80    use xios_output_mod, only: send_xios_field
[2934]81#endif
[3161]82
[3160]83use comsoil_h,         only: nsoilmx
84use writediagsoil_mod, only: writediagsoil
[2932]85
[3160]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
[2932]97    call writediagsoil(ngrid,field_name,title,units,3,field)
[3160]98else
99    call writediagfi(ngrid,field_name,title,units,3,field)
100endif
101
[3092]102#ifdef CPP_XIOS
[3055]103    ! only send the field to xios if the user asked for it
[3160]104    if (xios_is_active_field(field_name)) call send_xios_field(field_name,field)
[2932]105#endif
[3092]106
[3160]107END SUBROUTINE write_output_d2
[2932]108
[3160]109!-----------------------------------------------------------------------
[3092]110
[3160]111SUBROUTINE write_output_i0(field_name,title,units,field)
[3161]112! For a scalar
[3160]113
[2970]114#ifdef CPP_XIOS
[3160]115    use xios_output_mod, only: xios_is_active_field
116    use xios_output_mod, only: send_xios_field
[2970]117#endif
[3092]118
[3160]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)/))
[3092]130#ifdef CPP_XIOS
[3055]131    ! only send the field to xios if the user asked for it
[3160]132    if (xios_is_active_field(field_name)) call send_xios_field(field_name,(/real(field)/))
[2970]133#endif
[3092]134
[3160]135END SUBROUTINE write_output_i0
[2970]136
[3160]137!-----------------------------------------------------------------------
[3092]138
[3160]139SUBROUTINE write_output_i1(field_name,title,units,field)
140! For a surface field
141
[2970]142#ifdef CPP_XIOS
[3160]143    use xios_output_mod, only: xios_is_active_field
144    use xios_output_mod, only: send_xios_field
[2970]145#endif
[3092]146
[3160]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))
[3092]158#ifdef CPP_XIOS
[3055]159    ! only send the field to xios if the user asked for it
[3160]160    if (xios_is_active_field(field_name)) call send_xios_field(field_name,real(field))
[2970]161#endif
[3092]162
[3160]163END SUBROUTINE write_output_i1
[2970]164
[3160]165!-----------------------------------------------------------------------
[3092]166
[3160]167SUBROUTINE write_output_i2(field_name,title,units,field)
168! For a "3D" horizontal-vertical field
169
[2970]170#ifdef CPP_XIOS
[3160]171    use xios_output_mod, only: xios_is_active_field
172    use xios_output_mod, only: send_xios_field
[2970]173#endif
174
[3160]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
[2970]189    call writediagsoil(ngrid,field_name,title,units,3,real(field))
[3160]190else
191    call writediagfi(ngrid,field_name,title,units,3,real(field))
192endif
[3092]193#ifdef CPP_XIOS
[3055]194    ! only send the field to xios if the user asked for it
[3160]195    if (xios_is_active_field(field_name)) call send_xios_field(field_name,real(field))
[2970]196#endif
[3092]197
[3160]198END SUBROUTINE write_output_i2
[2970]199
[3160]200!-----------------------------------------------------------------------
[3092]201
[3160]202SUBROUTINE write_output_l0(field_name,title,units,field)
[3161]203! For a scalar
[3160]204 
[2970]205#ifdef CPP_XIOS
[3160]206    use xios_output_mod, only: xios_is_active_field
207    use xios_output_mod, only: send_xios_field
[2970]208#endif
209
[3160]210implicit none
[3092]211
[3160]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)
[3092]226#ifdef CPP_XIOS
[3055]227    ! only send the field to xios if the user asked for it
[3160]228    if (xios_is_active_field(field_name)) call send_xios_field(field_name,field_real)
[2970]229#endif
[3092]230
[3160]231END SUBROUTINE write_output_l0
[2970]232
[3160]233!-----------------------------------------------------------------------
[3092]234
[3160]235SUBROUTINE write_output_l1(field_name,title,units,field)
236! For a surface field
237
[2970]238#ifdef CPP_XIOS
[3160]239    use xios_output_mod, only: xios_is_active_field
240    use xios_output_mod, only: send_xios_field
[2970]241#endif
[3160]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
[2970]252! Local argument used to convert logical to real
[3160]253real, dimension(ngrid) :: field_real
[2970]254
[3160]255field_real = 0.
256where (field) field_real = 1.
[3092]257
[3160]258call writediagfi(ngrid,field_name,title,units,2,field_real)
[3092]259#ifdef CPP_XIOS
[3055]260    ! only send the field to xios if the user asked for it
[3160]261    if (xios_is_active_field(field_name)) call send_xios_field(field_name,field_real)
[2970]262#endif
[3092]263
[3160]264END SUBROUTINE write_output_l1
[2970]265
[3160]266!-----------------------------------------------------------------------
[3092]267
[3160]268SUBROUTINE write_output_l2(field_name,title,units,field)
269! For a "3D" horizontal-vertical field
270
[2970]271#ifdef CPP_XIOS
[3160]272    use xios_output_mod, only: xios_is_active_field
273    use xios_output_mod, only: send_xios_field
[2970]274#endif
[3160]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
[2970]288! Local argument used to convert logical to real
[3160]289real, allocatable, dimension(:,:) :: field_real
[2970]290
[3160]291allocate(field_real(size(field,1),size(field,2)))
292field_real = 0.
293where (field) field_real = 1.
[2970]294
[3160]295if (size(field,2) == nsoilmx) then
[2970]296    call writediagsoil(ngrid,field_name,title,units,3,field_real)
[3160]297else
298    call writediagfi(ngrid,field_name,title,units,3,field_real)
299endif
[2970]300
[3092]301#ifdef CPP_XIOS
[3055]302    ! only send the field to xios if the user asked for it
[3160]303    if (xios_is_active_field(field_name)) call send_xios_field(field_name,field_real)
[2970]304#endif
[3092]305
[3160]306deallocate(field_real)
[2976]307
[3160]308END SUBROUTINE write_output_l2
[2970]309
[2932]310END MODULE write_output_mod
Note: See TracBrowser for help on using the repository browser.