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

Last change on this file since 3118 was 3092, checked in by jbclement, 14 months ago

Marc PCM:

  • Correction of a bug in "writediagfi.F": the case of using the 1D model with parallelization was not anticipated so that the "diagfi.nc" file was filled with NaNf?;
  • Addition of the file "start1D.txt" as an example in the directory deftank/;
  • Some "cosmetic" modifications in "improvedclouds_mod.F", "write_output_mod.F90" and "testphys1d.F90".

JBC

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