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

Last change on this file was 4058, checked in by emillour, 5 days ago

Mars PCM:
Add a "output_diagfi" (default .true.) flag so that the user can decide if
a diagfi.nc file should be outputted.
Also do some cleaning around the few remaining calls to writediagfi
or reference to it througout the code; write_output() should be used.
EM

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