FMS  2024.03
Flexible Modeling System
mpp_io_unstructured_write.inc
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
20 !> @file
21 !> @brief Parallel file writes for unstructured grids, used in @ref mpp_io_mod
22 
23 !> @addtogroup mpp_io_mod
24 !> @{
25 
26 !>Write data for a 1D field associated with an unstructured mpp domain to a
27 !!restart file.
29  field, &
30  domain, &
31  fdata, &
32  nelems_io, &
33  tstamp, &
34  default_data)
35 
36  !Inputs/outputs
37  integer(i4_kind),intent(in) :: funit !<A file unit for the to which the
38  !! data will be written
39  type(fieldtype),intent(inout) :: field !<A field whose data will be written
40  type(domainug),intent(inout) :: domain !<An unstructured mpp domain associatd
41  !! with the inputted file
42  real(KIND=r8_kind),dimension(:),intent(inout) :: fdata !<The data that will be written to the file
43  integer,dimension(:),intent(in) :: nelems_io !<Number of grid points in the compressed
44  !! dimension for each rank (correct
45  !!sizes only exist for the
46  !!root rank of I/O domain pelist)
47  real(KIND=r8_kind),intent(in),optional :: tstamp !<A time value
48  real(KIND=r8_kind),intent(in), optional :: default_data !<Fill value for the inputted field
49 
50  !Local variables
51  real(KIND=r8_kind) :: fill !<Fill value for the inputted field (default: zero)
52  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
53  integer(i4_kind) :: io_domain_npes !<The total number of ranks in an I/O domain pelist
54  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
55  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
56  !! to be written
57  !!by the root rank of the pelist
58  real(kind=r8_kind),dimension(:),allocatable :: rbuff !<Buffer used to gather the data onto
59  !! the root rank of the pelist
60  real(kind=r8_kind),dimension(:),allocatable :: cdata !<Array used to write the data to the
61  !! file after the gather is performed
62  integer(i4_kind) :: i !<Loop variable
63 
64  !Start the mpp timer.
65  !mpp_write_clock is a module variable.
66  call mpp_clock_begin(mpp_write_clock)
67 
68  !Make sure that the module is initialized.
69  if (.not. module_is_initialized) then
70  call mpp_error(fatal, &
71  "mpp_io_unstructured_write_r_1D:" &
72  //" you must first call mpp_io_init.")
73  endif
74 
75  !Make sure that another NetCDF file is not currently using the inputted
76  !file unit.
77  if (.not. mpp_file(funit)%valid) then
78  call mpp_error(fatal, &
79  "mpp_io_unstructured_write_r_1D:" &
80  //" the inputted file unit is already in use.")
81  endif
82 
83  !Set the fill value for the field.
84  fill = 0.0
85  if (present(default_data)) then
86  fill = default_data
87  endif
88 
89  !Point to the I/O domain associated with the inputted unstructured mpp
90  !domain.
91  io_domain => null()
92  io_domain => mpp_get_ug_io_domain(domain)
93 
94  !Get the pelist associated with the I/O domain.
95  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
96  allocate(pelist(io_domain_npes))
97  call mpp_get_ug_domain_pelist(io_domain, &
98  pelist)
99 
100  !Make sure that only the root rank of the pelist will write to the file.
101  !This check is needed because data is only gathered on the lowest rank
102  !of the pelist.
103  if (mpp_pe() .eq. pelist(1) .and. .not. &
104  mpp_file(funit)%write_on_this_pe) then
105  call mpp_error(fatal, &
106  "mpp_io_unstructured_write_r_1D:" &
107  //" the root rank of the pelist must be allowed" &
108  //" to perform the write.")
109  endif
110  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
111  call mpp_error(fatal, &
112  "mpp_io_unstructured_write_r_1D:" &
113  //" the non-root ranks of the pelist are not" &
114  //" allowed to perform the write.")
115  endif
116 
117  !Allocate an array which will be used to gather the data to be written
118  !onto the root rank of the pelist.
119  if (mpp_pe() .eq. pelist(1)) then
120  nelems = sum(nelems_io)
121  allocate(rbuff(nelems))
122  else
123  allocate(rbuff(1))
124  endif
125 
126  !Perform the gather of data onto the root rank (pelist(1)).
127  call mpp_gather(fdata, &
128  size(fdata), &
129  rbuff, &
130  nelems_io, &
131  pelist)
132 
133  !Write out the data to the file. This is only done by the root rank
134  !of the pelist.
135  if (mpp_pe() .eq. pelist(1)) then
136  allocate(cdata(nelems))
137  cdata = fill
138  do i = 1,nelems
139  cdata(i) = rbuff(i)
140  enddo
141  field%size(1) = nelems
142  call write_record_r8(funit, &
143  field, &
144  nelems, &
145  cdata, &
146  time_in=tstamp)
147  deallocate(cdata)
148  endif
149 
150  !Deallocate local allocatables.
151  deallocate(rbuff)
152  deallocate(pelist)
153 
154  !Stop the mpp timer.
155  call mpp_clock_end(mpp_write_clock)
156 
157  return
158 end subroutine mpp_io_unstructured_write_r8_1d
159 
160 !------------------------------------------------------------------------------
161 !>Write data for a 2D field associated with an unstructured mpp domain to a
162 !!restart file.
164  field, &
165  domain, &
166  fdata, &
167  nelems_io, &
168  tstamp, &
169  default_data)
170 
171  !Inputs/outputs
172  integer(i4_kind),intent(in) :: funit !<A file unit for the to which the
173  !! data will be written
174  type(fieldtype),intent(inout) :: field !<A field whose data will be written
175  type(domainug),intent(inout) :: domain !<An unstructured mpp domain associatd
176  !! with the inputted file
177  real(KIND=r8_kind),dimension(:,:),intent(inout) :: fdata !<The data that will be written to the file
178  integer,dimension(:),intent(in) :: nelems_io !<Number of grid points in the compressed
179  !! dimension for each rank
180  !!(correct sizes only exist
181  !!for the root rank of I/O domain pelist)
182  real(KIND=r8_kind),intent(in),optional :: tstamp !<A time value
183  real(KIND=r8_kind),intent(in), optional :: default_data !<Fill value for the inputted field
184 
185  !Local variables
186  real(KIND=r8_kind) :: fill !<Fill value for the inputted field (default: zero)
187  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
188  integer(i4_kind) :: io_domain_npes !<The total number of ranks in an I/O domain pelist
189  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
190  integer(i4_kind) :: dim_size_1 !<Number of data points in the first
191  !! dimension (size(fdata,1))
192  integer(i4_kind) :: dim_size_2 !<Number of data points in the second
193  !! dimension (size(fdata,2))
194  real(kind=r8_kind),dimension(:),allocatable :: sbuff !<Buffer used to gather the data
195  !! onto the root rank of the pelist
196  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
197  !! to be written
198  !!by the root rank of the pelist
199  real(kind=r8_kind),dimension(:),allocatable :: rbuff !<Buffer used to gather the data
200  !! onto the root rank of the pelist
201  real(kind=r8_kind),dimension(:,:),allocatable :: cdata !<Array used to write the data to
202  !! the file after the gather is performed
203  integer(i4_kind) :: offset_r !<Offset for rbuff used to reorder
204  !! data before netCDF write
205  integer(i4_kind) :: offset_c !<Offset for cdata used to reorder
206  !! data before netCDF write
207  integer(i4_kind) :: i !<Loop variable
208  integer(i4_kind) :: j !<Loop variable
209  integer(i4_kind) :: k !<Loop variable
210 
211  !Start the mpp timer.
212  !mpp_write_clock is a module variable.
213  call mpp_clock_begin(mpp_write_clock)
214 
215  !Make sure that the module is initialized.
216  if (.not. module_is_initialized) then
217  call mpp_error(fatal, &
218  "mpp_io_unstructured_write_r_2D:" &
219  //" you must first call mpp_io_init.")
220  endif
221 
222  !Make sure that another NetCDF file is not currently using the inputted
223  !file unit.
224  if (.not. mpp_file(funit)%valid) then
225  call mpp_error(fatal, &
226  "mpp_io_unstructured_write_r_2D:" &
227  //" the inputted file unit is already in use.")
228  endif
229 
230  !Set the fill value for the field.
231  fill = 0.0
232  if (present(default_data)) then
233  fill = default_data
234  endif
235 
236  !Point to the I/O domain associated with the inputted unstructured mpp
237  !domain.
238  io_domain => null()
239  io_domain => mpp_get_ug_io_domain(domain)
240 
241  !Get the pelist associated with the I/O domain.
242  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
243  allocate(pelist(io_domain_npes))
244  call mpp_get_ug_domain_pelist(io_domain, &
245  pelist)
246 
247  !Make sure that only the root rank of the pelist will write to the file.
248  !This check is needed because data is only gathered on the lowest rank
249  !of the pelist.
250  if (mpp_pe() .eq. pelist(1) .and. .not. &
251  mpp_file(funit)%write_on_this_pe) then
252  call mpp_error(fatal, &
253  "mpp_io_unstructured_write_r_2D:" &
254  //" the root rank of the pelist must be allowed" &
255  //" to perform the write.")
256  endif
257  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
258  call mpp_error(fatal, &
259  "mpp_io_unstructured_write_r_2D:" &
260  //" the non-root ranks of the pelist are not" &
261  //" allowed to perform the write.")
262  endif
263 
264  !Load the data elements for each rank into a one dimensional array, which
265  !will be used to gather the data onto the root rank of the pelist.
266  allocate(sbuff(size(fdata)))
267  dim_size_1 = size(fdata,1)
268  dim_size_2 = size(fdata,2)
269  do j = 1,dim_size_2
270  do i = 1,dim_size_1
271  sbuff((j-1)*dim_size_1+i) = fdata(i,j)
272  enddo
273  enddo
274 
275  !Allocate an array which will be used to gather the data to be written
276  !onto the root rank of the pelist.
277  if (mpp_pe() .eq. pelist(1)) then
278  nelems = sum(nelems_io)
279  allocate(rbuff(nelems*dim_size_2))
280  else
281  allocate(rbuff(1))
282  endif
283 
284  !Perform the gather of data onto the root rank (pelist(1)).
285  call mpp_gather(sbuff, &
286  size(sbuff), &
287  rbuff, &
288  nelems_io*dim_size_2, &
289  pelist)
290 
291  !Reorder the gather data so that is of the form (nelems,dim_size_2). Write
292  !out the data to the file. This is only done by the root rank of the
293  !pelist.
294  if (mpp_pe() .eq. pelist(1)) then
295  allocate(cdata(nelems,dim_size_2))
296  cdata = fill
297  do j = 1,dim_size_2
298  offset_c = 0
299  do k = 1,io_domain_npes
300  if (k .gt. 1) then
301  offset_r = (j-1)*nelems_io(k) + dim_size_2*(sum(nelems_io(1:k-1)))
302  else
303  offset_r = (j-1)*nelems_io(k)
304  endif
305  do i = 1,nelems_io(k)
306  cdata(i+offset_c,j) = rbuff(i+offset_r)
307  enddo
308  offset_c = offset_c + nelems_io(k)
309  enddo
310  enddo
311  field%size(1) = nelems
312  call write_record_r8(funit, &
313  field, &
314  nelems*dim_size_2, &
315  cdata, &
316  time_in=tstamp)
317  deallocate(cdata)
318  endif
319 
320  !Deallocate local allocatables.
321  deallocate(sbuff)
322  deallocate(rbuff)
323  deallocate(pelist)
324 
325  !Stop the mpp timer.
326  call mpp_clock_end(mpp_write_clock)
327 
328  return
329 end subroutine mpp_io_unstructured_write_r8_2d
330 
331 !------------------------------------------------------------------------------
332 !>Write data for a 3D field associated with an unstructured mpp domain to a
333 !!restart file.
335  field, &
336  domain, &
337  fdata, &
338  nelems_io, &
339  tstamp, &
340  default_data)
341 
342  !Inputs/outputs
343  integer(i4_kind),intent(in) :: funit !<A file unit for the to which
344  !! the data will be written
345  type(fieldtype),intent(inout) :: field !<A field whose data will be written
346  type(domainug),intent(inout) :: domain !<An unstructured mpp domain associatd
347  !! with the inputted file
348  real(KIND=r8_kind),dimension(:,:,:),intent(inout) :: fdata !<The data that will be written to the file
349  integer,dimension(:),intent(in) :: nelems_io !<Number of grid points in the
350  !! compressed dimension for each rank
351  !!(correct sizes only
352  !!exist for the root rank of I/O domain pelist)
353  real(KIND=r8_kind),intent(in),optional :: tstamp !<A time value
354  real(KIND=r8_kind),intent(in), optional :: default_data !<Fill value for the inputted field
355 
356  !Local variables
357  real(KIND=r8_kind) :: fill !<Fill value for the inputted field
358  !! (defaults: zero)
359  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
360  integer(i4_kind) :: io_domain_npes !<The total number of ranks in
361  !! an I/O domain pelist
362  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
363  integer(i4_kind) :: dim_size_1 !<Number of data points in the
364  !! first dimension (size(fdata,1))
365  integer(i4_kind) :: dim_size_2 !<Number of data points in the
366  !! second dimension (size(fdata,2))
367  integer(i4_kind) :: dim_size_3 !<Number of data points in the
368  !! second dimension (size(fdata,3))
369  real(kind=r8_kind),dimension(:),allocatable :: sbuff !<Buffer used to gather the data
370  !! onto the root rank of the pelist
371  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
372  !! to be written
373  !!by the root rank of the pelist
374  real(kind=r8_kind),dimension(:),allocatable :: rbuff !<Buffer used to gather the data
375  !! onto the root rank of the pelist
376  real(kind=r8_kind),dimension(:,:,:),allocatable :: cdata !<Array used to write the data
377  !! to the file after the gather is performed
378  integer(i4_kind) :: offset_r !<Offset for rbuff used to reorder
379  !! data before netCDF write
380  integer(i4_kind) :: offset_c !<Offset for cdata used to reorder
381  !! data before netCDF write
382  integer(i4_kind) :: i !<Loop variable
383  integer(i4_kind) :: j !<Loop variable
384  integer(i4_kind) :: k !<Loop variable
385  integer(i4_kind) :: m !<Loop variable
386 
387  !Start the mpp timer.
388  !mpp_write_clock is a module variable.
389  call mpp_clock_begin(mpp_write_clock)
390 
391  !Make sure that the module is initialized.
392  if (.not. module_is_initialized) then
393  call mpp_error(fatal, &
394  "mpp_io_unstructured_write_r_3D:" &
395  //" you must first call mpp_io_init.")
396  endif
397 
398  !Make sure that another NetCDF file is not currently using the inputted
399  !file unit.
400  if (.not. mpp_file(funit)%valid) then
401  call mpp_error(fatal, &
402  "mpp_io_unstructured_write_r_3D:" &
403  //" the inputted file unit is already in use.")
404  endif
405 
406  !Set the fill value for the field.
407  fill = 0.0
408  if (present(default_data)) then
409  fill = default_data
410  endif
411 
412  !Point to the I/O domain associated with the inputted unstructured mpp
413  !domain.
414  io_domain => null()
415  io_domain => mpp_get_ug_io_domain(domain)
416 
417  !Get the pelist associated with the I/O domain.
418  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
419  allocate(pelist(io_domain_npes))
420  call mpp_get_ug_domain_pelist(io_domain, &
421  pelist)
422 
423  !Make sure that only the root rank of the pelist will write to the file.
424  !This check is needed because data is only gathered on the lowest rank
425  !of the pelist.
426  if (mpp_pe() .eq. pelist(1) .and. .not. &
427  mpp_file(funit)%write_on_this_pe) then
428  call mpp_error(fatal, &
429  "mpp_io_unstructured_write_r_3D:" &
430  //" the root rank of the pelist must be allowed" &
431  //" to perform the write.")
432  endif
433  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
434  call mpp_error(fatal, &
435  "mpp_io_unstructured_write_r_3D:" &
436  //" the non-root ranks of the pelist are not" &
437  //" allowed to perform the write.")
438  endif
439 
440  !Load the data elements for each rank into a one dimensional array, which
441  !will be used to gather the data onto the root rank of the pelist.
442  allocate(sbuff(size(fdata)))
443  dim_size_1 = size(fdata,1)
444  dim_size_2 = size(fdata,2)
445  dim_size_3 = size(fdata,3)
446  do k = 1,dim_size_3
447  do j = 1,dim_size_2
448  do i = 1,dim_size_1
449  sbuff((k-1)*dim_size_2*dim_size_1+(j-1)*dim_size_1+i) = fdata(i,j,k)
450  enddo
451  enddo
452  enddo
453 
454  !Allocate an array which will be used to gather the data to be written
455  !onto the root rank of the pelist.
456  if (mpp_pe() .eq. pelist(1)) then
457  nelems = sum(nelems_io)
458  allocate(rbuff(nelems*dim_size_2*dim_size_3))
459  else
460  allocate(rbuff(1))
461  endif
462 
463  !Perform the gather of data onto the root rank (pelist(1)).
464  call mpp_gather(sbuff, &
465  size(sbuff), &
466  rbuff, &
467  nelems_io*dim_size_2*dim_size_3, &
468  pelist)
469 
470  !Reorder the gather data so that is of the form (nelems,dim_size_2). Write
471  !out the data to the file. This is only done by the root rank of the
472  !pelist.
473  if (mpp_pe() .eq. pelist(1)) then
474  allocate(cdata(nelems,dim_size_2,dim_size_3))
475  cdata = fill
476  do m = 1,dim_size_3
477  do j = 1,dim_size_2
478  offset_c = 0
479  do k = 1,io_domain_npes
480  if (k .gt. 1) then
481  offset_r = (m-1)*dim_size_2*nelems_io(k) + &
482  (j-1)*nelems_io(k) + &
483  dim_size_2*dim_size_3*(sum(nelems_io(1:k-1)))
484  else
485  offset_r = (m-1)*dim_size_2*nelems_io(k) + &
486  (j-1)*nelems_io(k)
487  endif
488  do i = 1,nelems_io(k)
489  cdata(i+offset_c,j,m) = rbuff(i+offset_r)
490  enddo
491  offset_c = offset_c + nelems_io(k)
492  enddo
493  enddo
494  enddo
495  field%size(1) = nelems
496  call write_record_r8(funit, &
497  field, &
498  nelems*dim_size_2*dim_size_3, &
499  cdata, &
500  time_in=tstamp)
501  deallocate(cdata)
502  endif
503 
504  !Deallocate local allocatables.
505  deallocate(sbuff)
506  deallocate(rbuff)
507  deallocate(pelist)
508 
509  !Stop the mpp timer.
510  call mpp_clock_end(mpp_write_clock)
511 
512  return
513 end subroutine mpp_io_unstructured_write_r8_3d
514 
515 !------------------------------------------------------------------------------
516 !>Write data for a 4D field associated with an unstructured mpp domain to a
517 !!restart file.
519  field, &
520  domain, &
521  fdata, &
522  nelems_io_in, &
523  tstamp, &
524  default_data)
525 
526  !Inputs/outputs
527  integer(i4_kind),intent(in) :: funit !<A file unit for the to which
528  !! the data will be written
529  type(fieldtype),intent(inout) :: field !<A field whose data will be written
530  type(domainug),intent(inout) :: domain !<An unstructured mpp domain
531  !! associatd with the inputted file
532  real(KIND=r8_kind),dimension(:,:,:,:),intent(inout) :: fdata !<The data that will be written to the file
533  integer,dimension(:),intent(in),optional :: nelems_io_in !<Number of grid points in the
534  !! unstructured dimension for each rank
535  !!(correct sizes only
536  !!(cexist for the root
537  !!rank of I/O domain pelist
538  real(KIND=r8_kind),intent(in),optional :: tstamp !<A time value
539  real(KIND=r8_kind),intent(in), optional :: default_data !<Fill value for the inputted field
540 
541  !Local variables
542  real(KIND=r8_kind) :: fill !<Fill value for the inputted field (default: zero)
543  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
544  integer(i4_kind) :: io_domain_npes !<The total number of ranks in an I/O domain pelist
545  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
546  integer(i4_kind),dimension(:),allocatable :: nelems_io !<Number of grid points in the unstructured
547  !! dimension for each rank
548  integer(i4_kind) :: compute_size !<Size of the unstructured compute
549  !! domain for the current rank
550  integer(i4_kind) :: size_fdata_dim_2 !<Number of data points in a non-unstructured
551  !! dimension (size(fdata,2))
552  integer(i4_kind) :: size_fdata_dim_3 !<Number of data points in a non-unstructured
553  !! dimension (size(fdata,3))
554  integer(i4_kind) :: size_fdata_dim_4 !<Number of data points in a non-unstructured
555  !! dimension (size(fdata,3))
556  integer(i4_kind) :: mynelems !<Number of data points in the unstructured
557  !! dimension (size(fdata,1))
558  real(kind=r8_kind),dimension(:),allocatable :: sbuff !<Buffer used to gather the data
559  !! onto the root rank of the pelist
560  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
561  !! to be written
562  !!by the root rank of the pelist
563  real(kind=r8_kind),dimension(:),allocatable :: rbuff !<Buffer to gather the data onto
564  !! root rank of pelist
565  real(kind=r8_kind),dimension(:,:,:,:),allocatable :: cdata !<Array to write the data to file
566  !! after gather is performed
567  integer(i4_kind) :: i !<Loop variable
568  integer(i4_kind) :: j !<Loop variable
569  integer(i4_kind) :: k !<Loop variable
570  integer(i4_kind) :: n !<Loop variable
571 
572  !Start the mpp timer.
573  !mpp_write_clock is a module variable.
574  call mpp_clock_begin(mpp_write_clock)
575 
576  !Make sure that the module is initialized.
577  if (.not. module_is_initialized) then
578  call mpp_error(fatal, &
579  "mpp_io_unstructured_write_compressed_r_4D:" &
580  //" you must first call mpp_io_init.")
581  endif
582 
583  !Make sure that another NetCDF file is not currently using the inputted
584  !file unit.
585  if (.not. mpp_file(funit)%valid) then
586  call mpp_error(fatal, &
587  "mpp_io_unstructured_write_compressed_r_4D:" &
588  //" the inputted file unit is already in use.")
589  endif
590 
591  !Set the fill value for the field.
592  fill = 0.0
593  if (present(default_data)) then
594  fill = default_data
595  endif
596 
597  !Point to the I/O domain associated with the inputted unstructured mpp
598  !domain.
599  io_domain => null()
600  io_domain => mpp_get_ug_io_domain(domain)
601 
602  !Get the pelist associated with the I/O domain.
603  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
604  allocate(pelist(io_domain_npes))
605  call mpp_get_ug_domain_pelist(io_domain, &
606  pelist)
607 
608  !Make sure that only the root rank of the pelist will write to the file.
609  !This check is needed because data is only gathered on the lowest rank
610  !of the pelist.
611  if (mpp_pe() .eq. pelist(1) .and. .not. &
612  mpp_file(funit)%write_on_this_pe) then
613  call mpp_error(fatal, &
614  "mpp_io_unstructured_write_compressed_r_4D:" &
615  //" the root rank of the pelist must be allowed" &
616  //" to perform the write.")
617  endif
618  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
619  call mpp_error(fatal, &
620  "mpp_io_unstructured_write_compressed_r_4D:" &
621  //" the non-root ranks of the pelist are not" &
622  //" allowed to perform the write.")
623  endif
624 
625  !For the 3D unstructured case, data is assumed to be of the form
626  !fdata = fdata(unstructured,z,cc). The number of data elements in the
627  !unstructured dimension (size(fdata,1)) may differ between ranks.
628  !If not passed in, the number of data elements in the unstructured
629  !dimension must be gathered on the root rank of the pelist. The number
630  !data elements in the unstructured dimension should be equal to the size
631  !of the unstructured computed domain.
632  if (present(nelems_io_in)) then
633  allocate(nelems_io(size(nelems_io_in)))
634  nelems_io = nelems_io_in
635  else
636  allocate(nelems_io(io_domain_npes))
637  nelems_io = 0
638  call mpp_get_ug_compute_domains(io_domain, &
639  size=nelems_io)
640  endif
641 
642  !The number of data elements in the non-unstructured dimensions are
643  !required to be the same for all ranks. Perform gathers to check this.
644  size_fdata_dim_2 = size(fdata,2)
645  size_fdata_dim_3 = size(fdata,3)
646  size_fdata_dim_4 = size(fdata,4)
647 
648  !Allocate arrays which will be used to gather the data to be written
649  !onto the root rank of the pelist.
650  mynelems = size(fdata,1)
651  allocate(sbuff(mynelems*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4))
652  if (mpp_pe() .eq. pelist(1)) then
653  nelems = sum(nelems_io)
654  allocate(rbuff(nelems*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4))
655  else
656  allocate(rbuff(1))
657  endif
658 
659  !Load the data into the sbuff array. The data is transposed so that the
660  !gather may be performed more easily.
661  do k = 1,mynelems
662  do j = 1,size_fdata_dim_2
663  do i = 1,size_fdata_dim_3
664  do n = 1,size_fdata_dim_4
665  sbuff((k-1)*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4 &
666  + (j-1)*size_fdata_dim_3*size_fdata_dim_4 &
667  + (i-1)*size_fdata_dim_4 + n) = fdata(k,j,i,n)
668  enddo
669  enddo
670  enddo
671  enddo
672 
673  !Perform the gather of data onto the root rank (pelist(1)).
674  call mpp_gather(sbuff, &
675  size(sbuff), &
676  rbuff, &
677  nelems_io*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4, &
678  pelist)
679 
680  !Write out the data to the file. This is only done by the root rank
681  !of the pelist.
682  if (mpp_pe() .eq. pelist(1)) then
683  allocate(cdata(nelems,size_fdata_dim_2,size_fdata_dim_3,size_fdata_dim_4))
684  cdata = fill
685  do n = 1,size_fdata_dim_4
686  do k = 1,size_fdata_dim_3
687  do j = 1,size_fdata_dim_2
688  do i = 1,nelems
689  cdata(i,j,k,n) = rbuff((i-1)*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4 &
690  + (j-1)*size_fdata_dim_3*size_fdata_dim_4 &
691  + (k-1)*size_fdata_dim_4 + n)
692  enddo
693  enddo
694  enddo
695  enddo
696  field%size(1) = nelems
697  call write_record_r8(funit, &
698  field, &
699  nelems*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4, &
700  cdata, &
701  time_in=tstamp)
702  deallocate(cdata)
703  endif
704 
705  !Deallocate local allocatables.
706  deallocate(sbuff)
707  deallocate(rbuff)
708  deallocate(pelist)
709  deallocate(nelems_io)
710 
711  !Stop the mpp timer.
712  call mpp_clock_end(mpp_write_clock)
713 
714  return
715 end subroutine mpp_io_unstructured_write_r8_4d
716 
717 !------------------------------------------------------------------------------
718 
719 !----------
720 
721 !------------------------------------------------------------------------------
722 !>Write data for a 1D field associated with an unstructured mpp domain to a
723 !!restart file.
725  field, &
726  domain, &
727  fdata, &
728  nelems_io, &
729  tstamp, &
730  default_data)
731 
732  !Inputs/outputs
733  integer(i4_kind),intent(in) :: funit !<A file unit for the to which the
734  !! data will be written
735  type(fieldtype),intent(inout) :: field !<A field whose data will be written
736  type(domainug),intent(inout) :: domain !<An unstructured mpp domain associatd
737  !! with the inputted file
738  real(KIND=r4_kind),dimension(:),intent(inout) :: fdata !<The data that will be written to the file
739  integer,dimension(:),intent(in) :: nelems_io !<Number of grid points in the compressed
740  !! dimension for each rank (correct
741  !!sizes only exist for the
742  !!root rank of I/O domain pelist)
743  real(KIND=r4_kind),intent(in),optional :: tstamp !<A time value
744  real(KIND=r4_kind),intent(in), optional :: default_data !<Fill value for the inputted field
745 
746  !Local variables
747  real(KIND=r4_kind) :: fill !<Fill value for the inputted field (default: zero)
748  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
749  integer(i4_kind) :: io_domain_npes !<The total number of ranks in an I/O domain pelist
750  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
751  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
752  !! to be written
753  !!by the root rank of the pelist
754  real(kind=r4_kind),dimension(:),allocatable :: rbuff !<Buffer used to gather the data onto
755  !! the root rank of the pelist
756  real(kind=r4_kind),dimension(:),allocatable :: cdata !<Array used to write the data to the
757  !! file after the gather is performed
758  integer(i4_kind) :: i !<Loop variable
759 
760  !Start the mpp timer.
761  !mpp_write_clock is a module variable.
762  call mpp_clock_begin(mpp_write_clock)
763 
764  !Make sure that the module is initialized.
765  if (.not. module_is_initialized) then
766  call mpp_error(fatal, &
767  "mpp_io_unstructured_write_r_1D:" &
768  //" you must first call mpp_io_init.")
769  endif
770 
771  !Make sure that another NetCDF file is not currently using the inputted
772  !file unit.
773  if (.not. mpp_file(funit)%valid) then
774  call mpp_error(fatal, &
775  "mpp_io_unstructured_write_r_1D:" &
776  //" the inputted file unit is already in use.")
777  endif
778 
779  !Set the fill value for the field.
780  fill = 0.0
781  if (present(default_data)) then
782  fill = default_data
783  endif
784 
785  !Point to the I/O domain associated with the inputted unstructured mpp
786  !domain.
787  io_domain => null()
788  io_domain => mpp_get_ug_io_domain(domain)
789 
790  !Get the pelist associated with the I/O domain.
791  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
792  allocate(pelist(io_domain_npes))
793  call mpp_get_ug_domain_pelist(io_domain, &
794  pelist)
795 
796  !Make sure that only the root rank of the pelist will write to the file.
797  !This check is needed because data is only gathered on the lowest rank
798  !of the pelist.
799  if (mpp_pe() .eq. pelist(1) .and. .not. &
800  mpp_file(funit)%write_on_this_pe) then
801  call mpp_error(fatal, &
802  "mpp_io_unstructured_write_r_1D:" &
803  //" the root rank of the pelist must be allowed" &
804  //" to perform the write.")
805  endif
806  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
807  call mpp_error(fatal, &
808  "mpp_io_unstructured_write_r_1D:" &
809  //" the non-root ranks of the pelist are not" &
810  //" allowed to perform the write.")
811  endif
812 
813  !Allocate an array which will be used to gather the data to be written
814  !onto the root rank of the pelist.
815  if (mpp_pe() .eq. pelist(1)) then
816  nelems = sum(nelems_io)
817  allocate(rbuff(nelems))
818  else
819  allocate(rbuff(1))
820  endif
821 
822  !Perform the gather of data onto the root rank (pelist(1)).
823  call mpp_gather(fdata, &
824  size(fdata), &
825  rbuff, &
826  nelems_io, &
827  pelist)
828 
829  !Write out the data to the file. This is only done by the root rank
830  !of the pelist.
831  if (mpp_pe() .eq. pelist(1)) then
832  allocate(cdata(nelems))
833  cdata = fill
834  do i = 1,nelems
835  cdata(i) = rbuff(i)
836  enddo
837  field%size(1) = nelems
838  call write_record_r4(funit, &
839  field, &
840  nelems, &
841  cdata, &
842  time_in=tstamp)
843  deallocate(cdata)
844  endif
845 
846  !Deallocate local allocatables.
847  deallocate(rbuff)
848  deallocate(pelist)
849 
850  !Stop the mpp timer.
851  call mpp_clock_end(mpp_write_clock)
852 
853  return
854 end subroutine mpp_io_unstructured_write_r4_1d
855 
856 !------------------------------------------------------------------------------
857 !>Write data for a 2D field associated with an unstructured mpp domain to a
858 !!restart file.
860  field, &
861  domain, &
862  fdata, &
863  nelems_io, &
864  tstamp, &
865  default_data)
866 
867  !Inputs/outputs
868  integer(i4_kind),intent(in) :: funit !<A file unit for the to which the
869  !! data will be written
870  type(fieldtype),intent(inout) :: field !<A field whose data will be written
871  type(domainug),intent(inout) :: domain !<An unstructured mpp domain associatd
872  !! with the inputted file
873  real(KIND=r4_kind),dimension(:,:),intent(inout) :: fdata !<The data that will be written to the file
874  integer,dimension(:),intent(in) :: nelems_io !<Number of grid points in the compressed
875  !! dimension for each rank
876  !!(correct sizes only exist
877  !!for the root rank of I/O domain pelist)
878  real(KIND=r4_kind),intent(in),optional :: tstamp !<A time value
879  real(KIND=r4_kind),intent(in), optional :: default_data !<Fill value for the inputted field
880 
881  !Local variables
882  real(KIND=r4_kind) :: fill !<Fill value for the inputted field (default: zero)
883  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
884  integer(i4_kind) :: io_domain_npes !<The total number of ranks in an I/O domain pelist
885  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
886  integer(i4_kind) :: dim_size_1 !<Number of data points in the first
887  !! dimension (size(fdata,1))
888  integer(i4_kind) :: dim_size_2 !<Number of data points in the second
889  !! dimension (size(fdata,2))
890  real(kind=r4_kind),dimension(:),allocatable :: sbuff !<Buffer used to gather the data
891  !! onto the root rank of the pelist
892  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
893  !! to be written
894  !!by the root rank of the pelist
895  real(kind=r4_kind),dimension(:),allocatable :: rbuff !<Buffer used to gather the data
896  !! onto the root rank of the pelist
897  real(kind=r4_kind),dimension(:,:),allocatable :: cdata !<Array used to write the data to
898  !! the file after the gather is performed
899  integer(i4_kind) :: offset_r !<Offset for rbuff used to reorder
900  !! data before netCDF write
901  integer(i4_kind) :: offset_c !<Offset for cdata used to reorder
902  !! data before netCDF write
903  integer(i4_kind) :: i !<Loop variable
904  integer(i4_kind) :: j !<Loop variable
905  integer(i4_kind) :: k !<Loop variable
906 
907  !Start the mpp timer.
908  !mpp_write_clock is a module variable.
909  call mpp_clock_begin(mpp_write_clock)
910 
911  !Make sure that the module is initialized.
912  if (.not. module_is_initialized) then
913  call mpp_error(fatal, &
914  "mpp_io_unstructured_write_r_2D:" &
915  //" you must first call mpp_io_init.")
916  endif
917 
918  !Make sure that another NetCDF file is not currently using the inputted
919  !file unit.
920  if (.not. mpp_file(funit)%valid) then
921  call mpp_error(fatal, &
922  "mpp_io_unstructured_write_r_2D:" &
923  //" the inputted file unit is already in use.")
924  endif
925 
926  !Set the fill value for the field.
927  fill = 0.0
928  if (present(default_data)) then
929  fill = default_data
930  endif
931 
932  !Point to the I/O domain associated with the inputted unstructured mpp
933  !domain.
934  io_domain => null()
935  io_domain => mpp_get_ug_io_domain(domain)
936 
937  !Get the pelist associated with the I/O domain.
938  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
939  allocate(pelist(io_domain_npes))
940  call mpp_get_ug_domain_pelist(io_domain, &
941  pelist)
942 
943  !Make sure that only the root rank of the pelist will write to the file.
944  !This check is needed because data is only gathered on the lowest rank
945  !of the pelist.
946  if (mpp_pe() .eq. pelist(1) .and. .not. &
947  mpp_file(funit)%write_on_this_pe) then
948  call mpp_error(fatal, &
949  "mpp_io_unstructured_write_r_2D:" &
950  //" the root rank of the pelist must be allowed" &
951  //" to perform the write.")
952  endif
953  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
954  call mpp_error(fatal, &
955  "mpp_io_unstructured_write_r_2D:" &
956  //" the non-root ranks of the pelist are not" &
957  //" allowed to perform the write.")
958  endif
959 
960  !Load the data elements for each rank into a one dimensional array, which
961  !will be used to gather the data onto the root rank of the pelist.
962  allocate(sbuff(size(fdata)))
963  dim_size_1 = size(fdata,1)
964  dim_size_2 = size(fdata,2)
965  do j = 1,dim_size_2
966  do i = 1,dim_size_1
967  sbuff((j-1)*dim_size_1+i) = fdata(i,j)
968  enddo
969  enddo
970 
971  !Allocate an array which will be used to gather the data to be written
972  !onto the root rank of the pelist.
973  if (mpp_pe() .eq. pelist(1)) then
974  nelems = sum(nelems_io)
975  allocate(rbuff(nelems*dim_size_2))
976  else
977  allocate(rbuff(1))
978  endif
979 
980  !Perform the gather of data onto the root rank (pelist(1)).
981  call mpp_gather(sbuff, &
982  size(sbuff), &
983  rbuff, &
984  nelems_io*dim_size_2, &
985  pelist)
986 
987  !Reorder the gather data so that is of the form (nelems,dim_size_2). Write
988  !out the data to the file. This is only done by the root rank of the
989  !pelist.
990  if (mpp_pe() .eq. pelist(1)) then
991  allocate(cdata(nelems,dim_size_2))
992  cdata = fill
993  do j = 1,dim_size_2
994  offset_c = 0
995  do k = 1,io_domain_npes
996  if (k .gt. 1) then
997  offset_r = (j-1)*nelems_io(k) + dim_size_2*(sum(nelems_io(1:k-1)))
998  else
999  offset_r = (j-1)*nelems_io(k)
1000  endif
1001  do i = 1,nelems_io(k)
1002  cdata(i+offset_c,j) = rbuff(i+offset_r)
1003  enddo
1004  offset_c = offset_c + nelems_io(k)
1005  enddo
1006  enddo
1007  field%size(1) = nelems
1008  call write_record_r4(funit, &
1009  field, &
1010  nelems*dim_size_2, &
1011  cdata, &
1012  time_in=tstamp)
1013  deallocate(cdata)
1014  endif
1015 
1016  !Deallocate local allocatables.
1017  deallocate(sbuff)
1018  deallocate(rbuff)
1019  deallocate(pelist)
1020 
1021  !Stop the mpp timer.
1022  call mpp_clock_end(mpp_write_clock)
1023 
1024  return
1025 end subroutine mpp_io_unstructured_write_r4_2d
1026 
1027 !------------------------------------------------------------------------------
1028 !>Write data for a 3D field associated with an unstructured mpp domain to a
1029 !!restart file.
1031  field, &
1032  domain, &
1033  fdata, &
1034  nelems_io, &
1035  tstamp, &
1036  default_data)
1037 
1038  !Inputs/outputs
1039  integer(i4_kind),intent(in) :: funit !<A file unit for the to which
1040  !! the data will be written
1041  type(fieldtype),intent(inout) :: field !<A field whose data will be written
1042  type(domainug),intent(inout) :: domain !<An unstructured mpp domain associatd
1043  !! with the inputted file
1044  real(KIND=r4_kind),dimension(:,:,:),intent(inout) :: fdata !<The data that will be written to the file
1045  integer,dimension(:),intent(in) :: nelems_io !<Number of grid points in the
1046  !! compressed dimension for each rank
1047  !!(correct sizes only
1048  !!exist for the root rank of I/O domain pelist)
1049  real(KIND=r4_kind),intent(in),optional :: tstamp !<A time value
1050  real(KIND=r4_kind),intent(in), optional :: default_data !<Fill value for the inputted field
1051 
1052  !Local variables
1053  real(KIND=r4_kind) :: fill !<Fill value for the inputted field
1054  !! (defaults: zero)
1055  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
1056  integer(i4_kind) :: io_domain_npes !<The total number of ranks in
1057  !! an I/O domain pelist
1058  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
1059  integer(i4_kind) :: dim_size_1 !<Number of data points in the
1060  !! first dimension (size(fdata,1))
1061  integer(i4_kind) :: dim_size_2 !<Number of data points in the
1062  !! second dimension (size(fdata,2))
1063  integer(i4_kind) :: dim_size_3 !<Number of data points in the
1064  !! second dimension (size(fdata,3))
1065  real(kind=r4_kind),dimension(:),allocatable :: sbuff !<Buffer used to gather the data
1066  !! onto the root rank of the pelist
1067  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
1068  !! to be written
1069  !!by the root rank of the pelist
1070  real(kind=r4_kind),dimension(:),allocatable :: rbuff !<Buffer used to gather the data
1071  !! onto the root rank of the pelist
1072  real(kind=r4_kind),dimension(:,:,:),allocatable :: cdata !<Array used to write the data
1073  !! to the file after the gather is performed
1074  integer(i4_kind) :: offset_r !<Offset for rbuff used to reorder
1075  !! data before netCDF write
1076  integer(i4_kind) :: offset_c !<Offset for cdata used to reorder
1077  !! data before netCDF write
1078  integer(i4_kind) :: i !<Loop variable
1079  integer(i4_kind) :: j !<Loop variable
1080  integer(i4_kind) :: k !<Loop variable
1081  integer(i4_kind) :: m !<Loop variable
1082 
1083  !Start the mpp timer.
1084  !mpp_write_clock is a module variable.
1085  call mpp_clock_begin(mpp_write_clock)
1086 
1087  !Make sure that the module is initialized.
1088  if (.not. module_is_initialized) then
1089  call mpp_error(fatal, &
1090  "mpp_io_unstructured_write_r_3D:" &
1091  //" you must first call mpp_io_init.")
1092  endif
1093 
1094  !Make sure that another NetCDF file is not currently using the inputted
1095  !file unit.
1096  if (.not. mpp_file(funit)%valid) then
1097  call mpp_error(fatal, &
1098  "mpp_io_unstructured_write_r_3D:" &
1099  //" the inputted file unit is already in use.")
1100  endif
1101 
1102  !Set the fill value for the field.
1103  fill = 0.0
1104  if (present(default_data)) then
1105  fill = default_data
1106  endif
1107 
1108  !Point to the I/O domain associated with the inputted unstructured mpp
1109  !domain.
1110  io_domain => null()
1111  io_domain => mpp_get_ug_io_domain(domain)
1112 
1113  !Get the pelist associated with the I/O domain.
1114  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
1115  allocate(pelist(io_domain_npes))
1116  call mpp_get_ug_domain_pelist(io_domain, &
1117  pelist)
1118 
1119  !Make sure that only the root rank of the pelist will write to the file.
1120  !This check is needed because data is only gathered on the lowest rank
1121  !of the pelist.
1122  if (mpp_pe() .eq. pelist(1) .and. .not. &
1123  mpp_file(funit)%write_on_this_pe) then
1124  call mpp_error(fatal, &
1125  "mpp_io_unstructured_write_r_3D:" &
1126  //" the root rank of the pelist must be allowed" &
1127  //" to perform the write.")
1128  endif
1129  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
1130  call mpp_error(fatal, &
1131  "mpp_io_unstructured_write_r_3D:" &
1132  //" the non-root ranks of the pelist are not" &
1133  //" allowed to perform the write.")
1134  endif
1135 
1136  !Load the data elements for each rank into a one dimensional array, which
1137  !will be used to gather the data onto the root rank of the pelist.
1138  allocate(sbuff(size(fdata)))
1139  dim_size_1 = size(fdata,1)
1140  dim_size_2 = size(fdata,2)
1141  dim_size_3 = size(fdata,3)
1142  do k = 1,dim_size_3
1143  do j = 1,dim_size_2
1144  do i = 1,dim_size_1
1145  sbuff((k-1)*dim_size_2*dim_size_1+(j-1)*dim_size_1+i) = fdata(i,j,k)
1146  enddo
1147  enddo
1148  enddo
1149 
1150  !Allocate an array which will be used to gather the data to be written
1151  !onto the root rank of the pelist.
1152  if (mpp_pe() .eq. pelist(1)) then
1153  nelems = sum(nelems_io)
1154  allocate(rbuff(nelems*dim_size_2*dim_size_3))
1155  else
1156  allocate(rbuff(1))
1157  endif
1158 
1159  !Perform the gather of data onto the root rank (pelist(1)).
1160  call mpp_gather(sbuff, &
1161  size(sbuff), &
1162  rbuff, &
1163  nelems_io*dim_size_2*dim_size_3, &
1164  pelist)
1165 
1166  !Reorder the gather data so that is of the form (nelems,dim_size_2). Write
1167  !out the data to the file. This is only done by the root rank of the
1168  !pelist.
1169  if (mpp_pe() .eq. pelist(1)) then
1170  allocate(cdata(nelems,dim_size_2,dim_size_3))
1171  cdata = fill
1172  do m = 1,dim_size_3
1173  do j = 1,dim_size_2
1174  offset_c = 0
1175  do k = 1,io_domain_npes
1176  if (k .gt. 1) then
1177  offset_r = (m-1)*dim_size_2*nelems_io(k) + &
1178  (j-1)*nelems_io(k) + &
1179  dim_size_2*dim_size_3*(sum(nelems_io(1:k-1)))
1180  else
1181  offset_r = (m-1)*dim_size_2*nelems_io(k) + &
1182  (j-1)*nelems_io(k)
1183  endif
1184  do i = 1,nelems_io(k)
1185  cdata(i+offset_c,j,m) = rbuff(i+offset_r)
1186  enddo
1187  offset_c = offset_c + nelems_io(k)
1188  enddo
1189  enddo
1190  enddo
1191  field%size(1) = nelems
1192  call write_record_r4(funit, &
1193  field, &
1194  nelems*dim_size_2*dim_size_3, &
1195  cdata, &
1196  time_in=tstamp)
1197  deallocate(cdata)
1198  endif
1199 
1200  !Deallocate local allocatables.
1201  deallocate(sbuff)
1202  deallocate(rbuff)
1203  deallocate(pelist)
1204 
1205  !Stop the mpp timer.
1206  call mpp_clock_end(mpp_write_clock)
1207 
1208  return
1209 end subroutine mpp_io_unstructured_write_r4_3d
1210 
1211 !------------------------------------------------------------------------------
1212 !>Write data for a 4D field associated with an unstructured mpp domain to a
1213 !!restart file.
1215  field, &
1216  domain, &
1217  fdata, &
1218  nelems_io_in, &
1219  tstamp, &
1220  default_data)
1221 
1222  !Inputs/outputs
1223  integer(i4_kind),intent(in) :: funit !<A file unit for the to which
1224  !! the data will be written
1225  type(fieldtype),intent(inout) :: field !<A field whose data will be written
1226  type(domainug),intent(inout) :: domain !<An unstructured mpp domain
1227  !! associatd with the inputted file
1228  real(KIND=r4_kind),dimension(:,:,:,:),intent(inout) :: fdata !<The data that will be written to the file
1229  integer,dimension(:),intent(in),optional :: nelems_io_in !<Number of grid points in the
1230  !! unstructured dimension for each rank
1231  !!(correct sizes only
1232  !!(cexist for the root
1233  !!rank of I/O domain pelist
1234  real(KIND=r4_kind),intent(in),optional :: tstamp !<A time value
1235  real(KIND=r4_kind),intent(in), optional :: default_data !<Fill value for the inputted field
1236 
1237  !Local variables
1238  real(KIND=r4_kind) :: fill !<Fill value for the inputted field (default: zero)
1239  type(domainug),pointer :: io_domain !<Pointer to the unstructured I/O domain
1240  integer(i4_kind) :: io_domain_npes !<The total number of ranks in an I/O domain pelist
1241  integer(i4_kind),dimension(:),allocatable :: pelist !<A pelist
1242  integer(i4_kind),dimension(:),allocatable :: nelems_io !<Number of grid points in the unstructured
1243  !! dimension for each rank
1244  integer(i4_kind) :: compute_size !<Size of the unstructured compute
1245  !! domain for the current rank
1246  integer(i4_kind) :: size_fdata_dim_2 !<Number of data points in a non-unstructured
1247  !! dimension (size(fdata,2))
1248  integer(i4_kind) :: size_fdata_dim_3 !<Number of data points in a non-unstructured
1249  !! dimension (size(fdata,3))
1250  integer(i4_kind) :: size_fdata_dim_4 !<Number of data points in a non-unstructured
1251  !! dimension (size(fdata,3))
1252  integer(i4_kind) :: mynelems !<Number of data points in the unstructured
1253  !! dimension (size(fdata,1))
1254  real(kind=r4_kind),dimension(:),allocatable :: sbuff !<Buffer used to gather the data
1255  !! onto the root rank of the pelist
1256  integer(i4_kind) :: nelems !<Total number of data points (sum(nelems_io))
1257  !! to be written
1258  !!by the root rank of the pelist
1259  real(kind=r4_kind),dimension(:),allocatable :: rbuff !<Buffer to gather the data onto
1260  !! root rank of pelist
1261  real(kind=r4_kind),dimension(:,:,:,:),allocatable :: cdata !<Array to write the data to file
1262  !! after gather is performed
1263  integer(i4_kind) :: i !<Loop variable
1264  integer(i4_kind) :: j !<Loop variable
1265  integer(i4_kind) :: k !<Loop variable
1266  integer(i4_kind) :: n !<Loop variable
1267 
1268  !Start the mpp timer.
1269  !mpp_write_clock is a module variable.
1270  call mpp_clock_begin(mpp_write_clock)
1271 
1272  !Make sure that the module is initialized.
1273  if (.not. module_is_initialized) then
1274  call mpp_error(fatal, &
1275  "mpp_io_unstructured_write_compressed_r_4D:" &
1276  //" you must first call mpp_io_init.")
1277  endif
1278 
1279  !Make sure that another NetCDF file is not currently using the inputted
1280  !file unit.
1281  if (.not. mpp_file(funit)%valid) then
1282  call mpp_error(fatal, &
1283  "mpp_io_unstructured_write_compressed_r_4D:" &
1284  //" the inputted file unit is already in use.")
1285  endif
1286 
1287  !Set the fill value for the field.
1288  fill = 0.0
1289  if (present(default_data)) then
1290  fill = default_data
1291  endif
1292 
1293  !Point to the I/O domain associated with the inputted unstructured mpp
1294  !domain.
1295  io_domain => null()
1296  io_domain => mpp_get_ug_io_domain(domain)
1297 
1298  !Get the pelist associated with the I/O domain.
1299  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
1300  allocate(pelist(io_domain_npes))
1301  call mpp_get_ug_domain_pelist(io_domain, &
1302  pelist)
1303 
1304  !Make sure that only the root rank of the pelist will write to the file.
1305  !This check is needed because data is only gathered on the lowest rank
1306  !of the pelist.
1307  if (mpp_pe() .eq. pelist(1) .and. .not. &
1308  mpp_file(funit)%write_on_this_pe) then
1309  call mpp_error(fatal, &
1310  "mpp_io_unstructured_write_compressed_r_4D:" &
1311  //" the root rank of the pelist must be allowed" &
1312  //" to perform the write.")
1313  endif
1314  if (mpp_pe() .ne. pelist(1) .and. mpp_file(funit)%write_on_this_pe) then
1315  call mpp_error(fatal, &
1316  "mpp_io_unstructured_write_compressed_r_4D:" &
1317  //" the non-root ranks of the pelist are not" &
1318  //" allowed to perform the write.")
1319  endif
1320 
1321  !For the 3D unstructured case, data is assumed to be of the form
1322  !fdata = fdata(unstructured,z,cc). The number of data elements in the
1323  !unstructured dimension (size(fdata,1)) may differ between ranks.
1324  !If not passed in, the number of data elements in the unstructured
1325  !dimension must be gathered on the root rank of the pelist. The number
1326  !data elements in the unstructured dimension should be equal to the size
1327  !of the unstructured computed domain.
1328  if (present(nelems_io_in)) then
1329  allocate(nelems_io(size(nelems_io_in)))
1330  nelems_io = nelems_io_in
1331  else
1332  allocate(nelems_io(io_domain_npes))
1333  nelems_io = 0
1334  call mpp_get_ug_compute_domains(io_domain, &
1335  size=nelems_io)
1336  endif
1337 
1338  !The number of data elements in the non-unstructured dimensions are
1339  !required to be the same for all ranks. Perform gathers to check this.
1340  size_fdata_dim_2 = size(fdata,2)
1341  size_fdata_dim_3 = size(fdata,3)
1342  size_fdata_dim_4 = size(fdata,4)
1343 
1344  !Allocate arrays which will be used to gather the data to be written
1345  !onto the root rank of the pelist.
1346  mynelems = size(fdata,1)
1347  allocate(sbuff(mynelems*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4))
1348  if (mpp_pe() .eq. pelist(1)) then
1349  nelems = sum(nelems_io)
1350  allocate(rbuff(nelems*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4))
1351  else
1352  allocate(rbuff(1))
1353  endif
1354 
1355  !Load the data into the sbuff array. The data is transposed so that the
1356  !gather may be performed more easily.
1357  do k = 1,mynelems
1358  do j = 1,size_fdata_dim_2
1359  do i = 1,size_fdata_dim_3
1360  do n = 1,size_fdata_dim_4
1361  sbuff((k-1)*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4 &
1362  + (j-1)*size_fdata_dim_3*size_fdata_dim_4 &
1363  + (i-1)*size_fdata_dim_4 + n) = fdata(k,j,i,n)
1364  enddo
1365  enddo
1366  enddo
1367  enddo
1368 
1369  !Perform the gather of data onto the root rank (pelist(1)).
1370  call mpp_gather(sbuff, &
1371  size(sbuff), &
1372  rbuff, &
1373  nelems_io*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4, &
1374  pelist)
1375 
1376  !Write out the data to the file. This is only done by the root rank
1377  !of the pelist.
1378  if (mpp_pe() .eq. pelist(1)) then
1379  allocate(cdata(nelems,size_fdata_dim_2,size_fdata_dim_3,size_fdata_dim_4))
1380  cdata = fill
1381  do n = 1,size_fdata_dim_4
1382  do k = 1,size_fdata_dim_3
1383  do j = 1,size_fdata_dim_2
1384  do i = 1,nelems
1385  cdata(i,j,k,n) = rbuff((i-1)*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4 &
1386  + (j-1)*size_fdata_dim_3*size_fdata_dim_4 &
1387  + (k-1)*size_fdata_dim_4 + n)
1388  enddo
1389  enddo
1390  enddo
1391  enddo
1392  field%size(1) = nelems
1393  call write_record_r4(funit, &
1394  field, &
1395  nelems*size_fdata_dim_2*size_fdata_dim_3*size_fdata_dim_4, &
1396  cdata, &
1397  time_in=tstamp)
1398  deallocate(cdata)
1399  endif
1400 
1401  !Deallocate local allocatables.
1402  deallocate(sbuff)
1403  deallocate(rbuff)
1404  deallocate(pelist)
1405  deallocate(nelems_io)
1406 
1407  !Stop the mpp timer.
1408  call mpp_clock_end(mpp_write_clock)
1409 
1410  return
1411 end subroutine mpp_io_unstructured_write_r4_4d
1412 !> @}
subroutine mpp_io_unstructured_write_r8_2d(funit, field, domain, fdata, nelems_io, tstamp, default_data)
Write data for a 2D field associated with an unstructured mpp domain to a restart file.
subroutine mpp_io_unstructured_write_r8_4d(funit, field, domain, fdata, nelems_io_in, tstamp, default_data)
Write data for a 4D field associated with an unstructured mpp domain to a restart file.
subroutine mpp_io_unstructured_write_r8_3d(funit, field, domain, fdata, nelems_io, tstamp, default_data)
Write data for a 3D field associated with an unstructured mpp domain to a restart file.
subroutine mpp_io_unstructured_write_r4_4d(funit, field, domain, fdata, nelems_io_in, tstamp, default_data)
Write data for a 4D field associated with an unstructured mpp domain to a restart file.
subroutine mpp_io_unstructured_write_r4_3d(funit, field, domain, fdata, nelems_io, tstamp, default_data)
Write data for a 3D field associated with an unstructured mpp domain to a restart file.
subroutine mpp_io_unstructured_write_r4_2d(funit, field, domain, fdata, nelems_io, tstamp, default_data)
Write data for a 2D field associated with an unstructured mpp domain to a restart file.
subroutine mpp_io_unstructured_write_r4_1d(funit, field, domain, fdata, nelems_io, tstamp, default_data)
Write data for a 1D field associated with an unstructured mpp domain to a restart file.
subroutine mpp_io_unstructured_write_r8_1d(funit, field, domain, fdata, nelems_io, tstamp, default_data)
Write data for a 1D field associated with an unstructured mpp domain to a restart file.
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407