source: LMDZ6/trunk/libf/phylmd/ecrad/ifs/easy_netcdf_read_mpi.F90 @ 5229

Last change on this file since 5229 was 4773, checked in by idelkadi, 11 months ago
  • Update of Ecrad in LMDZ The same organization of the Ecrad offline version is retained in order to facilitate the updating of Ecrad in LMDZ and the comparison between online and offline results. version 1.6.1 of Ecrad (https://github.com/lguez/ecrad.git)
  • Implementation of the double call of Ecrad in LMDZ


File size: 9.1 KB
Line 
1! easy_netcdf_read_mpi.f90 - Read netcdf file on one task and share with other tasks
2!
3! (C) Copyright 2017- ECMWF.
4!
5! This software is licensed under the terms of the Apache Licence Version 2.0
6! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7!
8! In applying this licence, ECMWF does not waive the privileges and immunities
9! granted to it by virtue of its status as an intergovernmental organisation
10! nor does it submit to any jurisdiction.
11!
12! Author:  Robin Hogan
13! Email:   r.j.hogan@ecmwf.int
14!
15
16module easy_netcdf_read_mpi
17
18  use easy_netcdf,   only : netcdf_file_raw => netcdf_file
19  use parkind1,      only : jpim, jprb
20  use radiation_io,  only : nulout, nulerr, my_abort => radiation_abort
21
22  implicit none
23
24  ! MPI tag for radiation and physics communication
25  integer(kind=jpim), parameter :: mtagrad = 2800
26
27  !---------------------------------------------------------------------
28  ! An object of this type provides convenient read or write access to
29  ! a NetCDF file
30  type netcdf_file
31    type(netcdf_file_raw) :: file
32    logical :: is_master_task = .true.
33  contains
34    procedure :: open => open_netcdf_file
35    procedure :: close => close_netcdf_file
36    procedure :: get_real_scalar
37    procedure :: get_real_vector
38    procedure :: get_real_matrix
39    procedure :: get_real_array3
40    generic   :: get => get_real_scalar, get_real_vector, &
41         &              get_real_matrix, get_real_array3
42    procedure :: get_global_attribute
43
44    procedure :: set_verbose
45    procedure :: transpose_matrices
46    procedure :: exists
47  end type netcdf_file
48
49contains
50
51  ! --- GENERIC SUBROUTINES ---
52
53  !---------------------------------------------------------------------
54  ! Open a NetCDF file with name "file_name", optionally specifying the
55  ! verbosity level (0-5)
56  subroutine open_netcdf_file(this, file_name, iverbose)
57
58    USE MPL_MODULE, ONLY : MPL_BROADCAST, MPL_MYRANK
59
60    class(netcdf_file)            :: this
61    character(len=*), intent(in)  :: file_name
62    integer, intent(in), optional :: iverbose
63
64    integer                       :: istatus
65
66    ! Store verbosity level in object
67    if (present(iverbose)) then
68      this%file%iverbose = iverbose
69    else
70      ! By default announce files being opened and closed, but not
71      ! variables read/written
72      this%file%iverbose = 2
73    end if
74
75    ! By default we don't transpose 2D arrays on read
76    this%file%do_transpose_2d = .false.
77
78    if (MPL_MYRANK() == 1) then
79      this%is_master_task = .true.
80      call this%file%open(file_name, iverbose)
81    else
82      this%is_master_task = .false.
83    end if
84
85  end subroutine open_netcdf_file
86
87
88  !---------------------------------------------------------------------
89  ! Close the NetCDF file
90  subroutine close_netcdf_file(this)
91    class(netcdf_file) :: this
92    integer            :: istatus
93
94    if (this%is_master_task) then
95      call this%file%close()
96    end if
97
98  end subroutine close_netcdf_file
99
100
101  !---------------------------------------------------------------------
102  ! Set the verbosity level from 0 to 5, where the codes have the
103  ! following meaning: 0=errors only, 1=warning, 2=info, 3=progress,
104  ! 4=detailed, 5=debug
105  subroutine set_verbose(this, ival)
106    class(netcdf_file) :: this
107    integer, optional  :: ival
108
109    if (present(ival)) then
110      this%file%iverbose = ival
111    else
112      this%file%iverbose = 2
113    end if
114
115  end subroutine set_verbose
116
117
118
119  !---------------------------------------------------------------------
120  ! Specify whether 2D arrays should be transposed on read
121  subroutine transpose_matrices(this, do_transpose)
122    class(netcdf_file) :: this
123    logical, optional  :: do_transpose
124
125    if (present(do_transpose)) then
126      this%file%do_transpose_2d = do_transpose
127    else
128      this%file%do_transpose_2d = .true.
129    end if
130
131  end subroutine transpose_matrices
132
133
134
135  ! --- READING SUBROUTINES ---
136
137  !---------------------------------------------------------------------
138  ! Return true if the variable is present, false otherwise
139  function exists(this, var_name) result(is_present)
140
141    USE MPL_MODULE, ONLY : MPL_BROADCAST, MPL_NPROC
142
143    class(netcdf_file)           :: this
144    character(len=*), intent(in) :: var_name
145
146    logical :: is_present
147
148    if (this%is_master_task) then
149      is_present = this%file%exists(var_name)
150    end if
151
152    if (MPL_NPROC() > 1) then
153      CALL MPL_BROADCAST(is_present, mtagrad, 1, &
154           &  CDSTRING='EASY_NETCDF_READ_MPI:EXISTS')
155    end if
156
157  end function exists
158
159
160  !---------------------------------------------------------------------
161  ! The method "get" will read either a scalar, vector or matrix
162  ! depending on the rank of the output argument. This version reads a
163  ! scalar.
164  subroutine get_real_scalar(this, var_name, scalar)
165
166    USE MPL_MODULE, ONLY : MPL_BROADCAST, MPL_NPROC
167
168    class(netcdf_file)           :: this
169    character(len=*), intent(in) :: var_name
170    real(jprb), intent(out)      :: scalar
171
172    if (this%is_master_task) then
173      call this%file%get(var_name, scalar)
174    end if
175
176    if (MPL_NPROC() > 1) then
177      CALL MPL_BROADCAST(scalar, mtagrad, 1, &
178           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_REAL_SCALAR')
179    end if
180
181  end subroutine get_real_scalar
182
183
184  !---------------------------------------------------------------------
185  ! Read a 1D array into "vector", which must be allocatable and will
186  ! be reallocated if necessary
187  subroutine get_real_vector(this, var_name, vector)
188
189    USE MPL_MODULE, ONLY : MPL_BROADCAST, MPL_NPROC
190
191    class(netcdf_file)           :: this
192    character(len=*), intent(in) :: var_name
193    real(jprb), allocatable, intent(out) :: vector(:)
194
195    integer                      :: n  ! Length of vector
196
197    n = 0
198
199    if (this%is_master_task) then
200      call this%file%get(var_name, vector)
201      n = size(vector)
202    end if
203
204    if (MPL_NPROC() > 1) then
205      CALL MPL_BROADCAST(n, mtagrad, 1, &
206           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_REAL_VECTOR:SIZE')
207
208      if (.not. this%is_master_task) then
209        allocate(vector(n))
210      end if
211
212      CALL MPL_BROADCAST(vector, mtagrad, 1, &
213           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_REAL_VECTOR')
214    end if
215
216  end subroutine get_real_vector
217
218
219  !---------------------------------------------------------------------
220  ! Read 2D array into "matrix", which must be allocatable and will be
221  ! reallocated if necessary.  Whether to transpose is specifed by the
222  ! final optional argument, but can also be specified by the
223  ! do_transpose_2d class data member.
224  subroutine get_real_matrix(this, var_name, matrix, do_transp)
225
226    USE MPL_MODULE, ONLY : MPL_BROADCAST, MPL_NPROC
227
228    class(netcdf_file)           :: this
229    character(len=*), intent(in) :: var_name
230    real(jprb), allocatable, intent(out) :: matrix(:,:)
231    logical, optional, intent(in):: do_transp ! Transpose data?
232
233    integer                      :: n(2)
234
235    n = 0
236
237    if (this%is_master_task) then
238      call this%file%get(var_name, matrix, do_transp)
239      n = shape(matrix)
240    end if
241
242    if (MPL_NPROC() > 1) then
243      CALL MPL_BROADCAST(n, mtagrad, 1, &
244           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_REAL_MATRIX:SIZE')
245
246      if (.not. this%is_master_task) then
247        allocate(matrix(n(1),n(2)))
248      end if
249
250      CALL MPL_BROADCAST(matrix, mtagrad, 1, &
251           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_REAL_MATRIX')
252    end if
253
254  end subroutine get_real_matrix
255
256
257  !---------------------------------------------------------------------
258  ! Read 3D array into "var", which must be allocatable and will be
259  ! reallocated if necessary.  Whether to pemute is specifed by the
260  ! final optional argument
261  subroutine get_real_array3(this, var_name, var, ipermute)
262
263    USE MPL_MODULE, ONLY : MPL_BROADCAST, MPL_NPROC
264
265    class(netcdf_file)                   :: this
266    character(len=*), intent(in)         :: var_name
267    real(jprb), allocatable, intent(out) :: var(:,:,:)
268    integer, optional, intent(in)        :: ipermute(3)
269
270    integer                              :: n(3)
271
272    n = 0
273
274    if (this%is_master_task) then
275      call this%file%get(var_name, var, ipermute)
276      n = shape(var)
277    end if
278
279    if (MPL_NPROC() > 1) then
280      CALL MPL_BROADCAST(n, mtagrad, 1, &
281           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_REAL_ARRAY3:SIZE')
282
283      if (.not. this%is_master_task) then
284        allocate(var(n(1),n(2),n(3)))
285      end if
286
287      CALL MPL_BROADCAST(var, mtagrad, 1, &
288           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_REAL_ARRAY3')
289    end if
290
291  end subroutine get_real_array3
292
293
294  !---------------------------------------------------------------------
295  ! Get a global attribute as a character string
296  subroutine get_global_attribute(this, attr_name, attr_str)
297
298    USE MPL_MODULE, ONLY : MPL_BROADCAST, MPL_NPROC
299
300    class(netcdf_file) :: this
301
302    character(len=*), intent(in)    :: attr_name
303    character(len=*), intent(inout) :: attr_str
304
305    if (this%is_master_task) then
306      call this%file%get_global_attribute(attr_name, attr_str)
307    end if
308
309    if (MPL_NPROC() > 1) then
310      CALL MPL_BROADCAST(attr_str, mtagrad, 1, &
311           &  CDSTRING='EASY_NETCDF_READ_MPI:GET_GLOBAL_ATTRIBUTE')
312    end if
313
314  end subroutine get_global_attribute
315
316end module easy_netcdf_read_mpi
Note: See TracBrowser for help on using the repository browser.