Line data Source code
1 : ! Utility functions in support of PIO io interface
2 : module cam_pio_utils
3 :
4 : use pio, only: io_desc_t, iosystem_desc_t, file_desc_t, var_desc_t
5 : use pio, only: pio_freedecomp, pio_rearr_subset, pio_rearr_box
6 : use shr_kind_mod, only: r4 => shr_kind_r4, r8 => shr_kind_r8
7 : use cam_logfile, only: iulog
8 : use perf_mod, only: t_startf, t_stopf
9 : use spmd_utils, only: masterproc
10 :
11 : implicit none
12 : private
13 : save
14 :
15 : public :: cam_pio_createfile ! Create a new NetCDF file for PIO writing
16 : public :: cam_pio_openfile ! Open an existing NetCDF file
17 : public :: cam_pio_closefile ! Close an open PIO file handle
18 : public :: cam_pio_fileexists ! Check if file exists
19 : public :: cam_pio_newdecomp ! Create a new PIO decompsition (mapping)
20 : public :: init_pio_subsystem ! called from cam_comp
21 : public :: cam_pio_get_decomp ! Find an existing decomp or create a new one
22 : public :: cam_pio_handle_error ! If error, print a custom error message
23 : public :: cam_pio_set_fill ! Set the PIO fill value to PIO_FILL
24 : public :: cam_pio_inq_var_fill ! Return the buffer fill value
25 :
26 : public :: cam_permute_array
27 : public :: calc_permutation
28 :
29 : ! Convenience interfaces
30 : public :: cam_pio_def_dim
31 : public :: cam_pio_def_var
32 : public :: cam_pio_get_var
33 :
34 : ! General utility
35 : public :: cam_pio_var_info
36 : public :: cam_pio_find_var
37 : public :: cam_pio_check_var
38 :
39 : public :: clean_iodesc_list
40 :
41 : ! For help debugging code
42 : public :: cam_pio_dump_field
43 :
44 : integer :: pio_iotype
45 : integer :: pio_rearranger
46 : integer :: pio_ioformat
47 :
48 : ! This variable should be private ?
49 : type(iosystem_desc_t), pointer, public :: pio_subsystem => null()
50 :
51 : ! Some private string length parameters
52 : integer, parameter :: errormsg_str_len = 128
53 :
54 : ! The iodesc_list allows us to cache existing PIO decompositions
55 : ! The tag needs the dim lengths, the dtype and map id (+ optional permutation)
56 : integer, parameter :: tag_len = 48
57 : type iodesc_list
58 : character(tag_len) :: tag
59 : type(io_desc_t), pointer :: iodesc => NULL()
60 : type(iodesc_list), pointer :: next => NULL()
61 : end type iodesc_list
62 :
63 : type(iodesc_list), target :: iodesc_list_top
64 :
65 : ! Create a special type to hold a var_desc_t pointer so we can have an
66 : ! array of them
67 : type, public :: vdesc_ptr
68 : type(var_desc_t), pointer :: vd => NULL()
69 : end type vdesc_ptr
70 :
71 : interface cam_pio_def_var
72 : module procedure cam_pio_def_var_0d
73 : module procedure cam_pio_def_var_md
74 : end interface
75 :
76 : interface cam_pio_get_var
77 : module procedure cam_pio_get_var_2d_r8
78 : module procedure cam_pio_get_var_2d_r8_perm
79 : module procedure cam_pio_get_var_3d_r8
80 : module procedure cam_pio_get_var_3d_r8_perm
81 : end interface
82 :
83 : interface cam_pio_inq_var_fill
84 : module procedure inq_var_fill_i4
85 : module procedure inq_var_fill_r4
86 : module procedure inq_var_fill_r8
87 : end interface cam_pio_inq_var_fill
88 :
89 : interface calc_permutation
90 : module procedure calc_permutation_int
91 : module procedure calc_permutation_char
92 : end interface
93 :
94 : interface cam_permute_array
95 : module procedure permute_array_int
96 : module procedure permute_array_r8
97 : end interface
98 :
99 : interface cam_pio_dump_field
100 : module procedure dump_field_2d_d
101 : module procedure dump_field_3d_d
102 : module procedure dump_field_4d_d
103 : module procedure dump_field_6d_d
104 : end interface
105 :
106 : contains
107 :
108 : ! use_scam_limits is a private interface used to gather information about
109 : ! single-column usage and limits for use by the cam_pio_get_var interfaces
110 : ! This still only works for lat/lon dycores
111 0 : logical function use_scam_limits(File, start, kount, dimnames)
112 : use shr_scam_mod, only: shr_scam_getCloseLatLon
113 : use scamMod, only: scmlat, scmlon, single_column
114 : use cam_abortutils, only: endrun
115 :
116 : ! Dummy arguments
117 : type(file_desc_t), intent(inout) :: File
118 : integer, intent(inout) :: start(:)
119 : integer, intent(inout) :: kount(:)
120 : character(len=*), optional, intent(in) :: dimnames(:)
121 :
122 : ! Local variables
123 : character(len=*), parameter :: subname='USE_SCAM_LIMITS'
124 : real(r8) :: closelat, closelon
125 : integer :: latidx, lonidx
126 : integer :: i
127 : logical :: latfound
128 :
129 0 : use_scam_limits = single_column
130 0 : if (use_scam_limits) then
131 : call shr_scam_getCloseLatLon(File, scmlat, scmlon, closelat, closelon, &
132 0 : latidx, lonidx)
133 0 : if (present(dimnames)) then
134 0 : if (trim(dimnames(1)) == 'lon') then
135 0 : start(1) = lonidx ! First dim always lon for Eulerian dycore
136 : ! This could be generalized -- for now, stick with single column
137 0 : kount(1) = 1
138 : else
139 0 : call endrun(trim(subname)//': lon should be first dimension')
140 : end if
141 0 : latfound = .false.
142 0 : do i = 2, size(dimnames)
143 0 : if (size(start) < i) then
144 0 : call endrun(trim(subname)//': start too small')
145 : end if
146 0 : if (trim(dimnames(i)) == 'lat') then
147 0 : start(i) = latidx
148 : ! This could be generalized -- for now, stick with single column
149 0 : kount(i) = 1
150 0 : latfound = .true.
151 : end if
152 : end do
153 0 : if (.not. latfound) then
154 0 : call endrun(trim(subname)//': lat dimension not found')
155 : end if
156 : else
157 : ! No dimnames, assume standard positions (lon,lat)
158 0 : start(1) = lonidx
159 0 : start(2) = latidx
160 : ! This could be generalized -- for now, stick with single column
161 0 : kount(1:2) = 1
162 : end if
163 : end if
164 :
165 0 : end function use_scam_limits
166 :
167 : ! calc_permutation: Calculate a permutation array if filedims and arraydims
168 : ! are in a different order
169 : ! E.g.: If filedims is (lon, lat, lev, time) and
170 : ! arraydims is (lon, lev, lat), then
171 : ! perm is (1, 3, 2) and isperm is set to .true.
172 0 : subroutine calc_permutation_int(filedims, arraydims, perm, isperm)
173 0 : use cam_abortutils, only: endrun
174 :
175 : ! Dummy variables
176 : integer, intent(in) :: filedims(:)
177 : integer, intent(in) :: arraydims(:)
178 : integer, intent(out) :: perm(:)
179 : logical, intent(out) :: isperm
180 :
181 : ! Local variables
182 : character(len=*), parameter :: subname='CALC_PERMUTATION_INT'
183 : integer :: i, j
184 : integer :: adims, fdims
185 :
186 0 : perm = 0
187 0 : isperm = .false.
188 0 : adims = size(arraydims)
189 0 : fdims = size(filedims)
190 :
191 0 : if (size(perm) < adims) then
192 0 : call endrun(trim(subname)//': perm smaller than arraydims')
193 : end if
194 :
195 0 : if (fdims < adims) then
196 0 : call endrun(trim(subname)//': filedims smaller than arraydims')
197 : end if
198 :
199 0 : do i = 1, adims
200 0 : if (arraydims(i) == filedims(i)) then
201 0 : perm(i) = i
202 : else
203 0 : isperm = .true.
204 0 : do j = 1, fdims
205 0 : if (arraydims(i) == filedims(j)) then
206 0 : perm(i) = j
207 0 : exit
208 0 : else if (j == fdims) then
209 0 : call endrun(trim(subname)//': No match for array dimension')
210 : ! No else, just try the next j index
211 : end if
212 : end do
213 : end if
214 : end do
215 :
216 0 : end subroutine calc_permutation_int
217 :
218 27648 : subroutine calc_permutation_char(filedims, arraydims, perm, isperm)
219 : use cam_abortutils, only: endrun
220 :
221 : ! Dummy variables
222 : character(len=*), intent(in) :: filedims(:)
223 : character(len=*), intent(in) :: arraydims(:)
224 : integer, intent(out) :: perm(:)
225 : logical, intent(out) :: isperm
226 :
227 : ! Local variables
228 : character(len=*), parameter :: subname='CALC_PERMUTATION_CHAR'
229 : integer :: i, j
230 : integer :: adims, fdims
231 :
232 82944 : perm = 0
233 27648 : isperm = .false.
234 27648 : adims = size(arraydims)
235 27648 : fdims = size(filedims)
236 :
237 27648 : if (size(perm) < adims) then
238 0 : call endrun(trim(subname)//': perm smaller than arraydims')
239 : end if
240 :
241 27648 : if (fdims < adims) then
242 0 : call endrun(trim(subname)//': filedims smaller than arraydims')
243 : end if
244 :
245 82944 : ILOOP : do i = 1, adims
246 82944 : if (trim(arraydims(i)) == trim(filedims(i))) then
247 55296 : perm(i) = i
248 : else
249 0 : isperm = .true.
250 0 : do j = 1, fdims
251 0 : if (trim(arraydims(i)) == trim(filedims(j))) then
252 0 : perm(i) = j
253 0 : exit
254 0 : else if (j == fdims) then
255 : ! We have no match but for character strings, just say no perm
256 0 : isperm = .false.
257 0 : exit ILOOP
258 : ! No else, just try the next j index
259 : end if
260 : end do
261 : end if
262 : end do ILOOP
263 :
264 27648 : end subroutine calc_permutation_char
265 :
266 0 : subroutine permute_array_int(array, perm)
267 :
268 : ! Dummy arguments
269 : integer, intent(inout) :: array(:)
270 : integer, intent(in) :: perm(:)
271 :
272 : ! Local variables
273 0 : integer, allocatable :: temp(:)
274 : integer :: nelem, i
275 :
276 0 : nelem = size(array)
277 0 : allocate(temp(nelem))
278 0 : temp = array
279 0 : do i = 1, nelem
280 0 : array(i) = temp(perm(i))
281 : end do
282 :
283 0 : deallocate(temp)
284 0 : end subroutine permute_array_int
285 :
286 0 : subroutine permute_array_r8(array, perm)
287 :
288 : ! Dummy arguments
289 : real(r8), intent(inout) :: array(:)
290 : integer, intent(in) :: perm(:)
291 :
292 : ! Local variables
293 0 : real(r8), allocatable :: temp(:)
294 : integer :: nelem, i
295 :
296 0 : nelem = size(array)
297 0 : allocate(temp(nelem))
298 0 : temp = array
299 0 : do i = 1, nelem
300 0 : array(i) = temp(perm(i))
301 : end do
302 :
303 0 : deallocate(temp)
304 0 : end subroutine permute_array_r8
305 :
306 125393664 : subroutine cam_pio_handle_error(ierr, errorstr)
307 : use shr_kind_mod, only: SHR_KIND_CL
308 : use cam_abortutils, only: endrun
309 : use pio, only: pio_noerr
310 :
311 : ! Dummy arguments
312 : integer, intent(in) :: ierr
313 : character(len=*), intent(in) :: errorstr
314 :
315 : ! Local variables
316 : character(len=SHR_KIND_CL) :: errormsg
317 :
318 125393664 : if (ierr /= PIO_NOERR) then
319 0 : write(errormsg, '(a,i6,2a)') '(PIO:', ierr, ') ', trim(errorstr)
320 0 : call endrun(errormsg)
321 : end if
322 :
323 125393664 : end subroutine cam_pio_handle_error
324 :
325 : !-----------------------------------------------------------------------
326 : !
327 : ! cam_pio_var_info: Retrieve variable properties
328 : !
329 : !-----------------------------------------------------------------------
330 37632 : subroutine cam_pio_var_info(ncid, varid, ndims, dimids, dimlens, dimnames, varname, unlimDimID)
331 : use pio, only: PIO_inq_varndims, PIO_inq_vardimid, PIO_inq_dimlen
332 : use pio, only: PIO_inquire, PIO_inq_dimname
333 : use pio, only: PIO_seterrorhandling, PIO_BCAST_ERROR
334 : use cam_abortutils, only: endrun
335 :
336 :
337 : ! Dummy arguments
338 : type(file_desc_t), intent(inout) :: ncid
339 : type(var_desc_t), intent(in) :: varid
340 : integer, intent(out) :: ndims
341 : integer, intent(out) :: dimids(:)
342 : integer, intent(out) :: dimlens(:)
343 : character(len=*), optional, intent(out) :: dimnames(:)
344 : integer, optional, intent(out) :: unlimDimID
345 : character(len=*), optional, intent(in) :: varname
346 :
347 : ! Local variables
348 : integer :: ret ! PIO return value
349 : integer :: i
350 : integer :: err_handling
351 : character(len=128) :: errsuff
352 : !-----------------------------------------------------------------------
353 : ! We will handle errors for this routine
354 :
355 37632 : call PIO_seterrorhandling(ncid, PIO_BCAST_ERROR, err_handling)
356 :
357 264960 : dimids = -1
358 37632 : ndims = 0
359 264192 : dimlens = 0
360 :
361 37632 : if (present(varname)) then
362 36864 : errsuff = ' for '//trim(varname)
363 : else
364 768 : errsuff = ''
365 : end if
366 : ! Check dimensions
367 37632 : ret = PIO_inq_varndims(ncid, varid, ndims)
368 37632 : call cam_pio_handle_error(ret, 'CAM_PIO_VAR_INFO: Error with num dimensions')
369 37632 : if (size(dimids) < ndims) then
370 0 : call endrun('CAM_PIO_VAR_INFO: dimids too small'//trim(errsuff))
371 : end if
372 37632 : ret = PIO_inq_vardimid(ncid, varid, dimids(1:ndims))
373 37632 : call cam_pio_handle_error(ret, 'CAM_PIO_VAR_INFO: Error with inq dim ids'//trim(errsuff))
374 37632 : if (size(dimlens) < ndims) then
375 0 : call endrun('CAM_PIO_VAR_INFO: dimlens too small'//trim(errsuff))
376 : end if
377 124416 : do i = 1, ndims
378 86784 : ret = PIO_inq_dimlen(ncid, dimids(i), dimlens(i))
379 86784 : call cam_pio_handle_error(ret, 'CAM_PIO_VAR_INFO: Error with inq dimlens')
380 124416 : if (present(dimnames)) then
381 75264 : ret = PIO_inq_dimname(ncid, dimids(i), dimnames(i))
382 75264 : call cam_pio_handle_error(ret, 'CAM_PIO_VAR_INFO: Error with inq dimnames')
383 : end if
384 : end do
385 37632 : if (present(unlimDimID)) then
386 0 : ret = PIO_inquire(ncid, unlimitedDimID=unlimDimID)
387 0 : call cam_pio_handle_error(ret, 'CAM_PIO_VAR_INFO: Error with inquire')
388 : end if
389 37632 : call PIO_seterrorhandling(ncid, err_handling)
390 :
391 37632 : end subroutine cam_pio_var_info
392 :
393 0 : subroutine cam_pio_find_var(ncid, varname, varid, found)
394 : use pio, only: pio_inq_varid, pio_noerr
395 : use pio, only: PIO_seterrorhandling, PIO_BCAST_ERROR
396 :
397 : ! Dummy arguments
398 : type(file_desc_t), intent(inout) :: ncid
399 : character(len=*), intent(in) :: varname
400 : type(var_desc_t), intent(out) :: varid
401 : logical, intent(out) :: found
402 :
403 : ! Local variables
404 : integer :: ret ! PIO return value
405 : integer :: err_handling
406 :
407 : !-----------------------------------------------------------------------
408 : ! We will handle errors for this routine
409 :
410 0 : call PIO_seterrorhandling(ncid, PIO_BCAST_ERROR, err_handling)
411 0 : ret = PIO_inq_varid(ncid, trim(varname), varid)
412 0 : found = (ret == PIO_NOERR)
413 0 : call PIO_seterrorhandling(ncid, err_handling)
414 :
415 0 : end subroutine cam_pio_find_var
416 :
417 :
418 : !-----------------------------------------------------------------------
419 : !
420 : ! cam_pio_check_var: Make sure var exists and retrieve properties
421 : !
422 : !-----------------------------------------------------------------------
423 43008 : subroutine cam_pio_check_var(ncid, varname, varid, ndims, dimids, dimlens, &
424 43008 : readvar, dimnames)
425 : use pio, only: PIO_inq_varid, PIO_NOERR
426 : use pio, only: PIO_seterrorhandling, PIO_BCAST_ERROR
427 : use shr_sys_mod, only: shr_sys_flush ! Standardized system subroutines
428 :
429 : ! Dummy arguments
430 : type(file_desc_t), intent(inout) :: ncid
431 : character(len=*), intent(in) :: varname
432 : type(var_desc_t), intent(out) :: varid
433 : integer, intent(out) :: ndims
434 : integer, intent(out) :: dimids(:)
435 : integer, intent(out) :: dimlens(:)
436 : logical, intent(out) :: readvar
437 : character(len=*), optional, intent(out) :: dimnames(:)
438 :
439 : ! Local variables
440 : integer :: ret ! PIO return value
441 : integer :: err_handling
442 :
443 : !-----------------------------------------------------------------------
444 : ! We will handle errors for this routine
445 43008 : call pio_seterrorhandling(ncid, PIO_BCAST_ERROR, err_handling)
446 :
447 301056 : dimids = -1
448 43008 : ndims = 0
449 301056 : dimlens = 0
450 43008 : ret = PIO_inq_varid(ncid, trim(varname), varid)
451 43008 : if (ret /= PIO_NOERR) then
452 6144 : readvar = .false.
453 6144 : if (masterproc) then
454 8 : write(iulog,*)'CAM_PIO_CHECK_VAR INFO: variable ',trim(varname),' is not on file'
455 8 : call shr_sys_flush(iulog)
456 : end if
457 : else
458 36864 : readvar = .true.
459 : call cam_pio_var_info(ncid, varid, ndims, dimids, dimlens, &
460 46080 : dimnames=dimnames, varname=varname)
461 : end if
462 43008 : call pio_seterrorhandling(ncid, err_handling)
463 :
464 43008 : end subroutine cam_pio_check_var
465 :
466 1536 : subroutine init_pio_subsystem()
467 : use shr_pio_mod, only: shr_pio_getiosys, shr_pio_getiotype
468 : use shr_pio_mod, only: shr_pio_getioformat
469 : use cam_instance, only: atm_id
470 :
471 1536 : pio_subsystem => shr_pio_getiosys(atm_id)
472 1536 : pio_iotype = shr_pio_getiotype(atm_id)
473 1536 : pio_ioformat = shr_pio_getioformat(atm_id)
474 :
475 1536 : if (masterproc) then
476 2 : write(iulog,*)' '
477 2 : write(iulog,*)'Initialize PIO subsystem:'
478 2 : write(iulog,*)' iotype = ', pio_iotype
479 2 : write(iulog,*)' ioformat = ', pio_ioformat
480 : end if
481 :
482 1536 : end subroutine init_pio_subsystem
483 :
484 : ! cam_pio_get_decomp: retrieve or create a PIO decomposition for the field
485 : ! described by ldims and dtype where dims is the field's
486 : ! local shape.
487 : ! fdims is the shape of the field in a NetCDF file.
488 : ! map describes the mapping of the distributed dimensions
489 : ! field_dist_in is used if the dimensions of the
490 : ! field array are not in map order
491 : ! file_dist_in is used if the dimensions of the
492 : ! field on file are not in map order
493 : !
494 18239232 : subroutine cam_pio_get_decomp(iodesc, ldims, fdims, dtype, map, &
495 18239232 : field_dist_in, file_dist_in, permute)
496 : use pio, only: pio_offset_kind
497 : use cam_abortutils, only: endrun
498 : use cam_map_utils, only: cam_filemap_t
499 :
500 : ! Dummy arguments
501 : type(io_desc_t), pointer :: iodesc ! intent(out)
502 : integer, intent(in) :: ldims(:) ! Local array
503 : integer, intent(in) :: fdims(:) ! File dims
504 : integer, intent(in) :: dtype
505 : type(cam_filemap_t), target, intent(in) :: map
506 : integer, optional, intent(in) :: field_dist_in(:)
507 : integer, optional, intent(in) :: file_dist_in(:)
508 : integer, optional, intent(in) :: permute(:)
509 :
510 : ! Local variables
511 : logical :: found
512 18239232 : integer(PIO_OFFSET_KIND), pointer :: dof(:)
513 : type(iodesc_list), pointer :: iodesc_p
514 : character(len=errormsg_str_len) :: errormsg
515 :
516 18239232 : call t_startf('get_decomp')
517 :
518 18239232 : nullify(iodesc_p)
519 18239232 : nullify(dof)
520 36478464 : call find_iodesc(ldims, fdims, dtype, map, iodesc_p, found, perm=permute)
521 :
522 18239232 : if (.not. found) then
523 : ! Create a new iodesc
524 22272 : if(masterproc) then
525 29 : write(iulog,*) 'Creating new decomp: ', iodesc_p%tag
526 : end if
527 :
528 22272 : call t_startf('get_filemap')
529 : call map%get_filemap(ldims, fdims, dof, &
530 82944 : src_in=field_dist_in, dest_in=file_dist_in, permutation_in=permute)
531 22272 : call t_stopf('get_filemap')
532 59136 : if (any(fdims == 0)) then
533 : ! Quick sanity check
534 0 : write(errormsg, *) 'bad fdims, ',fdims
535 0 : call endrun('cam_pio_get_decomp: '//errormsg)
536 : end if
537 22272 : if (associated(iodesc_p%iodesc)) then
538 : ! Quick sanity check
539 0 : call endrun('cam_pio_get_decomp: iodesc already allocated')
540 : end if
541 22272 : allocate(iodesc_p%iodesc)
542 22272 : call t_startf('newdecomp')
543 22272 : call cam_pio_newdecomp(iodesc_p%iodesc, fdims, dof, dtype)
544 22272 : call t_stopf('newdecomp')
545 :
546 22272 : deallocate(dof)
547 : nullify(dof)
548 : end if
549 : ! At this point, we should have a decomp, assign iodesc
550 18239232 : iodesc => iodesc_p%iodesc
551 18239232 : nullify(iodesc_p)
552 :
553 18239232 : call t_stopf('get_decomp')
554 :
555 18239232 : end subroutine cam_pio_get_decomp
556 :
557 519936 : subroutine cam_pio_newdecomp(iodesc, dims, dof, dtype)
558 : use pio, only: pio_initdecomp, pio_offset_kind, pio_iotype_pnetcdf
559 :
560 : type(io_desc_t), pointer :: iodesc
561 : integer, intent(in) :: dims(:)
562 : integer(kind=PIO_OFFSET_KIND), intent(in) :: dof(:)
563 : integer, intent(in) :: dtype
564 :
565 519936 : call pio_initdecomp(pio_subsystem, dtype, dims, dof, iodesc)
566 :
567 519936 : end subroutine cam_pio_newdecomp
568 :
569 18239232 : subroutine find_iodesc(ldimlens, fdimlens, dtype, map, iodesc_p, found, perm)
570 : use cam_abortutils, only: endrun
571 : use cam_map_utils, only: cam_filemap_t
572 :
573 : ! Dummy arguments
574 : integer, intent(in) :: ldimlens(:)
575 : integer, intent(in) :: fdimlens(:)
576 : integer, intent(in) :: dtype
577 : type(cam_filemap_t), intent(in) :: map
578 : type(iodesc_list), pointer :: iodesc_p
579 : logical, intent(out) :: found
580 : integer, optional, intent(in) :: perm(:)
581 :
582 : ! Local variables
583 : type(iodesc_list), pointer :: curr, prev
584 : integer :: i
585 : integer :: lcnt
586 : integer :: fcnt
587 : integer :: mapind
588 : integer :: nperm
589 : character(len=128) :: form
590 : character(len=tag_len) :: tag
591 : character(len=*), parameter :: formc = 'i0,"(i0,""!""),""!"",",'
592 : character(len=*), parameter :: forme = '"""d"",i0,""!i"",i0,""!"""'
593 : character(len=*), parameter :: form2 = '("(",'//formc//formc//forme//',")")'
594 : character(len=*), parameter :: form3 = '("(",'//formc//formc//formc//forme//',")")'
595 :
596 18239232 : found = .false.
597 18239232 : curr => iodesc_list_top
598 :
599 : ! Retrieve the (hopefully) unique tag for this iodesc
600 : ! If a decomp was created using an earlier version of the map (hey, that
601 : ! might happen), we won't find it using this search because the current
602 : ! index is part of the search tag
603 18239232 : mapind = map%get_index()
604 18239232 : lcnt = size(ldimlens)
605 18239232 : fcnt = size(fdimlens)
606 18239232 : if (present(perm)) then
607 0 : if (size(perm) /= lcnt) then
608 0 : write(form, '(i0,a,i0)') size(perm), ', should be ', lcnt
609 0 : call endrun('FIND_IODESC: perm has wrong size, '//form)
610 : end if
611 0 : nperm = lcnt
612 : else
613 18239232 : nperm = 0
614 : end if
615 : if (present(perm)) then
616 0 : write(form, form3) lcnt, fcnt, nperm
617 0 : write(tag, form) (ldimlens(i),i=1,lcnt), (fdimlens(i),i=1,fcnt), (perm(i),i=1,lcnt), dtype, mapind
618 : else
619 18239232 : write(form, form2) lcnt, fcnt
620 95043072 : write(tag, form) (ldimlens(i),i=1,lcnt), (fdimlens(i),i=1,fcnt), dtype, mapind
621 : end if
622 :
623 86192640 : do while(associated(curr) .and. (.not. found))
624 86192640 : if(trim(tag) == trim(curr%tag)) then
625 18216960 : found = .true.
626 18216960 : iodesc_p => curr
627 : else
628 49736448 : prev => curr
629 49736448 : curr => curr%next
630 : end if
631 : end do
632 18239232 : if(.not. found) then
633 : ! We didn't find a match, make sure there is an unused iodesc_list
634 : ! object at the end of the list for the new decomp to be stored
635 22272 : curr => prev
636 22272 : if(associated(curr%iodesc)) then
637 19200 : allocate(curr%next)
638 19200 : curr => curr%next
639 19200 : nullify(curr%iodesc) ! Should already be null but . . .
640 19200 : nullify(curr%next) ! Should already be null but . . .
641 : end if
642 : ! This should be an unused object at the end of the list
643 22272 : curr%tag = tag
644 22272 : iodesc_p => curr
645 : end if
646 : ! if(masterproc) write(iulog,*) 'Using decomp: ',curr%tag
647 :
648 18239232 : end subroutine find_iodesc
649 :
650 :
651 : ! cam_pio_def_dim: Define a NetCDF dimension using the PIO interface
652 3597312 : subroutine cam_pio_def_dim(File, name, size, dimid, existOK)
653 : use cam_abortutils, only: endrun
654 : use pio, only: pio_inq_dimid, pio_def_dim, pio_inq_dimlen, PIO_NOERR
655 : use pio, only: PIO_seterrorhandling, PIO_BCAST_ERROR
656 :
657 : ! Dummy arguments
658 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
659 : character(len=*), intent(in) :: name ! Dimension name
660 : integer, intent(in) :: size ! Dimension length
661 : integer, intent(out) :: dimid ! NetCDF dimension ID
662 : logical, optional, intent(in) :: existOK ! OK if dim defined
663 :
664 : ! Local variables
665 : logical :: ok_if_dim_exists
666 : integer :: ierr
667 : integer :: err_handling
668 : integer :: dimlen
669 : character(len=errormsg_str_len) :: errormsg
670 : character(len=*), parameter :: subname = 'cam_pio_def_dim'
671 :
672 3597312 : if (present(existOK)) then
673 3105792 : ok_if_dim_exists = existOK
674 : else
675 : ok_if_dim_exists = .false.
676 : end if
677 :
678 : ! We will handle errors for this routine
679 3597312 : call pio_seterrorhandling(File, PIO_BCAST_ERROR, err_handling)
680 :
681 3597312 : ierr = pio_inq_dimid(File, trim(name), dimid)
682 3597312 : if (ierr == PIO_NOERR) then
683 1362432 : if (.not. ok_if_dim_exists) then
684 0 : write(errormsg, *) ': A dimension already exists for ', trim(name)
685 0 : call endrun(trim(subname)//errormsg)
686 : else
687 : ! It is OK for the dimension to exist but it better have the same size
688 1362432 : ierr = pio_inq_dimlen(File, dimid, dimlen)
689 1362432 : if (ierr /= PIO_NOERR) then
690 0 : write(errormsg, '(2a,i0,2a)') trim(subname), ': Error ', ierr, &
691 0 : ' finding dimension length for ', trim(name)
692 0 : call endrun(errormsg)
693 1362432 : else if (dimlen /= size) then
694 0 : write(errormsg, '(3a,2(i0,a))') ': Size mismatch for dimension, ', &
695 0 : trim(name), ': ', dimlen, ' (current), ', size, ' (desired)'
696 0 : call endrun(trim(subname)//errormsg)
697 : ! No else, existing dimension is OK
698 : end if
699 : end if
700 : else
701 : ! inq_dimid returned an error, define the dimension
702 2234880 : ierr = pio_def_dim(File, trim(name), size, dimid)
703 2234880 : call cam_pio_handle_error(ierr, trim(subname)//': Unable to define dimension '//trim(name))
704 : end if
705 :
706 : ! Back to whatever error handling was running before this routine
707 3597312 : call pio_seterrorhandling(File, err_handling)
708 :
709 3597312 : end subroutine cam_pio_def_dim
710 :
711 : ! cam_pio_def_var_0d: Define a NetCDF variable using the PIO interface
712 1536 : subroutine cam_pio_def_var_0d(File, name, dtype, vardesc, existOK)
713 :
714 : ! Dummy arguments
715 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
716 : character(len=*), intent(in) :: name ! Variable name
717 : integer, intent(in) :: dtype ! e.g., pio_int
718 : type(var_desc_t), intent(inout) :: vardesc ! Variable descriptor
719 : logical, optional, intent(in) :: existOK ! OK if var defined
720 :
721 : ! Local variables
722 : integer :: dimids(0)
723 :
724 1536 : call cam_pio_def_var(File, trim(name), dtype, dimids, vardesc, existOK)
725 1536 : end subroutine cam_pio_def_var_0d
726 :
727 : ! cam_pio_def_var_md: Define a NetCDF variable using the PIO interface
728 21321216 : subroutine cam_pio_def_var_md(File, name, dtype, dimids, vardesc, existOK)
729 : use cam_abortutils, only: endrun
730 : use pio, only: pio_inq_varid, pio_def_var, PIO_NOERR
731 : use pio, only: PIO_seterrorhandling, PIO_BCAST_ERROR
732 :
733 : ! Dummy arguments
734 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
735 : character(len=*), intent(in) :: name ! Variable name
736 : integer, intent(in) :: dtype ! e.g., pio_int
737 : integer, intent(in) :: dimids(:) ! NetCDF dim IDs
738 : type(var_desc_t), intent(inout) :: vardesc ! Var descriptor
739 : logical, optional, intent(in) :: existOK ! OK if var defined
740 :
741 : ! Local variables
742 : integer :: ierr
743 : integer :: err_handling
744 : logical :: ok_if_var_exists
745 : character(len=errormsg_str_len) :: errormsg
746 : character(len=*), parameter :: subname = 'cam_pio_def_var'
747 :
748 21321216 : if (present(existOK)) then
749 3380736 : ok_if_var_exists = existOK
750 : else
751 : ok_if_var_exists = .false.
752 : end if
753 :
754 : ! We will handle errors for this routine
755 21321216 : call pio_seterrorhandling(File, PIO_BCAST_ERROR, err_handling)
756 :
757 : ! Check to see if the variable already exists in the file
758 21321216 : ierr = pio_inq_varid(File, name, vardesc)
759 21321216 : if (ierr == PIO_NOERR) then
760 0 : if (.not. ok_if_var_exists) then
761 0 : write(errormsg, *) ': A variable already exists for ', trim(name)
762 0 : call endrun(trim(subname)//errormsg)
763 : end if
764 : else
765 : ! OK to define the variable
766 21321216 : if (size(dimids) > 0) then
767 21319680 : ierr = pio_def_var(File, trim(name), dtype, dimids, vardesc)
768 : else
769 1536 : ierr = pio_def_var(File, trim(name), dtype, vardesc)
770 : end if
771 21321216 : call cam_pio_handle_error(ierr, trim(subname)//': Unable to define variable '//trim(name))
772 : end if
773 :
774 : ! Back to whatever error handling was running before this routine
775 21321216 : call pio_seterrorhandling(File, err_handling)
776 :
777 21321216 : end subroutine cam_pio_def_var_md
778 :
779 0 : subroutine cam_pio_get_var_2d_r8(varname, File, field, start, kount, found)
780 : use cam_abortutils, only: endrun
781 : use pio, only: file_desc_t, var_desc_t, pio_get_var, PIO_MAX_NAME
782 : use pio, only: pio_inq_dimname
783 :
784 : ! Dummy arguments
785 : character(len=*), intent(in) :: varname
786 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
787 : real(r8), intent(inout) :: field(:,:)
788 : integer, optional, intent(in) :: start(2)
789 : integer, optional, intent(in) :: kount(2)
790 : logical, optional, intent(out) :: found
791 :
792 : ! Local variables
793 : character(len=*), parameter :: subname = 'cam_pio_get_var_2d_r8'
794 : character(len=PIO_MAX_NAME) :: tmpname
795 : type(var_desc_t) :: varid ! Var descriptor
796 : integer :: ierr
797 : integer :: strt(3)
798 : integer :: cnt(3)
799 : integer :: ndims
800 : integer :: dimids(3)
801 : logical :: exists
802 : character(len=PIO_MAX_NAME) :: filedims(4)
803 :
804 0 : if ( (present(start) .and. (.not. present(kount))) .or. &
805 : (present(kount) .and. (.not. present(start)))) then
806 0 : call endrun(trim(subname)//': start and kount must both be present')
807 : end if
808 :
809 0 : call cam_pio_find_var(File, trim(varname), varid, exists)
810 0 : if (present(found)) then
811 0 : found = exists
812 0 : else if (.not. exists) then
813 0 : call endrun(trim(subname)//': '//trim(varname)//' not found')
814 : end if
815 0 : if (exists) then
816 0 : call cam_pio_var_info(File, varid, ndims, dimids, cnt, dimnames=filedims, varname=varname)
817 :
818 0 : if (present(start)) then
819 : ! start and kount override other options and are not error checked
820 0 : strt(1:2) = start(1:2)
821 0 : strt(3) = 1
822 0 : cnt(1:2) = kount(1:2)
823 0 : cnt(3) = 1
824 : else
825 0 : strt = 1 ! cnt set by cam_pio_var_info
826 0 : exists = use_scam_limits(File, strt, cnt,filedims)
827 : end if
828 0 : if (ndims == 3) then
829 0 : ierr = pio_inq_dimname(File, dimids(3), tmpname)
830 0 : if (trim(tmpname) /= 'time') then
831 0 : call endrun(trim(subname)//': dimension mismatch for '//trim(varname))
832 : else
833 0 : ierr = pio_get_var(File, varid, strt, cnt, field)
834 : end if
835 0 : else if (ndims == 2) then
836 0 : ierr = pio_get_var(File, varid, strt, cnt, field)
837 0 : else if (ndims == 1) then
838 0 : ierr = pio_get_var(File, varid, strt(1:1), cnt(1:1), field(:,1))
839 : else
840 0 : call endrun(trim(subname)//': Incorrect variable rank')
841 : end if
842 : end if
843 :
844 0 : end subroutine cam_pio_get_var_2d_r8
845 :
846 0 : subroutine cam_pio_get_var_2d_r8_perm(varname, File, arraydims, field, &
847 : start, kount, found)
848 : use cam_abortutils, only: endrun
849 : use pio, only: file_desc_t, var_desc_t, pio_get_var, PIO_MAX_NAME
850 :
851 : ! Dummy arguments
852 : character(len=*), intent(in) :: varname
853 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
854 : character(len=*), intent(in) :: arraydims(2)
855 : real(r8), intent(inout) :: field(:,:)
856 : integer, optional, intent(in) :: start(2)
857 : integer, optional, intent(in) :: kount(2)
858 : logical, optional, intent(out) :: found
859 :
860 : ! Local variables
861 : character(len=*), parameter :: subname = 'cam_pio_get_var_2d_r8_perm'
862 : type(var_desc_t) :: varid ! Var descriptor
863 : integer :: ierr
864 : integer :: i, j, ind(2)
865 : integer :: strt(3)
866 : integer :: cnt(3)
867 : integer :: ndims
868 : integer :: dimids(3)
869 : integer :: perm(2)
870 : logical :: isperm
871 : logical :: exists
872 0 : real(r8), allocatable :: tmp_fld(:,:)
873 : character(len=PIO_MAX_NAME) :: filedims(3)
874 :
875 0 : if ( (present(start) .and. (.not. present(kount))) .or. &
876 : (present(kount) .and. (.not. present(start)))) then
877 0 : call endrun(trim(subname)//': start and kount must both be present')
878 : end if
879 :
880 0 : call cam_pio_find_var(File, trim(varname), varid, exists)
881 :
882 0 : if (present(found)) then
883 0 : found = exists
884 0 : else if (.not. exists) then
885 0 : call endrun(trim(subname)//': '//trim(varname)//' not found')
886 : end if
887 0 : if (exists) then
888 : call cam_pio_var_info(File, varid, ndims, dimids, cnt, &
889 0 : dimnames=filedims, varname=varname)
890 :
891 0 : if (present(start)) then
892 : ! start and kount override other options and are not error checked
893 0 : strt(1:2) = start
894 0 : strt(3) = 1
895 0 : cnt(1:2) = kount
896 : else
897 0 : strt = 1 ! cnt set by cam_pio_var_info
898 0 : exists = use_scam_limits(File, strt, cnt,filedims)
899 : end if
900 0 : if ( ((ndims == 2) .and. (trim(filedims(2)) /= 'time')) .or. &
901 : ((ndims == 3) .and. (trim(filedims(3)) == 'time'))) then
902 0 : call calc_permutation(filedims(1:2), arraydims, perm, isperm)
903 0 : if (isperm) then
904 0 : allocate(tmp_fld(cnt(1), cnt(2)))
905 0 : ierr = pio_get_var(File, varid, strt(1:ndims), cnt(1:ndims), tmp_fld)
906 0 : do j = 1, cnt(2)
907 0 : ind(2) = j
908 0 : do i = 1, cnt(1)
909 0 : ind(1) = i
910 0 : field(ind(perm(1)), ind(perm(2))) = tmp_fld(i, j)
911 : end do
912 : end do
913 : else
914 0 : ierr = pio_get_var(File, varid, strt(1:ndims), cnt(1:ndims), field)
915 : end if
916 : else
917 0 : call endrun(trim(subname)//': Incorrect variable rank')
918 : end if
919 : end if
920 :
921 0 : end subroutine cam_pio_get_var_2d_r8_perm
922 :
923 0 : subroutine cam_pio_get_var_3d_r8(varname, File, field, start, kount, found)
924 : use cam_abortutils, only: endrun
925 : use pio, only: file_desc_t, var_desc_t, pio_get_var, PIO_MAX_NAME
926 : use pio, only: pio_inq_dimname
927 :
928 : ! Dummy arguments
929 : character(len=*), intent(in) :: varname
930 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
931 : real(r8), intent(inout) :: field(:,:,:)
932 : integer, optional, intent(in) :: start(3)
933 : integer, optional, intent(in) :: kount(3)
934 : logical, optional, intent(out) :: found
935 :
936 : ! Local variables
937 : character(len=*), parameter :: subname = 'cam_pio_get_var_3d_r8'
938 : character(len=PIO_MAX_NAME) :: tmpname
939 : type(var_desc_t) :: varid ! Var descriptor
940 : integer :: ierr
941 : integer :: strt(4)
942 : integer :: cnt(4)
943 : integer :: ndims
944 : integer :: dimids(4)
945 : logical :: exists
946 : character(len=PIO_MAX_NAME) :: filedims(4)
947 :
948 0 : if ( (present(start) .and. (.not. present(kount))) .or. &
949 : (present(kount) .and. (.not. present(start)))) then
950 0 : call endrun(trim(subname)//': start and kount must both be present')
951 : end if
952 :
953 0 : call cam_pio_find_var(File, trim(varname), varid, exists)
954 :
955 0 : if (present(found)) then
956 0 : found = exists
957 0 : else if (.not. exists) then
958 0 : call endrun(trim(subname)//': '//trim(varname)//' not found')
959 : end if
960 0 : if (exists) then
961 0 : call cam_pio_var_info(File, varid, ndims, dimids, cnt,dimnames=filedims, varname=varname)
962 :
963 0 : if (present(start)) then
964 : ! start and kount override other options and are not error checked
965 0 : strt(1:3) = start(1:3)
966 0 : strt(4) = 1
967 0 : cnt(1:3) = kount(1:3)
968 0 : cnt(4) = 1
969 : else
970 0 : strt = 1 ! cnt set by cam_pio_var_info
971 0 : exists = use_scam_limits(File, strt, cnt,filedims)
972 : end if
973 :
974 0 : if (ndims == 4) then
975 0 : ierr = pio_inq_dimname(File, dimids(4), tmpname)
976 0 : if (trim(tmpname) /= 'time') then
977 0 : call endrun(trim(subname)//': dimension mismatch for '//trim(varname))
978 : else
979 0 : ierr = pio_get_var(File, varid, strt, cnt, field)
980 : end if
981 0 : else if (ndims == 3) then
982 0 : ierr = pio_get_var(File, varid, strt, cnt, field)
983 0 : else if (ndims == 2) then
984 0 : ierr = pio_get_var(File, varid, strt(1:ndims), cnt(1:ndims), field(:,:,1))
985 : else
986 0 : call endrun(trim(subname)//': Incorrect variable rank')
987 : end if
988 : end if
989 :
990 0 : end subroutine cam_pio_get_var_3d_r8
991 :
992 0 : subroutine cam_pio_get_var_3d_r8_perm(varname, File, arraydims, field, &
993 : start, kount, found)
994 : use cam_abortutils, only: endrun
995 : use pio, only: file_desc_t, var_desc_t, pio_get_var, PIO_MAX_NAME
996 :
997 : ! Dummy arguments
998 : character(len=*), intent(in) :: varname
999 : type(file_desc_t), intent(inout) :: File ! PIO file Handle
1000 : character(len=*), intent(in) :: arraydims(3)
1001 : real(r8), intent(inout) :: field(:,:,:)
1002 : integer, optional, intent(in) :: start(3)
1003 : integer, optional, intent(in) :: kount(3)
1004 : logical, optional, intent(out) :: found
1005 :
1006 : ! Local variables
1007 : character(len=*), parameter :: subname = 'cam_pio_get_var_3d_r8_perm'
1008 : type(var_desc_t) :: varid ! Var descriptor
1009 : integer :: ierr
1010 : integer :: i, j, k, ind(3)
1011 : integer :: strt(4)
1012 : integer :: cnt(4)
1013 : integer :: ndims
1014 : integer :: dimids(4)
1015 : integer :: perm(3)
1016 : logical :: exists
1017 : logical :: isperm
1018 0 : real(r8), allocatable :: tmp_fld(:,:,:)
1019 : character(len=PIO_MAX_NAME) :: filedims(4)
1020 :
1021 0 : if ( (present(start) .and. (.not. present(kount))) .or. &
1022 : (present(kount) .and. (.not. present(start)))) then
1023 0 : call endrun(trim(subname)//': start and kount must both be present')
1024 : end if
1025 :
1026 0 : call cam_pio_find_var(File, trim(varname), varid, exists)
1027 :
1028 0 : if (present(found)) then
1029 0 : found = exists
1030 0 : else if (.not. exists) then
1031 0 : call endrun(trim(subname)//': '//trim(varname)//' not found')
1032 : end if
1033 0 : if (exists) then
1034 : call cam_pio_var_info(File, varid, ndims, dimids, cnt, &
1035 0 : dimnames=filedims, varname=varname)
1036 :
1037 0 : if (present(start)) then
1038 : ! start and kount override other options and are not error checked
1039 0 : strt(1:3) = start
1040 0 : strt(4) = 1
1041 0 : cnt(1:3) = kount
1042 : else
1043 0 : strt = 1 ! cnt set by cam_pio_var_info
1044 0 : exists = use_scam_limits(File, strt, cnt,filedims)
1045 : end if
1046 :
1047 0 : if ( ((ndims == 3) .and. (trim(filedims(3)) /= 'time')) .or. &
1048 : ((ndims == 4) .and. (trim(filedims(4)) == 'time'))) then
1049 0 : call calc_permutation(filedims(1:3), arraydims, perm, isperm)
1050 0 : if (isperm) then
1051 0 : allocate(tmp_fld(cnt(1), cnt(2), cnt(3)))
1052 0 : ierr = pio_get_var(File, varid, strt(1:ndims), cnt(1:ndims), tmp_fld)
1053 0 : do k = 1, cnt(3)
1054 0 : ind(3) = k
1055 0 : do j = 1, cnt(2)
1056 0 : ind(2) = j
1057 0 : do i = 1, cnt(1)
1058 0 : ind(1) = i
1059 0 : field(ind(perm(1)), ind(perm(2)), ind(perm(3))) = tmp_fld(i, j, k)
1060 : end do
1061 : end do
1062 : end do
1063 : else
1064 0 : ierr = pio_get_var(File, varid, strt(1:ndims), cnt(1:ndims), field)
1065 : end if
1066 : else
1067 0 : call endrun(trim(subname)//': Incorrect variable rank')
1068 : end if
1069 : end if
1070 :
1071 0 : end subroutine cam_pio_get_var_3d_r8_perm
1072 :
1073 : ! clean_iodesc_list: Deallocate all entries in the iodesc list
1074 1536 : subroutine clean_iodesc_list()
1075 : type(iodesc_list), pointer :: this, prev
1076 :
1077 1536 : if(associated(iodesc_list_top%iodesc)) then
1078 : ! iodesc_list_top is not allocated so leave it (just empty)
1079 1536 : this => iodesc_list_top
1080 1536 : iodesc_list_top%tag = ''
1081 1536 : call pio_freedecomp(pio_subsystem, this%iodesc)
1082 1536 : deallocate(this%iodesc)
1083 : nullify(this%iodesc)
1084 1536 : this => this%next
1085 1536 : nullify(iodesc_list_top%next)
1086 :
1087 : ! All the other list items were allocated, blow them away
1088 6912 : do while(associated(this))
1089 5376 : call pio_freedecomp(pio_subsystem, this%iodesc)
1090 5376 : deallocate(this%iodesc)
1091 5376 : prev => this
1092 5376 : this => this%next
1093 5376 : deallocate(prev)
1094 : end do
1095 : end if
1096 1536 : end subroutine clean_iodesc_list
1097 :
1098 248832 : subroutine cam_pio_createfile(file, fname, mode_in)
1099 : use pio, only : pio_createfile, file_desc_t, pio_noerr, pio_clobber, pio_iotask_rank
1100 : use cam_abortutils, only : endrun
1101 :
1102 : ! Dummy arguments
1103 : type(file_desc_t), intent(inout) :: file
1104 : character(len=*), intent(in) :: fname
1105 : integer, optional, intent(in) :: mode_in
1106 :
1107 : ! Local variables
1108 : integer :: ierr
1109 : integer :: mode
1110 :
1111 248832 : mode = ior(PIO_CLOBBER, pio_ioformat)
1112 248832 : if (present(mode_in)) then
1113 248832 : mode = ior(mode, mode_in)
1114 : end if
1115 :
1116 248832 : ierr = pio_createfile(pio_subsystem, file, pio_iotype, fname, mode)
1117 :
1118 248832 : if(ierr /= PIO_NOERR) then
1119 0 : call endrun('Failed to open file,'//trim(fname)//', to write')
1120 248832 : else if(pio_iotask_rank(pio_subsystem) == 0) then
1121 324 : write(iulog, *) 'Opened file ', trim(fname), ' to write', file%fh
1122 : end if
1123 :
1124 248832 : end subroutine cam_pio_createfile
1125 :
1126 46848 : subroutine cam_pio_openfile(file, fname, mode)
1127 : use pio, only: pio_openfile, file_desc_t, pio_noerr, pio_iotask_rank
1128 : use cam_abortutils, only: endrun
1129 :
1130 : type(file_desc_t), intent(inout), target :: file
1131 : character(len=*), intent(in) :: fname
1132 : integer, intent(in) :: mode
1133 :
1134 : integer :: ierr
1135 :
1136 46848 : ierr = pio_openfile(pio_subsystem, file, pio_iotype, fname, mode)
1137 :
1138 46848 : if(ierr/= PIO_NOERR) then
1139 0 : call endrun('Failed to open '//trim(fname)//' to read')
1140 46848 : else if(pio_iotask_rank(pio_subsystem) == 0) then
1141 61 : write(iulog,*) 'Opened existing file ', trim(fname), file%fh
1142 : end if
1143 :
1144 46848 : end subroutine cam_pio_openfile
1145 :
1146 252672 : subroutine cam_pio_closefile(file)
1147 :
1148 : use pio, only : pio_closefile, file_desc_t
1149 :
1150 : type(file_desc_t), intent(inout), target :: file
1151 :
1152 252672 : call pio_closefile(file)
1153 :
1154 252672 : end subroutine cam_pio_closefile
1155 :
1156 0 : logical function cam_pio_fileexists(fname)
1157 : use pio, only: pio_openfile, file_desc_t, pio_noerr, PIO_NOWRITE
1158 : use pio, only: pio_seterrorhandling, PIO_BCAST_ERROR
1159 : use pio, only : pio_closefile
1160 :
1161 : character(len=*), intent(in) :: fname
1162 :
1163 : type(file_desc_t) :: file
1164 : integer :: ierr
1165 : integer :: err_handling
1166 :
1167 : ! We will handle errors for this routine
1168 :
1169 0 : call pio_seterrorhandling(pio_subsystem, PIO_BCAST_ERROR, err_handling)
1170 :
1171 0 : ierr = pio_openfile(pio_subsystem, file, pio_iotype, fname, PIO_NOWRITE)
1172 0 : cam_pio_fileexists = (ierr == PIO_NOERR)
1173 0 : if (cam_pio_fileexists) then
1174 0 : call pio_closefile(file)
1175 : end if
1176 :
1177 : ! Back to whatever error handling was running before this routine
1178 0 : call pio_seterrorhandling(pio_subsystem, err_handling)
1179 :
1180 0 : end function cam_pio_fileexists
1181 :
1182 1536 : integer function cam_pio_set_fill(File, fillmode, old_mode) result(ierr)
1183 : #ifdef PIO2
1184 : use pio, only: PIO_FILL, pio_set_fill
1185 : #endif
1186 : ! Dummy arguments
1187 : type(File_desc_t), intent(in) :: File
1188 : integer, optional, intent(in) :: fillmode
1189 : integer, optional, intent(out) :: old_mode
1190 : ! Local variables
1191 : integer :: oldfill
1192 : integer :: fillval
1193 :
1194 : #ifdef PIO2
1195 1536 : if (present(fillmode)) then
1196 0 : fillval = fillmode
1197 : else
1198 1536 : fillval = PIO_FILL
1199 : end if
1200 1536 : ierr = pio_set_fill(File, fillval, oldfill)
1201 1536 : if (present(old_mode)) then
1202 0 : old_mode = oldfill
1203 : end if
1204 : #else
1205 : ierr = 0
1206 : if (present(old_mode)) then
1207 : old_mode = 0
1208 : end if
1209 : #endif
1210 1536 : end function cam_pio_set_fill
1211 :
1212 0 : integer function inq_var_fill_i4(File, vdesc, fillvalue, no_fill) result(ierr)
1213 : #ifdef PIO2
1214 : use pio, only: pio_inq_var_fill
1215 : #endif
1216 : use pio, only: PIO_NOERR
1217 :
1218 : ! Dummy arguments
1219 : type(File_desc_t), intent(in) :: File
1220 : type(var_desc_t), intent(in) :: vdesc
1221 : ! fillvalue needs to not be optional to avoid ambiguity
1222 : integer, target, intent(out) :: fillvalue
1223 : integer, optional, intent(out) :: no_fill
1224 : ! Local variable
1225 : integer :: no_fill_use
1226 :
1227 : #ifdef PIO2
1228 0 : ierr = pio_inq_var_fill(File, vdesc, no_fill_use, fillvalue)
1229 0 : if (present(no_fill)) then
1230 0 : no_fill = no_fill_use
1231 : end if
1232 : #else
1233 : ierr = PIO_NOERR
1234 : fillvalue = 0
1235 : #endif
1236 :
1237 0 : end function inq_var_fill_i4
1238 :
1239 0 : integer function inq_var_fill_r4(File, vdesc, fillvalue, no_fill) result(ierr)
1240 : #ifdef PIO2
1241 : use pio, only: pio_inq_var_fill
1242 : #endif
1243 : use pio, only: PIO_NOERR
1244 :
1245 : ! Dummy arguments
1246 : type(File_desc_t), intent(in) :: File
1247 : type(var_desc_t), intent(in) :: vdesc
1248 : ! fillvalue needs to not be optional to avoid ambiguity
1249 : real(r4), target, intent(out) :: fillvalue
1250 : integer, optional, intent(out) :: no_fill
1251 : ! Local variable
1252 : integer :: no_fill_use
1253 :
1254 : #ifdef PIO2
1255 0 : ierr = pio_inq_var_fill(File, vdesc, no_fill_use, fillvalue)
1256 0 : if (present(no_fill)) then
1257 0 : no_fill = no_fill_use
1258 : end if
1259 : #else
1260 : ierr = PIO_NOERR
1261 : fillvalue = 0.0_R4
1262 : #endif
1263 :
1264 0 : end function inq_var_fill_r4
1265 :
1266 16128 : integer function inq_var_fill_r8(File, vdesc, fillvalue, no_fill) result(ierr)
1267 : #ifdef PIO2
1268 : use pio, only: pio_inq_var_fill
1269 : #endif
1270 : use pio, only: PIO_NOERR
1271 :
1272 : ! Dummy arguments
1273 : type(File_desc_t), intent(in) :: File
1274 : type(var_desc_t), intent(in) :: vdesc
1275 : ! fillvalue needs to not be optional to avoid ambiguity
1276 : real(r8), target, intent(out) :: fillvalue
1277 : integer, optional, intent(out) :: no_fill
1278 : ! Local variable
1279 : integer :: no_fill_use
1280 :
1281 : #ifdef PIO2
1282 16128 : ierr = pio_inq_var_fill(File, vdesc, no_fill_use, fillvalue)
1283 16128 : if (present(no_fill)) then
1284 0 : no_fill = no_fill_use
1285 : end if
1286 : #else
1287 : ierr = PIO_NOERR
1288 : fillvalue = 0.0_R8
1289 : #endif
1290 :
1291 16128 : end function inq_var_fill_r8
1292 :
1293 0 : subroutine find_dump_filename(fieldname, filename)
1294 :
1295 : ! Dummy arguments
1296 : character(len=*), intent(in) :: fieldname
1297 : character(len=*), intent(inout) :: filename
1298 :
1299 : ! Local variable
1300 : integer :: fnum
1301 :
1302 : ! Find an unused filename for this variable
1303 0 : filename = trim(fieldname)//'_dump_1.nc'
1304 0 : fnum = 1
1305 0 : do while (cam_pio_fileexists(trim(filename)))
1306 0 : fnum = fnum + 1
1307 0 : write(filename, '(2a,i0,a)') trim(fieldname), '_dump_', fnum, '.nc'
1308 : end do
1309 0 : end subroutine find_dump_filename
1310 :
1311 0 : subroutine dump_field_2d_d(fieldname, dim1b, dim1e, dim2b, dim2e, field, &
1312 : compute_maxdim_in, fill_value)
1313 : use pio, only: pio_offset_kind
1314 : use pio, only: pio_double, pio_int, pio_write_darray
1315 : use pio, only: pio_put_att, pio_initdecomp, pio_enddef
1316 : use spmd_utils, only: iam, npes, mpi_max, mpi_integer, mpicom
1317 :
1318 : ! Dummy arguments
1319 : character(len=*), intent(in) :: fieldname
1320 : integer, intent(in) :: dim1b
1321 : integer, intent(in) :: dim1e
1322 : integer, intent(in) :: dim2b
1323 : integer, intent(in) :: dim2e
1324 : real(r8), target, intent(in) :: field(dim1b:dim1e,dim2b:dim2e)
1325 : logical, optional, intent(in) :: compute_maxdim_in
1326 : real(r8), optional, intent(in) :: fill_value
1327 :
1328 : ! Local variables
1329 : type(file_desc_t) :: file
1330 : type(var_desc_t) :: vdesc
1331 : type(var_desc_t) :: bnddesc
1332 : type(io_desc_t) :: iodesc
1333 : character(len=64) :: filename
1334 : real(r8) :: fillval
1335 0 : integer(PIO_OFFSET_KIND), allocatable :: ldof(:)
1336 : integer :: dimids(3)
1337 : integer :: bnddimid
1338 : integer :: bounds(4)
1339 : integer :: dimsizes(3)
1340 : integer :: ierr
1341 : integer :: i, j, m, lsize
1342 : logical :: compute_maxdim
1343 :
1344 : ! Find an unused filename for this variable
1345 0 : call find_dump_filename(fieldname, filename)
1346 :
1347 : ! Should we compute max dim sizes or assume they are all the same?
1348 0 : if (present(compute_maxdim_in)) then
1349 0 : compute_maxdim = compute_maxdim_in
1350 : else
1351 : compute_maxdim = .true.
1352 : end if
1353 :
1354 0 : if (present(fill_value)) then
1355 0 : fillval = fill_value
1356 : else
1357 0 : fillval = -900._r8
1358 : end if
1359 :
1360 : ! Open the file for writing
1361 0 : call cam_pio_createfile(file, trim(filename))
1362 :
1363 : ! Define dimensions
1364 0 : if (compute_maxdim) then
1365 : call MPI_allreduce((dim1e - dim1b + 1), dimsizes(1), 1, MPI_integer, &
1366 0 : mpi_max, mpicom, ierr)
1367 : call MPI_allreduce((dim2e - dim2b + 1), dimsizes(2), 1, MPI_integer, &
1368 0 : mpi_max, mpicom, ierr)
1369 : else
1370 0 : dimsizes(1) = dim1e - dim1b + 1
1371 0 : dimsizes(2) = dim2e - dim2b + 1
1372 : end if
1373 0 : dimsizes(3) = npes
1374 0 : do i = 1, size(dimids, 1)
1375 0 : write(filename, '(a,i0)') 'dim', i
1376 0 : call cam_pio_def_dim(file, trim(filename), dimsizes(i), dimids(i))
1377 : end do
1378 0 : call cam_pio_def_dim(file, 'bounds', size(bounds, 1), bnddimid)
1379 : ! Define the variables
1380 0 : call cam_pio_def_var(file, trim(fieldname), pio_double, dimids, vdesc)
1381 : call cam_pio_def_var(file, 'field_bounds', pio_int, &
1382 0 : (/ bnddimid, dimids(size(dimids, 1)) /), bnddesc)
1383 0 : if (present(fill_value)) then
1384 0 : ierr = pio_put_att(file, vdesc, '_FillValue', fill_value)
1385 : end if
1386 0 : ierr = pio_enddef(file)
1387 :
1388 : ! Compute the variable decomposition and write field
1389 0 : lsize = product(dimsizes(1:2))
1390 0 : allocate(ldof((dim2e - dim2b + 1) * (dim1e - dim1b + 1)))
1391 0 : m = 0
1392 0 : do j = dim2b, dim2e
1393 0 : do i = dim1b, dim1e
1394 0 : m = m + 1
1395 0 : ldof(m) = (iam * lsize) + (dimsizes(1)*(j - dim2b)) + (i - dim1b + 1)
1396 : end do
1397 : end do
1398 0 : call pio_initdecomp(pio_subsystem, PIO_DOUBLE, dimsizes, ldof, iodesc)
1399 : call pio_write_darray(file, vdesc, iodesc, &
1400 0 : field(dim1b:dim1e,dim2b:dim2e), ierr, fillval)
1401 0 : call pio_freedecomp(file, iodesc)
1402 0 : deallocate(ldof)
1403 : ! Compute the bounds decomposition and write field bounds
1404 0 : bounds(1) = dim1b
1405 0 : bounds(2) = dim1e
1406 0 : bounds(3) = dim2b
1407 0 : bounds(4) = dim2e
1408 0 : dimsizes(1) = size(bounds, 1)
1409 0 : dimsizes(2) = npes
1410 0 : allocate(ldof(size(bounds, 1)))
1411 0 : do i = 1, size(bounds, 1)
1412 0 : ldof(i) = (iam * size(bounds, 1)) + i
1413 : end do
1414 0 : call pio_initdecomp(pio_subsystem, PIO_INT, dimsizes(1:2), ldof, iodesc)
1415 0 : call pio_write_darray(file, bnddesc, iodesc, bounds, ierr, -900)
1416 0 : call pio_freedecomp(file, iodesc)
1417 0 : deallocate(ldof)
1418 :
1419 : ! All done
1420 0 : call cam_pio_closefile(file)
1421 0 : end subroutine dump_field_2d_d
1422 :
1423 0 : subroutine dump_field_3d_d(fieldname, dim1b, dim1e, dim2b, dim2e, &
1424 0 : dim3b, dim3e, field, compute_maxdim_in, fill_value)
1425 : use pio, only: pio_offset_kind
1426 : use pio, only: pio_double, pio_int, pio_write_darray
1427 : use pio, only: pio_put_att, pio_initdecomp, pio_enddef
1428 : use spmd_utils, only: iam, npes, mpi_max, mpi_integer, mpicom
1429 :
1430 : ! Dummy arguments
1431 : character(len=*), intent(in) :: fieldname
1432 : integer, intent(in) :: dim1b
1433 : integer, intent(in) :: dim1e
1434 : integer, intent(in) :: dim2b
1435 : integer, intent(in) :: dim2e
1436 : integer, intent(in) :: dim3b
1437 : integer, intent(in) :: dim3e
1438 : real(r8), target, intent(in) :: field(dim1b:dim1e,dim2b:dim2e,dim3b:dim3e)
1439 : logical, optional, intent(in) :: compute_maxdim_in
1440 : real(r8), optional, intent(in) :: fill_value
1441 :
1442 : ! Local variables
1443 : type(file_desc_t) :: file
1444 : type(var_desc_t) :: vdesc
1445 : type(var_desc_t) :: bnddesc
1446 : type(io_desc_t) :: iodesc
1447 : character(len=64) :: filename
1448 : real(r8) :: fillval
1449 0 : integer(PIO_OFFSET_KIND), allocatable :: ldof(:)
1450 : integer :: dimids(4)
1451 : integer :: bnddimid
1452 : integer :: bounds(6)
1453 : integer :: dimsizes(4)
1454 : integer :: ierr
1455 : integer :: i, j, k, m, lsize
1456 : logical :: compute_maxdim
1457 :
1458 : ! Find an unused filename for this variable
1459 0 : call find_dump_filename(fieldname, filename)
1460 :
1461 : ! Should we compute max dim sizes or assume they are all the same?
1462 0 : if (present(compute_maxdim_in)) then
1463 0 : compute_maxdim = compute_maxdim_in
1464 : else
1465 : compute_maxdim = .true.
1466 : end if
1467 :
1468 0 : if (present(fill_value)) then
1469 0 : fillval = fill_value
1470 : else
1471 0 : fillval = -900._r8
1472 : end if
1473 :
1474 : ! Open the file for writing
1475 0 : call cam_pio_createfile(file, trim(filename))
1476 :
1477 : ! Define dimensions
1478 0 : if (compute_maxdim) then
1479 : call MPI_allreduce((dim1e - dim1b + 1), dimsizes(1), 1, MPI_integer, &
1480 0 : mpi_max, mpicom, ierr)
1481 : call MPI_allreduce((dim2e - dim2b + 1), dimsizes(2), 1, MPI_integer, &
1482 0 : mpi_max, mpicom, ierr)
1483 : call MPI_allreduce((dim3e - dim3b + 1), dimsizes(3), 1, MPI_integer, &
1484 0 : mpi_max, mpicom, ierr)
1485 : else
1486 0 : dimsizes(1) = dim1e - dim1b + 1
1487 0 : dimsizes(2) = dim2e - dim2b + 1
1488 0 : dimsizes(3) = dim3e - dim3b + 1
1489 : end if
1490 0 : dimsizes(4) = npes
1491 0 : do i = 1, size(dimids, 1)
1492 0 : write(filename, '(a,i0)') 'dim', i
1493 0 : call cam_pio_def_dim(file, trim(filename), dimsizes(i), dimids(i))
1494 : end do
1495 0 : call cam_pio_def_dim(file, 'bounds', size(bounds, 1), bnddimid)
1496 : ! Define the variables
1497 0 : call cam_pio_def_var(file, trim(fieldname), pio_double, dimids, vdesc)
1498 : call cam_pio_def_var(file, 'field_bounds', pio_int, &
1499 0 : (/ bnddimid, dimids(size(dimids, 1)) /), bnddesc)
1500 0 : ierr = pio_put_att(file, vdesc, '_FillValue', fillval)
1501 0 : ierr = pio_enddef(file)
1502 :
1503 : ! Compute the variable decomposition
1504 0 : lsize = product(dimsizes(1:3))
1505 0 : allocate(ldof((dim3e-dim3b+1) * (dim2e-dim2b+1) * (dim1e-dim1b+1)))
1506 0 : m = 0
1507 0 : do k = dim3b, dim3e
1508 0 : do j = dim2b, dim2e
1509 0 : do i = dim1b, dim1e
1510 0 : m = m + 1
1511 0 : ldof(m) = (iam * lsize) + (i - dim1b + 1) + &
1512 0 : (dimsizes(1)*((j - dim2b) + (dimsizes(2)*(k - dim3b))))
1513 : end do
1514 : end do
1515 : end do
1516 0 : call pio_initdecomp(pio_subsystem, PIO_DOUBLE, dimsizes, ldof, iodesc)
1517 : call pio_write_darray(file, vdesc, iodesc, &
1518 0 : field(dim1b:dim1e,dim2b:dim2e,dim3b:dim3e), ierr, fillval)
1519 0 : call pio_freedecomp(file, iodesc)
1520 0 : deallocate(ldof)
1521 : ! Compute the bounds decomposition and write field bounds
1522 0 : bounds(1) = dim1b
1523 0 : bounds(2) = dim1e
1524 0 : bounds(3) = dim2b
1525 0 : bounds(4) = dim2e
1526 0 : bounds(5) = dim3b
1527 0 : bounds(6) = dim3e
1528 0 : dimsizes(1) = size(bounds, 1)
1529 0 : dimsizes(2) = npes
1530 0 : allocate(ldof(size(bounds, 1)))
1531 0 : do i = 1, size(bounds, 1)
1532 0 : ldof(i) = (iam * size(bounds, 1)) + i
1533 : end do
1534 0 : call pio_initdecomp(pio_subsystem, PIO_INT, dimsizes(1:2), ldof, iodesc)
1535 0 : call pio_write_darray(file, bnddesc, iodesc, bounds, ierr, -900)
1536 0 : call pio_freedecomp(file, iodesc)
1537 0 : deallocate(ldof)
1538 :
1539 : ! All done
1540 0 : call cam_pio_closefile(file)
1541 0 : end subroutine dump_field_3d_d
1542 :
1543 0 : subroutine dump_field_4d_d(fieldname, dim1b, dim1e, dim2b, dim2e, &
1544 0 : dim3b, dim3e, dim4b, dim4e, field, compute_maxdim_in, fill_value)
1545 : use pio, only: pio_offset_kind
1546 : use pio, only: pio_double, pio_int, pio_write_darray
1547 : use pio, only: pio_put_att, pio_initdecomp, pio_enddef
1548 : use spmd_utils, only: iam, npes, mpi_max, mpi_integer, mpicom
1549 :
1550 : ! Dummy arguments
1551 : character(len=*), intent(in) :: fieldname
1552 : integer, intent(in) :: dim1b
1553 : integer, intent(in) :: dim1e
1554 : integer, intent(in) :: dim2b
1555 : integer, intent(in) :: dim2e
1556 : integer, intent(in) :: dim3b
1557 : integer, intent(in) :: dim3e
1558 : integer, intent(in) :: dim4b
1559 : integer, intent(in) :: dim4e
1560 : real(r8), target, intent(in) :: field(dim1b:dim1e,dim2b:dim2e,dim3b:dim3e,dim4b:dim4e)
1561 : logical, optional, intent(in) :: compute_maxdim_in
1562 : real(r8), optional, intent(in) :: fill_value
1563 :
1564 : ! Local variables
1565 : type(file_desc_t) :: file
1566 : type(var_desc_t) :: vdesc
1567 : type(var_desc_t) :: bnddesc
1568 : type(io_desc_t) :: iodesc
1569 : character(len=64) :: filename
1570 : real(r8) :: fillval
1571 0 : integer(PIO_OFFSET_KIND), allocatable :: ldof(:)
1572 : integer :: dimids(5)
1573 : integer :: bnddimid
1574 : integer :: bounds(8)
1575 : integer :: dimsizes(5)
1576 : integer :: ierr
1577 : integer :: i, j, k, l, m, lsize
1578 : logical :: compute_maxdim
1579 :
1580 : ! Find an unused filename for this variable
1581 0 : call find_dump_filename(fieldname, filename)
1582 :
1583 : ! Should we compute max dim sizes or assume they are all the same?
1584 0 : if (present(compute_maxdim_in)) then
1585 0 : compute_maxdim = compute_maxdim_in
1586 : else
1587 : compute_maxdim = .true.
1588 : end if
1589 :
1590 0 : if (present(fill_value)) then
1591 0 : fillval = fill_value
1592 : else
1593 0 : fillval = -900._r8
1594 : end if
1595 :
1596 : ! Open the file for writing
1597 0 : call cam_pio_createfile(file, trim(filename))
1598 :
1599 : ! Define dimensions
1600 0 : if (compute_maxdim) then
1601 : call MPI_allreduce((dim1e - dim1b + 1), dimsizes(1), 1, MPI_integer, &
1602 0 : mpi_max, mpicom, ierr)
1603 : call MPI_allreduce((dim2e - dim2b + 1), dimsizes(2), 1, MPI_integer, &
1604 0 : mpi_max, mpicom, ierr)
1605 : call MPI_allreduce((dim3e - dim3b + 1), dimsizes(3), 1, MPI_integer, &
1606 0 : mpi_max, mpicom, ierr)
1607 : call MPI_allreduce((dim4e - dim4b + 1), dimsizes(4), 1, MPI_integer, &
1608 0 : mpi_max, mpicom, ierr)
1609 : else
1610 0 : dimsizes(1) = dim1e - dim1b + 1
1611 0 : dimsizes(2) = dim2e - dim2b + 1
1612 0 : dimsizes(3) = dim3e - dim3b + 1
1613 0 : dimsizes(4) = dim4e - dim4b + 1
1614 : end if
1615 0 : dimsizes(5) = npes
1616 0 : do i = 1, size(dimids, 1)
1617 0 : write(filename, '(a,i0)') 'dim', i
1618 0 : call cam_pio_def_dim(file, trim(filename), dimsizes(i), dimids(i))
1619 : end do
1620 0 : call cam_pio_def_dim(file, 'bounds', size(bounds, 1), bnddimid)
1621 : ! Define the variables
1622 0 : call cam_pio_def_var(file, trim(fieldname), pio_double, dimids, vdesc)
1623 : call cam_pio_def_var(file, 'field_bounds', pio_int, &
1624 0 : (/ bnddimid, dimids(size(dimids, 1)) /), bnddesc)
1625 0 : ierr = pio_put_att(file, vdesc, '_FillValue', fillval)
1626 0 : ierr = pio_enddef(file)
1627 :
1628 : ! Compute the variable decomposition
1629 0 : lsize = product(dimsizes(1:4))
1630 0 : allocate(ldof((dim4e-dim4b+1) * (dim3e-dim3b+1) * (dim2e-dim2b+1) * (dim1e-dim1b+1)))
1631 0 : m = 0
1632 0 : do l = dim4b, dim4e
1633 0 : do k = dim3b, dim3e
1634 0 : do j = dim2b, dim2e
1635 0 : do i = dim1b, dim1e
1636 0 : m = m + 1
1637 0 : ldof(m) = (iam * lsize) + (i - dim1b + 1) + &
1638 : (dimsizes(1)*((j - dim2b) + &
1639 0 : (dimsizes(2)*((k - dim3b) + dimsizes(3)*(l - dim4b)))))
1640 : end do
1641 : end do
1642 : end do
1643 : end do
1644 0 : call pio_initdecomp(pio_subsystem, PIO_DOUBLE, dimsizes, ldof, iodesc)
1645 : call pio_write_darray(file, vdesc, iodesc, &
1646 0 : field(dim1b:dim1e,dim2b:dim2e,dim3b:dim3e,dim4b:dim4e), ierr, fillval)
1647 0 : call pio_freedecomp(file, iodesc)
1648 0 : deallocate(ldof)
1649 : ! Compute the bounds decomposition and write field bounds
1650 0 : bounds(1) = dim1b
1651 0 : bounds(2) = dim1e
1652 0 : bounds(3) = dim2b
1653 0 : bounds(4) = dim2e
1654 0 : bounds(5) = dim3b
1655 0 : bounds(6) = dim3e
1656 0 : bounds(7) = dim4b
1657 0 : bounds(8) = dim4e
1658 0 : dimsizes(1) = size(bounds, 1)
1659 0 : dimsizes(2) = npes
1660 0 : allocate(ldof(size(bounds, 1)))
1661 0 : do i = 1, size(bounds, 1)
1662 0 : ldof(i) = (iam * size(bounds, 1)) + i
1663 : end do
1664 0 : call pio_initdecomp(pio_subsystem, PIO_INT, dimsizes(1:2), ldof, iodesc)
1665 0 : call pio_write_darray(file, bnddesc, iodesc, bounds, ierr, -900)
1666 0 : call pio_freedecomp(file, iodesc)
1667 0 : deallocate(ldof)
1668 :
1669 : ! All done
1670 0 : call cam_pio_closefile(file)
1671 0 : end subroutine dump_field_4d_d
1672 :
1673 0 : subroutine dump_field_6d_d(fieldname, dimbs, dimes, field, &
1674 : compute_maxdim_in, fill_value)
1675 : use pio, only: pio_offset_kind
1676 : use pio, only: pio_double, pio_int, pio_write_darray
1677 : use pio, only: pio_put_att, pio_initdecomp, pio_enddef
1678 : use spmd_utils, only: iam, npes, mpi_max, mpi_integer, mpicom
1679 :
1680 : ! Dummy arguments
1681 : character(len=*), intent(in) :: fieldname
1682 : integer, intent(in) :: dimbs(:)
1683 : integer, intent(in) :: dimes(:)
1684 : real(r8), target, intent(in) :: field(:,:,:,:,:,:)
1685 : logical, optional, intent(in) :: compute_maxdim_in
1686 : real(r8), optional, intent(in) :: fill_value
1687 :
1688 : ! Local variables
1689 : type(file_desc_t) :: file
1690 : type(var_desc_t) :: vdesc
1691 : type(var_desc_t) :: bnddesc
1692 : type(io_desc_t) :: iodesc
1693 : character(len=64) :: filename
1694 : real(r8) :: fillval
1695 0 : integer(PIO_OFFSET_KIND), allocatable :: ldof(:)
1696 : integer :: dimids(7)
1697 : integer :: bnddimid
1698 : integer :: bounds(14)
1699 : integer :: dimsizes(7)
1700 : integer :: ierr
1701 : integer :: i1, i2, i3, i4, i5, i6, i(6)
1702 : integer :: ind, lsize, j
1703 : logical :: compute_maxdim
1704 :
1705 : ! Find an unused filename for this variable
1706 0 : call find_dump_filename(fieldname, filename)
1707 :
1708 : ! Should we compute max dim sizes or assume they are all the same?
1709 0 : if (present(compute_maxdim_in)) then
1710 0 : compute_maxdim = compute_maxdim_in
1711 : else
1712 : compute_maxdim = .true.
1713 : end if
1714 :
1715 0 : if (present(fill_value)) then
1716 0 : fillval = fill_value
1717 : else
1718 0 : fillval = -900._r8
1719 : end if
1720 :
1721 : ! Open the file for writing
1722 0 : call cam_pio_createfile(file, trim(filename))
1723 :
1724 : ! Define dimensions
1725 0 : if (compute_maxdim) then
1726 0 : do i1 = 1, 6
1727 0 : call MPI_allreduce((dimes(i1) - dimbs(i1) + 1), dimsizes(i1), 1, &
1728 0 : MPI_integer, mpi_max, mpicom, ierr)
1729 : end do
1730 : else
1731 0 : do i1 = 1, 6
1732 0 : dimsizes(i1) = dimes(i1) - dimbs(i1) + 1
1733 : end do
1734 : end if
1735 0 : dimsizes(7) = npes
1736 0 : do ind = 1, 7
1737 0 : write(filename, '(a,i0)') 'dim', ind
1738 0 : call cam_pio_def_dim(file, trim(filename), dimsizes(ind), dimids(ind))
1739 : end do
1740 0 : call cam_pio_def_dim(file, 'bounds', size(bounds, 1), bnddimid)
1741 : ! Define the variables
1742 : call cam_pio_def_var(file, 'field_bounds', pio_int, &
1743 0 : (/ bnddimid, dimids(size(dimids, 1)) /), bnddesc)
1744 0 : call cam_pio_def_var(file, trim(fieldname), pio_double, dimids, vdesc)
1745 0 : ierr = pio_put_att(file, vdesc, '_FillValue', fillval)
1746 0 : ierr = pio_enddef(file)
1747 :
1748 : ! Compute the variable decomposition
1749 0 : lsize = 1
1750 0 : do ind = 1, 6
1751 0 : lsize = lsize * (dimes(ind) - dimbs(ind) + 1)
1752 : end do
1753 0 : allocate(ldof(lsize))
1754 0 : ind = 0
1755 0 : do i6 = dimbs(6), dimes(6)
1756 0 : i(6) = i6
1757 0 : do i5 = dimbs(5), dimes(5)
1758 0 : i(5) = i5
1759 0 : do i4 = dimbs(4), dimes(4)
1760 0 : i(4) = i4
1761 0 : do i3 = dimbs(3), dimes(3)
1762 0 : i(3) = i3
1763 0 : do i2 = dimbs(2), dimes(2)
1764 0 : i(2) = i2
1765 0 : do i1 = dimbs(1), dimes(1)
1766 0 : i(1) = i1
1767 0 : ind = ind + 1
1768 0 : ldof(ind) = iam
1769 0 : do j = 6, 1, -1
1770 0 : ldof(ind) = (ldof(ind) * dimsizes(j)) + (i(j) - dimbs(j))
1771 : end do
1772 : end do
1773 : end do
1774 : end do
1775 : end do
1776 : end do
1777 : end do
1778 0 : call pio_initdecomp(pio_subsystem, PIO_DOUBLE, dimsizes, ldof, iodesc)
1779 : call pio_write_darray(file, vdesc, iodesc, &
1780 0 : field(dimbs(1):dimes(1),dimbs(2):dimes(2),dimbs(3):dimes(3), &
1781 0 : dimbs(4):dimes(4),dimbs(5):dimes(5),dimbs(6):dimes(6)), ierr, fillval)
1782 0 : call pio_freedecomp(file, iodesc)
1783 0 : deallocate(ldof)
1784 : ! Compute the bounds decomposition and write field bounds
1785 0 : do ind = 1, 6
1786 0 : bounds(2*ind - 1) = dimbs(ind)
1787 0 : bounds(2*ind) = dimes(ind)
1788 : end do
1789 0 : bounds(13) = 1
1790 0 : bounds(14) = npes
1791 0 : dimsizes(1) = size(bounds, 1)
1792 0 : dimsizes(2) = npes
1793 0 : allocate(ldof(size(bounds, 1)))
1794 0 : do ind = 1, size(bounds, 1)
1795 0 : ldof(ind) = (iam * size(bounds, 1)) + ind
1796 : end do
1797 0 : call pio_initdecomp(pio_subsystem, PIO_INT, dimsizes(1:2), ldof, iodesc)
1798 0 : call pio_write_darray(file, bnddesc, iodesc, bounds, ierr, -900)
1799 0 : call pio_freedecomp(file, iodesc)
1800 0 : deallocate(ldof)
1801 :
1802 : ! All done
1803 0 : call cam_pio_closefile(file)
1804 :
1805 0 : end subroutine dump_field_6d_d
1806 :
1807 0 : end module cam_pio_utils
|