Line data Source code
1 : !#define DEBUG 1
2 : module physics_buffer
3 :
4 : !-----------------------------------------------------------------------
5 : !
6 : ! Purpose:
7 : ! Buffer managment for persistant variables
8 : !
9 : ! Author: J. Edwards
10 : !
11 : ! This file is used with genf90.pl to generate buffer.F90
12 : !
13 : !-----------------------------------------------------------------------
14 :
15 : use shr_kind_mod, only: r8 => shr_kind_r8, r4=> shr_kind_r4, i4=> shr_kind_i4
16 : use ppgrid, only: pcols, begchunk, endchunk, psubcols
17 : use cam_logfile, only: iulog
18 : use pio, only: var_desc_t
19 : use dyn_grid, only: ptimelevels
20 : use cam_abortutils, only: endrun
21 : use buffer, only: buffer_field_allocate, buffer_field_deallocate, buffer_get_field_ptr, buffer_set_field, &
22 : dtype_i4, dtype_r4, dtype_r8, buffer_field_default_type, buffer_field_is_alloc
23 :
24 : implicit none
25 : private
26 :
27 : ! ngrid_types is a private parameter denoting how many types of a field
28 : ! (e.g., grid, subcol) a pbuf can hold
29 : integer, parameter :: ngrid_types = 2
30 :
31 : ! --col_type parameters -- see pbuf_add_field and pbuf_get_field
32 : ! -- These indices should start at 1 and increment by 1 so that they
33 : ! -- may be used in a loop from 1 to ngrid_types
34 :
35 : integer, parameter, public :: col_type_grid = 1
36 : integer, parameter, public :: col_type_subcol = 2
37 :
38 : !
39 : ! PBUF field suffix strings for different grid types for restart files
40 : ! NB: There should be ngrid_types entries
41 : !
42 : character(len=5), parameter :: field_grid_suff(ngrid_types) = (/ " ", "_SCOL" /)
43 :
44 : integer :: old_time_idx = 1
45 :
46 : ! The API has strings 'GLOBAL' and 'PHYSPKG' which correspond to these
47 : ! integer constants if global_allocate_all is false fields allocated with
48 : ! PHYSPKG persistence are deallocated at the end of each physics time step
49 : ! and reallocated at the beginning of the next.
50 : ! If global_allocate_all is true then all fields are allocated at model
51 : ! start and persist until model completion.
52 :
53 : integer, parameter :: persistence_global = 1
54 : integer, parameter :: persistence_physpkg = 2
55 : logical :: global_allocate_all = .true.
56 :
57 : !
58 : ! physics_buffer_hdr carries the description info for each buffer, only one header
59 : ! is allocated for each field and each chunk of the field points to it.
60 : ! It is carried as a linkedlist for initialization only.
61 : !
62 :
63 : type physics_buffer_hdr
64 : character(len=16) :: name = ''
65 : integer :: dtype = -1
66 : integer :: persistence = -1
67 : integer :: dimsizes(6,ngrid_types) = 0
68 : type(var_desc_t) :: vardesc(ngrid_types) ! var id for restart files
69 : logical :: is_copy(ngrid_types) = .false.
70 : logical :: added = .false.
71 : type(physics_buffer_hdr), pointer :: nexthdr => null()
72 :
73 : end type physics_buffer_hdr
74 :
75 : !
76 : ! The default type for a buffer field is buffer_field_double (r8) since that is the
77 : ! type most often used in the model. The F90 transfer function is used to move
78 : ! data of other types into and out of the pbuf2d
79 : !
80 :
81 : type physics_buffer_desc
82 : private
83 : integer :: lchnk
84 : type(physics_buffer_hdr), pointer :: hdr => null()
85 : type(buffer_field_default_type) :: bfg(ngrid_types)
86 : logical :: used = .false.
87 : end type physics_buffer_desc
88 :
89 : interface pbuf_get_field
90 : ! TYPE int,double,real
91 : ! DIMS 1,2,3,4,5
92 : module procedure get_pbuf1d_field_by_index_{DIMS}d_{TYPE}
93 : ! TYPE int,double,real
94 : ! DIMS 1,2,3,4,5
95 : module procedure get_pbuf2d_field_by_index_{DIMS}d_{TYPE}
96 : end interface
97 :
98 : interface pbuf_get_field_restart
99 : ! TYPE int,double,real
100 : module procedure get_pbuf2d_field_restart_{TYPE}
101 : end interface
102 :
103 : interface pbuf_set_field
104 : ! TYPE int,double,real
105 : ! DIMS 1,2,3,4,5
106 : module procedure set_pbuf2d_field_by_index_{DIMS}d_{TYPE}
107 : ! TYPE int,double,real
108 : ! DIMS 1,2,3,4,5
109 : module procedure set_pbuf1d_field_by_index_{DIMS}d_{TYPE}
110 : ! TYPE int,double,real
111 : module procedure set_pbuf1d_field_const_by_index_{TYPE}
112 : ! TYPE int,double,real
113 : module procedure set_pbuf2d_field_const_by_index_{TYPE}
114 : end interface
115 :
116 : interface pbuf_add_field
117 : ! TYPE int,double,real
118 : module procedure pbuf_add_field_{TYPE}
119 : end interface
120 :
121 : public :: pbuf_initialize, &
122 : pbuf_readnl, &! read namelist options
123 : pbuf_init_time, &! Initialize dyn_time_lvls
124 : pbuf_old_tim_idx, &! return the index for the oldest time
125 : pbuf_update_tim_idx, &! update the index for the oldest time
126 : pbuf_col_type_index, &
127 : pbuf_get_field_name, &
128 : pbuf_get_field, &
129 : pbuf_add_field, &
130 : pbuf_register_subcol, &
131 : physics_buffer_desc, &
132 : pbuf_get_index, &
133 : pbuf_get_chunk, &
134 : pbuf_allocate, &
135 : pbuf_deallocate, &
136 : pbuf_set_field, &
137 : pbuf_init_restart, &
138 : pbuf_write_restart, &
139 : pbuf_read_restart, &
140 : pbuf_is_used, &
141 : dtype_r8, dtype_r4, dtype_i4
142 :
143 : ! For help debugging code
144 : public :: pbuf_dump_pbuf
145 :
146 : ! For snapshot use
147 : public :: pbuf_get_dim_strings
148 : public :: pbuf_cam_snapshot_register
149 :
150 : integer, public :: dyn_time_lvls ! number of time levels in physics buffer (dycore dependent)
151 :
152 :
153 : !
154 : ! Currentpbufflds is incremented in calls to pbuf_add_field and determines the size
155 : ! of the allocated pbuf2d
156 : !
157 :
158 : integer :: currentpbufflds=0
159 : type(physics_buffer_hdr), pointer :: hdrbuffertop => null()
160 :
161 : !
162 : ! Insures that we do not attempt to allocate physics_buffer more than once
163 : !
164 :
165 : logical :: buffer_initialized =.false.
166 :
167 : !
168 : ! private pio descriptor for time
169 : !
170 : type(var_desc_t) :: timeidx_desc
171 :
172 : ! Set to .true. for more output
173 : logical :: debug = .false.
174 :
175 : !===============================================================================
176 : CONTAINS
177 : !===============================================================================
178 :
179 0 : subroutine pbuf_readnl(nlfile)
180 :
181 : use namelist_utils, only: find_group_name
182 : use units, only: getunit, freeunit
183 : use spmd_utils, only: masterproc, mpicom, mstrid=>masterprocid, mpi_logical
184 :
185 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
186 :
187 : ! Local variables
188 : integer :: unitn, ierr
189 : character(len=*), parameter :: sub = 'pbuf_readnl'
190 :
191 : logical :: pbuf_global_allocate
192 :
193 : namelist /pbuf_nl/ pbuf_global_allocate
194 : !-----------------------------------------------------------------------------
195 :
196 0 : pbuf_global_allocate = global_allocate_all
197 :
198 : ! Read namelist
199 0 : if (masterproc) then
200 0 : unitn = getunit()
201 0 : open( unitn, file=trim(nlfile), status='old' )
202 0 : call find_group_name(unitn, 'pbuf_nl', status=ierr)
203 0 : if (ierr == 0) then
204 0 : read(unitn, pbuf_nl, iostat=ierr)
205 0 : if (ierr /= 0) then
206 0 : call endrun(sub//': FATAL: reading namelist')
207 : end if
208 : end if
209 0 : close(unitn)
210 0 : call freeunit(unitn)
211 : end if
212 :
213 0 : call mpi_bcast(pbuf_global_allocate, 1, mpi_logical, mstrid, mpicom, ierr)
214 0 : if (ierr /= 0) call endrun(sub//": FATAL: mpi_bcast: pbuf_global_allocate")
215 :
216 0 : global_allocate_all = pbuf_global_allocate
217 :
218 0 : if (masterproc) then
219 0 : write(iulog,*) 'physics buffer options:'
220 0 : write(iulog,*) ' pbuf_global_allocate =', global_allocate_all
221 : end if
222 :
223 0 : end subroutine pbuf_readnl
224 :
225 : !===============================================================================
226 :
227 : !
228 : ! Initialize dyn_time_lvls
229 : !
230 1536 : subroutine pbuf_init_time()
231 1536 : dyn_time_lvls = ptimelevels - 1
232 1536 : end subroutine pbuf_init_time
233 :
234 : !
235 : ! Return index of oldest time sample in the physics buffer.
236 : !
237 :
238 11204424 : function pbuf_old_tim_idx()
239 : integer :: pbuf_old_tim_idx
240 11204424 : pbuf_old_tim_idx = old_time_idx
241 11204424 : end function pbuf_old_tim_idx
242 :
243 : !
244 : ! Update index of old time sample in the physics buffer.
245 : !
246 :
247 369408 : subroutine pbuf_update_tim_idx()
248 369408 : old_time_idx = mod(old_time_idx, dyn_time_lvls) + 1
249 369408 : end subroutine pbuf_update_tim_idx
250 :
251 : !
252 : ! pbuf_col_type_index returns an index for use with pbuf calls
253 : !
254 : ! * col_type: is set to col_type_grid (if use_subcol=.false.) or
255 : ! col_type_subcol (if (use_subcol=.true.).
256 : !
257 :
258 0 : subroutine pbuf_col_type_index(use_subcol, col_type)
259 :
260 : logical, intent(in) :: use_subcol
261 : integer, intent(out) :: col_type
262 :
263 0 : if (use_subcol) then
264 0 : col_type = col_type_subcol
265 : else
266 0 : col_type = col_type_grid
267 : end if
268 :
269 0 : end subroutine pbuf_col_type_index
270 :
271 : !
272 : ! Return a pointer to the current chunks physics_buffer.
273 : !
274 :
275 0 : function pbuf_get_field_name(index)
276 : integer, intent(in) :: index
277 : character(len=16) :: pbuf_get_field_name
278 : integer :: i
279 : type(physics_buffer_hdr), pointer :: hdrbuffer
280 :
281 0 : hdrbuffer => hdrbuffertop
282 0 : do i=2,index
283 0 : hdrbuffer=>hdrbuffer%nexthdr
284 : end do
285 0 : pbuf_get_field_name = hdrbuffer%name
286 :
287 0 : end function pbuf_get_field_name
288 :
289 :
290 22716096 : function pbuf_get_chunk(pbuf2d, lchnk)
291 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
292 : integer, intent(in) :: lchnk
293 :
294 : type(physics_buffer_desc), pointer :: pbuf_get_chunk(:)
295 :
296 22716096 : pbuf_get_chunk => pbuf2d(:,lchnk)
297 :
298 :
299 : end function pbuf_get_chunk
300 :
301 : !
302 : ! Return .true. iff pbuf has an allocated grid field
303 : !
304 :
305 : logical function pbuf_field_has_gridcols(pbuf, index) result(rval)
306 :
307 : type(physics_buffer_desc), pointer :: pbuf(:)
308 : integer, intent(in) :: index
309 :
310 : ! If the field is a copy, return .false. even if it is allocated
311 : if (pbuf(index)%hdr%is_copy(col_type_grid)) then
312 : rval = .false.
313 : else
314 : rval = buffer_field_is_alloc(pbuf(index)%bfg(col_type_grid))
315 : end if
316 :
317 : end function pbuf_field_has_gridcols
318 :
319 : !
320 : ! Return .true. iff pbuf has an allocated subcolumn field
321 : !
322 :
323 : logical function pbuf_field_has_subcols(pbuf, index) result(rval)
324 :
325 : type(physics_buffer_desc), pointer :: pbuf(:)
326 : integer, intent(in) :: index
327 :
328 : ! If the field is a copy, return .false. even if it is allocated
329 : if (pbuf(index)%hdr%is_copy(col_type_subcol)) then
330 : rval = .false.
331 : else
332 : rval = buffer_field_is_alloc(pbuf(index)%bfg(col_type_subcol))
333 : end if
334 :
335 : end function pbuf_field_has_subcols
336 :
337 : !
338 : ! Return .true. iff pbuf has an allocated col_type column field
339 : !
340 :
341 : logical function pbuf_field_has_col_type(pbuf, index, col_type) result(rval)
342 :
343 : type(physics_buffer_desc), pointer :: pbuf(:)
344 : integer, intent(in) :: index
345 : integer, intent(in) :: col_type
346 :
347 : if (col_type == col_type_grid) then
348 : rval = pbuf_field_has_gridcols(pbuf, index)
349 : else if(col_type == col_type_subcol) then
350 : rval = pbuf_field_has_subcols(pbuf, index)
351 : else
352 : call endrun('pbuf_field_has_col_type: Invalid col_type')
353 : end if
354 :
355 : end function pbuf_field_has_col_type
356 :
357 : !
358 : ! Initialize the buffer, should be called after all pbuf_add_field calls
359 : ! have been completed and should only be called once in a run
360 : !
361 2304 : subroutine pbuf_initialize(pbuf2d)
362 : use phys_grid, only: phys_grid_initialized
363 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
364 :
365 : integer :: i, c
366 : type(physics_buffer_hdr), pointer :: hdrbuffer
367 :
368 : !
369 : ! Allocate memory
370 : !
371 2304 : if (buffer_initialized) return
372 1536 : if (.not. phys_grid_initialized()) then
373 0 : call endrun('pbuf_initialize: Physics decomposition not initialized')
374 : end if
375 : ! Allocate at least 1 to avoid unallocated error in ideal physics
376 739872 : allocate(pbuf2d(max(1,currentpbufflds),begchunk:endchunk))
377 1536 : if(currentpbufflds<1) return
378 :
379 7728 : do c = begchunk, endchunk
380 6192 : hdrbuffer => hdrbuffertop
381 732192 : do i = 1, currentpbufflds
382 724464 : if(.not. hdrbuffer%added) then
383 0 : call endrun('pbuf_initialize: pbuf, '//trim(hdrbuffer%name)//', not added')
384 : end if
385 724464 : pbuf2d(i,c)%lchnk = c
386 724464 : pbuf2d(i,c)%hdr => hdrbuffer
387 730656 : hdrbuffer => hdrbuffer%nexthdr
388 : end do
389 : end do
390 :
391 1536 : buffer_initialized=.true.
392 1536 : call pbuf_allocate(pbuf2d, 'global')
393 : #ifdef DEBUG
394 1536 : call pbuf2d_print(pbuf2d)
395 : #endif
396 1536 : return
397 2304 : end subroutine pbuf_initialize
398 :
399 :
400 372480 : subroutine pbuf_allocate(pbuf2d, persistence)
401 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
402 : character(len=*), intent(in) :: persistence
403 : integer :: i
404 : logical :: allocate_all
405 :
406 : ! allocate_all is used to force allocation of all fields at same time as allocation
407 : ! for global scope.
408 372480 : allocate_all = .false.
409 372480 : if ( global_allocate_all ) then
410 372480 : if ( persistence == 'global' ) then
411 : allocate_all = .true.
412 : else
413 370944 : return
414 : end if
415 : end if
416 :
417 : if(allocate_all) then
418 181248 : do i=1,currentpbufflds
419 181248 : select case(pbuf2d(i,begchunk)%hdr%dtype)
420 : case(TYPEDOUBLE)
421 173568 : call pbuf_allocate_field_double(pbuf2d, i, dtype_r8)
422 : case(TYPEREAL)
423 0 : call pbuf_allocate_field_real(pbuf2d, i, dtype_r4)
424 : ! case(i8)
425 : ! call pbuf_allocate_field_long(pbuf2d, i)
426 : case(TYPEINT)
427 179712 : call pbuf_allocate_field_int(pbuf2d, i, dtype_i4)
428 : end select
429 : end do
430 0 : else if(persistence .eq. 'physpkg') then
431 0 : do i=1,currentpbufflds
432 0 : if(pbuf2d(i,begchunk)%hdr%persistence==persistence_physpkg) then
433 0 : select case(pbuf2d(i,begchunk)%hdr%dtype)
434 : case(TYPEDOUBLE)
435 0 : call pbuf_allocate_field_double(pbuf2d, i, dtype_r8)
436 : case(TYPEREAL)
437 0 : call pbuf_allocate_field_real(pbuf2d, i, dtype_r4)
438 : ! case(i8)
439 : ! call pbuf_allocate_field_long(pbuf2d, i)
440 : case(TYPEINT)
441 0 : call pbuf_allocate_field_int(pbuf2d, i, dtype_i4)
442 : end select
443 : end if
444 : end do
445 0 : else if(persistence .eq. 'global') then
446 0 : do i=1,currentpbufflds
447 0 : if(pbuf2d(i,begchunk)%hdr%persistence==persistence_global) then
448 :
449 0 : select case(pbuf2d(i,begchunk)%hdr%dtype)
450 : case(TYPEDOUBLE)
451 0 : call pbuf_allocate_field_double(pbuf2d, i, dtype_r8)
452 : case(TYPEREAL)
453 0 : call pbuf_allocate_field_real(pbuf2d, i, dtype_r4)
454 : ! case(i8)
455 : ! call pbuf_allocate_field_long(pbuf2d, i)
456 : case(TYPEINT)
457 0 : call pbuf_allocate_field_int(pbuf2d, i, dtype_i4)
458 : end select
459 : end if
460 : end do
461 : end if
462 2304 : end subroutine pbuf_allocate
463 :
464 : ! TYPE int,double,real
465 179712 : subroutine pbuf_allocate_field_{TYPE}(pbuf2d, index, dtype)
466 : {VTYPE}, intent(in) :: dtype
467 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
468 : integer, intent(in) :: index
469 :
470 179712 : integer, pointer :: dimsizes(:)
471 : integer :: c, i
472 : logical :: is_copy
473 :
474 539136 : do i = 1, ngrid_types
475 : ! Note - dimsizes(:)=0 is special case to indicate "do not allocate" and is not fatal
476 : ! Since this is called by a single thread, setting the dimsizes pointer to first chunk is okay
477 359424 : dimsizes => pbuf2d(index,begchunk)%hdr%dimsizes(:,i)
478 : ! Note - We may have dimsizes /= 0 but still don't allocate if copy
479 359424 : is_copy = pbuf2d(index,begchunk)%hdr%is_copy(i)
480 :
481 2515968 : if (any(dimsizes(:) < 0)) &
482 : call endrun('pbuf_allocate_field: dimsizes must be greater than 0 for pbuf field '&
483 0 : //trim(pbuf2d(index,begchunk)%hdr%name))
484 :
485 1617408 : if (all(dimsizes(:) /= 0) .and. (.not. is_copy)) then
486 904176 : do c = begchunk, endchunk
487 904176 : call buffer_field_allocate(pbuf2d(index,c)%bfg(i), dimsizes, dtype)
488 : end do
489 : end if
490 : end do
491 :
492 179712 : end subroutine pbuf_allocate_field_{TYPE}
493 :
494 86800896 : pure logical function pbuf_do_deallocate(hdr, persistence, col_type) result(rval)
495 : type(physics_buffer_hdr), pointer :: hdr
496 : character(len=*), intent(in) :: persistence
497 : integer, intent(in) :: col_type
498 :
499 86800896 : if (persistence .eq. 'physpkg') then
500 86441472 : if (hdr%is_copy(col_type)) then
501 : ! If this is a copied field, it is always deallocated
502 : rval = .true.
503 86441472 : else if (global_allocate_all) then
504 : ! This is an all-global run so no deallocation
505 : rval = .false.
506 0 : else if(hdr%persistence == persistence_physpkg) then
507 : ! This pbuf is a physpkg type, deallocate
508 : rval = .true.
509 : else
510 : ! This pbuf is global, do not deallocate
511 0 : rval = .false.
512 : end if
513 : else
514 : ! We are doing a global deallocate, everything must go!
515 : rval = .true.
516 : end if
517 86800896 : end function pbuf_do_deallocate
518 :
519 370944 : subroutine pbuf_deallocate(pbuf2d, persistence)
520 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
521 : character(len=*) :: persistence
522 :
523 : integer :: i, j, c
524 :
525 43771392 : do i = 1, currentpbufflds
526 130572288 : do j = 1, ngrid_types
527 130201344 : if (pbuf_do_deallocate(pbuf2d(i,begchunk)%hdr, persistence, j)) then
528 359424 : if (buffer_field_is_alloc(pbuf2d(i,begchunk)%bfg(j))) then
529 353280 : select case(pbuf2d(i,begchunk)%hdr%dtype)
530 : case(TYPEDOUBLE)
531 873264 : do c = begchunk,endchunk
532 873264 : call buffer_field_deallocate(pbuf2d(i,c)%bfg(j), dtype_r8)
533 : end do
534 : case(TYPEREAL)
535 0 : do c = begchunk,endchunk
536 0 : call buffer_field_deallocate(pbuf2d(i,c)%bfg(j), dtype_r4)
537 : end do
538 : case(TYPEINT)
539 204480 : do c = begchunk,endchunk
540 30912 : call buffer_field_deallocate(pbuf2d(i,c)%bfg(j), dtype_i4)
541 : end do
542 : end select
543 : end if
544 : end if
545 : end do ! ngrid_types
546 : end do ! currentpbufflds
547 370944 : end subroutine pbuf_deallocate
548 :
549 : ! Find a pbuf header pointer based on the name input
550 : ! Automatically tacks on an extra header if <name> is not found.
551 179712 : subroutine find_pbuf_header(name, index, bufptr)
552 : character(len=*), intent(in) :: name
553 : type(physics_buffer_hdr), pointer :: bufptr
554 : integer, intent(out) :: index
555 :
556 : ! Local Variables
557 : logical :: buf_found
558 :
559 179712 : if(.not. associated(hdrbuffertop)) then
560 : ! This is the very first pbuf, allocate and set
561 26112 : allocate(hdrbuffertop)
562 1536 : currentpbufflds = 1
563 1536 : hdrbuffertop%name = name
564 : end if
565 :
566 179712 : bufptr=>hdrbuffertop
567 179712 : buf_found = .false.
568 179712 : index = 1
569 10604544 : do while(.not. buf_found)
570 10604544 : if(trim(bufptr%name) == trim(name)) then
571 : buf_found = .true.
572 : else
573 10423296 : if (associated(bufptr%nexthdr)) then
574 10245120 : bufptr=>bufptr%nexthdr
575 10245120 : index = index + 1
576 : else
577 : ! We ran off the end of the buffers, make a new one for <name>
578 : ! Sanity check, we should have checked exactly currentpbufflds
579 178176 : if (index /= currentpbufflds) then
580 0 : call endrun("find_pbuf_header: currentpbufflds indexing off")
581 : end if
582 178176 : currentpbufflds = currentpbufflds + 1
583 178176 : index = currentpbufflds
584 3028992 : allocate(bufptr%nexthdr)
585 178176 : bufptr=>bufptr%nexthdr
586 178176 : bufptr%name = trim(name)
587 178176 : buf_found = .true.
588 : end if
589 : end if
590 : end do
591 179712 : end subroutine find_pbuf_header
592 : !
593 : ! Register a field in the pbuf
594 : ! This should be called from physics register routines.
595 : ! persistence must be 'global' or 'physpkg'
596 : ! dtype can be any of r8, r4, i4 as defined in shr_kinds_mod.F90
597 : ! col_type is either col_type_grid or col_type_subcol.
598 : ! If no col_type, then grid field is defined (i.e., dimsizes set)
599 :
600 179712 : subroutine pbuf_register_field_int(name, pname, index, persistence, &
601 179712 : dtype, dimsizes, col_type, pbuf_add)
602 :
603 : ! Dummy Arguments
604 : character(len=*), intent(in) :: name
605 : character(len=*), intent(in) :: pname ! name of calling parameterization
606 : integer, intent(out) :: index
607 : character(len=*), optional, intent(in) :: persistence
608 : integer, optional, intent(in) :: dtype ! used to differentiate specific calls
609 : integer, optional, intent(in) :: dimsizes(:) ! dimension sizes of grid field
610 : integer, optional, intent(in) :: col_type
611 : logical, optional, intent(in) :: pbuf_add
612 :
613 :
614 : ! Local Variables
615 : type(physics_buffer_hdr), pointer :: bufptr
616 : integer :: dimcnt, col_type_use
617 : character(len=128) :: errmsg
618 :
619 179712 : if(buffer_initialized) then
620 0 : call endrun('Attempt to register pbuf field after buffer initialized')
621 : end if
622 :
623 179712 : call find_pbuf_header(name, index, bufptr)
624 :
625 179712 : if (present(persistence)) then
626 179712 : if(persistence .eq. "global") then
627 58368 : bufptr%persistence = persistence_global
628 : else
629 121344 : bufptr%persistence = persistence_physpkg
630 : end if
631 : end if
632 :
633 179712 : if (present(dtype)) then
634 179712 : bufptr%dtype = dtype
635 : end if
636 :
637 179712 : if (present(dimsizes)) then
638 : ! Normally, we only set buffer dimsizes if dimsizes is passed
639 179712 : dimcnt = size(dimsizes)
640 :
641 179712 : if (present(col_type)) then
642 179712 : col_type_use = col_type
643 : else
644 : col_type_use = col_type_grid
645 : end if
646 :
647 : ! Only allow up to 5 dimensions. #6 is reserved for subcolumn umpacking
648 : ! and #7 is for the physics chunk index.
649 179712 : if (dimcnt > 5) then
650 0 : call endrun('pbuf_register_field: Attempt to exceed maximum of 5 dimensions for '//trim(name))
651 : end if
652 :
653 : ! Assign dimensions dependent on col_types input
654 : ! Note that dimensions are initialized to zero and set if being used
655 179712 : if (col_type_use == col_type_grid) then
656 : ! grid is requested
657 1257984 : bufptr%dimsizes(:,col_type_grid) = 1
658 505344 : bufptr%dimsizes(1:dimcnt,col_type_grid) = dimsizes
659 : ! If someone previously registered the subcol, reset those dims
660 179712 : if (bufptr%dimsizes(1,col_type_subcol) == pcols*psubcols) then
661 0 : bufptr%dimsizes(2:dimcnt,col_type_subcol) = dimsizes(2:dimcnt)
662 : end if
663 0 : else if (col_type_use == col_type_subcol) then
664 0 : bufptr%dimsizes(:,col_type_subcol) = 1
665 0 : bufptr%dimsizes(1,col_type_subcol) = pcols*psubcols
666 : ! This case should only be used for a pbuf_add_field call adding a subcolumn-only field
667 0 : if (dimcnt > 1) then
668 0 : bufptr%dimsizes(2:dimcnt,col_type_subcol) = dimsizes(2:dimcnt)
669 0 : else if (bufptr%dimsizes(1,col_type_grid) > 0) then
670 : ! If grid field previously registered or added, update dims
671 0 : bufptr%dimsizes(2:,col_type_subcol) = bufptr%dimsizes(2:,col_type_grid)
672 : end if
673 : end if
674 : end if
675 :
676 179712 : if (present(pbuf_add)) then
677 179712 : if (pbuf_add) then
678 : ! An add request has been made but we have to make sure we have all info
679 179712 : if (all(bufptr%dimsizes(:,:) == 0)) then
680 0 : write(errmsg, *) 'pbuf_add_field: trying to add field with no dimensions'
681 0 : call endrun(errmsg)
682 : end if
683 179712 : if (bufptr%dtype < 0) then
684 0 : write(errmsg, *) 'pbuf_add_field: trying to add field with no type'
685 0 : call endrun(errmsg)
686 : end if
687 179712 : if (bufptr%persistence < 0) then
688 0 : write(errmsg, *) 'pbuf_add_field: trying to add field with bad persistence'
689 0 : call endrun(errmsg)
690 : end if
691 179712 : bufptr%added = .true.
692 : end if
693 : end if
694 :
695 179712 : end subroutine pbuf_register_field_int
696 :
697 : !
698 : ! Add a field to the pbuf, this should be called from physics register routines
699 : ! name is required to be unique,
700 : ! persistence must be 'global' or 'physpkg'
701 : ! dtype can be any of r8, r4, i4 as defined in shr_kinds_mod.F90
702 : ! If present, col_type must be either col_type_grid or col_type_subcol
703 :
704 : ! TYPE int,double,real
705 179712 : subroutine pbuf_add_field_{TYPE}(name, persistence, dtype, dimsizes, index, col_type)
706 :
707 : character(len=*), intent(in) :: name, persistence
708 : {VTYPE}, intent(in) :: dtype ! used only to differentiate specific calls
709 : integer, intent(in) :: dimsizes(:) ! dimension sizes of grid field
710 : integer, intent(in), optional :: col_type
711 : integer, intent(out) :: index
712 :
713 : ! Local Variables
714 : integer :: col_type_use
715 :
716 179712 : if(buffer_initialized) then
717 0 : call endrun('Attempt to add pbuf field after buffer initialized')
718 : end if
719 :
720 179712 : if (present(col_type)) then
721 0 : col_type_use = col_type
722 : else
723 179712 : col_type_use = col_type_grid
724 : end if
725 :
726 : call pbuf_register_field_int(trim(name), '', index, &
727 : persistence=persistence, dtype={ITYPE}, &
728 179712 : dimsizes=dimsizes, col_type=col_type_use, pbuf_add=.true.)
729 :
730 179712 : end subroutine pbuf_add_field_{TYPE}
731 :
732 0 : subroutine pbuf_register_subcol(name, pname, index)
733 : use subcol_utils, only: is_subcol_on
734 :
735 : ! Dummy Arguments
736 : character(len=*), intent(in) :: name
737 : character(len=*), intent(in) :: pname ! name of calling parameterization
738 : integer, intent(out) :: index
739 :
740 : ! Local variables
741 : integer :: dimsizes(1)
742 :
743 : ! You really should not call this routine if subcolumns are not on
744 0 : if (.not. is_subcol_on()) then
745 0 : call endrun('pbuf_register_subcol: subcolumns are not active')
746 : end if
747 : ! Create and pass dimsizes so that subcolums are registered
748 0 : dimsizes(1) = pcols * psubcols
749 : call pbuf_register_field_int(trim(name), trim(pname), index, &
750 0 : dimsizes=dimsizes, col_type=col_type_subcol)
751 :
752 0 : end subroutine pbuf_register_subcol
753 :
754 1536 : subroutine pbuf2d_print(pbuf2d)
755 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
756 :
757 1536 : call pbuf1d_print(pbuf_get_chunk(pbuf2d,begchunk))
758 :
759 :
760 0 : end subroutine pbuf2d_print
761 :
762 1536 : subroutine pbuf1d_print(pbuf)
763 : use spmd_utils, only: masterproc
764 :
765 : type(physics_buffer_desc), pointer :: pbuf(:)
766 :
767 : integer :: ind
768 : type(physics_buffer_desc), pointer :: pbufPtr
769 :
770 1536 : if (masterproc .and. debug) then
771 0 : write(iulog, *) __FILE__, __LINE__, currentpbufflds, size(pbuf)
772 0 : do ind = 1, currentpbufflds
773 0 : pbufPtr => pbuf(ind)
774 0 : write(iulog, *) __FILE__, __LINE__, ind, trim(pbufPtr%hdr%name), &
775 0 : pbufPtr%hdr%dtype, pbufPtr%hdr%persistence, pbufPtr%hdr%dimsizes
776 : end do
777 : end if
778 :
779 1536 : end subroutine pbuf1d_print
780 : !
781 : ! Given a pbuf field name return an integer index to the field.
782 : ! This index can be used to retrieve the field and is much faster
783 : ! than using the name in most cases
784 : !
785 :
786 22590936 : function pbuf_get_index(name, errcode) result(index)
787 : character(len=*), intent(in) :: name
788 : integer, intent(inout), optional :: errcode
789 : integer :: index
790 : integer :: i
791 : type(physics_buffer_hdr), pointer :: bufptr
792 :
793 :
794 22590936 : bufptr=>hdrbuffertop
795 22590936 : index = -1
796 1189803672 : do i=1,currentpbufflds
797 1188298320 : if(bufptr%name .eq. name) then
798 : index=i
799 : exit
800 : end if
801 1168718088 : bufptr=>bufptr%nexthdr
802 : end do
803 :
804 22590936 : if (present(errcode)) then
805 15015888 : errcode = index
806 7575048 : else if(index<0) then
807 0 : call endrun('Attempt to find undefined name in pbuf '//trim(name))
808 : end if
809 :
810 :
811 22590936 : end function pbuf_get_index
812 :
813 : !=========================================================================================
814 :
815 : !
816 : ! Given a pbuf2d chunk and an index return a pointer to a field chunk
817 : !
818 : !
819 :
820 :
821 : ! TYPE int,double,real
822 : ! DIMS 1,2,3,4,5
823 302850720 : subroutine get_pbuf1d_field_by_index_{DIMS}d_{TYPE}(pbuf, index, field, start,kount, col_type, copy_if_needed, errcode)
824 :
825 : ! Get the data based on the col_type which is specified. If no col_type, then grid field is returned
826 302850720 : use subcol_utils, only: subcol_field_copy
827 :
828 : type(physics_buffer_desc), pointer:: pbuf(:)
829 : integer, intent(in) :: index
830 : {VTYPE}, pointer :: field{DIMSTR}
831 : integer, intent(in), optional :: start(:),kount(:)
832 : integer, intent(in), optional :: col_type
833 : logical, intent(in), optional :: copy_if_needed
834 : integer, intent(out), optional :: errcode
835 :
836 302850720 : integer, pointer :: dimsizes(:)
837 302850720 : {VTYPE}, pointer :: gfield{DIMSTR}
838 302850720 : {VTYPE}, pointer :: sfield{DIMSTR}
839 : integer :: col_type_use
840 302850720 : integer, allocatable :: kount_grid(:) ! For use in copy_if_needed
841 : logical :: subset
842 : character(len=*), parameter :: subname='GET_PBUF1D_FIELD_BY_INDEX'
843 : character(len=128) :: errmsg
844 :
845 : ! Copy the generic type to one compatable with the request
846 302850720 : if((index < 1) .or. (index > size(pbuf))) then
847 0 : write(errmsg, '(2a,i0,a)') subname,': index (',index,') out of range'
848 0 : call endrun(errmsg)
849 : end if
850 :
851 : ! Default col_type is grid columns
852 302850720 : if (present(col_type)) then
853 235296 : col_type_use = col_type
854 : else
855 : col_type_use = col_type_grid
856 : end if
857 :
858 : ! If there is an errcode, start with an OK status (zero)
859 302850720 : if (present(errcode)) then
860 0 : errcode = 0
861 : end if
862 :
863 : ! Check whether subset of data requested (default is false)
864 302850720 : subset = .false.
865 302850720 : if (present(start) .and. present(kount)) subset = .true.
866 :
867 : ! Check for ill-formed request
868 302850720 : if ( (present(start) .and. .not. present(kount)) .or. &
869 : (.not. present(start) .and. present(kount)) ) then
870 0 : call endrun('pbuf_get_field: Both start and kount must be present for '//trim(pbuf(index)%hdr%name))
871 : end if
872 :
873 : ! See if we need to copy the grid field to subcolumns
874 302850720 : if (present(copy_if_needed)) then
875 0 : if (copy_if_needed) then
876 0 : if (col_type_use == col_type_subcol) then
877 : ! If a subcolumn field buffer does not exist, allocate one.
878 : ! Even if start and kount are being passed, allocate and copy the
879 : ! full grid field buffer so that a future access will succeed.
880 0 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_use))) then
881 0 : dimsizes => pbuf(index)%hdr%dimsizes(:,col_type_use)
882 0 : dimsizes(2:) = pbuf(index)%hdr%dimsizes(2:,col_type_grid)
883 0 : dimsizes(1) = pcols * psubcols
884 0 : select case(pbuf(index)%hdr%dtype)
885 : case(TYPEDOUBLE)
886 : call buffer_field_allocate(pbuf(index)%bfg(col_type_subcol), &
887 0 : dimsizes, dtype_r8)
888 : case(TYPEREAL)
889 : call buffer_field_allocate(pbuf(index)%bfg(col_type_subcol), &
890 0 : dimsizes, dtype_r4)
891 : case(TYPEINT)
892 : call buffer_field_allocate(pbuf(index)%bfg(col_type_subcol), &
893 0 : dimsizes, dtype_i4)
894 : end select
895 : ! Set copy only if we did the allocation after init time
896 0 : pbuf(index)%hdr%is_copy(col_type_subcol) = .true.
897 : end if
898 : else
899 0 : call endrun('pbuf_get_field: copy_if_needed only supported for subcolumns')
900 : end if
901 0 : if (pbuf(index)%hdr%is_copy(col_type_subcol)) then
902 : ! Only do the copy if we did the alloc (i.e., set the is_copy flag)
903 : ! Only copy the portion we are going to hand back (i.e., start, kount)
904 : ! Chances are that kount(1) = pcols*psubcols because we are looking
905 : ! for a subcolumn field (or we wouldn't be here). Now,
906 : ! subcol_field_copy requires kount(1) = pcols for the input and
907 : ! therefore, kount(1) = pcols*psubcols for the output. Check and
908 : ! make it work
909 :
910 0 : if (subset) then
911 :
912 : ! Create kount array for grid field
913 : ! Use input subcol kount array, replacing the first dimension with pcols
914 0 : if (size(kount) > size(pbuf(index)%hdr%dimsizes(:,col_type_subcol))) then
915 0 : call endrun('pbuf_get_field: kount input has too many dimensions')
916 : end if
917 0 : if (kount(1) /= pcols * psubcols) then
918 0 : call endrun('pbuf_get_field: kount(1) must be pcols*psubcols when using copy_if_needed=.true.')
919 : endif
920 :
921 0 : allocate(kount_grid(size(kount)))
922 0 : kount_grid(2:) = kount(2:)
923 0 : kount_grid(1) = pcols
924 :
925 : ! Don't need to create start array for grid field as start array for subcolumn field is identical
926 0 : if (size(start) > size(pbuf(index)%hdr%dimsizes(:,col_type_subcol))) then
927 0 : call endrun('pbuf_get_field: start input has too many dimensions')
928 : end if
929 0 : if (start(1) /= 1) then
930 0 : call endrun('pbuf_get_field: start(1) must be 1 when using copy_if_needed=.true.')
931 : end if
932 :
933 : ! Get the grid field
934 0 : pbuf(index)%used = .true.
935 0 : call buffer_get_field_ptr(pbuf(index)%bfg(col_type_grid), &
936 0 : gfield, start, kount_grid)
937 :
938 0 : deallocate(kount_grid)
939 :
940 : else
941 : ! Get the grid field
942 0 : pbuf(index)%used = .true.
943 0 : call buffer_get_field_ptr(pbuf(index)%bfg(col_type_grid), &
944 0 : gfield)
945 : end if
946 :
947 : ! Get the subcol field pointer (note optional start/kount retain their status in this call)
948 0 : pbuf(index)%used = .true.
949 0 : call buffer_get_field_ptr(pbuf(index)%bfg(col_type_subcol), &
950 0 : sfield, start, kount)
951 :
952 0 : call subcol_field_copy(gfield, pbuf(index)%lchnk, sfield)
953 :
954 : end if
955 : end if
956 : end if
957 :
958 : ! Copy or not, retrieve the requested field pointer
959 302850720 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_use))) then
960 0 : if (present(errcode)) then
961 0 : errcode = -1
962 : else
963 0 : if (col_type_use == col_type_grid) then
964 : call endrun('pbuf_get_field: probably missing a pbuf_add_field call. field not allocated for '&
965 0 : //trim(pbuf(index)%hdr%name))
966 0 : else if (col_type_use == col_type_subcol) then
967 : call endrun('pbuf_get_field: probably missing a pbuf_register_subcol. field not allocated for '&
968 0 : //trim(pbuf(index)%hdr%name))
969 : else
970 0 : call endrun('pbuf_get_field: field not allocated for '//trim(pbuf(index)%hdr%name))
971 : end if
972 : end if
973 : else
974 : ! Get the field pointer (note optional start/kount retain their status in this call)
975 302850720 : pbuf(index)%used = .true.
976 826396704 : call buffer_get_field_ptr(pbuf(index)%bfg(col_type_use),field,start,kount )
977 : end if
978 :
979 605701440 : end subroutine get_pbuf1d_field_by_index_{DIMS}d_{TYPE}
980 :
981 : ! TYPE int,double,real
982 : ! DIMS 1,2,3,4,5
983 1730664 : subroutine get_pbuf2d_field_by_index_{DIMS}d_{TYPE}(pbuf2d, lchnk, index, field, start, kount, col_type, errcode)
984 :
985 : ! Get the data based on the col_type which is specified. If no col_type, then grid field is returned
986 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
987 : integer, intent(in) :: lchnk
988 : integer, intent(in) :: index
989 : integer, intent(in), optional :: start(:),kount(:)
990 : integer, intent(in), optional :: col_type
991 : integer, intent(out), optional :: errcode
992 :
993 : {VTYPE}, pointer :: field{DIMSTR}
994 :
995 : ! Check for ill-formed request
996 1730664 : if ( (present(start) .and. .not. present(kount)) .or.&
997 : (.not. present(start) .and. present(kount)) ) then
998 0 : call endrun('pbuf_get_field: Both start and kount must be present for '//trim(pbuf2d(index,begchunk)%hdr%name))
999 : end if
1000 :
1001 : call pbuf_get_field(pbuf_get_chunk(pbuf2d,lchnk), index, field, &
1002 5191992 : start=start, kount=kount, col_type=col_type, errcode=errcode)
1003 1730664 : end subroutine get_pbuf2d_field_by_index_{DIMS}d_{TYPE}
1004 :
1005 : ! TYPE int,double,real
1006 235296 : subroutine get_pbuf2d_field_restart_{TYPE}(pbuf2d, lchnk, index, field, mdimsize, col_type)
1007 : use subcol_pack_mod, only: subcol_unpack
1008 :
1009 : ! NB: This routine is really only useful for write_restart_field
1010 : ! Get the data based on the col_type which is specified.
1011 : ! If no col_type, then grid field is assumed
1012 : ! For grid field, reference buffer and copy into chunk (field)
1013 : ! If col_type is a subcol, unpack subcolumn data
1014 : ! Field is then reshaped into (/pcols, psubcols*mdims/).
1015 :
1016 : ! Dummy variables
1017 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1018 : integer, intent(in) :: lchnk
1019 : integer, intent(in) :: index
1020 : integer, intent(in) :: mdimsize
1021 : integer, intent(in) :: col_type
1022 : {VTYPE}, intent(inout) :: field(:,:)
1023 :
1024 : ! Local variables
1025 235296 : {VTYPE}, allocatable :: exp_fld(:,:,:,:,:,:)
1026 235296 : {VTYPE}, pointer :: buf(:,:,:,:,:)
1027 : {VTYPE} :: fillvalue
1028 : character(len=128) :: errmsg
1029 235296 : integer, pointer :: dimsizes(:)
1030 :
1031 : #if ({ITYPE}==TYPEREAL)
1032 0 : fillvalue = 0._r4
1033 : #elif ({ITYPE}==TYPEDOUBLE)
1034 235296 : fillvalue = 0._r8
1035 : #else
1036 0 : fillvalue = 0
1037 : #endif
1038 :
1039 235296 : if (col_type == col_type_grid) then
1040 235296 : call pbuf_get_field(pbuf2d, lchnk, index, buf, col_type=col_type)
1041 705888 : field(:,:) = reshape(buf, (/pcols, mdimsize/))
1042 0 : else if (col_type == col_type_subcol) then
1043 : ! Don't initialize field, unpack will fill in unused slots
1044 0 : call pbuf_get_field(pbuf2d, lchnk, index, buf, col_type=col_type)
1045 0 : dimsizes => pbuf2d(index, lchnk)%hdr%dimsizes(:,col_type_subcol)
1046 0 : allocate(exp_fld(pcols, psubcols, dimsizes(2), dimsizes(3), dimsizes(4), dimsizes(5)))
1047 : ! unpack the subcolumns into their own dimension
1048 0 : call subcol_unpack(lchnk, buf, exp_fld, fillvalue)
1049 : ! Reshape back to pcols for outputting
1050 0 : field(:,:) = reshape(exp_fld, (/pcols, mdimsize/))
1051 0 : deallocate(exp_fld)
1052 : else
1053 0 : write(errmsg, *) "get_pbuf2d_field_restart_{TYPE}: Bad col_type:",col_type
1054 0 : call endrun(errmsg)
1055 : end if
1056 470592 : end subroutine get_pbuf2d_field_restart_{TYPE}
1057 :
1058 : ! TYPE int,double,real
1059 : ! DIMS 1,2,3,4,5
1060 13056 : subroutine set_pbuf2d_field_const_by_index_{TYPE}(pbuf2d,index,const, col_type)
1061 :
1062 : ! Set the field specified by the col_type
1063 :
1064 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1065 : integer, intent(in) :: index
1066 : {VTYPE},intent(in) :: const
1067 : integer,intent(in) ,optional :: col_type
1068 :
1069 : integer :: c
1070 :
1071 65688 : do c=begchunk,endchunk
1072 65688 : if(present(col_type)) then
1073 0 : call set_pbuf1d_field_const_by_index_{TYPE}(pbuf_get_chunk(pbuf2d,c),index,const, col_type=col_type)
1074 : else
1075 52632 : call set_pbuf1d_field_const_by_index_{TYPE}(pbuf_get_chunk(pbuf2d,c),index,const)
1076 : end if
1077 : end do
1078 :
1079 13056 : end subroutine set_pbuf2d_field_const_by_index_{TYPE}
1080 :
1081 : ! TYPE int,double,real
1082 : ! DIMS 1,2,3,4,5
1083 52632 : subroutine set_pbuf1d_field_const_by_index_{TYPE}(pbuf,index,const,start,kount, col_type)
1084 :
1085 : ! Set the field(s) specified by the col_type
1086 : type(physics_buffer_desc), pointer :: pbuf(:)
1087 : integer, intent(in) :: index
1088 : {VTYPE},intent(in) :: const
1089 : integer, intent(in), optional :: start(:),kount(:)
1090 : integer, intent(in), optional :: col_type
1091 :
1092 : integer :: col_type_use
1093 : logical :: subset
1094 :
1095 : ! Default col_type is grid
1096 52632 : if (present(col_type)) then
1097 0 : col_type_use = col_type
1098 : else
1099 : col_type_use = col_type_grid
1100 : end if
1101 :
1102 : ! Check whether subset of data requested (default is false)
1103 52632 : subset = .false.
1104 52632 : if (present(start) .and. present(kount)) subset = .true.
1105 :
1106 : ! Check for ill-formed request
1107 52632 : if ( (present(start) .and. .not. present(kount)) .or.&
1108 : (.not. present(start) .and. present(kount)) ) then
1109 0 : call endrun('pbuf_set_field: Both start and kount must be present for '//trim(pbuf(index)%hdr%name))
1110 : end if
1111 :
1112 : ! Set the appropriate grid or sub-column field. Check that the fields have been allocated.
1113 52632 : if(subset) then
1114 :
1115 0 : if (col_type_use == col_type_subcol) then
1116 : ! Set sub-column field
1117 0 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_subcol))) &
1118 0 : call endrun('pbuf_set_field: sub-column field not allocated for '//trim(pbuf(index)%hdr%name))
1119 0 : pbuf(index)%used = .true.
1120 0 : call buffer_set_field(pbuf(index)%bfg(col_type_subcol),const,start,kount)
1121 :
1122 0 : else if (col_type_use == col_type_grid) then
1123 : ! Set grid field
1124 0 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_grid))) &
1125 0 : call endrun('pbuf_set_field: grid field not allocated for '//trim(pbuf(index)%hdr%name))
1126 0 : pbuf(index)%used = .true.
1127 0 : call buffer_set_field(pbuf(index)%bfg(col_type_grid),const,start,kount)
1128 : else
1129 0 : call endrun('pbuf_set_field: Trying to set '//trim(pbuf(index)%hdr%name)//&
1130 0 : ' but col_type is neither col_type_grid nor col_type_subcol')
1131 : end if
1132 :
1133 : else
1134 :
1135 52632 : if (col_type_use == col_type_subcol) then
1136 : ! Set sub-column field
1137 0 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_subcol))) &
1138 0 : call endrun('pbuf_set_field: sub-column field not allocated for '//trim(pbuf(index)%hdr%name))
1139 0 : pbuf(index)%used = .true.
1140 0 : call buffer_set_field(pbuf(index)%bfg(col_type_subcol),const)
1141 :
1142 52632 : else if (col_type_use == col_type_grid) then
1143 : ! Set grid field
1144 52632 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_grid))) &
1145 0 : call endrun('pbuf_set_field: grid field not allocated for '//trim(pbuf(index)%hdr%name))
1146 52632 : pbuf(index)%used = .true.
1147 52632 : call buffer_set_field(pbuf(index)%bfg(col_type_grid),const)
1148 : else
1149 0 : call endrun('pbuf_set_field: Trying to set '//trim(pbuf(index)%hdr%name)//&
1150 0 : ' but col_type is neither col_type_grid nor col_type_subcol')
1151 : end if
1152 :
1153 : end if
1154 :
1155 52632 : end subroutine set_pbuf1d_field_const_by_index_{TYPE}
1156 :
1157 : ! TYPE int,double,real
1158 : ! DIMS 1,2,3,4,5
1159 4608 : subroutine set_pbuf2d_field_by_index_{DIMS}d_{TYPE}(pbuf2d,index,field, start, kount, col_type)
1160 :
1161 : ! Set the field(s) specified by the col_type
1162 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1163 : integer, intent(in) :: index
1164 : integer,intent(in),optional :: start(:), kount(:)
1165 : integer,intent(in),optional :: col_type
1166 :
1167 : logical :: subset
1168 :
1169 : integer :: c
1170 4608 : {VTYPE}, pointer :: fld{DIMSTR}
1171 :
1172 :
1173 : #if ({DIMS}==1)
1174 : {VTYPE},pointer :: field(:,:)
1175 : #elif ({DIMS}==2)
1176 : {VTYPE},pointer :: field(:,:,:)
1177 : #elif ({DIMS}==3)
1178 : {VTYPE},pointer :: field(:,:,:,:)
1179 : #elif ({DIMS}==4)
1180 : {VTYPE},pointer :: field(:,:,:,:,:)
1181 : #elif ({DIMS}==5)
1182 : {VTYPE},pointer :: field(:,:,:,:,:,:)
1183 : #endif
1184 :
1185 : ! Check whether subset of data requested (default is false)
1186 4608 : subset = .false.
1187 2304 : if (present(start) .and. present(kount)) subset = .true.
1188 :
1189 : ! Check for ill-formed request
1190 4608 : if ( (present(start) .and. .not. present(kount)) .or.&
1191 : (.not. present(start) .and. present(kount)) ) then
1192 0 : call endrun('pbuf_set_field: Both start and kount must be present for '//trim(pbuf2d(index,begchunk)%hdr%name))
1193 : end if
1194 :
1195 23184 : do c=begchunk,endchunk
1196 18576 : fld => get_field_chunk_{DIMS}d_{TYPE}(field,c)
1197 23184 : if(subset .and. present(col_type)) then
1198 0 : call pbuf_set_field(pbuf_get_chunk(pbuf2d,c),index,fld,start,kount, col_type)
1199 18576 : else if(subset) then
1200 9288 : call pbuf_set_field(pbuf_get_chunk(pbuf2d,c),index,fld,start,kount)
1201 9288 : else if(present(col_type)) then
1202 0 : call pbuf_set_field(pbuf_get_chunk(pbuf2d,c),index,fld,col_type=col_type)
1203 : else
1204 9288 : call pbuf_set_field(pbuf_get_chunk(pbuf2d,c),index,fld)
1205 : end if
1206 : end do
1207 4608 : end subroutine set_pbuf2d_field_by_index_{DIMS}d_{TYPE}
1208 :
1209 : ! TYPE int,double,real
1210 : ! DIMS 1,2,3,4,5
1211 18576 : function get_field_chunk_{DIMS}d_{TYPE}(field, c) result(fld)
1212 : ! module private helper function
1213 : {VTYPE}, pointer :: fld{DIMSTR}
1214 : integer, intent(in) :: c
1215 :
1216 : #if ({DIMS}==1)
1217 : {VTYPE},pointer :: field(:,:)
1218 9288 : fld => field(:,c)
1219 : #elif ({DIMS}==2)
1220 : {VTYPE},pointer :: field(:,:,:)
1221 9288 : fld => field(:,:,c)
1222 : #elif ({DIMS}==3)
1223 : {VTYPE},pointer :: field(:,:,:,:)
1224 0 : fld => field(:,:,:,c)
1225 : #elif ({DIMS}==4)
1226 : {VTYPE},pointer :: field(:,:,:,:,:)
1227 0 : fld => field(:,:,:,:,c)
1228 : #elif ({DIMS}==5)
1229 : {VTYPE},pointer :: field(:,:,:,:,:,:)
1230 0 : fld => field(:,:,:,:,:,c)
1231 : #endif
1232 :
1233 18576 : end function get_field_chunk_{DIMS}d_{TYPE}
1234 :
1235 :
1236 :
1237 : ! TYPE int,double,real
1238 : ! DIMS 1,2,3,4,5
1239 5984568 : subroutine set_pbuf1d_field_by_index_{DIMS}d_{TYPE}(pbuf,index,field, start, kount, col_type)
1240 :
1241 : type(physics_buffer_desc), pointer :: pbuf(:)
1242 : integer, intent(in) :: index
1243 : {VTYPE}, intent(in) :: field{DIMSTR}
1244 : integer,intent(in),optional :: start(:), kount(:)
1245 : integer,intent(in),optional :: col_type
1246 :
1247 : integer :: col_type_use
1248 : logical :: subset
1249 :
1250 : ! Default col_type is grid only
1251 5984568 : if (present(col_type)) then
1252 0 : col_type_use = col_type
1253 : else
1254 : col_type_use = col_type_grid
1255 : end if
1256 :
1257 : ! Check whether subset of data requested (default is false)
1258 5984568 : subset = .false.
1259 5984568 : if (present(start) .and. present(kount)) subset = .true.
1260 :
1261 : ! Check for ill-formed request
1262 5984568 : if ( (present(start) .and. .not. present(kount)) .or.&
1263 : (.not. present(start) .and. present(kount)) ) then
1264 0 : call endrun('pbuf_set_field: Both start and kount must be present for '//trim(pbuf(index)%hdr%name))
1265 : end if
1266 :
1267 5984568 : if(subset) then
1268 :
1269 : ! Set sub-column field
1270 2993832 : if (col_type_use == col_type_subcol) then
1271 0 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_subcol))) &
1272 0 : call endrun('pbuf_set_field: sub-column field not allocated for '//trim(pbuf(index)%hdr%name))
1273 0 : pbuf(index)%used = .true.
1274 0 : call buffer_set_field(pbuf(index)%bfg(col_type_subcol),field,start,kount)
1275 :
1276 : ! Set grid field
1277 2993832 : else if (col_type_use == col_type_grid) then
1278 2993832 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_grid))) &
1279 0 : call endrun('pbuf_set_field: grid field not allocated for '//trim(pbuf(index)%hdr%name))
1280 2993832 : pbuf(index)%used = .true.
1281 2993832 : call buffer_set_field(pbuf(index)%bfg(col_type_grid),field,start,kount)
1282 : else
1283 0 : call endrun('pbuf_set_field: Trying to set '//trim(pbuf(index)%hdr%name)//&
1284 0 : ' but col_type is neigher col_type_grid nor col_type_subcol ')
1285 : end if
1286 : else
1287 :
1288 : ! Set sub-column field
1289 2990736 : if (col_type_use == col_type_subcol) then
1290 0 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_subcol))) &
1291 0 : call endrun('pbuf_set_field: sub-column field not allocated for '//trim(pbuf(index)%hdr%name))
1292 0 : pbuf(index)%used = .true.
1293 0 : call buffer_set_field(pbuf(index)%bfg(col_type_subcol),field)
1294 :
1295 : ! Set grid field
1296 2990736 : else if (col_type_use == col_type_grid) then
1297 2990736 : if (.not. buffer_field_is_alloc(pbuf(index)%bfg(col_type_grid))) &
1298 0 : call endrun('pbuf_set_field: grid field not allocated for '//trim(pbuf(index)%hdr%name))
1299 2990736 : pbuf(index)%used = .true.
1300 2990736 : call buffer_set_field(pbuf(index)%bfg(col_type_grid),field)
1301 : else
1302 0 : call endrun('pbuf_set_field: Trying to set '//trim(pbuf(index)%hdr%name)//&
1303 0 : ' but col_type is neither col_type_grid nor col_type_subcol ')
1304 : end if
1305 : endif
1306 :
1307 5984568 : end subroutine set_pbuf1d_field_by_index_{DIMS}d_{TYPE}
1308 :
1309 :
1310 58368 : function pbuftype2piotype(pbuftype) result(piotype)
1311 : use pio, only : pio_double, pio_int, pio_real
1312 :
1313 : integer, intent(in) :: pbuftype
1314 : integer :: piotype
1315 :
1316 58368 : select case(pbuftype)
1317 : case (TYPEDOUBLE)
1318 0 : piotype = pio_double
1319 : case(TYPEINT)
1320 0 : piotype = pio_int
1321 : case(TYPEREAL)
1322 0 : piotype = pio_real
1323 : ! case(TYPELONG)
1324 : ! piotype = pio_int
1325 : case default
1326 0 : write(iulog, *) 'Dtype = ', pbuftype
1327 58368 : call endrun('No restart support for dtype')
1328 : end select
1329 58368 : end function pbuftype2piotype
1330 :
1331 :
1332 : !
1333 : ! Initialize a restart file to write - all additional dims in a field are
1334 : ! bundled into a single dimension for output and a dim pbuf_xxxxx is declared
1335 : ! in the file if it does not already exist.
1336 : !
1337 1536 : subroutine pbuf_init_restart(File, pbuf2d)
1338 : use pio, only: file_desc_t, pio_int
1339 : use cam_pio_utils, only: cam_pio_def_dim, cam_pio_def_var
1340 : use phys_grid, only: phys_decomp
1341 : use cam_grid_support, only: cam_grid_get_file_dimids, cam_grid_dimensions
1342 :
1343 :
1344 : ! Dummy Variables
1345 : type(file_desc_t), intent(inout) :: file
1346 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1347 :
1348 : ! Local Variables
1349 : type(physics_buffer_desc), pointer :: pbuf
1350 : integer :: i, grid_select
1351 : integer :: adimid(3) ! PIO IDs
1352 : integer :: hdimcnt ! # grid dims
1353 : integer :: dimcnt ! # array dims
1354 : integer :: mdimsize, piodtype
1355 :
1356 : character(len=10) :: dimname
1357 : character(len=24) :: varname
1358 :
1359 : ! Use adimid as a temp to find number of horizontal dims
1360 1536 : call cam_grid_dimensions(phys_decomp, adimid(1:2), hdimcnt)
1361 1536 : call cam_grid_get_file_dimids(phys_decomp, File, adimid)
1362 :
1363 181248 : do i = 1, currentpbufflds
1364 179712 : pbuf => pbuf2d(i,begchunk)
1365 :
1366 : ! Only save global pbufs for restart
1367 179712 : if(pbuf%hdr%persistence /= persistence_global) cycle
1368 :
1369 58368 : piodtype = pbuftype2piotype(pbuf%hdr%dtype)
1370 :
1371 176640 : do grid_select = 1, ngrid_types
1372 : ! For subcol fields, mdimsize includes psubcols in size
1373 817152 : mdimsize = product(pbuf%hdr%dimsizes(:,grid_select))/pcols
1374 116736 : if(mdimsize > 1) then
1375 32256 : write(dimname,'(a,i5.5)') 'pbuf_',mdimsize
1376 32256 : call cam_pio_def_dim(File, dimname, mdimsize, adimid(hdimcnt+1), &
1377 64512 : existOK=.true.)
1378 32256 : dimcnt = hdimcnt + 1
1379 : else
1380 84480 : dimcnt = hdimcnt
1381 : end if
1382 296448 : if (mdimsize > 0) then
1383 58368 : varname = trim(pbuf%hdr%name)//trim(field_grid_suff(grid_select))
1384 0 : call cam_pio_def_var(File, varname, piodtype, adimid(1:dimcnt), &
1385 58368 : pbuf%hdr%vardesc(grid_select), existOK=.false.)
1386 : end if
1387 : end do
1388 :
1389 : end do
1390 :
1391 : call cam_pio_def_var(File, 'pbuf_time_idx', pio_int, timeidx_desc, &
1392 1536 : existOK=.false.)
1393 :
1394 1536 : end subroutine pbuf_init_restart
1395 :
1396 :
1397 1536 : subroutine pbuf_write_restart(File, pbuf2d)
1398 1536 : use pio, only: file_desc_t, pio_put_var
1399 : use cam_grid_support, only: cam_grid_dimensions
1400 : use phys_grid, only: phys_decomp
1401 :
1402 : ! Dummy Variables
1403 : type(file_desc_t), intent(inout) :: file
1404 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1405 :
1406 : ! Local Variables
1407 1536 : type(physics_buffer_desc), pointer :: pbufhdr(:)
1408 :
1409 : integer :: index, dtype, ierr
1410 : integer :: gdimlens(2)
1411 : integer :: grank
1412 :
1413 1536 : pbufhdr => pbuf_get_chunk(pbuf2d, begchunk)
1414 4608 : gdimlens = 1
1415 1536 : call cam_grid_dimensions(phys_decomp, gdimlens(1:2), grank)
1416 :
1417 181248 : do index = 1, currentpbufflds
1418 181248 : if(pbufhdr(index)%hdr%persistence == persistence_global) then
1419 58368 : dtype = pbufhdr(index)%hdr%dtype
1420 58368 : select case(dtype)
1421 : case (TYPEDOUBLE)
1422 58368 : call write_restart_field_double(File, pbuf2d, index, gdimlens, grank)
1423 : case (TYPEREAL)
1424 0 : call write_restart_field_real(File, pbuf2d, index, gdimlens, grank)
1425 : case (TYPEINT)
1426 58368 : call write_restart_field_int(File, pbuf2d, index, gdimlens, grank)
1427 : end select
1428 : end if
1429 : end do
1430 3072 : ierr = pio_put_var(File, timeidx_desc, (/old_time_idx/))
1431 :
1432 3072 : end subroutine pbuf_write_restart
1433 :
1434 145920 : subroutine pbuf_restart_dimsizes(pbuf, gdimlens, grank, gridnum, adimlens, &
1435 : fdimlens, frank)
1436 :
1437 : ! Dummy arguments
1438 : type(physics_buffer_desc), pointer :: pbuf
1439 : integer, intent(in) :: gdimlens(2) ! Grid horiz. dim sizes
1440 : integer, intent(in) :: grank ! Array rank on file
1441 : integer, intent(in) :: gridnum ! pbuf grid selector
1442 : integer, intent(out) :: adimlens(3) ! Array dims
1443 : integer, intent(out) :: fdimlens(3) ! Array dims on file
1444 : integer, intent(out) :: frank ! array rank on file
1445 :
1446 : ! Local variable
1447 : integer :: mdimsize
1448 :
1449 1021440 : mdimsize = PRODUCT(pbuf%hdr%dimsizes(:,gridnum)) / pcols
1450 437760 : fdimlens(1:2) = gdimlens
1451 145920 : if (grank == 1) then
1452 : ! Unstructured grid
1453 145920 : fdimlens(2) = mdimsize
1454 145920 : fdimlens(3) = 1
1455 145920 : frank = 2
1456 : else
1457 0 : fdimlens(3) = mdimsize
1458 0 : frank = 3
1459 : end if
1460 : ! Restart does not write a dimension if mdimsize == 1
1461 145920 : if (mdimsize == 1) then
1462 39168 : frank = frank - 1
1463 : end if
1464 :
1465 145920 : adimlens(1) = pcols
1466 145920 : adimlens(2) = mdimsize
1467 145920 : adimlens(3) = endchunk - begchunk + 1
1468 :
1469 1536 : end subroutine pbuf_restart_dimsizes
1470 :
1471 : ! TYPE int,double,real
1472 58368 : subroutine write_restart_field_{TYPE}(File, pbuf2d, index, gdimlens, grank)
1473 58368 : use cam_grid_support, only: cam_grid_write_dist_array
1474 : use phys_grid, only: phys_decomp
1475 : use pio, only: file_desc_t
1476 :
1477 : ! Dummy Variables
1478 : type(file_desc_t), intent(inout) :: File
1479 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1480 : integer, intent(in) :: index
1481 : integer, intent(in) :: gdimlens(2)
1482 : integer, intent(in) :: grank
1483 :
1484 : ! Local Variables
1485 : type(physics_buffer_desc), pointer :: pbuf
1486 58368 : {VTYPE}, allocatable :: field(:,:,:)
1487 :
1488 : integer :: grid_select ! 1=grid, 2=subcol
1489 : integer :: c ! Chunk index
1490 : integer :: fdimlens(3) ! Array dimensions on NetCDF file
1491 : integer :: adimlens(3) ! Array dimensions
1492 : integer :: frank
1493 :
1494 58368 : pbuf => pbuf2d(index, begchunk)
1495 :
1496 175104 : do grid_select = 1, ngrid_types
1497 : call pbuf_restart_dimsizes(pbuf, gdimlens, grank, grid_select, &
1498 116736 : adimlens, fdimlens, frank)
1499 525312 : if ((PRODUCT(adimlens) > 0) .and. (.not. pbuf%hdr%is_copy(grid_select))) then
1500 291840 : allocate(field(adimlens(1), adimlens(2), begchunk:endchunk))
1501 293664 : do c = begchunk, endchunk
1502 : call pbuf_get_field_restart(pbuf2d, c, index, field(:,:,c), &
1503 293664 : adimlens(2), grid_select)
1504 : end do
1505 58368 : if (size(field,2) == 1) then
1506 : ! Special case for 2-D pbuf fields
1507 26112 : adimlens(2) = adimlens(3)
1508 : call cam_grid_write_dist_array(File, phys_decomp, adimlens(1:2), &
1509 26112 : fdimlens(1:frank), field(:,1,:), pbuf%hdr%vardesc(grid_select))
1510 : else
1511 : call cam_grid_write_dist_array(File, phys_decomp, adimlens, &
1512 32256 : fdimlens(1:frank), field, pbuf%hdr%vardesc(grid_select))
1513 : end if
1514 58368 : deallocate(field)
1515 : end if
1516 : end do
1517 :
1518 58368 : end subroutine write_restart_field_{TYPE}
1519 :
1520 :
1521 768 : subroutine pbuf_read_restart(File, pbuf2d)
1522 0 : use cam_grid_support, only: cam_grid_dimensions
1523 : use phys_grid, only: phys_decomp
1524 : use pio, only: file_desc_t, pio_inq_varid, pio_get_var
1525 :
1526 : ! Dummy Variables
1527 : type(File_desc_t), intent(inout) :: File
1528 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1529 :
1530 : ! Local Variables
1531 768 : type(physics_buffer_desc), pointer :: pbufhdr(:)
1532 :
1533 : integer :: index, dtype, ierr
1534 : integer :: gdimlens(2) ! Horizontal grid dimensions
1535 : integer :: grank
1536 :
1537 768 : call pbuf_initialize(pbuf2d)
1538 :
1539 768 : pbufhdr => pbuf_get_chunk(pbuf2d, begchunk)
1540 :
1541 2304 : gdimlens = 1
1542 768 : call cam_grid_dimensions(phys_decomp, gdimlens(1:2), grank)
1543 768 : ierr = pio_inq_varid(File, 'pbuf_time_idx', timeidx_desc)
1544 768 : ierr = pio_get_var(File, timeidx_desc, old_time_idx)
1545 :
1546 90624 : do index = 1, currentpbufflds
1547 90624 : if(pbufhdr(index)%hdr%persistence == persistence_global) then
1548 29184 : dtype = pbufhdr(index)%hdr%dtype
1549 29184 : select case(dtype)
1550 : case (TYPEDOUBLE)
1551 29184 : call read_restart_field_double(File, pbuf2d, index, gdimlens, grank)
1552 : case (TYPEREAL)
1553 0 : call read_restart_field_real(File, pbuf2d, index, gdimlens, grank)
1554 : case (TYPEINT)
1555 29184 : call read_restart_field_int(File, pbuf2d, index, gdimlens, grank)
1556 : end select
1557 : end if
1558 : end do
1559 :
1560 1536 : end subroutine pbuf_read_restart
1561 :
1562 : ! TYPE int,double,real
1563 29184 : subroutine read_restart_field_{TYPE} (File, pbuf2d, index, gdimlens, grank)
1564 29952 : use pio, only: file_desc_t
1565 : use pio, only: pio_inq_varid
1566 : use cam_pio_utils, only: cam_pio_handle_error
1567 : use cam_grid_support, only: cam_grid_read_dist_array
1568 : use subcol_pack_mod, only: subcol_pack
1569 : use phys_grid, only: phys_decomp
1570 :
1571 : ! Dummy Variables
1572 : type(file_desc_t), intent(inout) :: File
1573 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1574 : integer, intent(in) :: index
1575 : integer, intent(in) :: gdimlens(2)
1576 : integer, intent(in) :: grank
1577 :
1578 : ! Local Variables
1579 : type(physics_buffer_desc), pointer :: pbuf
1580 29184 : {VTYPE}, pointer :: fld6(:,:,:,:,:,:)
1581 29184 : {VTYPE}, allocatable :: fld5_pack(:,:,:,:,:)
1582 29184 : {VTYPE}, allocatable :: field(:,:,:)
1583 :
1584 : integer :: grid_select ! (1=grid, 2=subcol)
1585 : integer :: ierr, c
1586 : integer :: start(6)
1587 : integer :: dimsizes(6)
1588 : integer :: fdimlens(3) ! Array dimensions on NetCDF file
1589 : integer :: frank ! Array rank on file
1590 : integer :: adimlens(3) ! Array dimensions
1591 :
1592 : character(len=24) :: varname
1593 : character(len=*), parameter :: subname = 'read_restart_field_{TYPE}'
1594 :
1595 29184 : pbuf => pbuf2d(index, begchunk)
1596 204288 : start = 1
1597 :
1598 87552 : do grid_select = 1, ngrid_types
1599 :
1600 408576 : dimsizes(:) = pbuf%hdr%dimsizes(:, grid_select)
1601 233472 : if(all(dimsizes(:) == 0)) then
1602 : ! None of this grid type for this variable
1603 : cycle
1604 : end if
1605 :
1606 : call pbuf_restart_dimsizes(pbuf, gdimlens, grank, grid_select, &
1607 29184 : adimlens, fdimlens, frank)
1608 : ! Fix up dimensions for subcolumn field
1609 29184 : if (grid_select == col_type_subcol) then
1610 : ! Field stored as (pcols, psubcols*restDims, chunks)
1611 : ! Read then pack to (pcols*psubcols,restDims, chunks)
1612 0 : allocate(fld5_pack(dimsizes(1),dimsizes(2),dimsizes(3),dimsizes(4),dimsizes(5)))
1613 0 : do c = 5, 2, -1
1614 0 : dimsizes(c + 1) = dimsizes(c)
1615 : end do
1616 0 : dimsizes(1) = pcols
1617 0 : dimsizes(2) = psubcols
1618 : end if
1619 :
1620 29184 : varname = trim(pbuf%hdr%name)//trim(field_grid_suff(grid_select))
1621 29184 : ierr = pio_inq_varid(File, varname, pbuf%hdr%vardesc(grid_select))
1622 29184 : call cam_pio_handle_error(ierr, trim(subname)//': '//trim(varname)//' not found')
1623 :
1624 145920 : allocate(field(adimlens(1), adimlens(2), begchunk:endchunk))
1625 29184 : if (size(field,2) == 1) then
1626 : ! Special case for 2-D pbuf fields
1627 13056 : adimlens(2) = adimlens(3)
1628 : call cam_grid_read_dist_array(File, phys_decomp, adimlens(1:2), &
1629 13056 : fdimlens(1:frank), field(:,1,:), pbuf%hdr%vardesc(grid_select))
1630 : else
1631 : call cam_grid_read_dist_array(File, phys_decomp, adimlens, &
1632 16128 : fdimlens(1:frank), field, pbuf%hdr%vardesc(grid_select))
1633 : end if
1634 :
1635 146832 : do c = begchunk, endchunk
1636 146832 : if (grid_select == col_type_grid) then
1637 117648 : pbuf%used = .true.
1638 117648 : call buffer_get_field_ptr(pbuf2d(index,c)%bfg(grid_select), fld6, &
1639 235296 : start, dimsizes)
1640 823536 : fld6 = reshape(field(:,:,c), dimsizes)
1641 0 : else if (grid_select == col_type_subcol) then
1642 0 : call subcol_pack(c, reshape(field(:,:,c), dimsizes), fld5_pack)
1643 0 : pbuf%used = .true.
1644 0 : call buffer_set_field(pbuf2d(index,c)%bfg(grid_select), fld5_pack)
1645 : else
1646 0 : call endrun('read_restart_field_{TYPE}: invalid grid selector - must be either 1 or 2')
1647 : end if
1648 : end do
1649 29184 : nullify(fld6)
1650 29184 : deallocate(field)
1651 87552 : if (allocated(fld5_pack)) then
1652 0 : deallocate(fld5_pack)
1653 : end if
1654 : end do
1655 :
1656 58368 : end subroutine read_restart_field_{TYPE}
1657 :
1658 0 : logical function pbuf_is_used(pbuf)
1659 : type(physics_buffer_desc), intent(in) :: pbuf
1660 :
1661 0 : pbuf_is_used = pbuf%used
1662 0 : end function pbuf_is_used
1663 :
1664 0 : subroutine pbuf_dump_pbuf(pbuf2d, name, num)
1665 : use cam_pio_utils, only: cam_pio_dump_field
1666 : use spmd_utils, only: masterproc
1667 :
1668 : ! Dummy arguments
1669 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
1670 : character(len=*), optional, intent(in) :: name
1671 : integer, optional, intent(in) :: num
1672 :
1673 : ! Local variables
1674 : integer, parameter :: max_name = 64
1675 : character(len=max_name) :: field_name
1676 : integer :: index, grid_select, c
1677 : integer :: namelen
1678 : integer :: dimstart(6), dimend(6)
1679 : integer :: ierr
1680 : integer :: dtype
1681 0 : type(physics_buffer_desc), pointer :: pbufhdr(:)
1682 0 : real(r8), allocatable :: field(:,:,:,:,:,:)
1683 0 : real(r8), pointer :: fld5(:,:,:,:,:)
1684 :
1685 0 : pbufhdr => pbuf_get_chunk(pbuf2d, begchunk)
1686 0 : dimstart = 1
1687 :
1688 0 : do index = 1, currentpbufflds
1689 0 : dtype = pbufhdr(index)%hdr%dtype
1690 0 : do grid_select = 1, ngrid_types
1691 0 : namelen = len_trim(pbufhdr(index)%hdr%name)
1692 0 : namelen = namelen + len_trim(field_grid_suff(grid_select))
1693 0 : if (present(name)) then
1694 0 : namelen = namelen + len_trim(name)
1695 : end if
1696 0 : if (present(num)) then
1697 0 : namelen = namelen + int(log10(real(num))) + 2
1698 : end if
1699 0 : if (namelen > 64) then
1700 0 : call endrun("PBUF_DUMP_PBUF: Name string too long")
1701 : end if
1702 0 : if (present(name)) then
1703 0 : if (present(num)) then
1704 0 : write(field_name, '(4a,i0)') trim(pbufhdr(index)%hdr%name), &
1705 0 : trim(field_grid_suff(grid_select)), trim(name), '_', num
1706 : else
1707 0 : write(field_name, '(3a)') trim(pbufhdr(index)%hdr%name), &
1708 0 : trim(field_grid_suff(grid_select)), trim(name)
1709 : end if
1710 0 : else if (present(num)) then
1711 0 : write(field_name, '(3a,i0)') trim(pbufhdr(index)%hdr%name), &
1712 0 : trim(field_grid_suff(grid_select)), '_', num
1713 : else
1714 0 : write(field_name, '(2a)') trim(pbufhdr(index)%hdr%name), &
1715 0 : trim(field_grid_suff(grid_select))
1716 : end if
1717 0 : dimend = pbufhdr(index)%hdr%dimsizes(:, grid_select)
1718 0 : if (PRODUCT(dimend) > 0) then
1719 0 : if (dimend(6) /= 1) then
1720 0 : if (masterproc) then
1721 0 : write(iulog, *) 'PBUF_DUMP_PBUF: ', trim(field_name), dimend
1722 : end if
1723 0 : call endrun("PBUF_DUMP_PBUF: No support for 6D pbuf field")
1724 : end if
1725 0 : dimend(6) = endchunk - begchunk + 1
1726 0 : select case(dtype)
1727 : case (TYPEDOUBLE)
1728 0 : allocate(field(dimend(1), dimend(2), dimend(3), dimend(4), &
1729 0 : dimend(5), dimend(6)))
1730 0 : do c = begchunk, endchunk
1731 : call pbuf_get_field(pbuf2d, c, index, fld5, &
1732 0 : col_type=grid_select, errcode=ierr)
1733 0 : field(:,:,:,:,:,c-begchunk+1) = fld5(:,:,:,:,:)
1734 : end do
1735 : call cam_pio_dump_field(field_name, dimstart, dimend, &
1736 0 : field, compute_maxdim_in=.false.)
1737 0 : deallocate(field)
1738 : case (TYPEREAL)
1739 0 : if (masterproc) then
1740 0 : write(iulog, *) 'PBUF_DUMP_PBUF: No support for real fields'
1741 : end if
1742 : case (TYPEINT)
1743 0 : if (masterproc) then
1744 0 : write(iulog, *) 'PBUF_DUMP_PBUF: No support for integer fields'
1745 : end if
1746 : end select
1747 : end if
1748 : end do
1749 : end do
1750 0 : end subroutine pbuf_dump_pbuf
1751 :
1752 0 : subroutine pbuf_cam_snapshot_register()
1753 : !--------------------------------------------------------
1754 : ! This subroutine checks that all of the dimensions in pbuf fields have been registered
1755 : ! If they have not, then a new dimension is added with the "pbuf_" prepending the actual size
1756 : !--------------------------------------------------------
1757 :
1758 0 : use cam_grid_support, only: max_hcoordname_len
1759 : use cam_history_support, only: hist_dimension_name, add_hist_coord
1760 :
1761 : type(physics_buffer_hdr), pointer :: hdrbuffer
1762 : integer :: i, j
1763 : integer :: dim_size
1764 : character(len=max_hcoordname_len) :: dim_name
1765 : character(len=max_hcoordname_len) :: dim_string
1766 :
1767 0 : hdrbuffer => hdrbuffertop
1768 0 : do i=1, currentpbufflds
1769 0 : do j =1, 6 ! The max number of pbuf dimensions
1770 0 : dim_size=hdrbuffer%dimsizes(j,1)
1771 0 : if (dim_size /= 1) then ! Don't check subcolumn dimensions
1772 0 : dim_name = hist_dimension_name(dim_size)
1773 0 : if (dim_name == '') then ! found a dimension which needs a name
1774 0 : write(dim_string,'(a,i0)') 'pbuf_',dim_size
1775 0 : call add_hist_coord(dim_string, dim_size, dim_string)
1776 : end if
1777 : end if
1778 : end do
1779 0 : hdrbuffer => hdrbuffer%nexthdr
1780 : end do
1781 :
1782 0 : end subroutine pbuf_cam_snapshot_register
1783 :
1784 0 : subroutine pbuf_get_dim_strings(pbuf, string)
1785 :
1786 : !--------------------------------------------------------
1787 : ! This subroutine retrieves the pbuf dimension strings
1788 : !--------------------------------------------------------
1789 :
1790 0 : use cam_history_support, only: hist_dimension_name
1791 : use cam_grid_support, only: max_hcoordname_len
1792 :
1793 : type(physics_buffer_desc), intent(in) :: pbuf(:)
1794 : character(len=*), intent(inout) :: string(:,:)
1795 :
1796 : character(len=max_hcoordname_len) :: dim_name
1797 : integer :: i, j
1798 :
1799 0 : do i = 1, size(pbuf(:))
1800 0 : string(i,:)=''
1801 0 : do j = 1, 6 ! Hardwired number of pbuf dimensions
1802 0 : dim_name = hist_dimension_name(pbuf(i)%hdr%dimsizes(j,1))
1803 0 : if (dim_name /= '') then
1804 0 : string(i,j)=trim(dim_name)
1805 : end if
1806 : end do
1807 : end do
1808 :
1809 0 : end subroutine pbuf_get_dim_strings
1810 :
1811 :
1812 :
1813 0 : end module physics_buffer
|